Actual source code: ex17.c
petsc-3.7.3 2016-07-24
2: static char help[] = "Solves a linear system with KSP. This problem is\n\
3: intended to test the complex numbers version of various solvers.\n\n";
5: #include <petscksp.h>
7: typedef enum {TEST_1,TEST_2,TEST_3,HELMHOLTZ_1,HELMHOLTZ_2} TestType;
8: extern PetscErrorCode FormTestMatrix(Mat,PetscInt,TestType);
12: int main(int argc,char **args)
13: {
14: Vec x,b,u; /* approx solution, RHS, exact solution */
15: Mat A; /* linear system matrix */
16: KSP ksp; /* KSP context */
18: PetscInt n = 10,its, dim,p = 1,use_random;
19: PetscScalar none = -1.0,pfive = 0.5;
20: PetscReal norm;
21: PetscRandom rctx;
22: TestType type;
23: PetscBool flg;
25: PetscInitialize(&argc,&args,(char*)0,help);
26: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
27: PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL);
28: switch (p) {
29: case 1: type = TEST_1; dim = n; break;
30: case 2: type = TEST_2; dim = n; break;
31: case 3: type = TEST_3; dim = n; break;
32: case 4: type = HELMHOLTZ_1; dim = n*n; break;
33: case 5: type = HELMHOLTZ_2; dim = n*n; break;
34: default: type = TEST_1; dim = n;
35: }
37: /* Create vectors */
38: VecCreate(PETSC_COMM_WORLD,&x);
39: VecSetSizes(x,PETSC_DECIDE,dim);
40: VecSetFromOptions(x);
41: VecDuplicate(x,&b);
42: VecDuplicate(x,&u);
44: use_random = 1;
45: flg = PETSC_FALSE;
46: PetscOptionsGetBool(NULL,NULL,"-norandom",&flg,NULL);
47: if (flg) {
48: use_random = 0;
49: VecSet(u,pfive);
50: } else {
51: PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
52: PetscRandomSetFromOptions(rctx);
53: VecSetRandom(u,rctx);
54: }
56: /* Create and assemble matrix */
57: MatCreate(PETSC_COMM_WORLD,&A);
58: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,dim,dim);
59: MatSetFromOptions(A);
60: MatSetUp(A);
61: FormTestMatrix(A,n,type);
62: MatMult(A,u,b);
63: flg = PETSC_FALSE;
64: PetscOptionsGetBool(NULL,NULL,"-printout",&flg,NULL);
65: if (flg) {
66: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
67: VecView(u,PETSC_VIEWER_STDOUT_WORLD);
68: VecView(b,PETSC_VIEWER_STDOUT_WORLD);
69: }
71: /* Create KSP context; set operators and options; solve linear system */
72: KSPCreate(PETSC_COMM_WORLD,&ksp);
73: KSPSetOperators(ksp,A,A);
74: KSPSetFromOptions(ksp);
75: KSPSolve(ksp,b,x);
76: /* KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD); */
78: /* Check error */
79: VecAXPY(x,none,u);
80: VecNorm(x,NORM_2,&norm);
81: KSPGetIterationNumber(ksp,&its);
82: if (norm >= 1.e-12) {
83: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g, Iterations %D\n",(double)norm,its);
84: } else {
85: PetscPrintf(PETSC_COMM_WORLD,"Norm of error < 1.e-12, Iterations %D\n",its);
86: }
88: /* Free work space */
89: VecDestroy(&x); VecDestroy(&u);
90: VecDestroy(&b); MatDestroy(&A);
91: if (use_random) {PetscRandomDestroy(&rctx);}
92: KSPDestroy(&ksp);
93: PetscFinalize();
94: return 0;
95: }
99: PetscErrorCode FormTestMatrix(Mat A,PetscInt n,TestType type)
100: {
101: #if !defined(PETSC_USE_COMPLEX)
102: SETERRQ(PetscObjectComm((PetscObject)A),1,"FormTestMatrix: These problems require complex numbers.");
103: #else
105: PetscScalar val[5];
107: PetscInt i,j,Ii,J,col[5],Istart,Iend;
109: MatGetOwnershipRange(A,&Istart,&Iend);
110: if (type == TEST_1) {
111: val[0] = 1.0; val[1] = 4.0; val[2] = -2.0;
112: for (i=1; i<n-1; i++) {
113: col[0] = i-1; col[1] = i; col[2] = i+1;
114: MatSetValues(A,1,&i,3,col,val,INSERT_VALUES);
115: }
116: i = n-1; col[0] = n-2; col[1] = n-1;
117: MatSetValues(A,1,&i,2,col,val,INSERT_VALUES);
118: i = 0; col[0] = 0; col[1] = 1; val[0] = 4.0; val[1] = -2.0;
119: MatSetValues(A,1,&i,2,col,val,INSERT_VALUES);
120: } else if (type == TEST_2) {
121: val[0] = 1.0; val[1] = 0.0; val[2] = 2.0; val[3] = 1.0;
122: for (i=2; i<n-1; i++) {
123: col[0] = i-2; col[1] = i-1; col[2] = i; col[3] = i+1;
124: MatSetValues(A,1,&i,4,col,val,INSERT_VALUES);
125: }
126: i = n-1; col[0] = n-3; col[1] = n-2; col[2] = n-1;
127: MatSetValues(A,1,&i,3,col,val,INSERT_VALUES);
128: i = 1; col[0] = 0; col[1] = 1; col[2] = 2;
129: MatSetValues(A,1,&i,3,col,&val[1],INSERT_VALUES);
130: i = 0;
131: MatSetValues(A,1,&i,2,col,&val[2],INSERT_VALUES);
132: } else if (type == TEST_3) {
133: val[0] = PETSC_i * 2.0;
134: val[1] = 4.0; val[2] = 0.0; val[3] = 1.0; val[4] = 0.7;
135: for (i=1; i<n-3; i++) {
136: col[0] = i-1; col[1] = i; col[2] = i+1; col[3] = i+2; col[4] = i+3;
137: MatSetValues(A,1,&i,5,col,val,INSERT_VALUES);
138: }
139: i = n-3; col[0] = n-4; col[1] = n-3; col[2] = n-2; col[3] = n-1;
140: MatSetValues(A,1,&i,4,col,val,INSERT_VALUES);
141: i = n-2; col[0] = n-3; col[1] = n-2; col[2] = n-1;
142: MatSetValues(A,1,&i,3,col,val,INSERT_VALUES);
143: i = n-1; col[0] = n-2; col[1] = n-1;
144: MatSetValues(A,1,&i,2,col,val,INSERT_VALUES);
145: i = 0; col[0] = 0; col[1] = 1; col[2] = 2; col[3] = 3;
146: MatSetValues(A,1,&i,4,col,&val[1],INSERT_VALUES);
147: } else if (type == HELMHOLTZ_1) {
148: /* Problem domain: unit square: (0,1) x (0,1)
149: Solve Helmholtz equation:
150: -delta u - sigma1*u + i*sigma2*u = f,
151: where delta = Laplace operator
152: Dirichlet b.c.'s on all sides
153: */
154: PetscRandom rctx;
155: PetscReal h2,sigma1 = 5.0;
156: PetscScalar sigma2;
157: PetscOptionsGetReal(NULL,NULL,"-sigma1",&sigma1,NULL);
158: PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
159: PetscRandomSetFromOptions(rctx);
160: PetscRandomSetInterval(rctx,0.0,PETSC_i);
161: h2 = 1.0/((n+1)*(n+1));
162: for (Ii=Istart; Ii<Iend; Ii++) {
163: *val = -1.0; i = Ii/n; j = Ii - i*n;
164: if (i>0) {
165: J = Ii-n; MatSetValues(A,1,&Ii,1,&J,val,ADD_VALUES);
166: }
167: if (i<n-1) {
168: J = Ii+n; MatSetValues(A,1,&Ii,1,&J,val,ADD_VALUES);
169: }
170: if (j>0) {
171: J = Ii-1; MatSetValues(A,1,&Ii,1,&J,val,ADD_VALUES);
172: }
173: if (j<n-1) {
174: J = Ii+1; MatSetValues(A,1,&Ii,1,&J,val,ADD_VALUES);
175: }
176: PetscRandomGetValue(rctx,&sigma2);
177: *val = 4.0 - sigma1*h2 + sigma2*h2;
178: MatSetValues(A,1,&Ii,1,&Ii,val,ADD_VALUES);
179: }
180: PetscRandomDestroy(&rctx);
181: } else if (type == HELMHOLTZ_2) {
182: /* Problem domain: unit square: (0,1) x (0,1)
183: Solve Helmholtz equation:
184: -delta u - sigma1*u = f,
185: where delta = Laplace operator
186: Dirichlet b.c.'s on 3 sides
187: du/dn = i*alpha*u on (1,y), 0<y<1
188: */
189: PetscReal h2,sigma1 = 200.0;
190: PetscScalar alpha_h;
191: PetscOptionsGetReal(NULL,NULL,"-sigma1",&sigma1,NULL);
192: h2 = 1.0/((n+1)*(n+1));
193: alpha_h = (PETSC_i * 10.0) / (PetscReal)(n+1); /* alpha_h = alpha * h */
194: for (Ii=Istart; Ii<Iend; Ii++) {
195: *val = -1.0; i = Ii/n; j = Ii - i*n;
196: if (i>0) {
197: J = Ii-n; MatSetValues(A,1,&Ii,1,&J,val,ADD_VALUES);
198: }
199: if (i<n-1) {
200: J = Ii+n; MatSetValues(A,1,&Ii,1,&J,val,ADD_VALUES);
201: }
202: if (j>0) {
203: J = Ii-1; MatSetValues(A,1,&Ii,1,&J,val,ADD_VALUES);
204: }
205: if (j<n-1) {
206: J = Ii+1; MatSetValues(A,1,&Ii,1,&J,val,ADD_VALUES);
207: }
208: *val = 4.0 - sigma1*h2;
209: if (!((Ii+1)%n)) *val += alpha_h;
210: MatSetValues(A,1,&Ii,1,&Ii,val,ADD_VALUES);
211: }
212: } else SETERRQ(PetscObjectComm((PetscObject)A),1,"FormTestMatrix: unknown test matrix type");
214: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
215: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
216: #endif
218: return 0;
219: }