Actual source code: ex70.c
petsc-3.7.3 2016-08-01
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: /*---------------------------------------------------------------------------- */
55: #include <petscksp.h>
57: typedef struct {
58: PetscBool userPC, userKSP; /* user defined preconditioner and matrix for the Schur complement */
59: PetscInt nx, ny; /* nb of cells in x- and y-direction */
60: PetscReal hx, hy; /* mesh size in x- and y-direction */
61: Mat A; /* block matrix */
62: Mat subA[4]; /* the four blocks */
63: Mat myS; /* the approximation of the Schur complement */
64: Vec x, b, y; /* solution, rhs and temporary vector */
65: IS isg[2]; /* index sets of split "0" and "1" */
66: } Stokes;
68: PetscErrorCode StokesSetupMatBlock00(Stokes*); /* setup the block Q */
69: PetscErrorCode StokesSetupMatBlock01(Stokes*); /* setup the block G */
70: PetscErrorCode StokesSetupMatBlock10(Stokes*); /* setup the block D (equal to the transpose of G) */
71: PetscErrorCode StokesSetupMatBlock11(Stokes*); /* setup the block C (equal to zero) */
73: PetscErrorCode StokesGetPosition(Stokes*, PetscInt, PetscInt*, PetscInt*); /* row number j*nx+i corresponds to position (i,j) in grid */
75: PetscErrorCode StokesStencilLaplacian(Stokes*, PetscInt, PetscInt, PetscInt*, PetscInt*, PetscScalar*); /* stencil of the Laplacian operator */
76: PetscErrorCode StokesStencilGradientX(Stokes*, PetscInt, PetscInt, PetscInt*, PetscInt*, PetscScalar*); /* stencil of the Gradient operator (x-component) */
77: PetscErrorCode StokesStencilGradientY(Stokes*, PetscInt, PetscInt, PetscInt*, PetscInt*, PetscScalar*); /* stencil of the Gradient operator (y-component) */
79: PetscErrorCode StokesRhs(Stokes*); /* rhs vector */
80: PetscErrorCode StokesRhsMomX(Stokes*, PetscInt, PetscInt, PetscScalar*); /* right hand side of velocity (x-component) */
81: PetscErrorCode StokesRhsMomY(Stokes*, PetscInt, PetscInt, PetscScalar*); /* right hand side of velocity (y-component) */
82: PetscErrorCode StokesRhsMass(Stokes*, PetscInt, PetscInt, PetscScalar*); /* right hand side of pressure */
84: PetscErrorCode StokesSetupApproxSchur(Stokes*); /* approximation of the Schur complement */
86: PetscErrorCode StokesExactSolution(Stokes*); /* exact solution vector */
87: PetscErrorCode StokesWriteSolution(Stokes*); /* write solution to file */
89: /* exact solution for the velocity (x-component, y-component is zero) */
90: PetscScalar StokesExactVelocityX(const PetscScalar y)
91: {
92: return 4.0*y*(1.0-y);
93: }
95: /* exact solution for the pressure */
96: PetscScalar StokesExactPressure(const PetscScalar x)
97: {
98: return 8.0*(2.0-x);
99: }
101: PetscErrorCode StokesSetupPC(Stokes *s, KSP ksp)
102: {
103: KSP *subksp;
104: PC pc;
105: PetscInt n = 1;
109: KSPGetPC(ksp, &pc);
110: PCFieldSplitSetIS(pc, "0", s->isg[0]);
111: PCFieldSplitSetIS(pc, "1", s->isg[1]);
112: if (s->userPC) {
113: PCFieldSplitSetSchurPre(pc, PC_FIELDSPLIT_SCHUR_PRE_USER, s->myS);
114: }
115: if (s->userKSP) {
116: PCSetUp(pc);
117: PCFieldSplitGetSubKSP(pc, &n, &subksp);
118: KSPSetOperators(subksp[1], s->myS, s->myS);
119: PetscFree(subksp);
120: }
121: return(0);
122: }
124: PetscErrorCode StokesWriteSolution(Stokes *s)
125: {
126: PetscMPIInt size;
127: PetscInt n,i,j;
128: const PetscScalar *array;
129: PetscErrorCode ierr;
132: /* write data (*warning* only works sequential) */
133: MPI_Comm_size(MPI_COMM_WORLD,&size);
134: /*PetscPrintf(PETSC_COMM_WORLD," number of processors = %D\n",size);*/
135: if (size == 1) {
136: PetscViewer viewer;
137: VecGetArrayRead(s->x, &array);
138: PetscViewerASCIIOpen(PETSC_COMM_WORLD, "solution.dat", &viewer);
139: PetscViewerASCIIPrintf(viewer, "# x, y, u, v, p\n");
140: for (j = 0; j < s->ny; j++) {
141: for (i = 0; i < s->nx; i++) {
142: n = j*s->nx+i;
143: 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]));
144: }
145: }
146: VecRestoreArrayRead(s->x, &array);
147: PetscViewerDestroy(&viewer);
148: }
149: return(0);
150: }
152: PetscErrorCode StokesSetupIndexSets(Stokes *s)
153: {
157: /* the two index sets */
158: MatNestGetISs(s->A, s->isg, NULL);
159: /* ISView(isg[0],PETSC_VIEWER_STDOUT_WORLD); */
160: /* ISView(isg[1],PETSC_VIEWER_STDOUT_WORLD); */
161: return(0);
162: }
164: PetscErrorCode StokesSetupVectors(Stokes *s)
165: {
169: /* solution vector x */
170: VecCreate(PETSC_COMM_WORLD, &s->x);
171: VecSetSizes(s->x, PETSC_DECIDE, 3*s->nx*s->ny);
172: VecSetType(s->x, VECMPI);
173: /* VecSetRandom(s->x, NULL); */
174: /* VecView(s->x, (PetscViewer) PETSC_VIEWER_DEFAULT); */
176: /* exact solution y */
177: VecDuplicate(s->x, &s->y);
178: StokesExactSolution(s);
179: /* VecView(s->y, (PetscViewer) PETSC_VIEWER_DEFAULT); */
181: /* rhs vector b */
182: VecDuplicate(s->x, &s->b);
183: StokesRhs(s);
184: /*VecView(s->b, (PetscViewer) PETSC_VIEWER_DEFAULT);*/
185: return(0);
186: }
188: PetscErrorCode StokesGetPosition(Stokes *s, PetscInt row, PetscInt *i, PetscInt *j)
189: {
190: PetscInt n;
193: /* cell number n=j*nx+i has position (i,j) in grid */
194: n = row%(s->nx*s->ny);
195: *i = n%s->nx;
196: *j = (n-(*i))/s->nx;
197: return(0);
198: }
200: PetscErrorCode StokesExactSolution(Stokes *s)
201: {
202: PetscInt row, start, end, i, j;
203: PetscScalar val;
204: Vec y0,y1;
208: /* velocity part */
209: VecGetSubVector(s->y, s->isg[0], &y0);
210: VecGetOwnershipRange(y0, &start, &end);
211: for (row = start; row < end; row++) {
212: StokesGetPosition(s, row,&i,&j);
213: if (row < s->nx*s->ny) {
214: val = StokesExactVelocityX(j*s->hy+s->hy/2);
215: } else {
216: val = 0;
217: }
218: VecSetValue(y0, row, val, INSERT_VALUES);
219: }
220: VecRestoreSubVector(s->y, s->isg[0], &y0);
222: /* pressure part */
223: VecGetSubVector(s->y, s->isg[1], &y1);
224: VecGetOwnershipRange(y1, &start, &end);
225: for (row = start; row < end; row++) {
226: StokesGetPosition(s, row, &i, &j);
227: val = StokesExactPressure(i*s->hx+s->hx/2);
228: VecSetValue(y1, row, val, INSERT_VALUES);
229: }
230: VecRestoreSubVector(s->y, s->isg[1], &y1);
231: return(0);
232: }
234: PetscErrorCode StokesRhs(Stokes *s)
235: {
236: PetscInt row, start, end, i, j;
237: PetscScalar val;
238: Vec b0,b1;
242: /* velocity part */
243: VecGetSubVector(s->b, s->isg[0], &b0);
244: VecGetOwnershipRange(b0, &start, &end);
245: for (row = start; row < end; row++) {
246: StokesGetPosition(s, row, &i, &j);
247: if (row < s->nx*s->ny) {
248: StokesRhsMomX(s, i, j, &val);
249: } else if (row < 2*s->nx*s->ny) {
250: StokesRhsMomY(s, i, j, &val);
251: }
252: VecSetValue(b0, row, val, INSERT_VALUES);
253: }
254: VecRestoreSubVector(s->b, s->isg[0], &b0);
256: /* pressure part */
257: VecGetSubVector(s->b, s->isg[1], &b1);
258: VecGetOwnershipRange(b1, &start, &end);
259: for (row = start; row < end; row++) {
260: StokesGetPosition(s, row, &i, &j);
261: StokesRhsMass(s, i, j, &val);
262: VecSetValue(b1, row, val, INSERT_VALUES);
263: }
264: VecRestoreSubVector(s->b, s->isg[1], &b1);
265: return(0);
266: }
268: PetscErrorCode StokesSetupMatBlock00(Stokes *s)
269: {
270: PetscInt row, start, end, size, i, j;
271: PetscInt cols[5];
272: PetscScalar vals[5];
276: /* A[0] is 2N-by-2N */
277: MatCreate(PETSC_COMM_WORLD,&s->subA[0]);
278: MatSetOptionsPrefix(s->subA[0],"a00_");
279: MatSetSizes(s->subA[0],PETSC_DECIDE,PETSC_DECIDE,2*s->nx*s->ny,2*s->nx*s->ny);
280: MatSetType(s->subA[0],MATMPIAIJ);
281: MatMPIAIJSetPreallocation(s->subA[0],5,NULL,5,NULL);
282: MatGetOwnershipRange(s->subA[0], &start, &end);
284: for (row = start; row < end; row++) {
285: StokesGetPosition(s, row, &i, &j);
286: /* first part: rows 0 to (nx*ny-1) */
287: StokesStencilLaplacian(s, i, j, &size, cols, vals);
288: /* second part: rows (nx*ny) to (2*nx*ny-1) */
289: if (row >= s->nx*s->ny) {
290: for (i = 0; i < 5; i++) cols[i] = cols[i] + s->nx*s->ny;
291: }
292: for (i = 0; i < 5; i++) vals[i] = -1.0*vals[i]; /* dynamic viscosity coef mu=-1 */
293: MatSetValues(s->subA[0], 1, &row, size, cols, vals, INSERT_VALUES);
294: }
295: MatAssemblyBegin(s->subA[0], MAT_FINAL_ASSEMBLY);
296: MatAssemblyEnd(s->subA[0], MAT_FINAL_ASSEMBLY);
297: return(0);
298: }
300: PetscErrorCode StokesSetupMatBlock01(Stokes *s)
301: {
302: PetscInt row, start, end, size, i, j;
303: PetscInt cols[5];
304: PetscScalar vals[5];
308: /* A[1] is 2N-by-N */
309: MatCreate(PETSC_COMM_WORLD, &s->subA[1]);
310: MatSetOptionsPrefix(s->subA[1],"a01_");
311: MatSetSizes(s->subA[1],PETSC_DECIDE,PETSC_DECIDE,2*s->nx*s->ny,s->nx*s->ny);
312: MatSetType(s->subA[1],MATMPIAIJ);
313: MatMPIAIJSetPreallocation(s->subA[1],5,NULL,5,NULL);
314: MatGetOwnershipRange(s->subA[1],&start,&end);
316: MatSetOption(s->subA[1],MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);
318: for (row = start; row < end; row++) {
319: StokesGetPosition(s, row, &i, &j);
320: /* first part: rows 0 to (nx*ny-1) */
321: if (row < s->nx*s->ny) {
322: StokesStencilGradientX(s, i, j, &size, cols, vals);
323: }
324: /* second part: rows (nx*ny) to (2*nx*ny-1) */
325: else {
326: StokesStencilGradientY(s, i, j, &size, cols, vals);
327: }
328: MatSetValues(s->subA[1], 1, &row, size, cols, vals, INSERT_VALUES);
329: }
330: MatAssemblyBegin(s->subA[1], MAT_FINAL_ASSEMBLY);
331: MatAssemblyEnd(s->subA[1], MAT_FINAL_ASSEMBLY);
332: return(0);
333: }
335: PetscErrorCode StokesSetupMatBlock10(Stokes *s)
336: {
340: /* A[2] is minus transpose of A[1] */
341: MatTranspose(s->subA[1], MAT_INITIAL_MATRIX, &s->subA[2]);
342: MatScale(s->subA[2], -1.0);
343: MatSetOptionsPrefix(s->subA[2], "a10_");
344: return(0);
345: }
347: PetscErrorCode StokesSetupMatBlock11(Stokes *s)
348: {
352: /* A[3] is N-by-N null matrix */
353: MatCreate(PETSC_COMM_WORLD, &s->subA[3]);
354: MatSetOptionsPrefix(s->subA[3], "a11_");
355: MatSetSizes(s->subA[3], PETSC_DECIDE, PETSC_DECIDE, s->nx*s->ny, s->nx*s->ny);
356: MatSetType(s->subA[3], MATMPIAIJ);
357: MatMPIAIJSetPreallocation(s->subA[3], 0, NULL, 0, NULL);
358: MatAssemblyBegin(s->subA[3], MAT_FINAL_ASSEMBLY);
359: MatAssemblyEnd(s->subA[3], MAT_FINAL_ASSEMBLY);
360: return(0);
361: }
363: PetscErrorCode StokesSetupApproxSchur(Stokes *s)
364: {
365: Vec diag;
369: /* Schur complement approximation: myS = A11 - A10 diag(A00)^(-1) A01 */
370: /* note: A11 is zero */
371: /* note: in real life this matrix would be build directly, */
372: /* i.e. without MatMatMult */
374: /* inverse of diagonal of A00 */
375: VecCreate(PETSC_COMM_WORLD,&diag);
376: VecSetSizes(diag,PETSC_DECIDE,2*s->nx*s->ny);
377: VecSetType(diag,VECMPI);
378: MatGetDiagonal(s->subA[0],diag);
379: VecReciprocal(diag);
381: /* compute: - A10 diag(A00)^(-1) A01 */
382: MatDiagonalScale(s->subA[1],diag,NULL); /* (*warning* overwrites subA[1]) */
383: MatMatMult(s->subA[2],s->subA[1],MAT_INITIAL_MATRIX,PETSC_DEFAULT,&s->myS);
384: MatScale(s->myS,-1.0);
386: /* restore A10 */
387: MatGetDiagonal(s->subA[0],diag);
388: MatDiagonalScale(s->subA[1],diag,NULL);
389: VecDestroy(&diag);
390: return(0);
391: }
393: PetscErrorCode StokesSetupMatrix(Stokes *s)
394: {
398: StokesSetupMatBlock00(s);
399: StokesSetupMatBlock01(s);
400: StokesSetupMatBlock10(s);
401: StokesSetupMatBlock11(s);
402: MatCreateNest(PETSC_COMM_WORLD, 2, NULL, 2, NULL, s->subA, &s->A);
403: StokesSetupApproxSchur(s);
404: return(0);
405: }
407: PetscErrorCode StokesStencilLaplacian(Stokes *s, PetscInt i, PetscInt j, PetscInt *size, PetscInt *cols, PetscScalar *vals)
408: {
409: PetscInt p =j*s->nx+i, w=p-1, e=p+1, s2=p-s->nx, n=p+s->nx;
410: PetscScalar ae=s->hy/s->hx, aeb=0;
411: PetscScalar aw=s->hy/s->hx, awb=s->hy/(s->hx/2);
412: PetscScalar as=s->hx/s->hy, asb=s->hx/(s->hy/2);
413: PetscScalar an=s->hx/s->hy, anb=s->hx/(s->hy/2);
416: if (i==0 && j==0) { /* south-west corner */
417: *size =3;
418: cols[0]=p; vals[0]=-(ae+awb+asb+an);
419: cols[1]=e; vals[1]=ae;
420: cols[2]=n; vals[2]=an;
421: } else if (i==0 && j==s->ny-1) { /* north-west corner */
422: *size =3;
423: cols[0]=s2; vals[0]=as;
424: cols[1]=p; vals[1]=-(ae+awb+as+anb);
425: cols[2]=e; vals[2]=ae;
426: } else if (i==s->nx-1 && j==0) { /* south-east corner */
427: *size =3;
428: cols[0]=w; vals[0]=aw;
429: cols[1]=p; vals[1]=-(aeb+aw+asb+an);
430: cols[2]=n; vals[2]=an;
431: } else if (i==s->nx-1 && j==s->ny-1) { /* north-east corner */
432: *size =3;
433: cols[0]=s2; vals[0]=as;
434: cols[1]=w; vals[1]=aw;
435: cols[2]=p; vals[2]=-(aeb+aw+as+anb);
436: } else if (i==0) { /* west boundary */
437: *size =4;
438: cols[0]=s2; vals[0]=as;
439: cols[1]=p; vals[1]=-(ae+awb+as+an);
440: cols[2]=e; vals[2]=ae;
441: cols[3]=n; vals[3]=an;
442: } else if (i==s->nx-1) { /* east boundary */
443: *size =4;
444: cols[0]=s2; vals[0]=as;
445: cols[1]=w; vals[1]=aw;
446: cols[2]=p; vals[2]=-(aeb+aw+as+an);
447: cols[3]=n; vals[3]=an;
448: } else if (j==0) { /* south boundary */
449: *size =4;
450: cols[0]=w; vals[0]=aw;
451: cols[1]=p; vals[1]=-(ae+aw+asb+an);
452: cols[2]=e; vals[2]=ae;
453: cols[3]=n; vals[3]=an;
454: } else if (j==s->ny-1) { /* north boundary */
455: *size =4;
456: cols[0]=s2; vals[0]=as;
457: cols[1]=w; vals[1]=aw;
458: cols[2]=p; vals[2]=-(ae+aw+as+anb);
459: cols[3]=e; vals[3]=ae;
460: } else { /* interior */
461: *size =5;
462: cols[0]=s2; vals[0]=as;
463: cols[1]=w; vals[1]=aw;
464: cols[2]=p; vals[2]=-(ae+aw+as+an);
465: cols[3]=e; vals[3]=ae;
466: cols[4]=n; vals[4]=an;
467: }
468: return(0);
469: }
471: PetscErrorCode StokesStencilGradientX(Stokes *s, PetscInt i, PetscInt j, PetscInt *size, PetscInt *cols, PetscScalar *vals)
472: {
473: PetscInt p =j*s->nx+i, w=p-1, e=p+1;
474: PetscScalar ae= s->hy/2, aeb=s->hy;
475: PetscScalar aw=-s->hy/2, awb=0;
478: if (i==0 && j==0) { /* south-west corner */
479: *size =2;
480: cols[0]=p; vals[0]=-(ae+awb);
481: cols[1]=e; vals[1]=ae;
482: } else if (i==0 && j==s->ny-1) { /* north-west corner */
483: *size =2;
484: cols[0]=p; vals[0]=-(ae+awb);
485: cols[1]=e; vals[1]=ae;
486: } else if (i==s->nx-1 && j==0) { /* south-east corner */
487: *size =2;
488: cols[0]=w; vals[0]=aw;
489: cols[1]=p; vals[1]=-(aeb+aw);
490: } else if (i==s->nx-1 && j==s->ny-1) { /* north-east corner */
491: *size =2;
492: cols[0]=w; vals[0]=aw;
493: cols[1]=p; vals[1]=-(aeb+aw);
494: } else if (i==0) { /* west boundary */
495: *size =2;
496: cols[0]=p; vals[0]=-(ae+awb);
497: cols[1]=e; vals[1]=ae;
498: } else if (i==s->nx-1) { /* east boundary */
499: *size =2;
500: cols[0]=w; vals[0]=aw;
501: cols[1]=p; vals[1]=-(aeb+aw);
502: } else if (j==0) { /* south boundary */
503: *size =3;
504: cols[0]=w; vals[0]=aw;
505: cols[1]=p; vals[1]=-(ae+aw);
506: cols[2]=e; vals[2]=ae;
507: } else if (j==s->ny-1) { /* north boundary */
508: *size =3;
509: cols[0]=w; vals[0]=aw;
510: cols[1]=p; vals[1]=-(ae+aw);
511: cols[2]=e; vals[2]=ae;
512: } else { /* interior */
513: *size =3;
514: cols[0]=w; vals[0]=aw;
515: cols[1]=p; vals[1]=-(ae+aw);
516: cols[2]=e; vals[2]=ae;
517: }
518: return(0);
519: }
521: PetscErrorCode StokesStencilGradientY(Stokes *s, PetscInt i, PetscInt j, PetscInt *size, PetscInt *cols, PetscScalar *vals)
522: {
523: PetscInt p =j*s->nx+i, s2=p-s->nx, n=p+s->nx;
524: PetscScalar as=-s->hx/2, asb=0;
525: PetscScalar an= s->hx/2, anb=0;
528: if (i==0 && j==0) { /* south-west corner */
529: *size =2;
530: cols[0]=p; vals[0]=-(asb+an);
531: cols[1]=n; vals[1]=an;
532: } else if (i==0 && j==s->ny-1) { /* north-west corner */
533: *size =2;
534: cols[0]=s2; vals[0]=as;
535: cols[1]=p; vals[1]=-(as+anb);
536: } else if (i==s->nx-1 && j==0) { /* south-east corner */
537: *size =2;
538: cols[0]=p; vals[0]=-(asb+an);
539: cols[1]=n; vals[1]=an;
540: } else if (i==s->nx-1 && j==s->ny-1) { /* north-east corner */
541: *size =2;
542: cols[0]=s2; vals[0]=as;
543: cols[1]=p; vals[1]=-(as+anb);
544: } else if (i==0) { /* west boundary */
545: *size =3;
546: cols[0]=s2; vals[0]=as;
547: cols[1]=p; vals[1]=-(as+an);
548: cols[2]=n; vals[2]=an;
549: } else if (i==s->nx-1) { /* east boundary */
550: *size =3;
551: cols[0]=s2; vals[0]=as;
552: cols[1]=p; vals[1]=-(as+an);
553: cols[2]=n; vals[2]=an;
554: } else if (j==0) { /* south boundary */
555: *size =2;
556: cols[0]=p; vals[0]=-(asb+an);
557: cols[1]=n; vals[1]=an;
558: } else if (j==s->ny-1) { /* north boundary */
559: *size =2;
560: cols[0]=s2; vals[0]=as;
561: cols[1]=p; vals[1]=-(as+anb);
562: } else { /* interior */
563: *size =3;
564: cols[0]=s2; vals[0]=as;
565: cols[1]=p; vals[1]=-(as+an);
566: cols[2]=n; vals[2]=an;
567: }
568: return(0);
569: }
571: PetscErrorCode StokesRhsMomX(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val)
572: {
573: PetscScalar y = j*s->hy+s->hy/2;
574: PetscScalar awb = s->hy/(s->hx/2);
577: if (i == 0) { /* west boundary */
578: *val = awb*StokesExactVelocityX(y);
579: } else {
580: *val = 0.0;
581: }
582: return(0);
583: }
585: PetscErrorCode StokesRhsMomY(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val)
586: {
588: *val = 0.0;
589: return(0);
590: }
592: PetscErrorCode StokesRhsMass(Stokes *s, PetscInt i, PetscInt j, PetscScalar *val)
593: {
594: PetscScalar y = j*s->hy+s->hy/2;
595: PetscScalar aeb = s->hy;
598: if (i == 0) { /* west boundary */
599: *val = aeb*StokesExactVelocityX(y);
600: } else {
601: *val = 0.0;
602: }
603: return(0);
604: }
606: PetscErrorCode StokesCalcResidual(Stokes *s)
607: {
608: PetscReal val;
609: Vec b0, b1;
613: /* residual Ax-b (*warning* overwrites b) */
614: VecScale(s->b, -1.0);
615: MatMultAdd(s->A, s->x, s->b, s->b);
616: /* VecView(s->b, (PetscViewer)PETSC_VIEWER_DEFAULT); */
618: /* residual velocity */
619: VecGetSubVector(s->b, s->isg[0], &b0);
620: VecNorm(b0, NORM_2, &val);
621: PetscPrintf(PETSC_COMM_WORLD," residual u = %g\n",(double)val);
622: VecRestoreSubVector(s->b, s->isg[0], &b0);
624: /* residual pressure */
625: VecGetSubVector(s->b, s->isg[1], &b1);
626: VecNorm(b1, NORM_2, &val);
627: PetscPrintf(PETSC_COMM_WORLD," residual p = %g\n",(double)val);
628: VecRestoreSubVector(s->b, s->isg[1], &b1);
630: /* total residual */
631: VecNorm(s->b, NORM_2, &val);
632: PetscPrintf(PETSC_COMM_WORLD," residual [u,p] = %g\n", (double)val);
633: return(0);
634: }
636: PetscErrorCode StokesCalcError(Stokes *s)
637: {
638: PetscScalar scale = PetscSqrtReal((double)s->nx*s->ny);
639: PetscReal val;
640: Vec y0, y1;
644: /* error y-x */
645: VecAXPY(s->y, -1.0, s->x);
646: /* VecView(s->y, (PetscViewer)PETSC_VIEWER_DEFAULT); */
648: /* error in velocity */
649: VecGetSubVector(s->y, s->isg[0], &y0);
650: VecNorm(y0, NORM_2, &val);
651: PetscPrintf(PETSC_COMM_WORLD," discretization error u = %g\n",(double)(PetscRealPart(val/scale)));
652: VecRestoreSubVector(s->y, s->isg[0], &y0);
654: /* error in pressure */
655: VecGetSubVector(s->y, s->isg[1], &y1);
656: VecNorm(y1, NORM_2, &val);
657: PetscPrintf(PETSC_COMM_WORLD," discretization error p = %g\n",(double)(PetscRealPart(val/scale)));
658: VecRestoreSubVector(s->y, s->isg[1], &y1);
660: /* total error */
661: VecNorm(s->y, NORM_2, &val);
662: PetscPrintf(PETSC_COMM_WORLD," discretization error [u,p] = %g\n", (double)PetscRealPart((val/scale)));
663: return(0);
664: }
666: int main(int argc, char **argv)
667: {
668: Stokes s;
669: KSP ksp;
672: PetscInitialize(&argc, &argv, NULL, help);
673: s.nx = 4;
674: s.ny = 6;
675: PetscOptionsGetInt(NULL,NULL, "-nx", &s.nx, NULL);
676: PetscOptionsGetInt(NULL,NULL, "-ny", &s.ny, NULL);
677: s.hx = 2.0/s.nx;
678: s.hy = 1.0/s.ny;
679: s.userPC = s.userKSP = PETSC_FALSE;
680: PetscOptionsHasName(NULL,NULL, "-user_pc", &s.userPC);
681: PetscOptionsHasName(NULL,NULL, "-user_ksp", &s.userKSP);
683: StokesSetupMatrix(&s);
684: StokesSetupIndexSets(&s);
685: StokesSetupVectors(&s);
687: KSPCreate(PETSC_COMM_WORLD, &ksp);
688: KSPSetOperators(ksp, s.A, s.A);
689: KSPSetFromOptions(ksp);
690: StokesSetupPC(&s, ksp);
691: KSPSolve(ksp, s.b, s.x);
693: /* don't trust, verify! */
694: StokesCalcResidual(&s);
695: StokesCalcError(&s);
696: StokesWriteSolution(&s);
698: KSPDestroy(&ksp);
699: MatDestroy(&s.subA[0]);
700: MatDestroy(&s.subA[1]);
701: MatDestroy(&s.subA[2]);
702: MatDestroy(&s.subA[3]);
703: MatDestroy(&s.A);
704: VecDestroy(&s.x);
705: VecDestroy(&s.b);
706: VecDestroy(&s.y);
707: MatDestroy(&s.myS);
708: PetscFinalize();
709: return 0;
710: }