Actual source code: ex17.c
1: /*$Id: ex17.c,v 1.43 2001/08/07 21:30:50 bsmith Exp $*/
3: static char help[] = "Solves a linear system with SLES. This problem isn
4: intended to test the complex numbers version of various solvers.nn";
6: #include petscsles.h
8: typedef enum {TEST_1,TEST_2,TEST_3,HELMHOLTZ_1,HELMHOLTZ_2} TestType;
9: extern int FormTestMatrix(Mat,int,TestType);
11: int main(int argc,char **args)
12: {
13: Vec x,b,u; /* approx solution, RHS, exact solution */
14: Mat A; /* linear system matrix */
15: SLES sles; /* SLES context */
16: int ierr,n = 10,its, dim,p = 1,use_random;
17: PetscScalar none = -1.0,pfive = 0.5;
18: PetscReal norm;
19: PetscRandom rctx;
20: TestType type;
21: PetscTruth flg;
23: PetscInitialize(&argc,&args,(char *)0,help);
24: PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
25: PetscOptionsGetInt(PETSC_NULL,"-p",&p,PETSC_NULL);
26: switch (p) {
27: case 1: type = TEST_1; dim = n; break;
28: case 2: type = TEST_2; dim = n; break;
29: case 3: type = TEST_3; dim = n; break;
30: case 4: type = HELMHOLTZ_1; dim = n*n; break;
31: case 5: type = HELMHOLTZ_2; dim = n*n; break;
32: default: type = TEST_1; dim = n;
33: }
35: /* Create vectors */
36: VecCreate(PETSC_COMM_WORLD,&x);
37: VecSetSizes(x,PETSC_DECIDE,dim);
38: VecSetFromOptions(x);
39: VecDuplicate(x,&b);
40: VecDuplicate(x,&u);
42: use_random = 1;
43: PetscOptionsHasName(PETSC_NULL,"-norandom",&flg);
44: if (flg) {
45: use_random = 0;
46: VecSet(&pfive,u);
47: } else {
48: PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT,&rctx);
49: VecSetRandom(rctx,u);
50: }
52: /* Create and assemble matrix */
53: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,dim,dim,&A);
54: MatSetFromOptions(A);
55: FormTestMatrix(A,n,type);
56: MatMult(A,u,b);
57: PetscOptionsHasName(PETSC_NULL,"-printout",&flg);
58: if (flg) {
59: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
60: VecView(u,PETSC_VIEWER_STDOUT_WORLD);
61: VecView(b,PETSC_VIEWER_STDOUT_WORLD);
62: }
64: /* Create SLES context; set operators and options; solve linear system */
65: SLESCreate(PETSC_COMM_WORLD,&sles);
66: SLESSetOperators(sles,A,A,DIFFERENT_NONZERO_PATTERN);
67:
68: SLESSetFromOptions(sles);
69: SLESSolve(sles,b,x,&its);
70: SLESView(sles,PETSC_VIEWER_STDOUT_WORLD);
72: /* Check error */
73: VecAXPY(&none,u,x);
74: ierr = VecNorm(x,NORM_2,&norm);
75: PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A,Iterations %dn",norm,its);
77: /* Free work space */
78: VecDestroy(x); VecDestroy(u);
79: VecDestroy(b); MatDestroy(A);
80: if (use_random) {PetscRandomDestroy(rctx);}
81: SLESDestroy(sles);
82: PetscFinalize();
83: return 0;
84: }
86: int FormTestMatrix(Mat A,int n,TestType type)
87: {
88: #if !defined(PETSC_USE_COMPLEX)
89: SETERRQ(1,"FormTestMatrix: These problems require complex numbers.");
90: #else
92: PetscScalar val[5];
93: int i,j,I,J,ierr,col[5],Istart,Iend;
95: MatGetOwnershipRange(A,&Istart,&Iend);
96: if (type == TEST_1) {
97: val[0] = 1.0; val[1] = 4.0; val[2] = -2.0;
98: for (i=1; i<n-1; i++) {
99: col[0] = i-1; col[1] = i; col[2] = i+1;
100: MatSetValues(A,1,&i,3,col,val,INSERT_VALUES);
101: }
102: i = n-1; col[0] = n-2; col[1] = n-1;
103: MatSetValues(A,1,&i,2,col,val,INSERT_VALUES);
104: i = 0; col[0] = 0; col[1] = 1; val[0] = 4.0; val[1] = -2.0;
105: MatSetValues(A,1,&i,2,col,val,INSERT_VALUES);
106: }
107: else if (type == TEST_2) {
108: val[0] = 1.0; val[1] = 0.0; val[2] = 2.0; val[3] = 1.0;
109: for (i=2; i<n-1; i++) {
110: col[0] = i-2; col[1] = i-1; col[2] = i; col[3] = i+1;
111: MatSetValues(A,1,&i,4,col,val,INSERT_VALUES);
112: }
113: i = n-1; col[0] = n-3; col[1] = n-2; col[2] = n-1;
114: MatSetValues(A,1,&i,3,col,val,INSERT_VALUES);
115: i = 1; col[0] = 0; col[1] = 1; col[2] = 2;
116: MatSetValues(A,1,&i,3,col,&val[1],INSERT_VALUES);
117: i = 0;
118: MatSetValues(A,1,&i,2,col,&val[2],INSERT_VALUES);
119: }
120: else if (type == TEST_3) {
121: val[0] = PETSC_i * 2.0;
122: val[1] = 4.0; val[2] = 0.0; val[3] = 1.0; val[4] = 0.7;
123: for (i=1; i<n-3; i++) {
124: col[0] = i-1; col[1] = i; col[2] = i+1; col[3] = i+2; col[4] = i+3;
125: MatSetValues(A,1,&i,5,col,val,INSERT_VALUES);
126: }
127: i = n-3; col[0] = n-4; col[1] = n-3; col[2] = n-2; col[3] = n-1;
128: MatSetValues(A,1,&i,4,col,val,INSERT_VALUES);
129: i = n-2; col[0] = n-3; col[1] = n-2; col[2] = n-1;
130: MatSetValues(A,1,&i,3,col,val,INSERT_VALUES);
131: i = n-1; col[0] = n-2; col[1] = n-1;
132: MatSetValues(A,1,&i,2,col,val,INSERT_VALUES);
133: i = 0; col[0] = 0; col[1] = 1; col[2] = 2; col[3] = 3;
134: MatSetValues(A,1,&i,4,col,&val[1],INSERT_VALUES);
135: }
136: else if (type == HELMHOLTZ_1) {
137: /* Problem domain: unit square: (0,1) x (0,1)
138: Solve Helmholtz equation:
139: -delta u - sigma1*u + i*sigma2*u = f,
140: where delta = Laplace operator
141: Dirichlet b.c.'s on all sides
142: */
143: PetscRandom rctx;
144: PetscReal h2,sigma1 = 5.0;
145: PetscScalar sigma2;
146: PetscOptionsGetReal(PETSC_NULL,"-sigma1",&sigma1,PETSC_NULL);
147: PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT_IMAGINARY,&rctx);
148: h2 = 1.0/((n+1)*(n+1));
149: for (I=Istart; I<Iend; I++) {
150: *val = -1.0; i = I/n; j = I - i*n;
151: if (i>0) {
152: J = I-n; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
153: if (i<n-1) {
154: J = I+n; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
155: if (j>0) {
156: J = I-1; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
157: if (j<n-1) {
158: J = I+1; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
159: PetscRandomGetValue(rctx,&sigma2);
160: *val = 4.0 - sigma1*h2 + sigma2*h2;
161: MatSetValues(A,1,&I,1,&I,val,ADD_VALUES);
162: }
163: PetscRandomDestroy(rctx);
164: }
165: else if (type == HELMHOLTZ_2) {
166: /* Problem domain: unit square: (0,1) x (0,1)
167: Solve Helmholtz equation:
168: -delta u - sigma1*u = f,
169: where delta = Laplace operator
170: Dirichlet b.c.'s on 3 sides
171: du/dn = i*alpha*u on (1,y), 0<y<1
172: */
173: PetscReal h2,sigma1 = 200.0;
174: PetscScalar alpha_h;
175: PetscOptionsGetReal(PETSC_NULL,"-sigma1",&sigma1,PETSC_NULL);
176: h2 = 1.0/((n+1)*(n+1));
177: alpha_h = (PETSC_i * 10.0) / (PetscReal)(n+1); /* alpha_h = alpha * h */
178: for (I=Istart; I<Iend; I++) {
179: *val = -1.0; i = I/n; j = I - i*n;
180: if (i>0) {
181: J = I-n; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
182: if (i<n-1) {
183: J = I+n; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
184: if (j>0) {
185: J = I-1; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
186: if (j<n-1) {
187: J = I+1; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
188: *val = 4.0 - sigma1*h2;
189: if (!((I+1)%n)) *val += alpha_h;
190: MatSetValues(A,1,&I,1,&I,val,ADD_VALUES);
191: }
192: }
193: else SETERRQ(1,"FormTestMatrix: unknown test matrix type");
195: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
196: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
197: #endif
199: return 0;
200: }