Actual source code: ex70.c

petsc-3.10.5 2019-03-28
Report Typos and Errors
  1: static char help[] = "Poiseuille flow problem. Viscous, laminar flow in a 2D channel with parabolic velocity\n\
  2:                       profile and linear pressure drop, exact solution of the 2D Stokes\n";

  4: /*---------------------------------------------------------------------------- */
  5: /* M A R I T I M E  R E S E A R C H  I N S T I T U T E  N E T H E R L A N D S  */
  6: /*---------------------------------------------------------------------------- */
  7: /* author : Christiaan M. Klaij                                                */
  8: /*---------------------------------------------------------------------------- */
  9: /*                                                                             */
 10: /* Poiseuille flow problem.                                                    */
 11: /*                                                                             */
 12: /* Viscous, laminar flow in a 2D channel with parabolic velocity               */
 13: /* profile and linear pressure drop, exact solution of the 2D Stokes           */
 14: /* equations.                                                                  */
 15: /*                                                                             */
 16: /* Discretized with the cell-centered finite-volume method on a                */
 17: /* Cartesian grid with co-located variables. Variables ordered as              */
 18: /* [u1...uN v1...vN p1...pN]^T. Matrix [A00 A01; A10, A11] solved with         */
 19: /* PCFIELDSPLIT. Lower factorization is used to mimick the Semi-Implicit       */
 20: /* Method for Pressure Linked Equations (SIMPLE) used as preconditioner        */
 21: /* instead of solver.                                                          */
 22: /*                                                                             */
 23: /* Disclaimer: does not contain the pressure-weighed interpolation             */
 24: /* method needed to suppress spurious pressure modes in real-life              */
 25: /* problems.                                                                   */
 26: /*                                                                             */
 27: /* Usage:                                                                      */
 28: /*                                                                             */
 29: /* mpiexec -n 2 ./ex70 -nx 32 -ny 48 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -fieldsplit_1_pc_type none */
 30: /*                                                                             */
 31: /*   Runs with PCFIELDSPLIT on 32x48 grid, no PC for the Schur                 */
 32: /*   complement because A11 is zero. FGMRES is needed because                  */
 33: /*   PCFIELDSPLIT is a variable preconditioner.                                */
 34: /*                                                                             */
 35: /* mpiexec -n 2 ./ex70 -nx 32 -ny 48 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -user_pc */
 36: /*                                                                             */
 37: /*   Same as above but with user defined PC for the true Schur                 */
 38: /*   complement. PC based on the SIMPLE-type approximation (inverse of         */
 39: /*   A00 approximated by inverse of its diagonal).                             */
 40: /*                                                                             */
 41: /* mpiexec -n 2 ./ex70 -nx 32 -ny 48 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -user_ksp */
 42: /*                                                                             */
 43: /*   Replace the true Schur complement with a user defined Schur               */
 44: /*   complement based on the SIMPLE-type approximation. Same matrix is         */
 45: /*   used as PC.                                                               */
 46: /*                                                                             */
 47: /* mpiexec -n 2 ./ex70 -nx 32 -ny 48 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -fieldsplit_0_ksp_type gmres -fieldsplit_0_pc_type bjacobi -fieldsplit_1_pc_type jacobi -fieldsplit_1_inner_ksp_type preonly -fieldsplit_1_inner_pc_type jacobi -fieldsplit_1_upper_ksp_type preonly -fieldsplit_1_upper_pc_type jacobi */
 48: /*                                                                             */
 49: /*   Out-of-the-box SIMPLE-type preconditioning. The major advantage           */
 50: /*   is that the user neither needs to provide the approximation of            */
 51: /*   the Schur complement, nor the corresponding preconditioner.               */
 52: /*                                                                             */
 53: /*---------------------------------------------------------------------------- */


 56:  #include <petscksp.h>

 58: typedef struct {
 59:   PetscBool userPC, userKSP; /* user defined preconditioner and matrix for the Schur complement */
 60:   PetscInt  nx, ny;  /* nb of cells in x- and y-direction */
 61:   PetscReal hx, hy;  /* mesh size in x- and y-direction */
 62:   Mat       A;       /* block matrix */
 63:   Mat       subA[4]; /* the four blocks */
 64:   Mat       myS;     /* the approximation of the Schur complement */
 65:   Vec       x, b, y; /* solution, rhs and temporary vector */
 66:   IS        isg[2];  /* index sets of split "0" and "1" */
 67: } Stokes;

 69: PetscErrorCode StokesSetupMatBlock00(Stokes*);  /* setup the block Q */
 70: PetscErrorCode StokesSetupMatBlock01(Stokes*);  /* setup the block G */
 71: PetscErrorCode StokesSetupMatBlock10(Stokes*);  /* setup the block D (equal to the transpose of G) */
 72: PetscErrorCode StokesSetupMatBlock11(Stokes*);  /* setup the block C (equal to zero) */

 74: PetscErrorCode StokesGetPosition(Stokes*, PetscInt, PetscInt*, PetscInt*); /* row number j*nx+i corresponds to position (i,j) in grid */

 76: PetscErrorCode StokesStencilLaplacian(Stokes*, PetscInt, PetscInt, PetscInt*, PetscInt*, PetscScalar*);  /* stencil of the Laplacian operator */
 77: PetscErrorCode StokesStencilGradientX(Stokes*, PetscInt, PetscInt, PetscInt*, PetscInt*, PetscScalar*);  /* stencil of the Gradient operator (x-component) */
 78: PetscErrorCode StokesStencilGradientY(Stokes*, PetscInt, PetscInt, PetscInt*, PetscInt*, PetscScalar*);  /* stencil of the Gradient operator (y-component) */

 80: PetscErrorCode StokesRhs(Stokes*);                                               /* rhs vector */
 81: PetscErrorCode StokesRhsMomX(Stokes*, PetscInt, PetscInt, PetscScalar*);   /* right hand side of velocity (x-component) */
 82: PetscErrorCode StokesRhsMomY(Stokes*, PetscInt, PetscInt, PetscScalar*);   /* right hand side of velocity (y-component) */
 83: PetscErrorCode StokesRhsMass(Stokes*, PetscInt, PetscInt, PetscScalar*);   /* right hand side of pressure */

 85: PetscErrorCode StokesSetupApproxSchur(Stokes*);  /* approximation of the Schur complement */

 87: PetscErrorCode StokesExactSolution(Stokes*); /* exact solution vector */
 88: PetscErrorCode StokesWriteSolution(Stokes*); /* write solution to file */

 90: /* exact solution for the velocity (x-component, y-component is zero) */
 91: PetscScalar StokesExactVelocityX(const PetscScalar y)
 92: {
 93:   return 4.0*y*(1.0-y);
 94: }

 96: /* exact solution for the pressure */
 97: PetscScalar StokesExactPressure(const PetscScalar x)
 98: {
 99:   return 8.0*(2.0-x);
100: }

102: PetscErrorCode StokesSetupPC(Stokes *s, KSP ksp)
103: {
104:   KSP            *subksp;
105:   PC             pc;
106:   PetscInt       n = 1;

110:   KSPGetPC(ksp, &pc);
111:   PCFieldSplitSetIS(pc, "0", s->isg[0]);
112:   PCFieldSplitSetIS(pc, "1", s->isg[1]);
113:   if (s->userPC) {
114:     PCFieldSplitSetSchurPre(pc, PC_FIELDSPLIT_SCHUR_PRE_USER, s->myS);
115:   }
116:   if (s->userKSP) {
117:     PCSetUp(pc);
118:     PCFieldSplitGetSubKSP(pc, &n, &subksp);
119:     KSPSetOperators(subksp[1], s->myS, s->myS);
120:     PetscFree(subksp);
121:   }
122:   return(0);
123: }

125: PetscErrorCode StokesWriteSolution(Stokes *s)
126: {
127:   PetscMPIInt       size;
128:   PetscInt          n,i,j;
129:   const PetscScalar *array;
130:   PetscErrorCode    ierr;

133:   /* write data (*warning* only works sequential) */
134:   MPI_Comm_size(MPI_COMM_WORLD,&size);
135:   /*PetscPrintf(PETSC_COMM_WORLD," number of processors = %D\n",size);*/
136:   if (size == 1) {
137:     PetscViewer viewer;
138:     VecGetArrayRead(s->x, &array);
139:     PetscViewerASCIIOpen(PETSC_COMM_WORLD, "solution.dat", &viewer);
140:     PetscViewerASCIIPrintf(viewer, "# x, y, u, v, p\n");
141:     for (j = 0; j < s->ny; j++) {
142:       for (i = 0; i < s->nx; i++) {
143:         n    = j*s->nx+i;
144:         PetscViewerASCIIPrintf(viewer, "%.12g %.12g %.12g %.12g %.12g\n", (double)(i*s->hx+s->hx/2),(double)(j*s->hy+s->hy/2), (double)PetscRealPart(array[n]), (double)PetscRealPart(array[n+s->nx*s->ny]),(double)PetscRealPart(array[n+2*s->nx*s->ny]));
145:       }
146:     }
147:     VecRestoreArrayRead(s->x, &array);
148:     PetscViewerDestroy(&viewer);
149:   }
150:   return(0);
151: }

153: PetscErrorCode StokesSetupIndexSets(Stokes *s)
154: {

158:   /* the two index sets */
159:   MatNestGetISs(s->A, s->isg, NULL);
160:   /*  ISView(isg[0],PETSC_VIEWER_STDOUT_WORLD); */
161:   /*  ISView(isg[1],PETSC_VIEWER_STDOUT_WORLD); */
162:   return(0);
163: }

165: PetscErrorCode StokesSetupVectors(Stokes *s)
166: {

170:   /* solution vector x */
171:   VecCreate(PETSC_COMM_WORLD, &s->x);
172:   VecSetSizes(s->x, PETSC_DECIDE, 3*s->nx*s->ny);
173:   VecSetType(s->x, VECMPI);
174:   /*  VecSetRandom(s->x, NULL); */
175:   /*  VecView(s->x, (PetscViewer) PETSC_VIEWER_DEFAULT); */

177:   /* exact solution y */
178:   VecDuplicate(s->x, &s->y);
179:   StokesExactSolution(s);
180:   /*  VecView(s->y, (PetscViewer) PETSC_VIEWER_DEFAULT); */

182:   /* rhs vector b */
183:   VecDuplicate(s->x, &s->b);
184:   StokesRhs(s);
185:   /*VecView(s->b, (PetscViewer) PETSC_VIEWER_DEFAULT);*/
186:   return(0);
187: }

189: PetscErrorCode StokesGetPosition(Stokes *s, PetscInt row, PetscInt *i, PetscInt *j)
190: {
191:   PetscInt n;

194:   /* cell number n=j*nx+i has position (i,j) in grid */
195:   n  = row%(s->nx*s->ny);
196:   *i = n%s->nx;
197:   *j = (n-(*i))/s->nx;
198:   return(0);
199: }

201: PetscErrorCode StokesExactSolution(Stokes *s)
202: {
203:   PetscInt       row, start, end, i, j;
204:   PetscScalar    val;
205:   Vec            y0,y1;

209:   /* velocity part */
210:   VecGetSubVector(s->y, s->isg[0], &y0);
211:   VecGetOwnershipRange(y0, &start, &end);
212:   for (row = start; row < end; row++) {
213:     StokesGetPosition(s, row,&i,&j);
214:     if (row < s->nx*s->ny) {
215:       val = StokesExactVelocityX(j*s->hy+s->hy/2);
216:     } else {
217:       val = 0;
218:     }
219:     VecSetValue(y0, row, val, INSERT_VALUES);
220:   }
221:   VecRestoreSubVector(s->y, s->isg[0], &y0);

223:   /* pressure part */
224:   VecGetSubVector(s->y, s->isg[1], &y1);
225:   VecGetOwnershipRange(y1, &start, &end);
226:   for (row = start; row < end; row++) {
227:     StokesGetPosition(s, row, &i, &j);
228:     val  = StokesExactPressure(i*s->hx+s->hx/2);
229:     VecSetValue(y1, row, val, INSERT_VALUES);
230:   }
231:   VecRestoreSubVector(s->y, s->isg[1], &y1);
232:   return(0);
233: }

235: PetscErrorCode StokesRhs(Stokes *s)
236: {
237:   PetscInt       row, start, end, i, j;
238:   PetscScalar    val;
239:   Vec            b0,b1;

243:   /* velocity part */
244:   VecGetSubVector(s->b, s->isg[0], &b0);
245:   VecGetOwnershipRange(b0, &start, &end);
246:   for (row = start; row < end; row++) {
247:     StokesGetPosition(s, row, &i, &j);
248:     if (row < s->nx*s->ny) {
249:       StokesRhsMomX(s, i, j, &val);
250:     } else {
251:       StokesRhsMomY(s, i, j, &val);
252:     }
253:     VecSetValue(b0, row, val, INSERT_VALUES);
254:   }
255:   VecRestoreSubVector(s->b, s->isg[0], &b0);

257:   /* pressure part */
258:   VecGetSubVector(s->b, s->isg[1], &b1);
259:   VecGetOwnershipRange(b1, &start, &end);
260:   for (row = start; row < end; row++) {
261:     StokesGetPosition(s, row, &i, &j);
262:     StokesRhsMass(s, i, j, &val);
263:     VecSetValue(b1, row, val, INSERT_VALUES);
264:   }
265:   VecRestoreSubVector(s->b, s->isg[1], &b1);
266:   return(0);
267: }

269: PetscErrorCode StokesSetupMatBlock00(Stokes *s)
270: {
271:   PetscInt       row, start, end, sz, i, j;
272:   PetscInt       cols[5];
273:   PetscScalar    vals[5];

277:   /* A[0] is 2N-by-2N */
278:   MatCreate(PETSC_COMM_WORLD,&s->subA[0]);
279:   MatSetOptionsPrefix(s->subA[0],"a00_");
280:   MatSetSizes(s->subA[0],PETSC_DECIDE,PETSC_DECIDE,2*s->nx*s->ny,2*s->nx*s->ny);
281:   MatSetType(s->subA[0],MATMPIAIJ);
282:   MatMPIAIJSetPreallocation(s->subA[0],5,NULL,5,NULL);
283:   MatGetOwnershipRange(s->subA[0], &start, &end);

285:   for (row = start; row < end; row++) {
286:     StokesGetPosition(s, row, &i, &j);
287:     /* first part: rows 0 to (nx*ny-1) */
288:     StokesStencilLaplacian(s, i, j, &sz, cols, vals);
289:     /* second part: rows (nx*ny) to (2*nx*ny-1) */
290:     if (row >= s->nx*s->ny) {
291:       for (i = 0; i < sz; i++) cols[i] += s->nx*s->ny;
292:     }
293:     for (i = 0; i < sz; i++) vals[i] = -1.0*vals[i]; /* dynamic viscosity coef mu=-1 */
294:     MatSetValues(s->subA[0], 1, &row, sz, cols, vals, INSERT_VALUES);
295:   }
296:   MatAssemblyBegin(s->subA[0], MAT_FINAL_ASSEMBLY);
297:   MatAssemblyEnd(s->subA[0], MAT_FINAL_ASSEMBLY);
298:   return(0);
299: }

301: PetscErrorCode StokesSetupMatBlock01(Stokes *s)
302: {
303:   PetscInt       row, start, end, sz, i, j;
304:   PetscInt       cols[5];
305:   PetscScalar    vals[5];

309:   /* A[1] is 2N-by-N */
310:   MatCreate(PETSC_COMM_WORLD, &s->subA[1]);
311:   MatSetOptionsPrefix(s->subA[1],"a01_");
312:   MatSetSizes(s->subA[1],PETSC_DECIDE,PETSC_DECIDE,2*s->nx*s->ny,s->nx*s->ny);
313:   MatSetType(s->subA[1],MATMPIAIJ);
314:   MatMPIAIJSetPreallocation(s->subA[1],5,NULL,5,NULL);
315:   MatGetOwnershipRange(s->subA[1],&start,&end);

317:   MatSetOption(s->subA[1],MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);

319:   for (row = start; row < end; row++) {
320:     StokesGetPosition(s, row, &i, &j);
321:     /* first part: rows 0 to (nx*ny-1) */
322:     if (row < s->nx*s->ny) {
323:       StokesStencilGradientX(s, i, j, &sz, cols, vals);
324:     } else {    /* second part: rows (nx*ny) to (2*nx*ny-1) */
325:       StokesStencilGradientY(s, i, j, &sz, cols, vals);
326:     }
327:     MatSetValues(s->subA[1], 1, &row, sz, cols, vals, INSERT_VALUES);
328:   }
329:   MatAssemblyBegin(s->subA[1], MAT_FINAL_ASSEMBLY);
330:   MatAssemblyEnd(s->subA[1], MAT_FINAL_ASSEMBLY);
331:   return(0);
332: }

334: PetscErrorCode StokesSetupMatBlock10(Stokes *s)
335: {

339:   /* A[2] is minus transpose of A[1] */
340:   MatTranspose(s->subA[1], MAT_INITIAL_MATRIX, &s->subA[2]);
341:   MatScale(s->subA[2], -1.0);
342:   MatSetOptionsPrefix(s->subA[2], "a10_");
343:   return(0);
344: }

346: PetscErrorCode StokesSetupMatBlock11(Stokes *s)
347: {

351:   /* A[3] is N-by-N null matrix */
352:   MatCreate(PETSC_COMM_WORLD, &s->subA[3]);
353:   MatSetOptionsPrefix(s->subA[3], "a11_");
354:   MatSetSizes(s->subA[3], PETSC_DECIDE, PETSC_DECIDE, s->nx*s->ny, s->nx*s->ny);
355:   MatSetType(s->subA[3], MATMPIAIJ);
356:   MatMPIAIJSetPreallocation(s->subA[3], 0, NULL, 0, NULL);
357:   MatAssemblyBegin(s->subA[3], MAT_FINAL_ASSEMBLY);
358:   MatAssemblyEnd(s->subA[3], MAT_FINAL_ASSEMBLY);
359:   return(0);
360: }

362: PetscErrorCode StokesSetupApproxSchur(Stokes *s)
363: {
364:   Vec            diag;

368:   /* Schur complement approximation: myS = A11 - A10 diag(A00)^(-1) A01 */
369:   /* note: A11 is zero */
370:   /* note: in real life this matrix would be build directly, */
371:   /* i.e. without MatMatMult */

373:   /* inverse of diagonal of A00 */
374:   VecCreate(PETSC_COMM_WORLD,&diag);
375:   VecSetSizes(diag,PETSC_DECIDE,2*s->nx*s->ny);
376:   VecSetType(diag,VECMPI);
377:   MatGetDiagonal(s->subA[0],diag);
378:   VecReciprocal(diag);

380:   /* compute: - A10 diag(A00)^(-1) A01 */
381:   MatDiagonalScale(s->subA[1],diag,NULL); /* (*warning* overwrites subA[1]) */
382:   MatMatMult(s->subA[2],s->subA[1],MAT_INITIAL_MATRIX,PETSC_DEFAULT,&s->myS);
383:   MatScale(s->myS,-1.0);

385:   /* restore A10 */
386:   MatGetDiagonal(s->subA[0],diag);
387:   MatDiagonalScale(s->subA[1],diag,NULL);
388:   VecDestroy(&diag);
389:   return(0);
390: }

392: PetscErrorCode StokesSetupMatrix(Stokes *s)
393: {

397:   StokesSetupMatBlock00(s);
398:   StokesSetupMatBlock01(s);
399:   StokesSetupMatBlock10(s);
400:   StokesSetupMatBlock11(s);
401:   MatCreateNest(PETSC_COMM_WORLD, 2, NULL, 2, NULL, s->subA, &s->A);
402:   StokesSetupApproxSchur(s);
403:   return(0);
404: }

406: PetscErrorCode StokesStencilLaplacian(Stokes *s, PetscInt i, PetscInt j, PetscInt *sz, PetscInt *cols, PetscScalar *vals)
407: {
408:   PetscInt    p =j*s->nx+i, w=p-1, e=p+1, s2=p-s->nx, n=p+s->nx;
409:   PetscScalar ae=s->hy/s->hx, aeb=0;
410:   PetscScalar aw=s->hy/s->hx, awb=s->hy/(s->hx/2);
411:   PetscScalar as=s->hx/s->hy, asb=s->hx/(s->hy/2);
412:   PetscScalar an=s->hx/s->hy, anb=s->hx/(s->hy/2);

415:   if (i==0 && j==0) { /* south-west corner */
416:     *sz  =3;
417:     cols[0]=p; vals[0]=-(ae+awb+asb+an);
418:     cols[1]=e; vals[1]=ae;
419:     cols[2]=n; vals[2]=an;
420:   } else if (i==0 && j==s->ny-1) { /* north-west corner */
421:     *sz  =3;
422:     cols[0]=s2; vals[0]=as;
423:     cols[1]=p; vals[1]=-(ae+awb+as+anb);
424:     cols[2]=e; vals[2]=ae;
425:   } else if (i==s->nx-1 && j==0) { /* south-east corner */
426:     *sz  =3;
427:     cols[0]=w; vals[0]=aw;
428:     cols[1]=p; vals[1]=-(aeb+aw+asb+an);
429:     cols[2]=n; vals[2]=an;
430:   } else if (i==s->nx-1 && j==s->ny-1) { /* north-east corner */
431:     *sz  =3;
432:     cols[0]=s2; vals[0]=as;
433:     cols[1]=w; vals[1]=aw;
434:     cols[2]=p; vals[2]=-(aeb+aw+as+anb);
435:   } else if (i==0) { /* west boundary */
436:     *sz  =4;
437:     cols[0]=s2; vals[0]=as;
438:     cols[1]=p; vals[1]=-(ae+awb+as+an);
439:     cols[2]=e; vals[2]=ae;
440:     cols[3]=n; vals[3]=an;
441:   } else if (i==s->nx-1) { /* east boundary */
442:     *sz  =4;
443:     cols[0]=s2; vals[0]=as;
444:     cols[1]=w; vals[1]=aw;
445:     cols[2]=p; vals[2]=-(aeb+aw+as+an);
446:     cols[3]=n; vals[3]=an;
447:   } else if (j==0) { /* south boundary */
448:     *sz  =4;
449:     cols[0]=w; vals[0]=aw;
450:     cols[1]=p; vals[1]=-(ae+aw+asb+an);
451:     cols[2]=e; vals[2]=ae;
452:     cols[3]=n; vals[3]=an;
453:   } else if (j==s->ny-1) { /* north boundary */
454:     *sz  =4;
455:     cols[0]=s2; vals[0]=as;
456:     cols[1]=w; vals[1]=aw;
457:     cols[2]=p; vals[2]=-(ae+aw+as+anb);
458:     cols[3]=e; vals[3]=ae;
459:   } else { /* interior */
460:     *sz  =5;
461:     cols[0]=s2; vals[0]=as;
462:     cols[1]=w; vals[1]=aw;
463:     cols[2]=p; vals[2]=-(ae+aw+as+an);
464:     cols[3]=e; vals[3]=ae;
465:     cols[4]=n; vals[4]=an;
466:   }
467:   return(0);
468: }

470: PetscErrorCode StokesStencilGradientX(Stokes *s, PetscInt i, PetscInt j, PetscInt *sz, PetscInt *cols, PetscScalar *vals)
471: {
472:   PetscInt    p =j*s->nx+i, w=p-1, e=p+1;
473:   PetscScalar ae= s->hy/2, aeb=s->hy;
474:   PetscScalar aw=-s->hy/2, awb=0;

477:   if (i==0 && j==0) { /* south-west corner */
478:     *sz  =2;
479:     cols[0]=p; vals[0]=-(ae+awb);
480:     cols[1]=e; vals[1]=ae;
481:   } else if (i==0 && j==s->ny-1) { /* north-west corner */
482:     *sz  =2;
483:     cols[0]=p; vals[0]=-(ae+awb);
484:     cols[1]=e; vals[1]=ae;
485:   } else if (i==s->nx-1 && j==0) { /* south-east corner */
486:     *sz  =2;
487:     cols[0]=w; vals[0]=aw;
488:     cols[1]=p; vals[1]=-(aeb+aw);
489:   } else if (i==s->nx-1 && j==s->ny-1) { /* north-east corner */
490:     *sz  =2;
491:     cols[0]=w; vals[0]=aw;
492:     cols[1]=p; vals[1]=-(aeb+aw);
493:   } else if (i==0) { /* west boundary */
494:     *sz  =2;
495:     cols[0]=p; vals[0]=-(ae+awb);
496:     cols[1]=e; vals[1]=ae;
497:   } else if (i==s->nx-1) { /* east boundary */
498:     *sz  =2;
499:     cols[0]=w; vals[0]=aw;
500:     cols[1]=p; vals[1]=-(aeb+aw);
501:   } else if (j==0) { /* south boundary */
502:     *sz  =3;
503:     cols[0]=w; vals[0]=aw;
504:     cols[1]=p; vals[1]=-(ae+aw);
505:     cols[2]=e; vals[2]=ae;
506:   } else if (j==s->ny-1) { /* north boundary */
507:     *sz  =3;
508:     cols[0]=w; vals[0]=aw;
509:     cols[1]=p; vals[1]=-(ae+aw);
510:     cols[2]=e; vals[2]=ae;
511:   } else { /* interior */
512:     *sz  =3;
513:     cols[0]=w; vals[0]=aw;
514:     cols[1]=p; vals[1]=-(ae+aw);
515:     cols[2]=e; vals[2]=ae;
516:   }
517:   return(0);
518: }

520: PetscErrorCode StokesStencilGradientY(Stokes *s, PetscInt i, PetscInt j, PetscInt *sz, PetscInt *cols, PetscScalar *vals)
521: {
522:   PetscInt    p =j*s->nx+i, s2=p-s->nx, n=p+s->nx;
523:   PetscScalar as=-s->hx/2, asb=0;
524:   PetscScalar an= s->hx/2, anb=0;

527:   if (i==0 && j==0) { /* south-west corner */
528:     *sz  =2;
529:     cols[0]=p; vals[0]=-(asb+an);
530:     cols[1]=n; vals[1]=an;
531:   } else if (i==0 && j==s->ny-1) { /* north-west corner */
532:     *sz  =2;
533:     cols[0]=s2; vals[0]=as;
534:     cols[1]=p; vals[1]=-(as+anb);
535:   } else if (i==s->nx-1 && j==0) { /* south-east corner */
536:     *sz  =2;
537:     cols[0]=p; vals[0]=-(asb+an);
538:     cols[1]=n; vals[1]=an;
539:   } else if (i==s->nx-1 && j==s->ny-1) { /* north-east corner */
540:     *sz  =2;
541:     cols[0]=s2; vals[0]=as;
542:     cols[1]=p; vals[1]=-(as+anb);
543:   } else if (i==0) { /* west boundary */
544:     *sz  =3;
545:     cols[0]=s2; vals[0]=as;
546:     cols[1]=p; vals[1]=-(as+an);
547:     cols[2]=n; vals[2]=an;
548:   } else if (i==s->nx-1) { /* east boundary */
549:     *sz  =3;
550:     cols[0]=s2; vals[0]=as;
551:     cols[1]=p; vals[1]=-(as+an);
552:     cols[2]=n; vals[2]=an;
553:   } else if (j==0) { /* south boundary */
554:     *sz  =2;
555:     cols[0]=p; vals[0]=-(asb+an);
556:     cols[1]=n; vals[1]=an;
557:   } else if (j==s->ny-1) { /* north boundary */
558:     *sz  =2;
559:     cols[0]=s2; vals[0]=as;
560:     cols[1]=p; vals[1]=-(as+anb);
561:   } else { /* interior */
562:     *sz  =3;
563:     cols[0]=s2; vals[0]=as;
564:     cols[1]=p; vals[1]=-(as+an);
565:     cols[2]=n; vals[2]=an;
566:   }
567:   return(0);
568: }

570: PetscErrorCode StokesRhsMomX(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val)
571: {
572:   PetscScalar y   = j*s->hy+s->hy/2;
573:   PetscScalar awb = s->hy/(s->hx/2);

576:   if (i == 0) { /* west boundary */
577:     *val = awb*StokesExactVelocityX(y);
578:   } else {
579:     *val = 0.0;
580:   }
581:   return(0);
582: }

584: PetscErrorCode StokesRhsMomY(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val)
585: {
587:   *val = 0.0;
588:   return(0);
589: }

591: PetscErrorCode StokesRhsMass(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val)
592: {
593:   PetscScalar y   = j*s->hy+s->hy/2;
594:   PetscScalar aeb = s->hy;

597:   if (i == 0) { /* west boundary */
598:     *val = aeb*StokesExactVelocityX(y);
599:   } else {
600:     *val = 0.0;
601:   }
602:   return(0);
603: }

605: PetscErrorCode StokesCalcResidual(Stokes *s)
606: {
607:   PetscReal      val;
608:   Vec            b0, b1;

612:   /* residual Ax-b (*warning* overwrites b) */
613:   VecScale(s->b, -1.0);
614:   MatMultAdd(s->A, s->x, s->b, s->b);
615:   /*  VecView(s->b, (PetscViewer)PETSC_VIEWER_DEFAULT); */

617:   /* residual velocity */
618:   VecGetSubVector(s->b, s->isg[0], &b0);
619:   VecNorm(b0, NORM_2, &val);
620:   PetscPrintf(PETSC_COMM_WORLD," residual u = %g\n",(double)val);
621:   VecRestoreSubVector(s->b, s->isg[0], &b0);

623:   /* residual pressure */
624:   VecGetSubVector(s->b, s->isg[1], &b1);
625:   VecNorm(b1, NORM_2, &val);
626:   PetscPrintf(PETSC_COMM_WORLD," residual p = %g\n",(double)val);
627:   VecRestoreSubVector(s->b, s->isg[1], &b1);

629:   /* total residual */
630:   VecNorm(s->b, NORM_2, &val);
631:   PetscPrintf(PETSC_COMM_WORLD," residual [u,p] = %g\n", (double)val);
632:   return(0);
633: }

635: PetscErrorCode StokesCalcError(Stokes *s)
636: {
637:   PetscScalar    scale = PetscSqrtReal((double)s->nx*s->ny);
638:   PetscReal      val;
639:   Vec            y0, y1;

643:   /* error y-x */
644:   VecAXPY(s->y, -1.0, s->x);
645:   /* VecView(s->y, (PetscViewer)PETSC_VIEWER_DEFAULT); */

647:   /* error in velocity */
648:   VecGetSubVector(s->y, s->isg[0], &y0);
649:   VecNorm(y0, NORM_2, &val);
650:   PetscPrintf(PETSC_COMM_WORLD," discretization error u = %g\n",(double)(PetscRealPart(val/scale)));
651:   VecRestoreSubVector(s->y, s->isg[0], &y0);

653:   /* error in pressure */
654:   VecGetSubVector(s->y, s->isg[1], &y1);
655:   VecNorm(y1, NORM_2, &val);
656:   PetscPrintf(PETSC_COMM_WORLD," discretization error p = %g\n",(double)(PetscRealPart(val/scale)));
657:   VecRestoreSubVector(s->y, s->isg[1], &y1);

659:   /* total error */
660:   VecNorm(s->y, NORM_2, &val);
661:   PetscPrintf(PETSC_COMM_WORLD," discretization error [u,p] = %g\n", (double)PetscRealPart((val/scale)));
662:   return(0);
663: }

665: int main(int argc, char **argv)
666: {
667:   Stokes         s;
668:   KSP            ksp;

671:   PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
672:   s.nx     = 4;
673:   s.ny     = 6;
674:   PetscOptionsGetInt(NULL,NULL, "-nx", &s.nx, NULL);
675:   PetscOptionsGetInt(NULL,NULL, "-ny", &s.ny, NULL);
676:   s.hx     = 2.0/s.nx;
677:   s.hy     = 1.0/s.ny;
678:   s.userPC = s.userKSP = PETSC_FALSE;
679:   PetscOptionsHasName(NULL,NULL, "-user_pc", &s.userPC);
680:   PetscOptionsHasName(NULL,NULL, "-user_ksp", &s.userKSP);

682:   StokesSetupMatrix(&s);
683:   StokesSetupIndexSets(&s);
684:   StokesSetupVectors(&s);

686:   KSPCreate(PETSC_COMM_WORLD, &ksp);
687:   KSPSetOperators(ksp, s.A, s.A);
688:   KSPSetFromOptions(ksp);
689:   StokesSetupPC(&s, ksp);
690:   KSPSolve(ksp, s.b, s.x);

692:   /* don't trust, verify! */
693:   StokesCalcResidual(&s);
694:   StokesCalcError(&s);
695:   StokesWriteSolution(&s);

697:   KSPDestroy(&ksp);
698:   MatDestroy(&s.subA[0]);
699:   MatDestroy(&s.subA[1]);
700:   MatDestroy(&s.subA[2]);
701:   MatDestroy(&s.subA[3]);
702:   MatDestroy(&s.A);
703:   VecDestroy(&s.x);
704:   VecDestroy(&s.b);
705:   VecDestroy(&s.y);
706:   MatDestroy(&s.myS);
707:   PetscFinalize();
708:   return ierr;
709: }


712: /*TEST

714:    test:
715:       nsize: 2
716:       args: -nx 16 -ny 24 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -fieldsplit_1_pc_type none

718:    test:
719:       suffix: 2
720:       nsize: 2
721:       args: -nx 16 -ny 24 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -user_pc

723:    test:
724:       suffix: 3
725:       nsize: 2
726:       args: -nx 16 -ny 24 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -user_pc

728:    test:
729:       suffix: 4
730:       nsize: 2
731:       args: -nx 16 -ny 24 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -fieldsplit_0_ksp_type gmres -fieldsplit_0_pc_type bjacobi -fieldsplit_1_pc_type jacobi -fieldsplit_1_inner_ksp_type preonly -fieldsplit_1_inner_pc_type jacobi -fieldsplit_1_upper_ksp_type preonly -fieldsplit_1_upper_pc_type jacobi

733:    test:
734:       suffix: 4_pcksp
735:       nsize: 2
736:       args: -nx 16 -ny 24 -ksp_type fgmres -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type lower -fieldsplit_0_ksp_type gmres -fieldsplit_0_pc_type bjacobi -fieldsplit_1_pc_type jacobi -fieldsplit_1_inner_ksp_type preonly -fieldsplit_1_upper_ksp_type preonly -fieldsplit_1_upper_pc_type jacobi

738: TEST*/