Actual source code: ex17.c
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;
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: PetscTruth flg;
25: PetscInitialize(&argc,&args,(char *)0,help);
26: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
27: PetscOptionsGetInt(PETSC_NULL,"-p",&p,PETSC_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: PetscOptionsHasName(PETSC_NULL,"-norandom",&flg);
46: if (flg) {
47: use_random = 0;
48: VecSet(&pfive,u);
49: } else {
50: PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT,&rctx);
51: VecSetRandom(rctx,u);
52: }
54: /* Create and assemble matrix */
55: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,dim,dim,&A);
56: MatSetFromOptions(A);
57: FormTestMatrix(A,n,type);
58: MatMult(A,u,b);
59: PetscOptionsHasName(PETSC_NULL,"-printout",&flg);
60: if (flg) {
61: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
62: VecView(u,PETSC_VIEWER_STDOUT_WORLD);
63: VecView(b,PETSC_VIEWER_STDOUT_WORLD);
64: }
66: /* Create KSP context; set operators and options; solve linear system */
67: KSPCreate(PETSC_COMM_WORLD,&ksp);
68: KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN);
69: KSPSetFromOptions(ksp);
70: KSPSolve(ksp,b,x);
71: KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD);
73: /* Check error */
74: VecAXPY(&none,u,x);
75: VecNorm(x,NORM_2,&norm);
76: KSPGetIterationNumber(ksp,&its);
77: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A,Iterations %D\n",norm,its);
79: /* Free work space */
80: VecDestroy(x); VecDestroy(u);
81: VecDestroy(b); MatDestroy(A);
82: if (use_random) {PetscRandomDestroy(rctx);}
83: KSPDestroy(ksp);
84: PetscFinalize();
85: return 0;
86: }
90: PetscErrorCode FormTestMatrix(Mat A,PetscInt n,TestType type)
91: {
92: #if !defined(PETSC_USE_COMPLEX)
93: SETERRQ(1,"FormTestMatrix: These problems require complex numbers.");
94: #else
96: PetscScalar val[5];
98: PetscInt i,j,I,J,col[5],Istart,Iend;
100: MatGetOwnershipRange(A,&Istart,&Iend);
101: if (type == TEST_1) {
102: val[0] = 1.0; val[1] = 4.0; val[2] = -2.0;
103: for (i=1; i<n-1; i++) {
104: col[0] = i-1; col[1] = i; col[2] = i+1;
105: MatSetValues(A,1,&i,3,col,val,INSERT_VALUES);
106: }
107: i = n-1; col[0] = n-2; col[1] = n-1;
108: MatSetValues(A,1,&i,2,col,val,INSERT_VALUES);
109: i = 0; col[0] = 0; col[1] = 1; val[0] = 4.0; val[1] = -2.0;
110: MatSetValues(A,1,&i,2,col,val,INSERT_VALUES);
111: }
112: else if (type == TEST_2) {
113: val[0] = 1.0; val[1] = 0.0; val[2] = 2.0; val[3] = 1.0;
114: for (i=2; i<n-1; i++) {
115: col[0] = i-2; col[1] = i-1; col[2] = i; col[3] = i+1;
116: MatSetValues(A,1,&i,4,col,val,INSERT_VALUES);
117: }
118: i = n-1; col[0] = n-3; col[1] = n-2; col[2] = n-1;
119: MatSetValues(A,1,&i,3,col,val,INSERT_VALUES);
120: i = 1; col[0] = 0; col[1] = 1; col[2] = 2;
121: MatSetValues(A,1,&i,3,col,&val[1],INSERT_VALUES);
122: i = 0;
123: MatSetValues(A,1,&i,2,col,&val[2],INSERT_VALUES);
124: }
125: else if (type == TEST_3) {
126: val[0] = PETSC_i * 2.0;
127: val[1] = 4.0; val[2] = 0.0; val[3] = 1.0; val[4] = 0.7;
128: for (i=1; i<n-3; i++) {
129: col[0] = i-1; col[1] = i; col[2] = i+1; col[3] = i+2; col[4] = i+3;
130: MatSetValues(A,1,&i,5,col,val,INSERT_VALUES);
131: }
132: i = n-3; col[0] = n-4; col[1] = n-3; col[2] = n-2; col[3] = n-1;
133: MatSetValues(A,1,&i,4,col,val,INSERT_VALUES);
134: i = n-2; col[0] = n-3; col[1] = n-2; col[2] = n-1;
135: MatSetValues(A,1,&i,3,col,val,INSERT_VALUES);
136: i = n-1; col[0] = n-2; col[1] = n-1;
137: MatSetValues(A,1,&i,2,col,val,INSERT_VALUES);
138: i = 0; col[0] = 0; col[1] = 1; col[2] = 2; col[3] = 3;
139: MatSetValues(A,1,&i,4,col,&val[1],INSERT_VALUES);
140: }
141: else if (type == HELMHOLTZ_1) {
142: /* Problem domain: unit square: (0,1) x (0,1)
143: Solve Helmholtz equation:
144: -delta u - sigma1*u + i*sigma2*u = f,
145: where delta = Laplace operator
146: Dirichlet b.c.'s on all sides
147: */
148: PetscRandom rctx;
149: PetscReal h2,sigma1 = 5.0;
150: PetscScalar sigma2;
151: PetscOptionsGetReal(PETSC_NULL,"-sigma1",&sigma1,PETSC_NULL);
152: PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT_IMAGINARY,&rctx);
153: h2 = 1.0/((n+1)*(n+1));
154: for (I=Istart; I<Iend; I++) {
155: *val = -1.0; i = I/n; j = I - i*n;
156: if (i>0) {
157: J = I-n; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
158: if (i<n-1) {
159: J = I+n; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
160: if (j>0) {
161: J = I-1; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
162: if (j<n-1) {
163: J = I+1; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
164: PetscRandomGetValue(rctx,&sigma2);
165: *val = 4.0 - sigma1*h2 + sigma2*h2;
166: MatSetValues(A,1,&I,1,&I,val,ADD_VALUES);
167: }
168: PetscRandomDestroy(rctx);
169: }
170: else if (type == HELMHOLTZ_2) {
171: /* Problem domain: unit square: (0,1) x (0,1)
172: Solve Helmholtz equation:
173: -delta u - sigma1*u = f,
174: where delta = Laplace operator
175: Dirichlet b.c.'s on 3 sides
176: du/dn = i*alpha*u on (1,y), 0<y<1
177: */
178: PetscReal h2,sigma1 = 200.0;
179: PetscScalar alpha_h;
180: PetscOptionsGetReal(PETSC_NULL,"-sigma1",&sigma1,PETSC_NULL);
181: h2 = 1.0/((n+1)*(n+1));
182: alpha_h = (PETSC_i * 10.0) / (PetscReal)(n+1); /* alpha_h = alpha * h */
183: for (I=Istart; I<Iend; I++) {
184: *val = -1.0; i = I/n; j = I - i*n;
185: if (i>0) {
186: J = I-n; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
187: if (i<n-1) {
188: J = I+n; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
189: if (j>0) {
190: J = I-1; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
191: if (j<n-1) {
192: J = I+1; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
193: *val = 4.0 - sigma1*h2;
194: if (!((I+1)%n)) *val += alpha_h;
195: MatSetValues(A,1,&I,1,&I,val,ADD_VALUES);
196: }
197: }
198: else SETERRQ(1,"FormTestMatrix: unknown test matrix type");
200: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
201: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
202: #endif
204: return 0;
205: }