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:   flg        = PETSC_FALSE;
 46:   PetscOptionsGetTruth(PETSC_NULL,"-norandom",&flg,PETSC_NULL);
 47:   if (flg) {
 48:     use_random = 0;
 49:     VecSet(u,pfive);
 50:   } else {
 51:     PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
 52:     PetscRandomSetFromOptions(rctx);
 53:     VecSetRandom(u,rctx);
 54:   }

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

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

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

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

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

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

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

206:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
207:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
208: #endif

210:   return 0;
211: }