Actual source code: ex43.c

petsc-master 2016-07-25
Report Typos and Errors
  1: static char help[] = "Solves the incompressible, variable viscosity Stokes equation in 2d on the unit domain \n\
  2: using Q1Q1 elements, stabilized with Bochev's polynomial projection method. \n\
  3: The models defined utilise free slip boundary conditions on all sides. \n\
  4: Options: \n"
  5: "\
  6:      -mx : Number of elements in the x-direction \n\
  7:      -my : Number of elements in the y-direction \n\
  8:      -o : Specify output filename for solution (will be petsc binary format or paraview format if the extension is .vts) \n\
  9:      -gnuplot : Output Gauss point coordinates, coefficients and u,p solution in gnuplot format \n\
 10:      -c_str : Indicates the structure of the coefficients to use \n"
 11: "\
 12:           -c_str 0 => Coefficient definition for an analytic solution with a vertical jump in viscosity at x = xc \n\
 13:                       This problem is driven by the forcing function f(x,y) = (0, sin(nz pi y)cos(pi x) \n\
 14:                          Parameters: \n\
 15:                               -solcx_eta0  : Viscosity to the left of the interface \n\
 16:                               -solcx_eta1  : Viscosity to the right of the interface \n\
 17:                               -solcx_xc    : Location of the interface \n\
 18:                               -solcx_nz    : Wavenumber in the y direction \n"
 19: "\
 20:           -c_str 1 => Coefficient definition for a dense rectangular blob located at the center of the domain \n\
 21:                          Parameters: \n\
 22:                               -sinker_eta0 : Viscosity of the background fluid \n\
 23:                               -sinker_eta1 : Viscosity of the blob \n\
 24:                               -sinker_dx   : Width of the blob \n\
 25:                               -sinker_dy   : Height of the blob \n"
 26: "\
 27:           -c_str 2 => Coefficient definition for a dense circular blob located at the center of the domain \n\
 28:                          Parameters: \n\
 29:                               -sinker_eta0 : Viscosity of the background fluid \n\
 30:                               -sinker_eta1 : Viscosity of the blob \n\
 31:                               -sinker_r    : Radius of the blob \n"
 32: "\
 33:           -c_str 3 => Coefficient definition for a dense circular and rectangular inclusion (located at the center of the domain) \n\
 34:                               -sinker_eta0 : Viscosity of the background fluid \n\
 35:                               -sinker_eta1 : Viscosity of the two inclusions \n\
 36:                               -sinker_r    : Radius of the circular inclusion \n\
 37:                               -sinker_c0x  : Origin (x-coord) of the circular inclusion \n\
 38:                               -sinker_c0y  : Origin (y-coord) of the circular inclusion \n\
 39:                               -sinker_dx   : Width of the rectangular inclusion \n\
 40:                               -sinker_dy   : Height of the rectangular inclusion \n\
 41:                               -sinker_phi  : Rotation angle of the rectangular inclusion \n\
 42:      -use_gp_coords : Evaluate the viscosity and force term at the global coordinates of each quadrature point \n\
 43:                       By default, the viscosity and force term are evaulated at the element center and applied as a constant over the entire element \n";

 45: /* Contributed by Dave May */

 47: #include <petscksp.h>
 48: #include <petscdm.h>
 49: #include <petscdmda.h>

 51: /* A Maple-generated exact solution created by Mirko Velic (mirko.velic@sci.monash.edu.au) */
 52: #include "ex43-solcx.h"

 54: static PetscErrorCode DMDABCApplyFreeSlip(DM,Mat,Vec);


 57: #define NSD            2 /* number of spatial dimensions */
 58: #define NODES_PER_EL   4 /* nodes per element */
 59: #define U_DOFS         2 /* degrees of freedom per velocity node */
 60: #define P_DOFS         1 /* degrees of freedom per pressure node */
 61: #define GAUSS_POINTS   4

 63: /* Gauss point based evaluation 8+4+4+4 = 20 */
 64: typedef struct {
 65:   PetscScalar gp_coords[2*GAUSS_POINTS];
 66:   PetscScalar eta[GAUSS_POINTS];
 67:   PetscScalar fx[GAUSS_POINTS];
 68:   PetscScalar fy[GAUSS_POINTS];
 69: } GaussPointCoefficients;

 71: typedef struct {
 72:   PetscScalar u_dof;
 73:   PetscScalar v_dof;
 74:   PetscScalar p_dof;
 75: } StokesDOF;


 78: /*
 79:  Element: Local basis function ordering
 80:  1-----2
 81:  |     |
 82:  |     |
 83:  0-----3
 84:  */
 85: static void ConstructQ12D_Ni(PetscScalar _xi[],PetscScalar Ni[])
 86: {
 87:   PetscScalar xi  = _xi[0];
 88:   PetscScalar eta = _xi[1];

 90:   Ni[0] = 0.25*(1.0-xi)*(1.0-eta);
 91:   Ni[1] = 0.25*(1.0-xi)*(1.0+eta);
 92:   Ni[2] = 0.25*(1.0+xi)*(1.0+eta);
 93:   Ni[3] = 0.25*(1.0+xi)*(1.0-eta);
 94: }

 96: static void ConstructQ12D_GNi(PetscScalar _xi[],PetscScalar GNi[][NODES_PER_EL])
 97: {
 98:   PetscScalar xi  = _xi[0];
 99:   PetscScalar eta = _xi[1];

101:   GNi[0][0] = -0.25*(1.0-eta);
102:   GNi[0][1] = -0.25*(1.0+eta);
103:   GNi[0][2] =   0.25*(1.0+eta);
104:   GNi[0][3] =   0.25*(1.0-eta);

106:   GNi[1][0] = -0.25*(1.0-xi);
107:   GNi[1][1] =   0.25*(1.0-xi);
108:   GNi[1][2] =   0.25*(1.0+xi);
109:   GNi[1][3] = -0.25*(1.0+xi);
110: }

112: static void ConstructQ12D_GNx(PetscScalar GNi[][NODES_PER_EL],PetscScalar GNx[][NODES_PER_EL],PetscScalar coords[],PetscScalar *det_J)
113: {
114:   PetscScalar J00,J01,J10,J11,J;
115:   PetscScalar iJ00,iJ01,iJ10,iJ11;
116:   PetscInt    i;

118:   J00 = J01 = J10 = J11 = 0.0;
119:   for (i = 0; i < NODES_PER_EL; i++) {
120:     PetscScalar cx = coords[2*i];
121:     PetscScalar cy = coords[2*i+1];

123:     J00 = J00+GNi[0][i]*cx;      /* J_xx = dx/dxi */
124:     J01 = J01+GNi[0][i]*cy;      /* J_xy = dy/dxi */
125:     J10 = J10+GNi[1][i]*cx;      /* J_yx = dx/deta */
126:     J11 = J11+GNi[1][i]*cy;      /* J_yy = dy/deta */
127:   }
128:   J = (J00*J11)-(J01*J10);

130:   iJ00 =  J11/J;
131:   iJ01 = -J01/J;
132:   iJ10 = -J10/J;
133:   iJ11 =  J00/J;

135:   for (i = 0; i < NODES_PER_EL; i++) {
136:     GNx[0][i] = GNi[0][i]*iJ00+GNi[1][i]*iJ01;
137:     GNx[1][i] = GNi[0][i]*iJ10+GNi[1][i]*iJ11;
138:   }

140:   if (det_J) *det_J = J;
141: }

143: static void ConstructGaussQuadrature(PetscInt *ngp,PetscScalar gp_xi[][2],PetscScalar gp_weight[])
144: {
145:   *ngp         = 4;
146:   gp_xi[0][0]  = -0.57735026919; gp_xi[0][1] = -0.57735026919;
147:   gp_xi[1][0]  = -0.57735026919; gp_xi[1][1] =  0.57735026919;
148:   gp_xi[2][0]  =  0.57735026919; gp_xi[2][1] =  0.57735026919;
149:   gp_xi[3][0]  =  0.57735026919; gp_xi[3][1] = -0.57735026919;
150:   gp_weight[0] = 1.0;
151:   gp_weight[1] = 1.0;
152:   gp_weight[2] = 1.0;
153:   gp_weight[3] = 1.0;
154: }

156: /* procs to the left claim the ghost node as their element */
159: static PetscErrorCode DMDAGetLocalElementSize(DM da,PetscInt *mxl,PetscInt *myl,PetscInt *mzl)
160: {
161:   PetscInt m,n,p,M,N,P;
162:   PetscInt sx,sy,sz;

165:   DMDAGetInfo(da,0,&M,&N,&P,0,0,0,0,0,0,0,0,0);
166:   DMDAGetCorners(da,&sx,&sy,&sz,&m,&n,&p);

168:   if (mxl) {
169:     *mxl = m;
170:     if ((sx+m) == M) *mxl = m-1;  /* last proc */
171:   }
172:   if (myl) {
173:     *myl = n;
174:     if ((sy+n) == N) *myl = n-1;  /* last proc */
175:   }
176:   if (mzl) {
177:     *mzl = p;
178:     if ((sz+p) == P) *mzl = p-1;  /* last proc */
179:   }
180:   return(0);
181: }

185: static PetscErrorCode DMDAGetElementCorners(DM da,PetscInt *sx,PetscInt *sy,PetscInt *sz,PetscInt *mx,PetscInt *my,PetscInt *mz)
186: {
187:   PetscInt       si,sj,sk;

191:   DMDAGetGhostCorners(da,&si,&sj,&sk,0,0,0);

193:   *sx = si;
194:   if (si) *sx = si+1;

196:   *sy = sj;
197:   if (sj) *sy = sj+1;

199:   if (sz) {
200:     *sz = sk;
201:     if (sk) *sz = sk+1;
202:   }

204:   DMDAGetLocalElementSize(da,mx,my,mz);
205:   return(0);
206: }

208: /*
209: i,j are the element indices
210: The unknown is a vector quantity.
211: The s[].c is used to indicate the degree of freedom.
212: */
215: static PetscErrorCode DMDAGetElementEqnums_up(MatStencil s_u[],MatStencil s_p[],PetscInt i,PetscInt j)
216: {
218:   /* velocity */
219:   /* node 0 */
220:   s_u[0].i = i; s_u[0].j = j; s_u[0].c = 0;       /* Vx0 */
221:   s_u[1].i = i; s_u[1].j = j; s_u[1].c = 1;       /* Vy0 */

223:   /* node 1 */
224:   s_u[2].i = i; s_u[2].j = j+1; s_u[2].c = 0;     /* Vx1 */
225:   s_u[3].i = i; s_u[3].j = j+1; s_u[3].c = 1;     /* Vy1 */

227:   /* node 2 */
228:   s_u[4].i = i+1; s_u[4].j = j+1; s_u[4].c = 0;   /* Vx2 */
229:   s_u[5].i = i+1; s_u[5].j = j+1; s_u[5].c = 1;   /* Vy2 */

231:   /* node 3 */
232:   s_u[6].i = i+1; s_u[6].j = j; s_u[6].c = 0;     /* Vx3 */
233:   s_u[7].i = i+1; s_u[7].j = j; s_u[7].c = 1;     /* Vy3 */

235:   /* pressure */
236:   s_p[0].i = i;   s_p[0].j = j;   s_p[0].c = 2; /* P0 */
237:   s_p[1].i = i;   s_p[1].j = j+1; s_p[1].c = 2; /* P1 */
238:   s_p[2].i = i+1; s_p[2].j = j+1; s_p[2].c = 2; /* P2 */
239:   s_p[3].i = i+1; s_p[3].j = j;   s_p[3].c = 2; /* P3 */
240:   return(0);
241: }

245: static PetscErrorCode DMDAGetElementOwnershipRanges2d(DM da,PetscInt **_lx,PetscInt **_ly)
246: {
248:   PetscMPIInt    rank;
249:   PetscInt       proc_I,proc_J;
250:   PetscInt       cpu_x,cpu_y;
251:   PetscInt       local_mx,local_my;
252:   Vec            vlx,vly;
253:   PetscInt       *LX,*LY,i;
254:   PetscScalar    *_a;
255:   Vec            V_SEQ;
256:   VecScatter     ctx;

259:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

261:   DMDAGetInfo(da,0,0,0,0,&cpu_x,&cpu_y,0,0,0,0,0,0,0);

263:   proc_J = rank/cpu_x;
264:   proc_I = rank-cpu_x*proc_J;

266:   PetscMalloc1(cpu_x,&LX);
267:   PetscMalloc1(cpu_y,&LY);

269:   DMDAGetLocalElementSize(da,&local_mx,&local_my,NULL);
270:   VecCreate(PETSC_COMM_WORLD,&vlx);
271:   VecSetSizes(vlx,PETSC_DECIDE,cpu_x);
272:   VecSetFromOptions(vlx);

274:   VecCreate(PETSC_COMM_WORLD,&vly);
275:   VecSetSizes(vly,PETSC_DECIDE,cpu_y);
276:   VecSetFromOptions(vly);

278:   VecSetValue(vlx,proc_I,(PetscScalar)(local_mx+1.0e-9),INSERT_VALUES);
279:   VecSetValue(vly,proc_J,(PetscScalar)(local_my+1.0e-9),INSERT_VALUES);
280:   VecAssemblyBegin(vlx);VecAssemblyEnd(vlx);
281:   VecAssemblyBegin(vly);VecAssemblyEnd(vly);

283:   VecScatterCreateToAll(vlx,&ctx,&V_SEQ);
284:   VecScatterBegin(ctx,vlx,V_SEQ,INSERT_VALUES,SCATTER_FORWARD);
285:   VecScatterEnd(ctx,vlx,V_SEQ,INSERT_VALUES,SCATTER_FORWARD);
286:   VecGetArray(V_SEQ,&_a);
287:   for (i = 0; i < cpu_x; i++) LX[i] = (PetscInt)PetscRealPart(_a[i]);
288:   VecRestoreArray(V_SEQ,&_a);
289:   VecScatterDestroy(&ctx);
290:   VecDestroy(&V_SEQ);

292:   VecScatterCreateToAll(vly,&ctx,&V_SEQ);
293:   VecScatterBegin(ctx,vly,V_SEQ,INSERT_VALUES,SCATTER_FORWARD);
294:   VecScatterEnd(ctx,vly,V_SEQ,INSERT_VALUES,SCATTER_FORWARD);
295:   VecGetArray(V_SEQ,&_a);
296:   for (i = 0; i < cpu_y; i++) LY[i] = (PetscInt)PetscRealPart(_a[i]);
297:   VecRestoreArray(V_SEQ,&_a);
298:   VecScatterDestroy(&ctx);
299:   VecDestroy(&V_SEQ);

301:   *_lx = LX;
302:   *_ly = LY;

304:   VecDestroy(&vlx);
305:   VecDestroy(&vly);
306:   return(0);
307: }

311: static PetscErrorCode DMDACoordViewGnuplot2d(DM da,const char prefix[])
312: {
313:   DM             cda;
314:   Vec            coords;
315:   DMDACoor2d     **_coords;
316:   PetscInt       si,sj,nx,ny,i,j;
317:   FILE           *fp;
318:   char           fname[PETSC_MAX_PATH_LEN];
319:   PetscMPIInt    rank;

323:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
324:   PetscSNPrintf(fname,sizeof(fname),"%s-p%1.4d.dat",prefix,rank);
325:   PetscFOpen(PETSC_COMM_SELF,fname,"w",&fp);
326:   if (!fp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot open file");

328:   PetscFPrintf(PETSC_COMM_SELF,fp,"### Element geometry for processor %1.4d ### \n",rank);

330:   DMGetCoordinateDM(da,&cda);
331:   DMGetCoordinatesLocal(da,&coords);
332:   DMDAVecGetArrayRead(cda,coords,&_coords);
333:   DMDAGetGhostCorners(cda,&si,&sj,0,&nx,&ny,0);
334:   for (j = sj; j < sj+ny-1; j++) {
335:     for (i = si; i < si+nx-1; i++) {
336:       PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e \n",(double)PetscRealPart(_coords[j][i].x),(double)PetscRealPart(_coords[j][i].y));
337:       PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e \n",(double)PetscRealPart(_coords[j+1][i].x),(double)PetscRealPart(_coords[j+1][i].y));
338:       PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e \n",(double)PetscRealPart(_coords[j+1][i+1].x),(double)PetscRealPart(_coords[j+1][i+1].y));
339:       PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e \n",(double)PetscRealPart(_coords[j][i+1].x),(double)PetscRealPart(_coords[j][i+1].y));
340:       PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e \n\n",(double)PetscRealPart(_coords[j][i].x),(double)PetscRealPart(_coords[j][i].y));
341:     }
342:   }
343:   DMDAVecRestoreArrayRead(cda,coords,&_coords);

345:   PetscFClose(PETSC_COMM_SELF,fp);
346:   return(0);
347: }

351: static PetscErrorCode DMDAViewGnuplot2d(DM da,Vec fields,const char comment[],const char prefix[])
352: {
353:   DM             cda;
354:   Vec            coords,local_fields;
355:   DMDACoor2d     **_coords;
356:   FILE           *fp;
357:   char           fname[PETSC_MAX_PATH_LEN];
358:   PetscMPIInt    rank;
359:   PetscInt       si,sj,nx,ny,i,j;
360:   PetscInt       n_dofs,d;
361:   PetscScalar    *_fields;

365:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
366:   PetscSNPrintf(fname,sizeof(fname),"%s-p%1.4d.dat",prefix,rank);
367:   PetscFOpen(PETSC_COMM_SELF,fname,"w",&fp);
368:   if (!fp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot open file");

370:   PetscFPrintf(PETSC_COMM_SELF,fp,"### %s (processor %1.4d) ### \n",comment,rank);
371:   DMDAGetInfo(da,0,0,0,0,0,0,0,&n_dofs,0,0,0,0,0);
372:   PetscFPrintf(PETSC_COMM_SELF,fp,"### x y ");
373:   for (d = 0; d < n_dofs; d++) {
374:     const char *field_name;
375:     DMDAGetFieldName(da,d,&field_name);
376:     PetscFPrintf(PETSC_COMM_SELF,fp,"%s ",field_name);
377:   }
378:   PetscFPrintf(PETSC_COMM_SELF,fp,"###\n");

380:   DMGetCoordinateDM(da,&cda);
381:   DMGetCoordinatesLocal(da,&coords);
382:   DMDAVecGetArray(cda,coords,&_coords);
383:   DMDAGetGhostCorners(cda,&si,&sj,0,&nx,&ny,0);

385:   DMCreateLocalVector(da,&local_fields);
386:   DMGlobalToLocalBegin(da,fields,INSERT_VALUES,local_fields);
387:   DMGlobalToLocalEnd(da,fields,INSERT_VALUES,local_fields);
388:   VecGetArray(local_fields,&_fields);

390:   for (j = sj; j < sj+ny; j++) {
391:     for (i = si; i < si+nx; i++) {
392:       PetscScalar coord_x,coord_y;
393:       PetscScalar field_d;

395:       coord_x = _coords[j][i].x;
396:       coord_y = _coords[j][i].y;

398:       PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e ",(double)PetscRealPart(coord_x),(double)PetscRealPart(coord_y));
399:       for (d = 0; d < n_dofs; d++) {
400:         field_d = _fields[n_dofs*((i-si)+(j-sj)*(nx))+d];
401:         PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e ",(double)PetscRealPart(field_d));
402:       }
403:       PetscFPrintf(PETSC_COMM_SELF,fp,"\n");
404:     }
405:   }
406:   VecRestoreArray(local_fields,&_fields);
407:   VecDestroy(&local_fields);

409:   DMDAVecRestoreArray(cda,coords,&_coords);

411:   PetscFClose(PETSC_COMM_SELF,fp);
412:   return(0);
413: }

417: static PetscErrorCode DMDAViewCoefficientsGnuplot2d(DM da,Vec fields,const char comment[],const char prefix[])
418: {
419:   DM                     cda;
420:   Vec                    local_fields;
421:   FILE                   *fp;
422:   char                   fname[PETSC_MAX_PATH_LEN];
423:   PetscMPIInt            rank;
424:   PetscInt               si,sj,nx,ny,i,j,p;
425:   PetscInt               n_dofs,d;
426:   GaussPointCoefficients **_coefficients;
427:   PetscErrorCode         ierr;

430:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
431:   PetscSNPrintf(fname,sizeof(fname),"%s-p%1.4d.dat",prefix,rank);
432:   PetscFOpen(PETSC_COMM_SELF,fname,"w",&fp);
433:   if (!fp) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot open file");

435:   PetscFPrintf(PETSC_COMM_SELF,fp,"### %s (processor %1.4d) ### \n",comment,rank);
436:   DMDAGetInfo(da,0,0,0,0,0,0,0,&n_dofs,0,0,0,0,0);
437:   PetscFPrintf(PETSC_COMM_SELF,fp,"### x y ");
438:   for (d = 0; d < n_dofs; d++) {
439:     const char *field_name;
440:     DMDAGetFieldName(da,d,&field_name);
441:     PetscFPrintf(PETSC_COMM_SELF,fp,"%s ",field_name);
442:   }
443:   PetscFPrintf(PETSC_COMM_SELF,fp,"###\n");

445:   DMGetCoordinateDM(da,&cda);
446:   DMDAGetGhostCorners(cda,&si,&sj,0,&nx,&ny,0);

448:   DMCreateLocalVector(da,&local_fields);
449:   DMGlobalToLocalBegin(da,fields,INSERT_VALUES,local_fields);
450:   DMGlobalToLocalEnd(da,fields,INSERT_VALUES,local_fields);
451:   DMDAVecGetArray(da,local_fields,&_coefficients);

453:   for (j = sj; j < sj+ny; j++) {
454:     for (i = si; i < si+nx; i++) {
455:       PetscScalar coord_x,coord_y;

457:       for (p = 0; p < GAUSS_POINTS; p++) {
458:         coord_x = _coefficients[j][i].gp_coords[2*p];
459:         coord_y = _coefficients[j][i].gp_coords[2*p+1];

461:         PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e ",(double)PetscRealPart(coord_x),(double)PetscRealPart(coord_y));

463:         PetscFPrintf(PETSC_COMM_SELF,fp,"%1.6e %1.6e %1.6e",(double)PetscRealPart(_coefficients[j][i].eta[p]),(double)PetscRealPart(_coefficients[j][i].fx[p]),(double)PetscRealPart(_coefficients[j][i].fy[p]));
464:         PetscFPrintf(PETSC_COMM_SELF,fp,"\n");
465:       }
466:     }
467:   }
468:   DMDAVecRestoreArray(da,local_fields,&_coefficients);
469:   VecDestroy(&local_fields);

471:   PetscFClose(PETSC_COMM_SELF,fp);
472:   return(0);
473: }

475: static PetscInt ASS_MAP_wIwDI_uJuDJ(PetscInt wi,PetscInt wd,PetscInt w_NPE,PetscInt w_dof,PetscInt ui,PetscInt ud,PetscInt u_NPE,PetscInt u_dof)
476: {
477:   PetscInt ij;
478:   PetscInt r,c,nc;

480:   nc = u_NPE*u_dof;
481:   r = w_dof*wi+wd;
482:   c = u_dof*ui+ud;
483:   ij = r*nc+c;
484:   return ij;
485: }

487: /*
488:  D = [ 2.eta   0   0   ]
489:      [   0   2.eta 0   ]
490:      [   0     0   eta ]
491:  
492:  B = [ d_dx   0   ]
493:      [  0    d_dy ]
494:      [ d_dy  d_dx ]
495: */
496: static void FormStressOperatorQ1(PetscScalar Ke[],PetscScalar coords[],PetscScalar eta[])
497: {
498:   PetscInt       ngp;
499:   PetscScalar    gp_xi[GAUSS_POINTS][2];
500:   PetscScalar    gp_weight[GAUSS_POINTS];
501:   PetscInt       p,i,j,k;
502:   PetscScalar    GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
503:   PetscScalar    J_p,tildeD[3];
504:   PetscScalar    B[3][U_DOFS*NODES_PER_EL];

506:   /* define quadrature rule */
507:   ConstructGaussQuadrature(&ngp,gp_xi,gp_weight);

509:   /* evaluate integral */
510:   for (p = 0; p < ngp; p++) {
511:     ConstructQ12D_GNi(gp_xi[p],GNi_p);
512:     ConstructQ12D_GNx(GNi_p,GNx_p,coords,&J_p);

514:     for (i = 0; i < NODES_PER_EL; i++) {
515:       PetscScalar d_dx_i = GNx_p[0][i];
516:       PetscScalar d_dy_i = GNx_p[1][i];

518:       B[0][2*i] = d_dx_i;B[0][2*i+1] = 0.0;
519:       B[1][2*i] = 0.0;B[1][2*i+1] = d_dy_i;
520:       B[2][2*i] = d_dy_i;B[2][2*i+1] = d_dx_i;
521:     }

523:     tildeD[0] = 2.0*gp_weight[p]*J_p*eta[p];
524:     tildeD[1] = 2.0*gp_weight[p]*J_p*eta[p];
525:     tildeD[2] =       gp_weight[p]*J_p*eta[p];

527:     /* form Bt tildeD B */
528:     /*
529:     Ke_ij = Bt_ik . D_kl . B_lj
530:           = B_ki . D_kl . B_lj
531:           = B_ki . D_kk . B_kj
532:     */
533:     for (i = 0; i < 8; i++) {
534:       for (j = 0; j < 8; j++) {
535:         for (k = 0; k < 3; k++) { /* Note D is diagonal for stokes */
536:           Ke[i+8*j] = Ke[i+8*j]+B[k][i]*tildeD[k]*B[k][j];
537:         }
538:       }
539:     }
540:   }
541: }

543: static void FormGradientOperatorQ1(PetscScalar Ke[],PetscScalar coords[])
544: {
545:   PetscInt    ngp;
546:   PetscScalar gp_xi[GAUSS_POINTS][2];
547:   PetscScalar gp_weight[GAUSS_POINTS];
548:   PetscInt    p,i,j,di;
549:   PetscScalar Ni_p[NODES_PER_EL];
550:   PetscScalar GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
551:   PetscScalar J_p,fac;

553:   /* define quadrature rule */
554:   ConstructGaussQuadrature(&ngp,gp_xi,gp_weight);

556:   /* evaluate integral */
557:   for (p = 0; p < ngp; p++) {
558:     ConstructQ12D_Ni(gp_xi[p],Ni_p);
559:     ConstructQ12D_GNi(gp_xi[p],GNi_p);
560:     ConstructQ12D_GNx(GNi_p,GNx_p,coords,&J_p);
561:     fac = gp_weight[p]*J_p;

563:     for (i = 0; i < NODES_PER_EL; i++) { /* u nodes */
564:       for (di = 0; di < NSD; di++) { /* u dofs */
565:         for (j = 0; j < 4; j++) {  /* p nodes, p dofs = 1 (ie no loop) */
566:           PetscInt IJ;
567:           IJ = ASS_MAP_wIwDI_uJuDJ(i,di,NODES_PER_EL,2,j,0,NODES_PER_EL,1);

569:           Ke[IJ] = Ke[IJ]-GNx_p[di][i]*Ni_p[j]*fac;
570:         }
571:       }
572:     }
573:   }
574: }

576: static void FormDivergenceOperatorQ1(PetscScalar De[],PetscScalar coords[])
577: {
578:   PetscScalar Ge[U_DOFS*NODES_PER_EL*P_DOFS*NODES_PER_EL];
579:   PetscInt    i,j;
580:   PetscInt    nr_g,nc_g;

582:   PetscMemzero(Ge,sizeof(PetscScalar)*U_DOFS*NODES_PER_EL*P_DOFS*NODES_PER_EL);
583:   FormGradientOperatorQ1(Ge,coords);

585:   nr_g = U_DOFS*NODES_PER_EL;
586:   nc_g = P_DOFS*NODES_PER_EL;

588:   for (i = 0; i < nr_g; i++) {
589:     for (j = 0; j < nc_g; j++) {
590:       De[nr_g*j+i] = Ge[nc_g*i+j];
591:     }
592:   }
593: }

595: static void FormStabilisationOperatorQ1(PetscScalar Ke[],PetscScalar coords[],PetscScalar eta[])
596: {
597:   PetscInt    ngp;
598:   PetscScalar gp_xi[GAUSS_POINTS][2];
599:   PetscScalar gp_weight[GAUSS_POINTS];
600:   PetscInt    p,i,j;
601:   PetscScalar Ni_p[NODES_PER_EL];
602:   PetscScalar GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
603:   PetscScalar J_p,fac,eta_avg;

605:   /* define quadrature rule */
606:   ConstructGaussQuadrature(&ngp,gp_xi,gp_weight);

608:   /* evaluate integral */
609:   for (p = 0; p < ngp; p++) {
610:     ConstructQ12D_Ni(gp_xi[p],Ni_p);
611:     ConstructQ12D_GNi(gp_xi[p],GNi_p);
612:     ConstructQ12D_GNx(GNi_p,GNx_p,coords,&J_p);
613:     fac = gp_weight[p]*J_p;

615:     for (i = 0; i < NODES_PER_EL; i++) {
616:       for (j = 0; j < NODES_PER_EL; j++) {
617:         Ke[NODES_PER_EL*i+j] = Ke[NODES_PER_EL*i+j]-fac*(Ni_p[i]*Ni_p[j]-0.0625);
618:       }
619:     }
620:   }

622:   /* scale */
623:   eta_avg = 0.0;
624:   for (p = 0; p < ngp; p++) eta_avg += eta[p];
625:   eta_avg = (1.0/((PetscScalar)ngp))*eta_avg;
626:   fac     = 1.0/eta_avg;
627:   for (i = 0; i < NODES_PER_EL; i++) {
628:     for (j = 0; j < NODES_PER_EL; j++) {
629:       Ke[NODES_PER_EL*i+j] = fac*Ke[NODES_PER_EL*i+j];
630:     }
631:   }
632: }

634: static void FormScaledMassMatrixOperatorQ1(PetscScalar Ke[],PetscScalar coords[],PetscScalar eta[])
635: {
636:   PetscInt    ngp;
637:   PetscScalar gp_xi[GAUSS_POINTS][2];
638:   PetscScalar gp_weight[GAUSS_POINTS];
639:   PetscInt    p,i,j;
640:   PetscScalar Ni_p[NODES_PER_EL];
641:   PetscScalar GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
642:   PetscScalar J_p,fac,eta_avg;

644:   /* define quadrature rule */
645:   ConstructGaussQuadrature(&ngp,gp_xi,gp_weight);

647:   /* evaluate integral */
648:   for (p = 0; p < ngp; p++) {
649:     ConstructQ12D_Ni(gp_xi[p],Ni_p);
650:     ConstructQ12D_GNi(gp_xi[p],GNi_p);
651:     ConstructQ12D_GNx(GNi_p,GNx_p,coords,&J_p);
652:     fac = gp_weight[p]*J_p;

654:     for (i = 0; i < NODES_PER_EL; i++) {
655:       for (j = 0; j < NODES_PER_EL; j++) {
656:         Ke[NODES_PER_EL*i+j] = Ke[NODES_PER_EL*i+j]-fac*Ni_p[i]*Ni_p[j];
657:       }
658:     }
659:   }

661:   /* scale */
662:   eta_avg = 0.0;
663:   for (p = 0; p < ngp; p++) eta_avg += eta[p];
664:   eta_avg = (1.0/((PetscScalar)ngp))*eta_avg;
665:   fac     = 1.0/eta_avg;
666:   for (i = 0; i < NODES_PER_EL; i++) {
667:     for (j = 0; j < NODES_PER_EL; j++) {
668:       Ke[NODES_PER_EL*i+j] = fac*Ke[NODES_PER_EL*i+j];
669:     }
670:   }
671: }

673: static void FormMomentumRhsQ1(PetscScalar Fe[],PetscScalar coords[],PetscScalar fx[],PetscScalar fy[])
674: {
675:   PetscInt    ngp;
676:   PetscScalar gp_xi[GAUSS_POINTS][2];
677:   PetscScalar gp_weight[GAUSS_POINTS];
678:   PetscInt    p,i;
679:   PetscScalar Ni_p[NODES_PER_EL];
680:   PetscScalar GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
681:   PetscScalar J_p,fac;

683:   /* define quadrature rule */
684:   ConstructGaussQuadrature(&ngp,gp_xi,gp_weight);

686:   /* evaluate integral */
687:   for (p = 0; p < ngp; p++) {
688:     ConstructQ12D_Ni(gp_xi[p],Ni_p);
689:     ConstructQ12D_GNi(gp_xi[p],GNi_p);
690:     ConstructQ12D_GNx(GNi_p,GNx_p,coords,&J_p);
691:     fac = gp_weight[p]*J_p;

693:     for (i = 0; i < NODES_PER_EL; i++) {
694:       Fe[NSD*i]   += fac*Ni_p[i]*fx[p];
695:       Fe[NSD*i+1] += fac*Ni_p[i]*fy[p];
696:     }
697:   }
698: }

702: static PetscErrorCode GetElementCoords(DMDACoor2d **_coords,PetscInt ei,PetscInt ej,PetscScalar el_coords[])
703: {
705:   /* get coords for the element */
706:   el_coords[NSD*0] = _coords[ej][ei].x;     el_coords[NSD*0+1] = _coords[ej][ei].y;
707:   el_coords[NSD*1] = _coords[ej+1][ei].x;   el_coords[NSD*1+1] = _coords[ej+1][ei].y;
708:   el_coords[NSD*2] = _coords[ej+1][ei+1].x; el_coords[NSD*2+1] = _coords[ej+1][ei+1].y;
709:   el_coords[NSD*3] = _coords[ej][ei+1].x;   el_coords[NSD*3+1] = _coords[ej][ei+1].y;
710:   return(0);
711: }

715: static PetscErrorCode AssembleA_Stokes(Mat A,DM stokes_da,DM properties_da,Vec properties)
716: {
717:   DM                     cda;
718:   Vec                    coords;
719:   DMDACoor2d             **_coords;
720:   MatStencil             u_eqn[NODES_PER_EL*U_DOFS]; /* 2 degrees of freedom */
721:   MatStencil             p_eqn[NODES_PER_EL*P_DOFS]; /* 1 degrees of freedom */
722:   PetscInt               sex,sey,mx,my;
723:   PetscInt               ei,ej;
724:   PetscScalar            Ae[NODES_PER_EL*U_DOFS*NODES_PER_EL*U_DOFS];
725:   PetscScalar            Ge[NODES_PER_EL*U_DOFS*NODES_PER_EL*P_DOFS];
726:   PetscScalar            De[NODES_PER_EL*P_DOFS*NODES_PER_EL*U_DOFS];
727:   PetscScalar            Ce[NODES_PER_EL*P_DOFS*NODES_PER_EL*P_DOFS];
728:   PetscScalar            el_coords[NODES_PER_EL*NSD];
729:   Vec                    local_properties;
730:   GaussPointCoefficients **props;
731:   PetscScalar            *prop_eta;
732:   PetscErrorCode         ierr;

735:   /* setup for coords */
736:   DMGetCoordinateDM(stokes_da,&cda);
737:   DMGetCoordinatesLocal(stokes_da,&coords);
738:   DMDAVecGetArray(cda,coords,&_coords);

740:   /* setup for coefficients */
741:   DMCreateLocalVector(properties_da,&local_properties);
742:   DMGlobalToLocalBegin(properties_da,properties,INSERT_VALUES,local_properties);
743:   DMGlobalToLocalEnd(properties_da,properties,INSERT_VALUES,local_properties);
744:   DMDAVecGetArray(properties_da,local_properties,&props);

746:   DMDAGetElementCorners(stokes_da,&sex,&sey,NULL,&mx,&my,NULL);
747:   for (ej = sey; ej < sey+my; ej++) {
748:     for (ei = sex; ei < sex+mx; ei++) {
749:       /* get coords for the element */
750:       GetElementCoords(_coords,ei,ej,el_coords);

752:       /* get coefficients for the element */
753:       prop_eta = props[ej][ei].eta;

755:       /* initialise element stiffness matrix */
756:       PetscMemzero(Ae,sizeof(PetscScalar)*NODES_PER_EL*U_DOFS*NODES_PER_EL*U_DOFS);
757:       PetscMemzero(Ge,sizeof(PetscScalar)*NODES_PER_EL*U_DOFS*NODES_PER_EL*P_DOFS);
758:       PetscMemzero(De,sizeof(PetscScalar)*NODES_PER_EL*P_DOFS*NODES_PER_EL*U_DOFS);
759:       PetscMemzero(Ce,sizeof(PetscScalar)*NODES_PER_EL*P_DOFS*NODES_PER_EL*P_DOFS);

761:       /* form element stiffness matrix */
762:       FormStressOperatorQ1(Ae,el_coords,prop_eta);
763:       FormGradientOperatorQ1(Ge,el_coords);
764:       FormDivergenceOperatorQ1(De,el_coords);
765:       FormStabilisationOperatorQ1(Ce,el_coords,prop_eta);

767:       /* insert element matrix into global matrix */
768:       DMDAGetElementEqnums_up(u_eqn,p_eqn,ei,ej);
769:       MatSetValuesStencil(A,NODES_PER_EL*U_DOFS,u_eqn,NODES_PER_EL*U_DOFS,u_eqn,Ae,ADD_VALUES);
770:       MatSetValuesStencil(A,NODES_PER_EL*U_DOFS,u_eqn,NODES_PER_EL*P_DOFS,p_eqn,Ge,ADD_VALUES);
771:       MatSetValuesStencil(A,NODES_PER_EL*P_DOFS,p_eqn,NODES_PER_EL*U_DOFS,u_eqn,De,ADD_VALUES);
772:       MatSetValuesStencil(A,NODES_PER_EL*P_DOFS,p_eqn,NODES_PER_EL*P_DOFS,p_eqn,Ce,ADD_VALUES);
773:     }
774:   }
775:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
776:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

778:   DMDAVecRestoreArray(cda,coords,&_coords);

780:   DMDAVecRestoreArray(properties_da,local_properties,&props);
781:   VecDestroy(&local_properties);
782:   return(0);
783: }

787: static PetscErrorCode AssembleA_PCStokes(Mat A,DM stokes_da,DM properties_da,Vec properties)
788: {
789:   DM                     cda;
790:   Vec                    coords;
791:   DMDACoor2d             **_coords;
792:   MatStencil             u_eqn[NODES_PER_EL*U_DOFS]; /* 2 degrees of freedom */
793:   MatStencil             p_eqn[NODES_PER_EL*P_DOFS]; /* 1 degrees of freedom */
794:   PetscInt               sex,sey,mx,my;
795:   PetscInt               ei,ej;
796:   PetscScalar            Ae[NODES_PER_EL*U_DOFS*NODES_PER_EL*U_DOFS];
797:   PetscScalar            Ge[NODES_PER_EL*U_DOFS*NODES_PER_EL*P_DOFS];
798:   PetscScalar            De[NODES_PER_EL*P_DOFS*NODES_PER_EL*U_DOFS];
799:   PetscScalar            Ce[NODES_PER_EL*P_DOFS*NODES_PER_EL*P_DOFS];
800:   PetscScalar            el_coords[NODES_PER_EL*NSD];
801:   Vec                    local_properties;
802:   GaussPointCoefficients **props;
803:   PetscScalar            *prop_eta;
804:   PetscErrorCode         ierr;

807:   /* setup for coords */
808:   DMGetCoordinateDM(stokes_da,&cda);
809:   DMGetCoordinatesLocal(stokes_da,&coords);
810:   DMDAVecGetArray(cda,coords,&_coords);

812:   /* setup for coefficients */
813:   DMCreateLocalVector(properties_da,&local_properties);
814:   DMGlobalToLocalBegin(properties_da,properties,INSERT_VALUES,local_properties);
815:   DMGlobalToLocalEnd(properties_da,properties,INSERT_VALUES,local_properties);
816:   DMDAVecGetArray(properties_da,local_properties,&props);

818:   DMDAGetElementCorners(stokes_da,&sex,&sey,0,&mx,&my,0);
819:   for (ej = sey; ej < sey+my; ej++) {
820:     for (ei = sex; ei < sex+mx; ei++) {
821:       /* get coords for the element */
822:       GetElementCoords(_coords,ei,ej,el_coords);

824:       /* get coefficients for the element */
825:       prop_eta = props[ej][ei].eta;

827:       /* initialise element stiffness matrix */
828:       PetscMemzero(Ae,sizeof(PetscScalar)*NODES_PER_EL*U_DOFS*NODES_PER_EL*U_DOFS);
829:       PetscMemzero(Ge,sizeof(PetscScalar)*NODES_PER_EL*U_DOFS*NODES_PER_EL*P_DOFS);
830:       PetscMemzero(De,sizeof(PetscScalar)*NODES_PER_EL*P_DOFS*NODES_PER_EL*U_DOFS);
831:       PetscMemzero(Ce,sizeof(PetscScalar)*NODES_PER_EL*P_DOFS*NODES_PER_EL*P_DOFS);

833:       /* form element stiffness matrix */
834:       FormStressOperatorQ1(Ae,el_coords,prop_eta);
835:       FormGradientOperatorQ1(Ge,el_coords);
836:       FormScaledMassMatrixOperatorQ1(Ce,el_coords,prop_eta);

838:       /* insert element matrix into global matrix */
839:       DMDAGetElementEqnums_up(u_eqn,p_eqn,ei,ej);
840:       MatSetValuesStencil(A,NODES_PER_EL*U_DOFS,u_eqn,NODES_PER_EL*U_DOFS,u_eqn,Ae,ADD_VALUES);
841:       MatSetValuesStencil(A,NODES_PER_EL*U_DOFS,u_eqn,NODES_PER_EL*P_DOFS,p_eqn,Ge,ADD_VALUES);
842:       MatSetValuesStencil(A,NODES_PER_EL*P_DOFS,p_eqn,NODES_PER_EL*P_DOFS,p_eqn,Ce,ADD_VALUES);
843:     }
844:   }
845:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
846:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

848:   DMDAVecRestoreArray(cda,coords,&_coords);

850:   DMDAVecRestoreArray(properties_da,local_properties,&props);
851:   VecDestroy(&local_properties);
852:   return(0);
853: }

857: static PetscErrorCode DMDASetValuesLocalStencil_ADD_VALUES(StokesDOF **fields_F,MatStencil u_eqn[],MatStencil p_eqn[],PetscScalar Fe_u[],PetscScalar Fe_p[])
858: {
859:   PetscInt n;

862:   for (n = 0; n < 4; n++) {
863:     fields_F[u_eqn[2*n].j][u_eqn[2*n].i].u_dof     = fields_F[u_eqn[2*n].j][u_eqn[2*n].i].u_dof+Fe_u[2*n];
864:     fields_F[u_eqn[2*n+1].j][u_eqn[2*n+1].i].v_dof = fields_F[u_eqn[2*n+1].j][u_eqn[2*n+1].i].v_dof+Fe_u[2*n+1];
865:     fields_F[p_eqn[n].j][p_eqn[n].i].p_dof         = fields_F[p_eqn[n].j][p_eqn[n].i].p_dof+Fe_p[n];
866:   }
867:   return(0);
868: }

872: static PetscErrorCode AssembleF_Stokes(Vec F,DM stokes_da,DM properties_da,Vec properties)
873: {
874:   DM                     cda;
875:   Vec                    coords;
876:   DMDACoor2d             **_coords;
877:   MatStencil             u_eqn[NODES_PER_EL*U_DOFS]; /* 2 degrees of freedom */
878:   MatStencil             p_eqn[NODES_PER_EL*P_DOFS]; /* 1 degrees of freedom */
879:   PetscInt               sex,sey,mx,my;
880:   PetscInt               ei,ej;
881:   PetscScalar            Fe[NODES_PER_EL*U_DOFS];
882:   PetscScalar            He[NODES_PER_EL*P_DOFS];
883:   PetscScalar            el_coords[NODES_PER_EL*NSD];
884:   Vec                    local_properties;
885:   GaussPointCoefficients **props;
886:   PetscScalar            *prop_fx,*prop_fy;
887:   Vec                    local_F;
888:   StokesDOF              **ff;
889:   PetscErrorCode         ierr;

892:   /* setup for coords */
893:   DMGetCoordinateDM(stokes_da,&cda);
894:   DMGetCoordinatesLocal(stokes_da,&coords);
895:   DMDAVecGetArray(cda,coords,&_coords);

897:   /* setup for coefficients */
898:   DMGetLocalVector(properties_da,&local_properties);
899:   DMGlobalToLocalBegin(properties_da,properties,INSERT_VALUES,local_properties);
900:   DMGlobalToLocalEnd(properties_da,properties,INSERT_VALUES,local_properties);
901:   DMDAVecGetArray(properties_da,local_properties,&props);

903:   /* get acces to the vector */
904:   DMGetLocalVector(stokes_da,&local_F);
905:   VecZeroEntries(local_F);
906:   DMDAVecGetArray(stokes_da,local_F,&ff);

908:   DMDAGetElementCorners(stokes_da,&sex,&sey,NULL,&mx,&my,NULL);
909:   for (ej = sey; ej < sey+my; ej++) {
910:     for (ei = sex; ei < sex+mx; ei++) {
911:       /* get coords for the element */
912:       GetElementCoords(_coords,ei,ej,el_coords);

914:       /* get coefficients for the element */
915:       prop_fx = props[ej][ei].fx;
916:       prop_fy = props[ej][ei].fy;

918:       /* initialise element stiffness matrix */
919:       PetscMemzero(Fe,sizeof(PetscScalar)*NODES_PER_EL*U_DOFS);
920:       PetscMemzero(He,sizeof(PetscScalar)*NODES_PER_EL*P_DOFS);

922:       /* form element stiffness matrix */
923:       FormMomentumRhsQ1(Fe,el_coords,prop_fx,prop_fy);

925:       /* insert element matrix into global matrix */
926:       DMDAGetElementEqnums_up(u_eqn,p_eqn,ei,ej);

928:       DMDASetValuesLocalStencil_ADD_VALUES(ff,u_eqn,p_eqn,Fe,He);
929:     }
930:   }

932:   DMDAVecRestoreArray(stokes_da,local_F,&ff);
933:   DMLocalToGlobalBegin(stokes_da,local_F,ADD_VALUES,F);
934:   DMLocalToGlobalEnd(stokes_da,local_F,ADD_VALUES,F);
935:   DMRestoreLocalVector(stokes_da,&local_F);

937:   DMDAVecRestoreArray(cda,coords,&_coords);

939:   DMDAVecRestoreArray(properties_da,local_properties,&props);
940:   DMRestoreLocalVector(properties_da,&local_properties);
941:   return(0);
942: }

946: static PetscErrorCode DMDACreateSolCx(PetscReal eta0,PetscReal eta1,PetscReal xc,PetscInt nz,PetscInt mx,PetscInt my,DM *_da,Vec *_X)
947: {
948:   DM             da,cda;
949:   Vec            X;
950:   StokesDOF      **_stokes;
951:   Vec            coords;
952:   DMDACoor2d     **_coords;
953:   PetscInt       si,sj,ei,ej,i,j;

957:   DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,mx+1,my+1,PETSC_DECIDE,PETSC_DECIDE,3,1,NULL,NULL,&da);
958:   DMDASetFieldName(da,0,"anlytic_Vx");
959:   DMDASetFieldName(da,1,"anlytic_Vy");
960:   DMDASetFieldName(da,2,"analytic_P");

962:   DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.,0.);

964:   DMGetCoordinatesLocal(da,&coords);
965:   DMGetCoordinateDM(da,&cda);
966:   DMDAVecGetArray(cda,coords,&_coords);

968:   DMCreateGlobalVector(da,&X);
969:   DMDAVecGetArray(da,X,&_stokes);

971:   DMDAGetCorners(da,&si,&sj,0,&ei,&ej,0);
972:   for (j = sj; j < sj+ej; j++) {
973:     for (i = si; i < si+ei; i++) {
974:       PetscReal pos[2],pressure,vel[2],total_stress[3],strain_rate[3];

976:       pos[0] = PetscRealPart(_coords[j][i].x);
977:       pos[1] = PetscRealPart(_coords[j][i].y);

979:       evaluate_solCx(pos,eta0,eta1,xc,nz,vel,&pressure,total_stress,strain_rate);

981:       _stokes[j][i].u_dof = vel[0];
982:       _stokes[j][i].v_dof = vel[1];
983:       _stokes[j][i].p_dof = pressure;
984:     }
985:   }
986:   DMDAVecRestoreArray(da,X,&_stokes);
987:   DMDAVecRestoreArray(cda,coords,&_coords);

989:   *_da = da;
990:   *_X  = X;
991:   return(0);
992: }

996: static PetscErrorCode StokesDAGetNodalFields(StokesDOF **fields,PetscInt ei,PetscInt ej,StokesDOF nodal_fields[])
997: {
999:   /* get the nodal fields */
1000:   nodal_fields[0].u_dof = fields[ej][ei].u_dof;     nodal_fields[0].v_dof = fields[ej][ei].v_dof;     nodal_fields[0].p_dof = fields[ej][ei].p_dof;
1001:   nodal_fields[1].u_dof = fields[ej+1][ei].u_dof;   nodal_fields[1].v_dof = fields[ej+1][ei].v_dof;   nodal_fields[1].p_dof = fields[ej+1][ei].p_dof;
1002:   nodal_fields[2].u_dof = fields[ej+1][ei+1].u_dof; nodal_fields[2].v_dof = fields[ej+1][ei+1].v_dof; nodal_fields[2].p_dof = fields[ej+1][ei+1].p_dof;
1003:   nodal_fields[3].u_dof = fields[ej][ei+1].u_dof;   nodal_fields[3].v_dof = fields[ej][ei+1].v_dof;   nodal_fields[3].p_dof = fields[ej][ei+1].p_dof;
1004:   return(0);
1005: }

1009: static PetscErrorCode DMDAIntegrateErrors(DM stokes_da,Vec X,Vec X_analytic)
1010: {
1011:   DM             cda;
1012:   Vec            coords,X_analytic_local,X_local;
1013:   DMDACoor2d     **_coords;
1014:   PetscInt       sex,sey,mx,my;
1015:   PetscInt       ei,ej;
1016:   PetscScalar    el_coords[NODES_PER_EL*NSD];
1017:   StokesDOF      **stokes_analytic,**stokes;
1018:   StokesDOF      stokes_analytic_e[4],stokes_e[4];

1020:   PetscScalar    GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
1021:   PetscScalar    Ni_p[NODES_PER_EL];
1022:   PetscInt       ngp;
1023:   PetscScalar    gp_xi[GAUSS_POINTS][2];
1024:   PetscScalar    gp_weight[GAUSS_POINTS];
1025:   PetscInt       p,i;
1026:   PetscScalar    J_p,fac;
1027:   PetscScalar    h,p_e_L2,u_e_L2,u_e_H1,p_L2,u_L2,u_H1,tp_L2,tu_L2,tu_H1;
1028:   PetscInt       M;
1029:   PetscReal      xymin[2],xymax[2];

1033:   /* define quadrature rule */
1034:   ConstructGaussQuadrature(&ngp,gp_xi,gp_weight);

1036:   /* setup for coords */
1037:   DMGetCoordinateDM(stokes_da,&cda);
1038:   DMGetCoordinatesLocal(stokes_da,&coords);
1039:   DMDAVecGetArray(cda,coords,&_coords);

1041:   /* setup for analytic */
1042:   DMCreateLocalVector(stokes_da,&X_analytic_local);
1043:   DMGlobalToLocalBegin(stokes_da,X_analytic,INSERT_VALUES,X_analytic_local);
1044:   DMGlobalToLocalEnd(stokes_da,X_analytic,INSERT_VALUES,X_analytic_local);
1045:   DMDAVecGetArray(stokes_da,X_analytic_local,&stokes_analytic);

1047:   /* setup for solution */
1048:   DMCreateLocalVector(stokes_da,&X_local);
1049:   DMGlobalToLocalBegin(stokes_da,X,INSERT_VALUES,X_local);
1050:   DMGlobalToLocalEnd(stokes_da,X,INSERT_VALUES,X_local);
1051:   DMDAVecGetArray(stokes_da,X_local,&stokes);

1053:   DMDAGetInfo(stokes_da,0,&M,0,0,0,0,0,0,0,0,0,0,0);
1054:   DMDAGetBoundingBox(stokes_da,xymin,xymax);

1056:   h = (xymax[0]-xymin[0])/((PetscReal)M);

1058:   tp_L2 = tu_L2 = tu_H1 = 0.0;

1060:   DMDAGetElementCorners(stokes_da,&sex,&sey,NULL,&mx,&my,NULL);
1061:   for (ej = sey; ej < sey+my; ej++) {
1062:     for (ei = sex; ei < sex+mx; ei++) {
1063:       /* get coords for the element */
1064:       GetElementCoords(_coords,ei,ej,el_coords);
1065:       StokesDAGetNodalFields(stokes,ei,ej,stokes_e);
1066:       StokesDAGetNodalFields(stokes_analytic,ei,ej,stokes_analytic_e);

1068:       /* evaluate integral */
1069:       p_e_L2 = 0.0;
1070:       u_e_L2 = 0.0;
1071:       u_e_H1 = 0.0;
1072:       for (p = 0; p < ngp; p++) {
1073:         ConstructQ12D_Ni(gp_xi[p],Ni_p);
1074:         ConstructQ12D_GNi(gp_xi[p],GNi_p);
1075:         ConstructQ12D_GNx(GNi_p,GNx_p,el_coords,&J_p);
1076:         fac = gp_weight[p]*J_p;

1078:         for (i = 0; i < NODES_PER_EL; i++) {
1079:           PetscScalar u_error,v_error;

1081:           p_e_L2 = p_e_L2+fac*Ni_p[i]*(stokes_e[i].p_dof-stokes_analytic_e[i].p_dof)*(stokes_e[i].p_dof-stokes_analytic_e[i].p_dof);

1083:           u_error = stokes_e[i].u_dof-stokes_analytic_e[i].u_dof;
1084:           v_error = stokes_e[i].v_dof-stokes_analytic_e[i].v_dof;
1085:           u_e_L2 += fac*Ni_p[i]*(u_error*u_error+v_error*v_error);

1087:           u_e_H1 = u_e_H1+fac*(GNx_p[0][i]*u_error*GNx_p[0][i]*u_error              /* du/dx */
1088:                                +GNx_p[1][i]*u_error*GNx_p[1][i]*u_error               /* du/dy */
1089:                                +GNx_p[0][i]*v_error*GNx_p[0][i]*v_error               /* dv/dx */
1090:                                +GNx_p[1][i]*v_error*GNx_p[1][i]*v_error);             /* dv/dy */
1091:         }
1092:       }

1094:       tp_L2 += p_e_L2;
1095:       tu_L2 += u_e_L2;
1096:       tu_H1 += u_e_H1;
1097:     }
1098:   }
1099:   MPI_Allreduce(&tp_L2,&p_L2,1,MPIU_SCALAR,MPIU_SUM,PETSC_COMM_WORLD);
1100:   MPI_Allreduce(&tu_L2,&u_L2,1,MPIU_SCALAR,MPIU_SUM,PETSC_COMM_WORLD);
1101:   MPI_Allreduce(&tu_H1,&u_H1,1,MPIU_SCALAR,MPIU_SUM,PETSC_COMM_WORLD);
1102:   p_L2 = PetscSqrtScalar(p_L2);
1103:   u_L2 = PetscSqrtScalar(u_L2);
1104:   u_H1 = PetscSqrtScalar(u_H1);

1106:   PetscPrintf(PETSC_COMM_WORLD,"%1.4e   %1.4e   %1.4e   %1.4e \n",(double)PetscRealPart(h),(double)PetscRealPart(p_L2),(double)PetscRealPart(u_L2),(double)PetscRealPart(u_H1));

1108:   DMDAVecRestoreArray(cda,coords,&_coords);

1110:   DMDAVecRestoreArray(stokes_da,X_analytic_local,&stokes_analytic);
1111:   VecDestroy(&X_analytic_local);
1112:   DMDAVecRestoreArray(stokes_da,X_local,&stokes);
1113:   VecDestroy(&X_local);
1114:   return(0);
1115: }

1119: static PetscErrorCode solve_stokes_2d_coupled(PetscInt mx,PetscInt my)
1120: {
1121:   DM                     da_Stokes,da_prop;
1122:   PetscInt               u_dof,p_dof,dof,stencil_width;
1123:   Mat                    A,B;
1124:   PetscInt               mxl,myl;
1125:   DM                     prop_cda,vel_cda;
1126:   Vec                    prop_coords,vel_coords;
1127:   PetscInt               si,sj,nx,ny,i,j,p;
1128:   Vec                    f,X;
1129:   PetscInt               prop_dof,prop_stencil_width;
1130:   Vec                    properties,l_properties;
1131:   PetscReal              dx,dy;
1132:   PetscInt               M,N;
1133:   DMDACoor2d             **_prop_coords,**_vel_coords;
1134:   GaussPointCoefficients **element_props;
1135:   PetscInt               its;
1136:   KSP                    ksp_S;
1137:   PetscInt               coefficient_structure = 0;
1138:   PetscInt               cpu_x,cpu_y,*lx = NULL,*ly = NULL;
1139:   PetscBool              use_gp_coords = PETSC_FALSE,set,output_gnuplot = PETSC_FALSE;
1140:   char                   filename[PETSC_MAX_PATH_LEN];
1141:   PetscErrorCode         ierr;


1145:   PetscOptionsGetBool(NULL,NULL,"-gnuplot",&output_gnuplot,NULL);
1146: 
1147:   /* Generate the da for velocity and pressure */
1148:   /*
1149:   We use Q1 elements for the temperature.
1150:   FEM has a 9-point stencil (BOX) or connectivity pattern
1151:   Num nodes in each direction is mx+1, my+1
1152:   */
1153:   u_dof         = U_DOFS; /* Vx, Vy - velocities */
1154:   p_dof         = P_DOFS; /* p - pressure */
1155:   dof           = u_dof+p_dof;
1156:   stencil_width = 1;
1157:   DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,mx+1,my+1,PETSC_DECIDE,PETSC_DECIDE,dof,stencil_width,NULL,NULL,&da_Stokes);
1158:   DMDASetFieldName(da_Stokes,0,"Vx");
1159:   DMDASetFieldName(da_Stokes,1,"Vy");
1160:   DMDASetFieldName(da_Stokes,2,"P");

1162:   /* unit box [0,1] x [0,1] */
1163:   DMDASetUniformCoordinates(da_Stokes,0.0,1.0,0.0,1.0,0.,0.);

1165:   /* Generate element properties, we will assume all material properties are constant over the element */
1166:   /* local number of elements */
1167:   DMDAGetLocalElementSize(da_Stokes,&mxl,&myl,NULL);

1169:   /* !!! IN PARALLEL WE MUST MAKE SURE THE TWO DMDA's ALIGN !!!  */
1170:   DMDAGetInfo(da_Stokes,0,0,0,0,&cpu_x,&cpu_y,0,0,0,0,0,0,0);
1171:   DMDAGetElementOwnershipRanges2d(da_Stokes,&lx,&ly);

1173:   prop_dof           = (int)(sizeof(GaussPointCoefficients)/sizeof(PetscScalar)); /* gauss point setup */
1174:   prop_stencil_width = 0;
1175:   DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,mx,my,cpu_x,cpu_y,prop_dof,prop_stencil_width,lx,ly,&da_prop);
1176:   PetscFree(lx);
1177:   PetscFree(ly);

1179:   /* define centroid positions */
1180:   DMDAGetInfo(da_prop,0,&M,&N,0,0,0,0,0,0,0,0,0,0);
1181:   dx   = 1.0/((PetscReal)(M));
1182:   dy   = 1.0/((PetscReal)(N));

1184:   DMDASetUniformCoordinates(da_prop,0.0+0.5*dx,1.0-0.5*dx,0.0+0.5*dy,1.0-0.5*dy,0.,0);

1186:   /* define coefficients */
1187:   PetscOptionsGetInt(NULL,NULL,"-c_str",&coefficient_structure,NULL);

1189:   DMCreateGlobalVector(da_prop,&properties);
1190:   DMCreateLocalVector(da_prop,&l_properties);
1191:   DMDAVecGetArray(da_prop,l_properties,&element_props);

1193:   DMGetCoordinateDM(da_prop,&prop_cda);
1194:   DMGetCoordinatesLocal(da_prop,&prop_coords);
1195:   DMDAVecGetArray(prop_cda,prop_coords,&_prop_coords);

1197:   DMDAGetGhostCorners(prop_cda,&si,&sj,0,&nx,&ny,0);

1199:   DMGetCoordinateDM(da_Stokes,&vel_cda);
1200:   DMGetCoordinatesLocal(da_Stokes,&vel_coords);
1201:   DMDAVecGetArray(vel_cda,vel_coords,&_vel_coords);

1203:   /* interpolate the coordinates */
1204:   for (j = sj; j < sj+ny; j++) {
1205:     for (i = si; i < si+nx; i++) {
1206:       PetscInt    ngp;
1207:       PetscScalar gp_xi[GAUSS_POINTS][2],gp_weight[GAUSS_POINTS];
1208:       PetscScalar el_coords[8];

1210:       GetElementCoords(_vel_coords,i,j,el_coords);
1211:       ConstructGaussQuadrature(&ngp,gp_xi,gp_weight);

1213:       for (p = 0; p < GAUSS_POINTS; p++) {
1214:         PetscScalar gp_x,gp_y;
1215:         PetscInt    n;
1216:         PetscScalar xi_p[2],Ni_p[4];

1218:         xi_p[0] = gp_xi[p][0];
1219:         xi_p[1] = gp_xi[p][1];
1220:         ConstructQ12D_Ni(xi_p,Ni_p);

1222:         gp_x = 0.0;
1223:         gp_y = 0.0;
1224:         for (n = 0; n < NODES_PER_EL; n++) {
1225:           gp_x = gp_x+Ni_p[n]*el_coords[2*n];
1226:           gp_y = gp_y+Ni_p[n]*el_coords[2*n+1];
1227:         }
1228:         element_props[j][i].gp_coords[2*p]   = gp_x;
1229:         element_props[j][i].gp_coords[2*p+1] = gp_y;
1230:       }
1231:     }
1232:   }

1234:   /* define the coefficients */
1235:   PetscOptionsGetBool(NULL,NULL,"-use_gp_coords",&use_gp_coords,NULL);

1237:   for (j = sj; j < sj+ny; j++) {
1238:     for (i = si; i < si+nx; i++) {
1239:       PetscReal centroid_x = PetscRealPart(_prop_coords[j][i].x); /* centroids of cell */
1240:       PetscReal centroid_y = PetscRealPart(_prop_coords[j][i].y);
1241:       PetscReal coord_x,coord_y;

1243:       if (coefficient_structure == 0) {
1244:         PetscReal opts_eta0,opts_eta1,opts_xc;
1245:         PetscInt  opts_nz;

1247:         opts_eta0 = 1.0;
1248:         opts_eta1 = 1.0;
1249:         opts_xc   = 0.5;
1250:         opts_nz   = 1;

1252:         PetscOptionsGetReal(NULL,NULL,"-solcx_eta0",&opts_eta0,NULL);
1253:         PetscOptionsGetReal(NULL,NULL,"-solcx_eta1",&opts_eta1,NULL);
1254:         PetscOptionsGetReal(NULL,NULL,"-solcx_xc",&opts_xc,NULL);
1255:         PetscOptionsGetInt(NULL,NULL,"-solcx_nz",&opts_nz,NULL);

1257:         for (p = 0; p < GAUSS_POINTS; p++) {
1258:           coord_x = centroid_x;
1259:           coord_y = centroid_y;
1260:           if (use_gp_coords) {
1261:             coord_x = PetscRealPart(element_props[j][i].gp_coords[2*p]);
1262:             coord_y = PetscRealPart(element_props[j][i].gp_coords[2*p+1]);
1263:           }

1265:           element_props[j][i].eta[p] = opts_eta0;
1266:           if (coord_x > opts_xc) element_props[j][i].eta[p] = opts_eta1;

1268:           element_props[j][i].fx[p] = 0.0;
1269:           element_props[j][i].fy[p] = PetscSinReal(opts_nz*PETSC_PI*coord_y)*PetscCosReal(1.0*PETSC_PI*coord_x);
1270:         }
1271:       } else if (coefficient_structure == 1) { /* square sinker */
1272:         PetscReal opts_eta0,opts_eta1,opts_dx,opts_dy;

1274:         opts_eta0 = 1.0;
1275:         opts_eta1 = 1.0;
1276:         opts_dx   = 0.50;
1277:         opts_dy   = 0.50;

1279:         PetscOptionsGetReal(NULL,NULL,"-sinker_eta0",&opts_eta0,NULL);
1280:         PetscOptionsGetReal(NULL,NULL,"-sinker_eta1",&opts_eta1,NULL);
1281:         PetscOptionsGetReal(NULL,NULL,"-sinker_dx",&opts_dx,NULL);
1282:         PetscOptionsGetReal(NULL,NULL,"-sinker_dy",&opts_dy,NULL);

1284:         for (p = 0; p < GAUSS_POINTS; p++) {
1285:           coord_x = centroid_x;
1286:           coord_y = centroid_y;
1287:           if (use_gp_coords) {
1288:             coord_x = PetscRealPart(element_props[j][i].gp_coords[2*p]);
1289:             coord_y = PetscRealPart(element_props[j][i].gp_coords[2*p+1]);
1290:           }

1292:           element_props[j][i].eta[p] = opts_eta0;
1293:           element_props[j][i].fx[p]  = 0.0;
1294:           element_props[j][i].fy[p]  = 0.0;

1296:           if ((coord_x > -0.5*opts_dx+0.5) && (coord_x < 0.5*opts_dx+0.5)) {
1297:             if ((coord_y > -0.5*opts_dy+0.5) && (coord_y < 0.5*opts_dy+0.5)) {
1298:               element_props[j][i].eta[p] =  opts_eta1;
1299:               element_props[j][i].fx[p]  =  0.0;
1300:               element_props[j][i].fy[p]  = -1.0;
1301:             }
1302:           }
1303:         }
1304:       } else if (coefficient_structure == 2) { /* circular sinker */
1305:         PetscReal opts_eta0,opts_eta1,opts_r,radius2;

1307:         opts_eta0 = 1.0;
1308:         opts_eta1 = 1.0;
1309:         opts_r    = 0.25;

1311:         PetscOptionsGetReal(NULL,NULL,"-sinker_eta0",&opts_eta0,NULL);
1312:         PetscOptionsGetReal(NULL,NULL,"-sinker_eta1",&opts_eta1,NULL);
1313:         PetscOptionsGetReal(NULL,NULL,"-sinker_r",&opts_r,NULL);

1315:         for (p = 0; p < GAUSS_POINTS; p++) {
1316:           coord_x = centroid_x;
1317:           coord_y = centroid_y;
1318:           if (use_gp_coords) {
1319:             coord_x = PetscRealPart(element_props[j][i].gp_coords[2*p]);
1320:             coord_y = PetscRealPart(element_props[j][i].gp_coords[2*p+1]);
1321:           }

1323:           element_props[j][i].eta[p] = opts_eta0;
1324:           element_props[j][i].fx[p]  = 0.0;
1325:           element_props[j][i].fy[p]  = 0.0;

1327:           radius2 = (coord_x-0.5)*(coord_x-0.5)+(coord_y-0.5)*(coord_y-0.5);
1328:           if (radius2 < opts_r*opts_r) {
1329:             element_props[j][i].eta[p] =  opts_eta1;
1330:             element_props[j][i].fx[p]  =  0.0;
1331:             element_props[j][i].fy[p]  = -1.0;
1332:           }
1333:         }
1334:       } else if (coefficient_structure == 3) { /* circular and rectangular inclusion */
1335:         PetscReal opts_eta0,opts_eta1,opts_r,opts_dx,opts_dy,opts_c0x,opts_c0y,opts_s0x,opts_s0y,opts_phi,radius2;

1337:         opts_eta0 = 1.0;
1338:         opts_eta1 = 1.0;
1339:         opts_r    = 0.25;
1340:         opts_c0x  = 0.35;       /* circle center */
1341:         opts_c0y  = 0.35;
1342:         opts_s0x  = 0.7;       /* square center */
1343:         opts_s0y  = 0.7;
1344:         opts_dx   = 0.25;
1345:         opts_dy   = 0.25;
1346:         opts_phi  = 25;

1348:         PetscOptionsGetReal(NULL,NULL,"-sinker_eta0",&opts_eta0,NULL);
1349:         PetscOptionsGetReal(NULL,NULL,"-sinker_eta1",&opts_eta1,NULL);
1350:         PetscOptionsGetReal(NULL,NULL,"-sinker_r",&opts_r,NULL);
1351:         PetscOptionsGetReal(NULL,NULL,"-sinker_c0x",&opts_c0x,NULL);
1352:         PetscOptionsGetReal(NULL,NULL,"-sinker_c0y",&opts_c0y,NULL);
1353:         PetscOptionsGetReal(NULL,NULL,"-sinker_s0x",&opts_s0x,NULL);
1354:         PetscOptionsGetReal(NULL,NULL,"-sinker_s0y",&opts_s0y,NULL);
1355:         PetscOptionsGetReal(NULL,NULL,"-sinker_dx",&opts_dx,NULL);
1356:         PetscOptionsGetReal(NULL,NULL,"-sinker_dy",&opts_dy,NULL);
1357:         PetscOptionsGetReal(NULL,NULL,"-sinker_phi",&opts_phi,NULL);
1358:         opts_phi *= PETSC_PI / 180;

1360:         for (p = 0; p < GAUSS_POINTS; p++) {
1361:           coord_x = centroid_x;
1362:           coord_y = centroid_y;
1363:           if (use_gp_coords) {
1364:             coord_x = PetscRealPart(element_props[j][i].gp_coords[2*p]);
1365:             coord_y = PetscRealPart(element_props[j][i].gp_coords[2*p+1]);
1366:           }

1368:           element_props[j][i].eta[p] = opts_eta0;
1369:           element_props[j][i].fx[p]  = 0.0;
1370:           element_props[j][i].fy[p]  = -0.2;

1372:           radius2 = PetscSqr(coord_x - opts_c0x) + PetscSqr(coord_y - opts_c0y);
1373:           if (radius2 < opts_r*opts_r
1374:               || (PetscAbs(+(coord_x - opts_s0x)*PetscCosReal(opts_phi) + (coord_y - opts_s0y)*PetscSinReal(opts_phi)) < opts_dx/2
1375:                   && PetscAbs(-(coord_x - opts_s0x)*PetscSinReal(opts_phi) + (coord_y - opts_s0y)*PetscCosReal(opts_phi)) < opts_dy/2)) {
1376:             element_props[j][i].eta[p] =  opts_eta1;
1377:             element_props[j][i].fx[p]  =  0.0;
1378:             element_props[j][i].fy[p]  = -1.0;
1379:           }
1380:         }
1381:       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Unknown coefficient_structure");
1382:     }
1383:   }
1384:   DMDAVecRestoreArray(prop_cda,prop_coords,&_prop_coords);

1386:   DMDAVecRestoreArray(vel_cda,vel_coords,&_vel_coords);

1388:   DMDAVecRestoreArray(da_prop,l_properties,&element_props);
1389:   DMLocalToGlobalBegin(da_prop,l_properties,ADD_VALUES,properties);
1390:   DMLocalToGlobalEnd(da_prop,l_properties,ADD_VALUES,properties);

1392:   if (output_gnuplot) {
1393:     DMDACoordViewGnuplot2d(da_Stokes,"mesh");
1394:     DMDAViewCoefficientsGnuplot2d(da_prop,properties,"Coeffcients for Stokes eqn.","properties");
1395:   }

1397:   /* Generate a matrix with the correct non-zero pattern of type AIJ. This will work in parallel and serial */
1398:   DMSetMatType(da_Stokes,MATAIJ);
1399:   DMCreateMatrix(da_Stokes,&A);
1400:   DMCreateMatrix(da_Stokes,&B);
1401:   DMCreateGlobalVector(da_Stokes,&f);
1402:   DMCreateGlobalVector(da_Stokes,&X);

1404:   /* assemble A11 */
1405:   MatZeroEntries(A);
1406:   MatZeroEntries(B);
1407:   VecZeroEntries(f);

1409:   AssembleA_Stokes(A,da_Stokes,da_prop,properties);
1410:   AssembleA_PCStokes(B,da_Stokes,da_prop,properties);
1411:   /* build force vector */
1412:   AssembleF_Stokes(f,da_Stokes,da_prop,properties);

1414:   DMDABCApplyFreeSlip(da_Stokes,A,f);
1415:   DMDABCApplyFreeSlip(da_Stokes,B,NULL);

1417:   /* SOLVE */
1418:   KSPCreate(PETSC_COMM_WORLD,&ksp_S);
1419:   KSPSetOptionsPrefix(ksp_S,"stokes_");
1420:   KSPSetDM(ksp_S,da_Stokes);
1421:   KSPSetDMActive(ksp_S,PETSC_FALSE);
1422:   KSPSetOperators(ksp_S,A,B);
1423:   KSPSetFromOptions(ksp_S);
1424:   {
1425:     PC             pc;
1426:     const PetscInt ufields[] = {0,1},pfields[1] = {2};
1427:     KSPGetPC(ksp_S,&pc);
1428:     PCFieldSplitSetBlockSize(pc,3);
1429:     PCFieldSplitSetFields(pc,"u",2,ufields,ufields);
1430:     PCFieldSplitSetFields(pc,"p",1,pfields,pfields);
1431:   }

1433:   {
1434:     PC        pc_S;
1435:     KSP       *sub_ksp,ksp_U;
1436:     PetscInt  nsplits;
1437:     DM        da_U;
1438:     PetscBool is_pcfs;

1440:     KSPSetUp(ksp_S);
1441:     KSPGetPC(ksp_S,&pc_S);

1443:     is_pcfs = PETSC_FALSE;
1444:     PetscObjectTypeCompare((PetscObject)pc_S,PCFIELDSPLIT,&is_pcfs);

1446:     if (is_pcfs) {
1447:       PCFieldSplitGetSubKSP(pc_S,&nsplits,&sub_ksp);
1448:       ksp_U = sub_ksp[0];
1449:       PetscFree(sub_ksp);

1451:       if (nsplits == 2) {
1452:         DMDAGetReducedDMDA(da_Stokes,2,&da_U);

1454:         KSPSetDM(ksp_U,da_U);
1455:         KSPSetDMActive(ksp_U,PETSC_FALSE);
1456:         DMDestroy(&da_U);
1457:       }
1458:     }
1459:   }

1461:   KSPSolve(ksp_S,f,X);

1463:   PetscOptionsGetString(NULL,NULL,"-o",filename,sizeof(filename),&set);
1464:   if (set) {
1465:     char        *ext;
1466:     PetscViewer viewer;
1467:     PetscBool   flg;
1468:     PetscViewerCreate(PETSC_COMM_WORLD,&viewer);
1469:     PetscStrrchr(filename,'.',&ext);
1470:     PetscStrcmp("vts",ext,&flg);
1471:     if (flg) {
1472:       PetscViewerSetType(viewer,PETSCVIEWERVTK);
1473:     } else {
1474:       PetscViewerSetType(viewer,PETSCVIEWERBINARY);
1475:     }
1476:     PetscViewerFileSetMode(viewer,FILE_MODE_WRITE);
1477:     PetscViewerFileSetName(viewer,filename);
1478:     VecView(X,viewer);
1479:     PetscViewerDestroy(&viewer);
1480:   }
1481:   if (output_gnuplot) {
1482:     DMDAViewGnuplot2d(da_Stokes,X,"Velocity solution for Stokes eqn.","X");
1483:   }

1485:   KSPGetIterationNumber(ksp_S,&its);

1487:   if (coefficient_structure == 0) {
1488:     PetscReal opts_eta0,opts_eta1,opts_xc;
1489:     PetscInt  opts_nz,N;
1490:     DM        da_Stokes_analytic;
1491:     Vec       X_analytic;
1492:     PetscReal nrm1[3],nrm2[3],nrmI[3];

1494:     opts_eta0 = 1.0;
1495:     opts_eta1 = 1.0;
1496:     opts_xc   = 0.5;
1497:     opts_nz   = 1;

1499:     PetscOptionsGetReal(NULL,NULL,"-solcx_eta0",&opts_eta0,NULL);
1500:     PetscOptionsGetReal(NULL,NULL,"-solcx_eta1",&opts_eta1,NULL);
1501:     PetscOptionsGetReal(NULL,NULL,"-solcx_xc",&opts_xc,NULL);
1502:     PetscOptionsGetInt(NULL,NULL,"-solcx_nz",&opts_nz,NULL);

1504:     DMDACreateSolCx(opts_eta0,opts_eta1,opts_xc,opts_nz,mx,my,&da_Stokes_analytic,&X_analytic);
1505:     if (output_gnuplot) {
1506:       DMDAViewGnuplot2d(da_Stokes_analytic,X_analytic,"Analytic solution for Stokes eqn.","X_analytic");
1507:     }
1508:     DMDAIntegrateErrors(da_Stokes_analytic,X,X_analytic);

1510:     VecAXPY(X_analytic,-1.0,X);
1511:     VecGetSize(X_analytic,&N);
1512:     N    = N/3;

1514:     VecStrideNorm(X_analytic,0,NORM_1,&nrm1[0]);
1515:     VecStrideNorm(X_analytic,0,NORM_2,&nrm2[0]);
1516:     VecStrideNorm(X_analytic,0,NORM_INFINITY,&nrmI[0]);

1518:     VecStrideNorm(X_analytic,1,NORM_1,&nrm1[1]);
1519:     VecStrideNorm(X_analytic,1,NORM_2,&nrm2[1]);
1520:     VecStrideNorm(X_analytic,1,NORM_INFINITY,&nrmI[1]);

1522:     VecStrideNorm(X_analytic,2,NORM_1,&nrm1[2]);
1523:     VecStrideNorm(X_analytic,2,NORM_2,&nrm2[2]);
1524:     VecStrideNorm(X_analytic,2,NORM_INFINITY,&nrmI[2]);

1526:     DMDestroy(&da_Stokes_analytic);
1527:     VecDestroy(&X_analytic);
1528:   }

1530:   KSPDestroy(&ksp_S);
1531:   VecDestroy(&X);
1532:   VecDestroy(&f);
1533:   MatDestroy(&A);
1534:   MatDestroy(&B);

1536:   DMDestroy(&da_Stokes);
1537:   DMDestroy(&da_prop);

1539:   VecDestroy(&properties);
1540:   VecDestroy(&l_properties);
1541:   return(0);
1542: }

1546: int main(int argc,char **args)
1547: {
1549:   PetscInt       mx,my;

1551:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
1552:   mx   = my = 10;
1553:   PetscOptionsGetInt(NULL,NULL,"-mx",&mx,NULL);
1554:   PetscOptionsGetInt(NULL,NULL,"-my",&my,NULL);
1555:   solve_stokes_2d_coupled(mx,my);
1556:   PetscFinalize();
1557:   return ierr;
1558: }

1560: /* -------------------------- helpers for boundary conditions -------------------------------- */
1563: static PetscErrorCode BCApplyZero_EAST(DM da,PetscInt d_idx,Mat A,Vec b)
1564: {
1565:   DM                     cda;
1566:   Vec                    coords;
1567:   PetscInt               si,sj,nx,ny,i,j;
1568:   PetscInt               M,N;
1569:   DMDACoor2d             **_coords;
1570:   const PetscInt         *g_idx;
1571:   PetscInt               *bc_global_ids;
1572:   PetscScalar            *bc_vals;
1573:   PetscInt               nbcs;
1574:   PetscInt               n_dofs;
1575:   PetscErrorCode         ierr;
1576:   ISLocalToGlobalMapping ltogm;

1579:   DMGetLocalToGlobalMapping(da,&ltogm);
1580:   ISLocalToGlobalMappingGetIndices(ltogm,&g_idx);

1582:   DMGetCoordinateDM(da,&cda);
1583:   DMGetCoordinatesLocal(da,&coords);
1584:   DMDAVecGetArray(cda,coords,&_coords);
1585:   DMDAGetGhostCorners(cda,&si,&sj,0,&nx,&ny,0);
1586:   DMDAGetInfo(da,0,&M,&N,0,0,0,0,&n_dofs,0,0,0,0,0);

1588:   PetscMalloc1(ny*n_dofs,&bc_global_ids);
1589:   PetscMalloc1(ny*n_dofs,&bc_vals);

1591:   /* init the entries to -1 so VecSetValues will ignore them */
1592:   for (i = 0; i < ny*n_dofs; i++) bc_global_ids[i] = -1;

1594:   i = nx-1;
1595:   for (j = 0; j < ny; j++) {
1596:     PetscInt local_id;

1598:     local_id = i+j*nx;

1600:     bc_global_ids[j] = g_idx[n_dofs*local_id+d_idx];

1602:     bc_vals[j] =  0.0;
1603:   }
1604:   ISLocalToGlobalMappingRestoreIndices(ltogm,&g_idx);
1605:   nbcs = 0;
1606:   if ((si+nx) == (M)) nbcs = ny;

1608:   if (b) {
1609:     VecSetValues(b,nbcs,bc_global_ids,bc_vals,INSERT_VALUES);
1610:     VecAssemblyBegin(b);
1611:     VecAssemblyEnd(b);
1612:   }
1613:   if (A) {
1614:     MatZeroRowsColumns(A,nbcs,bc_global_ids,1.0,0,0);
1615:   }

1617:   PetscFree(bc_vals);
1618:   PetscFree(bc_global_ids);

1620:   DMDAVecRestoreArray(cda,coords,&_coords);
1621:   return(0);
1622: }

1626: static PetscErrorCode BCApplyZero_WEST(DM da,PetscInt d_idx,Mat A,Vec b)
1627: {
1628:   DM                     cda;
1629:   Vec                    coords;
1630:   PetscInt               si,sj,nx,ny,i,j;
1631:   PetscInt               M,N;
1632:   DMDACoor2d             **_coords;
1633:   const PetscInt         *g_idx;
1634:   PetscInt               *bc_global_ids;
1635:   PetscScalar            *bc_vals;
1636:   PetscInt               nbcs;
1637:   PetscInt               n_dofs;
1638:   PetscErrorCode         ierr;
1639:   ISLocalToGlobalMapping ltogm;

1642:   DMGetLocalToGlobalMapping(da,&ltogm);
1643:   ISLocalToGlobalMappingGetIndices(ltogm,&g_idx);

1645:   DMGetCoordinateDM(da,&cda);
1646:   DMGetCoordinatesLocal(da,&coords);
1647:   DMDAVecGetArray(cda,coords,&_coords);
1648:   DMDAGetGhostCorners(cda,&si,&sj,0,&nx,&ny,0);
1649:   DMDAGetInfo(da,0,&M,&N,0,0,0,0,&n_dofs,0,0,0,0,0);

1651:   PetscMalloc1(ny*n_dofs,&bc_global_ids);
1652:   PetscMalloc1(ny*n_dofs,&bc_vals);

1654:   /* init the entries to -1 so VecSetValues will ignore them */
1655:   for (i = 0; i < ny*n_dofs; i++) bc_global_ids[i] = -1;

1657:   i = 0;
1658:   for (j = 0; j < ny; j++) {
1659:     PetscInt local_id;

1661:     local_id = i+j*nx;

1663:     bc_global_ids[j] = g_idx[n_dofs*local_id+d_idx];

1665:     bc_vals[j] =  0.0;
1666:   }
1667:   ISLocalToGlobalMappingRestoreIndices(ltogm,&g_idx);
1668:   nbcs = 0;
1669:   if (si == 0) nbcs = ny;

1671:   if (b) {
1672:     VecSetValues(b,nbcs,bc_global_ids,bc_vals,INSERT_VALUES);
1673:     VecAssemblyBegin(b);
1674:     VecAssemblyEnd(b);
1675:   }

1677:   if (A) {
1678:     MatZeroRowsColumns(A,nbcs,bc_global_ids,1.0,0,0);
1679:   }

1681:   PetscFree(bc_vals);
1682:   PetscFree(bc_global_ids);

1684:   DMDAVecRestoreArray(cda,coords,&_coords);
1685:   return(0);
1686: }

1690: static PetscErrorCode BCApplyZero_NORTH(DM da,PetscInt d_idx,Mat A,Vec b)
1691: {
1692:   DM                     cda;
1693:   Vec                    coords;
1694:   PetscInt               si,sj,nx,ny,i,j;
1695:   PetscInt               M,N;
1696:   DMDACoor2d             **_coords;
1697:   const PetscInt         *g_idx;
1698:   PetscInt               *bc_global_ids;
1699:   PetscScalar            *bc_vals;
1700:   PetscInt               nbcs;
1701:   PetscInt               n_dofs;
1702:   PetscErrorCode         ierr;
1703:   ISLocalToGlobalMapping ltogm;

1706:   DMGetLocalToGlobalMapping(da,&ltogm);
1707:   ISLocalToGlobalMappingGetIndices(ltogm,&g_idx);

1709:   DMGetCoordinateDM(da,&cda);
1710:   DMGetCoordinatesLocal(da,&coords);
1711:   DMDAVecGetArray(cda,coords,&_coords);
1712:   DMDAGetGhostCorners(cda,&si,&sj,0,&nx,&ny,0);
1713:   DMDAGetInfo(da,0,&M,&N,0,0,0,0,&n_dofs,0,0,0,0,0);

1715:   PetscMalloc1(nx,&bc_global_ids);
1716:   PetscMalloc1(nx,&bc_vals);

1718:   /* init the entries to -1 so VecSetValues will ignore them */
1719:   for (i = 0; i < nx; i++) bc_global_ids[i] = -1;

1721:   j = ny-1;
1722:   for (i = 0; i < nx; i++) {
1723:     PetscInt local_id;

1725:     local_id = i+j*nx;

1727:     bc_global_ids[i] = g_idx[n_dofs*local_id+d_idx];

1729:     bc_vals[i] =  0.0;
1730:   }
1731:   ISLocalToGlobalMappingRestoreIndices(ltogm,&g_idx);
1732:   nbcs = 0;
1733:   if ((sj+ny) == (N)) nbcs = nx;

1735:   if (b) {
1736:     VecSetValues(b,nbcs,bc_global_ids,bc_vals,INSERT_VALUES);
1737:     VecAssemblyBegin(b);
1738:     VecAssemblyEnd(b);
1739:   }
1740:   if (A) {
1741:     MatZeroRowsColumns(A,nbcs,bc_global_ids,1.0,0,0);
1742:   }

1744:   PetscFree(bc_vals);
1745:   PetscFree(bc_global_ids);

1747:   DMDAVecRestoreArray(cda,coords,&_coords);
1748:   return(0);
1749: }

1753: static PetscErrorCode BCApplyZero_SOUTH(DM da,PetscInt d_idx,Mat A,Vec b)
1754: {
1755:   DM                     cda;
1756:   Vec                    coords;
1757:   PetscInt               si,sj,nx,ny,i,j;
1758:   PetscInt               M,N;
1759:   DMDACoor2d             **_coords;
1760:   const PetscInt         *g_idx;
1761:   PetscInt               *bc_global_ids;
1762:   PetscScalar            *bc_vals;
1763:   PetscInt               nbcs;
1764:   PetscInt               n_dofs;
1765:   PetscErrorCode         ierr;
1766:   ISLocalToGlobalMapping ltogm;

1769:   DMGetLocalToGlobalMapping(da,&ltogm);
1770:   ISLocalToGlobalMappingGetIndices(ltogm,&g_idx);

1772:   DMGetCoordinateDM(da,&cda);
1773:   DMGetCoordinatesLocal(da,&coords);
1774:   DMDAVecGetArray(cda,coords,&_coords);
1775:   DMDAGetGhostCorners(cda,&si,&sj,0,&nx,&ny,0);
1776:   DMDAGetInfo(da,0,&M,&N,0,0,0,0,&n_dofs,0,0,0,0,0);

1778:   PetscMalloc1(nx,&bc_global_ids);
1779:   PetscMalloc1(nx,&bc_vals);

1781:   /* init the entries to -1 so VecSetValues will ignore them */
1782:   for (i = 0; i < nx; i++) bc_global_ids[i] = -1;

1784:   j = 0;
1785:   for (i = 0; i < nx; i++) {
1786:     PetscInt local_id;

1788:     local_id = i+j*nx;

1790:     bc_global_ids[i] = g_idx[n_dofs*local_id+d_idx];

1792:     bc_vals[i] =  0.0;
1793:   }
1794:   ISLocalToGlobalMappingRestoreIndices(ltogm,&g_idx);
1795:   nbcs = 0;
1796:   if (sj == 0) nbcs = nx;

1798:   if (b) {
1799:     VecSetValues(b,nbcs,bc_global_ids,bc_vals,INSERT_VALUES);
1800:     VecAssemblyBegin(b);
1801:     VecAssemblyEnd(b);
1802:   }
1803:   if (A) {
1804:     MatZeroRowsColumns(A,nbcs,bc_global_ids,1.0,0,0);
1805:   }

1807:   PetscFree(bc_vals);
1808:   PetscFree(bc_global_ids);

1810:   DMDAVecRestoreArray(cda,coords,&_coords);
1811:   return(0);
1812: }

1814: /* Impose free slip boundary conditions; u_{i} n_{i} = 0, tau_{ij} t_j = 0 */
1817: static PetscErrorCode DMDABCApplyFreeSlip(DM da_Stokes,Mat A,Vec f)
1818: {

1822:   BCApplyZero_NORTH(da_Stokes,1,A,f);
1823:   BCApplyZero_EAST(da_Stokes,0,A,f);
1824:   BCApplyZero_SOUTH(da_Stokes,1,A,f);
1825:   BCApplyZero_WEST(da_Stokes,0,A,f);
1826:   return(0);
1827: }