Actual source code: ex17.c

  1: /*$Id: ex17.c,v 1.40 2001/03/23 23:23:50 balay 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:   Scalar      none = -1.0,pfive = 0.5;
 18:   double      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,PETSC_DECIDE,dim,&x);
 37:   VecSetFromOptions(x);
 38:   VecDuplicate(x,&b);
 39:   VecDuplicate(x,&u);

 41:   use_random = 1;
 42:   PetscOptionsHasName(PETSC_NULL,"-norandom",&flg);
 43:   if (flg) {
 44:     use_random = 0;
 45:     VecSet(&pfive,u);
 46:   } else {
 47:     PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT,&rctx);
 48:     VecSetRandom(rctx,u);
 49:   }

 51:   /* Create and assemble matrix */
 52:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,dim,dim,&A);
 53:   MatSetFromOptions(A);
 54:   FormTestMatrix(A,n,type);
 55:   MatMult(A,u,b);
 56:   PetscOptionsHasName(PETSC_NULL,"-printout",&flg);
 57:   if (flg) {
 58:     MatView(A,PETSC_VIEWER_STDOUT_WORLD);
 59:     VecView(u,PETSC_VIEWER_STDOUT_WORLD);
 60:     VecView(b,PETSC_VIEWER_STDOUT_WORLD);
 61:   }

 63:   /* Create SLES context; set operators and options; solve linear system */
 64:   SLESCreate(PETSC_COMM_WORLD,&sles);
 65:   SLESSetOperators(sles,A,A,DIFFERENT_NONZERO_PATTERN);
 66: 
 67:   SLESSetFromOptions(sles);
 68:   SLESSolve(sles,b,x,&its);
 69:   SLESView(sles,PETSC_VIEWER_STDOUT_WORLD);

 71:   /* Check error */
 72:   VecAXPY(&none,u,x);
 73:   ierr  = VecNorm(x,NORM_2,&norm);
 74:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A,Iterations %dn",norm,its);

 76:   /* Free work space */
 77:   VecDestroy(x); VecDestroy(u);
 78:   VecDestroy(b); MatDestroy(A);
 79:   if (use_random) {PetscRandomDestroy(rctx);}
 80:   SLESDestroy(sles);
 81:   PetscFinalize();
 82:   return 0;
 83: }

 85: int FormTestMatrix(Mat A,int n,TestType type)
 86: {
 87: #if !defined(PETSC_USE_COMPLEX)
 88:   SETERRQ(1,"FormTestMatrix: These problems require complex numbers.");
 89: #else

 91:   Scalar val[5],h;
 92:   int    i,j,I,J,ierr,col[5],Istart,Iend;

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

194:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
195:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
196: #endif

198:   return 0;
199: }