Actual source code: ex2.c

  1: static char help[] = "Solve a toy 2D problem on a staggered grid\n\n";
  2: /*

  4:   To demonstrate the basic functionality of DMStag, solves an isoviscous
  5:   incompressible Stokes problem on a rectangular 2D domain, using a manufactured
  6:   solution.

  8:   u_xx + u_yy - p_x = f^x
  9:   v_xx + v_yy - p_y = f^y
 10:   u_x + v_y         = g

 12:   g is 0 in the physical case.

 14:   Boundary conditions give prescribed flow perpendicular to the boundaries,
 15:   and zero derivative perpendicular to them (free slip).

 17:   Use the -pinpressure option to fix a pressure node, instead of providing
 18:   a constant-pressure nullspace. This allows use of direct solvers, e.g. to
 19:   use UMFPACK,

 21:      ./ex2 -pinpressure 1 -pc_type lu -pc_factor_mat_solver_type umfpack

 23:   This example demonstrates the use of DMProduct to efficiently store coordinates
 24:   on an orthogonal grid.

 26: */
 27: #include <petscdm.h>
 28: #include <petscksp.h>
 29: #include <petscdmstag.h>

 31: /* Shorter, more convenient names for DMStagStencilLocation entries */
 32: #define DOWN_LEFT  DMSTAG_DOWN_LEFT
 33: #define DOWN       DMSTAG_DOWN
 34: #define DOWN_RIGHT DMSTAG_DOWN_RIGHT
 35: #define LEFT       DMSTAG_LEFT
 36: #define ELEMENT    DMSTAG_ELEMENT
 37: #define RIGHT      DMSTAG_RIGHT
 38: #define UP_LEFT    DMSTAG_UP_LEFT
 39: #define UP         DMSTAG_UP
 40: #define UP_RIGHT   DMSTAG_UP_RIGHT

 42: static PetscErrorCode CreateSystem(DM,Mat*,Vec*,PetscBool);
 43: static PetscErrorCode CreateReferenceSolution(DM,Vec*);
 44: static PetscErrorCode AttachNullspace(DM,Mat);
 45: static PetscErrorCode CheckSolution(Vec,Vec);

 47: /* Manufactured solution. Chosen to be higher order than can be solved exactly,
 48: and to have a zero derivative for flow parallel to the boundaries. That is,
 49: d(ux)/dy = 0 at the top and bottom boundaries, and d(uy)/dx = 0 at the right
 50: and left boundaries. */
 51: static PetscScalar uxRef(PetscScalar x,PetscScalar y){return 0.0*x + y*y - 2.0*y*y*y + y*y*y*y;}      /* no x-dependence  */
 52: static PetscScalar uyRef(PetscScalar x,PetscScalar y) {return x*x - 2.0*x*x*x + x*x*x*x + 0.0*y;}      /* no y-dependence  */
 53: static PetscScalar pRef (PetscScalar x,PetscScalar y) {return -1.0*(x-0.5) + -3.0/2.0*y*y + 0.5;}    /* zero integral    */
 54: static PetscScalar fx   (PetscScalar x,PetscScalar y) {return 0.0*x + 2.0 -12.0*y + 12.0*y*y + 1.0;} /* no x-dependence  */
 55: static PetscScalar fy   (PetscScalar x,PetscScalar y) {return 2.0 -12.0*x + 12.0*x*x + 3.0*y;}
 56: static PetscScalar g    (PetscScalar x,PetscScalar y) {return 0.0*x*y;}                              /* identically zero */

 58: int main(int argc,char **argv)
 59: {
 60:   DM             dmSol;
 61:   Vec            sol,solRef,rhs;
 62:   Mat            A;
 63:   KSP            ksp;
 64:   PC             pc;
 65:   PetscBool      pinPressure;

 67:   /* Initialize PETSc and process command line arguments */
 68:   PetscInitialize(&argc,&argv,(char*)0,help);
 69:   pinPressure = PETSC_TRUE;
 70:   PetscOptionsGetBool(NULL,NULL,"-pinpressure",&pinPressure,NULL);

 72:   /* Create 2D DMStag for the solution, and set up. */
 73:   {
 74:     const PetscInt dof0 = 0, dof1 = 1,dof2 = 1; /* 1 dof on each edge and element center */
 75:     const PetscInt stencilWidth = 1;
 76:     DMStagCreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,7,9,PETSC_DECIDE,PETSC_DECIDE,dof0,dof1,dof2,DMSTAG_STENCIL_BOX,stencilWidth,NULL,NULL,&dmSol);
 77:     DMSetFromOptions(dmSol);
 78:     DMSetUp(dmSol);
 79:   }

 81:   /* Define uniform coordinates as a product of 1D arrays */
 82:   DMStagSetUniformCoordinatesProduct(dmSol,0.0,1.0,0.0,1.0,0.0,0.0);

 84:   /* Compute (manufactured) reference solution */
 85:   CreateReferenceSolution(dmSol,&solRef);

 87:   /* Assemble system */
 88:   CreateSystem(dmSol,&A,&rhs,pinPressure);

 90:   /* Attach a constant-pressure nullspace to the operator
 91:   (logically, this should be in CreateSystem, but we separate it here to highlight
 92:    the features of DMStag exposed, in particular DMStagMigrateVec()) */
 93:   if (!pinPressure) {
 94:     AttachNullspace(dmSol,A);
 95:   }

 97:   /* Solve, using the default FieldSplit (Approximate Block Factorization) Preconditioner
 98:      This is not intended to be an example of a good solver!  */
 99:   DMCreateGlobalVector(dmSol,&sol);
100:   KSPCreate(PETSC_COMM_WORLD,&ksp);
101:   KSPSetType(ksp,KSPFGMRES);
102:   KSPSetOperators(ksp,A,A);
103:   KSPGetPC(ksp,&pc);
104:   PCSetType(pc,PCFIELDSPLIT);
105:   PCFieldSplitSetDetectSaddlePoint(pc,PETSC_TRUE);
106:   KSPSetFromOptions(ksp);
107:   KSPSolve(ksp,rhs,sol);

109:   /* Check Solution */
110:   CheckSolution(sol,solRef);

112:   /* Clean up and finalize PETSc */
113:   KSPDestroy(&ksp);
114:   VecDestroy(&sol);
115:   VecDestroy(&solRef);
116:   VecDestroy(&rhs);
117:   MatDestroy(&A);
118:   DMDestroy(&dmSol);
119:   PetscFinalize();
120:   return 0;
121: }

123: /*
124: Note: this system is not well-scaled! Generally one would adjust the equations
125:  to try to get matrix entries to be of comparable order, regardless of grid spacing
126:  or choice of coefficients.
127: */
128: static PetscErrorCode CreateSystem(DM dmSol,Mat *pA,Vec *pRhs, PetscBool pinPressure)
129: {
130:   PetscInt       N[2];
131:   PetscInt       ex,ey,startx,starty,nx,ny;
132:   PetscInt       iprev,icenter,inext;
133:   Mat            A;
134:   Vec            rhs;
135:   PetscReal      hx,hy;
136:   PetscScalar    **cArrX,**cArrY;

138:   /* Here, we showcase two different methods for manipulating local vector entries.
139:      One can use DMStagStencil objects with DMStagVecSetValuesStencil(),
140:      making sure to call VecAssemble[Begin/End]() after all values are set.
141:      Alternately, one can use DMStagVecGetArray[Read]() and DMStagVecRestoreArray[Read]().
142:      The first approach is used to build the rhs, and the second is used to
143:      obtain coordinate values. Working with the array is almost certainly more efficient,
144:      but only allows setting local entries, requires understanding which "slot" to use,
145:      and doesn't correspond as precisely to the matrix assembly process using DMStagStencil objects */

148:   DMCreateMatrix(dmSol,pA);
149:   A = *pA;
150:   DMCreateGlobalVector(dmSol,pRhs);
151:   rhs = *pRhs;
152:   DMStagGetCorners(dmSol,&startx,&starty,NULL,&nx,&ny,NULL,NULL,NULL,NULL);
153:   DMStagGetGlobalSizes(dmSol,&N[0],&N[1],NULL);
154:   hx = 1.0/N[0]; hy = 1.0/N[1];
155:   DMStagGetProductCoordinateArraysRead(dmSol,&cArrX,&cArrY,NULL);
156:   DMStagGetProductCoordinateLocationSlot(dmSol,ELEMENT,&icenter);
157:   DMStagGetProductCoordinateLocationSlot(dmSol,LEFT,&iprev);
158:   DMStagGetProductCoordinateLocationSlot(dmSol,RIGHT,&inext);

160:   /* Loop over all local elements. Note that it may be more efficient in real
161:      applications to loop over each boundary separately */
162:   for (ey = starty; ey<starty+ny; ++ey) { /* With DMStag, always iterate x fastest, y second fastest, z slowest */
163:     for (ex = startx; ex<startx+nx; ++ex) {

165:       if (ex == N[0]-1) {
166:         /* Right Boundary velocity Dirichlet */
167:         DMStagStencil row;
168:         PetscScalar   valRhs;
169:         const PetscScalar valA = 1.0;
170:         row.i = ex; row.j = ey; row.loc = RIGHT; row.c = 0;
171:         DMStagMatSetValuesStencil(dmSol,A,1,&row,1,&row,&valA,INSERT_VALUES);
172:         valRhs = uxRef(cArrX[ex][inext],cArrY[ey][icenter]);
173:         DMStagVecSetValuesStencil(dmSol,rhs,1,&row,&valRhs,INSERT_VALUES);
174:       }
175:       if (ey == N[1]-1) {
176:         /* Top boundary velocity Dirichlet */
177:         DMStagStencil row;
178:         PetscScalar   valRhs;
179:         const PetscScalar valA = 1.0;
180:         row.i = ex; row.j = ey; row.loc = UP; row.c = 0;
181:         DMStagMatSetValuesStencil(dmSol,A,1,&row,1,&row,&valA,INSERT_VALUES);
182:         valRhs = uyRef(cArrX[ex][icenter],cArrY[ey][inext]);
183:         DMStagVecSetValuesStencil(dmSol,rhs,1,&row,&valRhs,INSERT_VALUES);
184:       }

186:       if (ey == 0) {
187:         /* Bottom boundary velocity Dirichlet */
188:         DMStagStencil row;
189:         PetscScalar   valRhs;
190:         const PetscScalar valA = 1.0;
191:         row.i = ex; row.j = ey; row.loc = DOWN; row.c = 0;
192:         DMStagMatSetValuesStencil(dmSol,A,1,&row,1,&row,&valA,INSERT_VALUES);
193:         valRhs = uyRef(cArrX[ex][icenter],cArrY[ey][iprev]);
194:         DMStagVecSetValuesStencil(dmSol,rhs,1,&row,&valRhs,INSERT_VALUES);
195:       } else {
196:         /* Y-momentum equation : (u_xx + u_yy) - p_y = f^y */
197:         DMStagStencil row,col[7];
198:         PetscScalar   valA[7],valRhs;
199:         PetscInt      nEntries;

201:         row.i    = ex  ; row.j    = ey  ; row.loc    = DOWN;    row.c     = 0;
202:         if (ex == 0) {
203:           nEntries = 6;
204:           col[0].i = ex  ; col[0].j = ey  ; col[0].loc = DOWN;    col[0].c  = 0; valA[0] = -1.0 / (hx*hx) -2.0 / (hy*hy);
205:           col[1].i = ex  ; col[1].j = ey-1; col[1].loc = DOWN;    col[1].c  = 0; valA[1] =  1.0 / (hy*hy);
206:           col[2].i = ex  ; col[2].j = ey+1; col[2].loc = DOWN;    col[2].c  = 0; valA[2] =  1.0 / (hy*hy);
207:           /* Missing left element */
208:           col[3].i = ex+1; col[3].j = ey  ; col[3].loc = DOWN;    col[3].c  = 0; valA[3] =  1.0 / (hx*hx);
209:           col[4].i = ex  ; col[4].j = ey-1; col[4].loc = ELEMENT; col[4].c  = 0; valA[4] =  1.0 / hy;
210:           col[5].i = ex  ; col[5].j = ey  ; col[5].loc = ELEMENT; col[5].c  = 0; valA[5] = -1.0 / hy;
211:         } else if (ex == N[0]-1) {
212:           /* Right boundary y velocity stencil */
213:           nEntries = 6;
214:           col[0].i = ex  ; col[0].j = ey  ; col[0].loc = DOWN;    col[0].c  = 0; valA[0] = -1.0 / (hx*hx) -2.0 / (hy*hy);
215:           col[1].i = ex  ; col[1].j = ey-1; col[1].loc = DOWN;    col[1].c  = 0; valA[1] =  1.0 / (hy*hy);
216:           col[2].i = ex  ; col[2].j = ey+1; col[2].loc = DOWN;    col[2].c  = 0; valA[2] =  1.0 / (hy*hy);
217:           col[3].i = ex-1; col[3].j = ey  ; col[3].loc = DOWN;    col[3].c  = 0; valA[3] =  1.0 / (hx*hx);
218:           /* Missing right element */
219:           col[4].i = ex  ; col[4].j = ey-1; col[4].loc = ELEMENT; col[4].c  = 0; valA[4] =  1.0 / hy;
220:           col[5].i = ex  ; col[5].j = ey  ; col[5].loc = ELEMENT; col[5].c  = 0; valA[5] = -1.0 / hy;
221:         } else {
222:           nEntries = 7;
223:           col[0].i = ex  ; col[0].j = ey  ; col[0].loc = DOWN;    col[0].c  = 0; valA[0] = -2.0 / (hx*hx) -2.0 / (hy*hy);
224:           col[1].i = ex  ; col[1].j = ey-1; col[1].loc = DOWN;    col[1].c  = 0; valA[1] =  1.0 / (hy*hy);
225:           col[2].i = ex  ; col[2].j = ey+1; col[2].loc = DOWN;    col[2].c  = 0; valA[2] =  1.0 / (hy*hy);
226:           col[3].i = ex-1; col[3].j = ey  ; col[3].loc = DOWN;    col[3].c  = 0; valA[3] =  1.0 / (hx*hx);
227:           col[4].i = ex+1; col[4].j = ey  ; col[4].loc = DOWN;    col[4].c  = 0; valA[4] =  1.0 / (hx*hx);
228:           col[5].i = ex  ; col[5].j = ey-1; col[5].loc = ELEMENT; col[5].c  = 0; valA[5] =  1.0 / hy;
229:           col[6].i = ex  ; col[6].j = ey  ; col[6].loc = ELEMENT; col[6].c  = 0; valA[6] = -1.0 / hy;
230:         }
231:         DMStagMatSetValuesStencil(dmSol,A,1,&row,nEntries,col,valA,INSERT_VALUES);
232:         valRhs = fy(cArrX[ex][icenter],cArrY[ey][iprev]);
233:         DMStagVecSetValuesStencil(dmSol,rhs,1,&row,&valRhs,INSERT_VALUES);
234:       }

236:       if (ex == 0) {
237:         /* Left velocity Dirichlet */
238:         DMStagStencil row;
239:         PetscScalar   valRhs;
240:         const PetscScalar valA = 1.0;
241:         row.i = ex; row.j = ey; row.loc = LEFT; row.c = 0;
242:         DMStagMatSetValuesStencil(dmSol,A,1,&row,1,&row,&valA,INSERT_VALUES);
243:         valRhs = uxRef(cArrX[ex][iprev],cArrY[ey][icenter]);
244:         DMStagVecSetValuesStencil(dmSol,rhs,1,&row,&valRhs,INSERT_VALUES);
245:       } else {
246:         /* X-momentum equation : (u_xx + u_yy) - p_x = f^x */
247:         DMStagStencil row,col[7];
248:         PetscScalar   valA[7],valRhs;
249:         PetscInt      nEntries;
250:         row.i    = ex  ; row.j    = ey  ; row.loc    = LEFT;    row.c     = 0;

252:         if (ey == 0) {
253:           nEntries = 6;
254:           col[0].i = ex  ; col[0].j = ey  ; col[0].loc = LEFT;    col[0].c  = 0; valA[0] = -2.0 /(hx*hx) -1.0 /(hy*hy);
255:           /* missing term from element below */
256:           col[1].i = ex  ; col[1].j = ey+1; col[1].loc = LEFT;    col[1].c  = 0; valA[1] =  1.0 / (hy*hy);
257:           col[2].i = ex-1; col[2].j = ey  ; col[2].loc = LEFT;    col[2].c  = 0; valA[2] =  1.0 / (hx*hx);
258:           col[3].i = ex+1; col[3].j = ey  ; col[3].loc = LEFT;    col[3].c  = 0; valA[3] =  1.0 / (hx*hx);
259:           col[4].i = ex-1; col[4].j = ey  ; col[4].loc = ELEMENT; col[4].c  = 0; valA[4] =  1.0 / hx;
260:           col[5].i = ex  ; col[5].j = ey  ; col[5].loc = ELEMENT; col[5].c  = 0; valA[5] = -1.0 / hx;
261:         } else if (ey == N[1]-1) {
262:           /* Top boundary x velocity stencil */
263:           nEntries = 6;
264:           row.i    = ex  ; row.j    = ey  ; row.loc    = LEFT;    row.c     = 0;
265:           col[0].i = ex  ; col[0].j = ey  ; col[0].loc = LEFT;    col[0].c  = 0; valA[0] = -2.0 / (hx*hx) -1.0 / (hy*hy);
266:           col[1].i = ex  ; col[1].j = ey-1; col[1].loc = LEFT;    col[1].c  = 0; valA[1] =  1.0 / (hy*hy);
267:           /* Missing element above term */
268:           col[2].i = ex-1; col[2].j = ey  ; col[2].loc = LEFT;    col[2].c  = 0; valA[2] =  1.0 / (hx*hx);
269:           col[3].i = ex+1; col[3].j = ey  ; col[3].loc = LEFT;    col[3].c  = 0; valA[3] =  1.0 / (hx*hx);
270:           col[4].i = ex-1; col[4].j = ey  ; col[4].loc = ELEMENT; col[4].c  = 0; valA[4] =  1.0 / hx;
271:           col[5].i = ex  ; col[5].j = ey  ; col[5].loc = ELEMENT; col[5].c  = 0; valA[5] = -1.0 / hx;
272:         } else {
273:           /* Note how this is identical to the stencil for U_y, with "DOWN" replaced by "LEFT" and the pressure derivative in the other direction */
274:           nEntries = 7;
275:           col[0].i = ex  ; col[0].j = ey  ; col[0].loc = LEFT;    col[0].c  = 0; valA[0] = -2.0 / (hx*hx) + -2.0 / (hy*hy);
276:           col[1].i = ex  ; col[1].j = ey-1; col[1].loc = LEFT;    col[1].c  = 0; valA[1] =  1.0 / (hy*hy);
277:           col[2].i = ex  ; col[2].j = ey+1; col[2].loc = LEFT;    col[2].c  = 0; valA[2] =  1.0 / (hy*hy);
278:           col[3].i = ex-1; col[3].j = ey  ; col[3].loc = LEFT;    col[3].c  = 0; valA[3] =  1.0 / (hx*hx);
279:           col[4].i = ex+1; col[4].j = ey  ; col[4].loc = LEFT;    col[4].c  = 0; valA[4] =  1.0 / (hx*hx);
280:           col[5].i = ex-1; col[5].j = ey  ; col[5].loc = ELEMENT; col[5].c  = 0; valA[5] =  1.0 / hx;
281:           col[6].i = ex  ; col[6].j = ey  ; col[6].loc = ELEMENT; col[6].c  = 0; valA[6] = -1.0 / hx;

283:         }
284:         DMStagMatSetValuesStencil(dmSol,A,1,&row,nEntries,col,valA,INSERT_VALUES);
285:         valRhs = fx(cArrX[ex][iprev],cArrY[ey][icenter]);
286:         DMStagVecSetValuesStencil(dmSol,rhs,1,&row,&valRhs,INSERT_VALUES);
287:       }

289:       /* P equation : u_x + v_y = g
290:          Note that this includes an explicit zero on the diagonal. This is only needed for
291:          direct solvers (not required if using an iterative solver and setting the constant-pressure nullspace) */
292:       if (pinPressure && ex == 0 && ey == 0) { /* Pin the first pressure node, if requested */
293:         DMStagStencil row;
294:         PetscScalar valA,valRhs;
295:         row.i = ex; row.j = ey; row.loc  = ELEMENT; row.c = 0;
296:         valA = 1.0;
297:         DMStagMatSetValuesStencil(dmSol,A,1,&row,1,&row,&valA,INSERT_VALUES);
298:         valRhs = pRef(cArrX[ex][icenter],cArrY[ey][icenter]);
299:         DMStagVecSetValuesStencil(dmSol,rhs,1,&row,&valRhs,INSERT_VALUES);
300:       } else {
301:         DMStagStencil row,col[5];
302:         PetscScalar   valA[5],valRhs;

304:         row.i    = ex; row.j    = ey; row.loc    = ELEMENT; row.c    = 0;
305:         col[0].i = ex; col[0].j = ey; col[0].loc = LEFT;    col[0].c = 0; valA[0] = -1.0 / hx;
306:         col[1].i = ex; col[1].j = ey; col[1].loc = RIGHT;   col[1].c = 0; valA[1] =  1.0 / hx;
307:         col[2].i = ex; col[2].j = ey; col[2].loc = DOWN;    col[2].c = 0; valA[2] = -1.0 / hy;
308:         col[3].i = ex; col[3].j = ey; col[3].loc = UP;      col[3].c = 0; valA[3] =  1.0 / hy;
309:         col[4] = row;                                                     valA[4] = 0.0;
310:         DMStagMatSetValuesStencil(dmSol,A,1,&row,5,col,valA,INSERT_VALUES);
311:         valRhs = g(cArrX[ex][icenter],cArrY[ey][icenter]);
312:         DMStagVecSetValuesStencil(dmSol,rhs,1,&row,&valRhs,INSERT_VALUES);
313:       }
314:     }
315:   }
316:   DMStagRestoreProductCoordinateArraysRead(dmSol,&cArrX,&cArrY,NULL);
317:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
318:   VecAssemblyBegin(rhs);
319:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
320:   VecAssemblyEnd(rhs);
321:   return 0;
322: }

324: /* Create a pressure-only DMStag and use it to generate a nullspace vector
325:    - Create a compatible DMStag with one dof per element (and nothing else).
326:    - Create a constant vector and normalize it
327:    - Migrate it to a vector on the original dmSol (making use of the fact
328:    that this will fill in zeros for "extra" dof)
329:    - Set the nullspace for the operator
330:    - Destroy everything (the operator keeps the references it needs) */
331: static PetscErrorCode AttachNullspace(DM dmSol,Mat A)
332: {
333:   DM             dmPressure;
334:   Vec            constantPressure,basis;
335:   PetscReal      nrm;
336:   MatNullSpace   matNullSpace;

339:   DMStagCreateCompatibleDMStag(dmSol,0,0,1,0,&dmPressure);
340:   DMGetGlobalVector(dmPressure,&constantPressure);
341:   VecSet(constantPressure,1.0);
342:   VecNorm(constantPressure,NORM_2,&nrm);
343:   VecScale(constantPressure,1.0/nrm);
344:   DMCreateGlobalVector(dmSol,&basis);
345:   DMStagMigrateVec(dmPressure,constantPressure,dmSol,basis);
346:   MatNullSpaceCreate(PetscObjectComm((PetscObject)dmSol),PETSC_FALSE,1,&basis,&matNullSpace);
347:   VecDestroy(&basis);
348:   VecDestroy(&constantPressure);
349:   MatSetNullSpace(A,matNullSpace);
350:   MatNullSpaceDestroy(&matNullSpace);
351:   return 0;
352: }

354: /* Create a reference solution.
355:    Here, we use the more direct method of iterating over arrays.  */
356: static PetscErrorCode CreateReferenceSolution(DM dmSol,Vec *pSolRef)
357: {
358:   PetscInt       startx,starty,nx,ny,nExtra[2],ex,ey;
359:   PetscInt       iuy,iux,ip,iprev,icenter;
360:   PetscScalar    ***arrSol,**cArrX,**cArrY;
361:   Vec            solRefLocal;

364:   DMCreateGlobalVector(dmSol,pSolRef);
365:   DMGetLocalVector(dmSol,&solRefLocal);

367:   /* Obtain indices to use in the raw arrays */
368:   DMStagGetLocationSlot(dmSol,DOWN,0,&iuy);
369:   DMStagGetLocationSlot(dmSol,LEFT,0,&iux);
370:   DMStagGetLocationSlot(dmSol,ELEMENT,0,&ip);

372:   /* Use high-level convenience functions to get raw arrays and indices for 1d coordinates */
373:   DMStagGetProductCoordinateArraysRead(dmSol,&cArrX,&cArrY,NULL);
374:   DMStagGetProductCoordinateLocationSlot(dmSol,ELEMENT,&icenter);
375:   DMStagGetProductCoordinateLocationSlot(dmSol,LEFT,&iprev);

377:   DMStagVecGetArray(dmSol,solRefLocal,&arrSol);
378:   DMStagGetCorners(dmSol,&startx,&starty,NULL,&nx,&ny,NULL,&nExtra[0],&nExtra[1],NULL);
379:   for (ey=starty; ey<starty + ny + nExtra[1]; ++ey) {
380:     for (ex=startx; ex<startx + nx + nExtra[0]; ++ex) {
381:       arrSol[ey][ex][iuy] = uyRef(cArrX[ex][icenter],cArrY[ey][iprev]);
382:       arrSol[ey][ex][iux] = uxRef(cArrX[ex][iprev],cArrY[ey][icenter]);
383:       if (ey < starty+ny && ex < startx+nx) { /* Don't fill on the dummy elements (though you could, and these values would just be ignored) */
384:         arrSol[ey][ex][ip]  = pRef(cArrX[ex][icenter],cArrY[ey][icenter]);
385:       }
386:     }
387:   }
388:   DMStagVecRestoreArray(dmSol,solRefLocal,&arrSol);
389:   DMStagRestoreProductCoordinateArraysRead(dmSol,&cArrX,&cArrY,NULL);
390:   DMLocalToGlobal(dmSol,solRefLocal,INSERT_VALUES,*pSolRef);
391:   DMRestoreLocalVector(dmSol,&solRefLocal);
392:   return 0;
393: }

395: static PetscErrorCode CheckSolution(Vec sol,Vec solRef)
396: {
397:   Vec            diff;
398:   PetscReal      normsolRef,errAbs,errRel;

401:   VecDuplicate(sol,&diff);
402:   VecCopy(sol,diff);
403:   VecAXPY(diff,-1.0,solRef);
404:   VecNorm(diff,NORM_2,&errAbs);
405:   VecNorm(solRef,NORM_2,&normsolRef);
406:   errRel = errAbs/normsolRef;
407:   PetscPrintf(PETSC_COMM_WORLD,"Error (abs): %g\nError (rel): %g\n",(double)errAbs,(double)errRel);
408:   VecDestroy(&diff);
409:   return 0;
410: }

412: /*TEST

414:    test:
415:       suffix: 1
416:       nsize: 4
417:       args: -ksp_monitor_short -ksp_converged_reason

419:    test:
420:       suffix: direct_umfpack
421:       requires: suitesparse
422:       nsize: 1
423:       args: -pinpressure 1 -stag_grid_x 8 -stag_grid_y 6 -ksp_monitor_short -pc_type lu -pc_factor_mat_solver_type umfpack

425: TEST*/