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 is\n\
  4: intended to test the complex numbers version of various solvers.\n\n";

 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);

 13: int main(int argc,char **args)
 14: {
 15:   Vec         x,b,u;      /* approx solution, RHS, exact solution */
 16:   Mat         A;            /* linear system matrix */
 17:   SLES        sles;         /* SLES context */
 18:   int         ierr,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 SLES context; set operators and options; solve linear system */
 67:   SLESCreate(PETSC_COMM_WORLD,&sles);
 68:   SLESSetOperators(sles,A,A,DIFFERENT_NONZERO_PATTERN);
 69: 
 70:   SLESSetFromOptions(sles);
 71:   SLESSolve(sles,b,x,&its);
 72:   SLESView(sles,PETSC_VIEWER_STDOUT_WORLD);

 74:   /* Check error */
 75:   VecAXPY(&none,u,x);
 76:   VecNorm(x,NORM_2,&norm);
 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:   SLESDestroy(sles);
 84:   PetscFinalize();
 85:   return 0;
 86: }

 90: int FormTestMatrix(Mat A,int 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];
 97:   int    i,j,I,J,ierr,col[5],Istart,Iend;

 99:   MatGetOwnershipRange(A,&Istart,&Iend);
100:   if (type == TEST_1) {
101:     val[0] = 1.0; val[1] = 4.0; val[2] = -2.0;
102:     for (i=1; i<n-1; i++) {
103:       col[0] = i-1; col[1] = i; col[2] = i+1;
104:       MatSetValues(A,1,&i,3,col,val,INSERT_VALUES);
105:     }
106:     i = n-1; col[0] = n-2; col[1] = n-1;
107:     MatSetValues(A,1,&i,2,col,val,INSERT_VALUES);
108:     i = 0; col[0] = 0; col[1] = 1; val[0] = 4.0; val[1] = -2.0;
109:     MatSetValues(A,1,&i,2,col,val,INSERT_VALUES);
110:   }
111:   else if (type == TEST_2) {
112:     val[0] = 1.0; val[1] = 0.0; val[2] = 2.0; val[3] = 1.0;
113:     for (i=2; i<n-1; i++) {
114:       col[0] = i-2; col[1] = i-1; col[2] = i; col[3] = i+1;
115:       MatSetValues(A,1,&i,4,col,val,INSERT_VALUES);
116:     }
117:     i = n-1; col[0] = n-3; col[1] = n-2; col[2] = n-1;
118:     MatSetValues(A,1,&i,3,col,val,INSERT_VALUES);
119:     i = 1; col[0] = 0; col[1] = 1; col[2] = 2;
120:     MatSetValues(A,1,&i,3,col,&val[1],INSERT_VALUES);
121:     i = 0;
122:     MatSetValues(A,1,&i,2,col,&val[2],INSERT_VALUES);
123:   }
124:   else if (type == TEST_3) {
125:     val[0] = PETSC_i * 2.0;
126:     val[1] = 4.0; val[2] = 0.0; val[3] = 1.0; val[4] = 0.7;
127:     for (i=1; i<n-3; i++) {
128:       col[0] = i-1; col[1] = i; col[2] = i+1; col[3] = i+2; col[4] = i+3;
129:       MatSetValues(A,1,&i,5,col,val,INSERT_VALUES);
130:     }
131:     i = n-3; col[0] = n-4; col[1] = n-3; col[2] = n-2; col[3] = n-1;
132:     MatSetValues(A,1,&i,4,col,val,INSERT_VALUES);
133:     i = n-2; col[0] = n-3; col[1] = n-2; col[2] = n-1;
134:     MatSetValues(A,1,&i,3,col,val,INSERT_VALUES);
135:     i = n-1; col[0] = n-2; col[1] = n-1;
136:     MatSetValues(A,1,&i,2,col,val,INSERT_VALUES);
137:     i = 0; col[0] = 0; col[1] = 1; col[2] = 2; col[3] = 3;
138:     MatSetValues(A,1,&i,4,col,&val[1],INSERT_VALUES);
139:   }
140:   else if (type == HELMHOLTZ_1) {
141:     /* Problem domain: unit square: (0,1) x (0,1)
142:        Solve Helmholtz equation:
143:           -delta u - sigma1*u + i*sigma2*u = f, 
144:            where delta = Laplace operator
145:        Dirichlet b.c.'s on all sides
146:      */
147:     PetscRandom rctx;
148:     PetscReal   h2,sigma1 = 5.0;
149:     PetscScalar sigma2;
150:     PetscOptionsGetReal(PETSC_NULL,"-sigma1",&sigma1,PETSC_NULL);
151:     PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT_IMAGINARY,&rctx);
152:     h2 = 1.0/((n+1)*(n+1));
153:     for (I=Istart; I<Iend; I++) {
154:       *val = -1.0; i = I/n; j = I - i*n;
155:       if (i>0) {
156:         J = I-n; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
157:       if (i<n-1) {
158:         J = I+n; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
159:       if (j>0) {
160:         J = I-1; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
161:       if (j<n-1) {
162:         J = I+1; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
163:       PetscRandomGetValue(rctx,&sigma2);
164:       *val = 4.0 - sigma1*h2 + sigma2*h2;
165:       MatSetValues(A,1,&I,1,&I,val,ADD_VALUES);
166:     }
167:     PetscRandomDestroy(rctx);
168:   }
169:   else if (type == HELMHOLTZ_2) {
170:     /* Problem domain: unit square: (0,1) x (0,1)
171:        Solve Helmholtz equation:
172:           -delta u - sigma1*u = f, 
173:            where delta = Laplace operator
174:        Dirichlet b.c.'s on 3 sides
175:        du/dn = i*alpha*u on (1,y), 0<y<1
176:      */
177:     PetscReal   h2,sigma1 = 200.0;
178:     PetscScalar alpha_h;
179:     PetscOptionsGetReal(PETSC_NULL,"-sigma1",&sigma1,PETSC_NULL);
180:     h2 = 1.0/((n+1)*(n+1));
181:     alpha_h = (PETSC_i * 10.0) / (PetscReal)(n+1);  /* alpha_h = alpha * h */
182:     for (I=Istart; I<Iend; I++) {
183:       *val = -1.0; i = I/n; j = I - i*n;
184:       if (i>0) {
185:         J = I-n; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
186:       if (i<n-1) {
187:         J = I+n; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
188:       if (j>0) {
189:         J = I-1; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
190:       if (j<n-1) {
191:         J = I+1; MatSetValues(A,1,&I,1,&J,val,ADD_VALUES);}
192:       *val = 4.0 - sigma1*h2;
193:       if (!((I+1)%n)) *val += alpha_h;
194:       MatSetValues(A,1,&I,1,&I,val,ADD_VALUES);
195:     }
196:   }
197:   else SETERRQ(1,"FormTestMatrix: unknown test matrix type");

199:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
200:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
201: #endif

203:   return 0;
204: }