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(u,pfive);
 49:   } else {
 50:     PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
 51:     PetscRandomSetFromOptions(rctx);
 52:     VecSetRandom(u,rctx);
 53:   }

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

 68:   /* Create KSP context; set operators and options; solve linear system */
 69:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 70:   KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN);
 71:   KSPSetFromOptions(ksp);
 72:   KSPSolve(ksp,b,x);
 73:   KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD);

 75:   /* Check error */
 76:   VecAXPY(x,none,u);
 77:   VecNorm(x,NORM_2,&norm);
 78:   KSPGetIterationNumber(ksp,&its);
 79:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %A,Iterations %D\n",norm,its);

 81:   /* Free work space */
 82:   VecDestroy(x); VecDestroy(u);
 83:   VecDestroy(b); MatDestroy(A);
 84:   if (use_random) {PetscRandomDestroy(rctx);}
 85:   KSPDestroy(ksp);
 86:   PetscFinalize();
 87:   return 0;
 88: }

 92: PetscErrorCode FormTestMatrix(Mat A,PetscInt n,TestType type)
 93: {
 94: #if !defined(PETSC_USE_COMPLEX)
 95:   SETERRQ(1,"FormTestMatrix: These problems require complex numbers.");
 96: #else

 98:   PetscScalar    val[5];
100:   PetscInt       i,j,Ii,J,col[5],Istart,Iend;

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

203:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
204:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
205: #endif

207:   return 0;
208: }