Actual source code: ex42.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  1: static char help[] = "Solves the incompressible, variable viscosity stokes equation in 3d using Q1Q1 elements, \n\
  2: stabilized with Bochev's polynomial projection method. Note that implementation here assumes \n\
  3: all boundaries are free-slip, i.e. zero normal flow and zero tangential stress \n\
  4:      -mx : number elements in x-direction \n\
  5:      -my : number elements in y-direction \n\
  6:      -mz : number elements in z-direction \n\
  7:      -stokes_ksp_monitor_blocks : active monitor for each component u,v,w,p \n\
  8:      -model : defines viscosity and forcing function. Choose either: 0 (isoviscous), 1 (manufactured-broken), 2 (solcx-2d), 3 (sinker) \n";

 10: /* Contributed by Dave May */

 12: #include <petscksp.h>
 13: #include <petscdmda.h>

 15: #define PROFILE_TIMING
 16: #define ASSEMBLE_LOWER_TRIANGULAR

 18: #define NSD            3 /* number of spatial dimensions */
 19: #define NODES_PER_EL   8 /* nodes per element */
 20: #define U_DOFS         3 /* degrees of freedom per velocity node */
 21: #define P_DOFS         1 /* degrees of freedom per pressure node */
 22: #define GAUSS_POINTS   8

 24: /* Gauss point based evaluation */
 25: typedef struct {
 26:   PetscScalar gp_coords[NSD*GAUSS_POINTS];
 27:   PetscScalar eta[GAUSS_POINTS];
 28:   PetscScalar fx[GAUSS_POINTS];
 29:   PetscScalar fy[GAUSS_POINTS];
 30:   PetscScalar fz[GAUSS_POINTS];
 31:   PetscScalar hc[GAUSS_POINTS];
 32: } GaussPointCoefficients;

 34: typedef struct {
 35:   PetscScalar u_dof;
 36:   PetscScalar v_dof;
 37:   PetscScalar w_dof;
 38:   PetscScalar p_dof;
 39: } StokesDOF;

 41: typedef struct _p_CellProperties *CellProperties;
 42: struct _p_CellProperties {
 43:   PetscInt               ncells;
 44:   PetscInt               mx,my,mz;
 45:   PetscInt               sex,sey,sez;
 46:   GaussPointCoefficients *gpc;
 47: };

 49: static PetscErrorCode DMDAGetElementCorners(DM da,PetscInt *sx,PetscInt *sy,PetscInt *sz,PetscInt *mx,PetscInt *my,PetscInt *mz);

 51: /* elements */
 54: PetscErrorCode CellPropertiesCreate(DM da_stokes,CellProperties *C)
 55: {
 57:   CellProperties cells;
 58:   PetscInt       mx,my,mz,sex,sey,sez;

 61:   PetscMalloc(sizeof(struct _p_CellProperties),&cells);

 63:   DMDAGetElementCorners(da_stokes,&sex,&sey,&sez,&mx,&my,&mz);

 65:   cells->mx     = mx;
 66:   cells->my     = my;
 67:   cells->mz     = mz;
 68:   cells->ncells = mx * my * mz;
 69:   cells->sex    = sex;
 70:   cells->sey    = sey;
 71:   cells->sez    = sez;

 73:   PetscMalloc(sizeof(GaussPointCoefficients)*mx*my*mz,&cells->gpc);

 75:   *C = cells;
 76:   return(0);
 77: }

 81: PetscErrorCode CellPropertiesDestroy(CellProperties *C)
 82: {
 84:   CellProperties cells;

 87:   if (!C) return(0);
 88:   cells = *C;
 89:   PetscFree(cells->gpc);
 90:   PetscFree(cells);
 91:   *C = NULL;
 92:   return(0);
 93: }

 97: PetscErrorCode CellPropertiesGetCell(CellProperties C,PetscInt II,PetscInt J,PetscInt K,GaussPointCoefficients **G)
 98: {
100:   *G = &C->gpc[(II-C->sex) + (J-C->sey)*C->mx + (K-C->sez)*C->mx*C->my];
101:   return(0);
102: }

104: /* FEM routines */
105: /*
106:  Element: Local basis function ordering
107:  1-----2
108:  |     |
109:  |     |
110:  0-----3
111:  */
112: static void ShapeFunctionQ13D_Evaluate(PetscScalar _xi[],PetscScalar Ni[])
113: {
114:   PetscReal xi   = PetscRealPart(_xi[0]);
115:   PetscReal eta  = PetscRealPart(_xi[1]);
116:   PetscReal zeta = PetscRealPart(_xi[2]);

118:   Ni[0] = 0.125 * (1.0 - xi) * (1.0 - eta) * (1.0 - zeta);
119:   Ni[1] = 0.125 * (1.0 - xi) * (1.0 + eta) * (1.0 - zeta);
120:   Ni[2] = 0.125 * (1.0 + xi) * (1.0 + eta) * (1.0 - zeta);
121:   Ni[3] = 0.125 * (1.0 + xi) * (1.0 - eta) * (1.0 - zeta);

123:   Ni[4] = 0.125 * (1.0 - xi) * (1.0 - eta) * (1.0 + zeta);
124:   Ni[5] = 0.125 * (1.0 - xi) * (1.0 + eta) * (1.0 + zeta);
125:   Ni[6] = 0.125 * (1.0 + xi) * (1.0 + eta) * (1.0 + zeta);
126:   Ni[7] = 0.125 * (1.0 + xi) * (1.0 - eta) * (1.0 + zeta);
127: }

129: static void ShapeFunctionQ13D_Evaluate_dxi(PetscScalar _xi[],PetscScalar GNi[][NODES_PER_EL])
130: {
131:   PetscReal xi   = PetscRealPart(_xi[0]);
132:   PetscReal eta  = PetscRealPart(_xi[1]);
133:   PetscReal zeta = PetscRealPart(_xi[2]);
134:   /* xi deriv */
135:   GNi[0][0] = -0.125 * (1.0 - eta) * (1.0 - zeta);
136:   GNi[0][1] = -0.125 * (1.0 + eta) * (1.0 - zeta);
137:   GNi[0][2] =  0.125 * (1.0 + eta) * (1.0 - zeta);
138:   GNi[0][3] =  0.125 * (1.0 - eta) * (1.0 - zeta);

140:   GNi[0][4] = -0.125 * (1.0 - eta) * (1.0 + zeta);
141:   GNi[0][5] = -0.125 * (1.0 + eta) * (1.0 + zeta);
142:   GNi[0][6] =  0.125 * (1.0 + eta) * (1.0 + zeta);
143:   GNi[0][7] =  0.125 * (1.0 - eta) * (1.0 + zeta);
144:   /* eta deriv */
145:   GNi[1][0] = -0.125 * (1.0 - xi) * (1.0 - zeta);
146:   GNi[1][1] =  0.125 * (1.0 - xi) * (1.0 - zeta);
147:   GNi[1][2] =  0.125 * (1.0 + xi) * (1.0 - zeta);
148:   GNi[1][3] = -0.125 * (1.0 + xi) * (1.0 - zeta);

150:   GNi[1][4] = -0.125 * (1.0 - xi) * (1.0 + zeta);
151:   GNi[1][5] =  0.125 * (1.0 - xi) * (1.0 + zeta);
152:   GNi[1][6] =  0.125 * (1.0 + xi) * (1.0 + zeta);
153:   GNi[1][7] = -0.125 * (1.0 + xi) * (1.0 + zeta);
154:   /* zeta deriv */
155:   GNi[2][0] = -0.125 * (1.0 - xi) * (1.0 - eta);
156:   GNi[2][1] = -0.125 * (1.0 - xi) * (1.0 + eta);
157:   GNi[2][2] = -0.125 * (1.0 + xi) * (1.0 + eta);
158:   GNi[2][3] = -0.125 * (1.0 + xi) * (1.0 - eta);

160:   GNi[2][4] = 0.125 * (1.0 - xi) * (1.0 - eta);
161:   GNi[2][5] = 0.125 * (1.0 - xi) * (1.0 + eta);
162:   GNi[2][6] = 0.125 * (1.0 + xi) * (1.0 + eta);
163:   GNi[2][7] = 0.125 * (1.0 + xi) * (1.0 - eta);
164: }

166: static void matrix_inverse_3x3(PetscScalar A[3][3],PetscScalar B[3][3])
167: {
168:   PetscScalar t4, t6, t8, t10, t12, t14, t17;

170:   t4  = A[2][0] * A[0][1];
171:   t6  = A[2][0] * A[0][2];
172:   t8  = A[1][0] * A[0][1];
173:   t10 = A[1][0] * A[0][2];
174:   t12 = A[0][0] * A[1][1];
175:   t14 = A[0][0] * A[1][2];
176:   t17 = 0.1e1 / (t4 * A[1][2] - t6 * A[1][1] - t8 * A[2][2] + t10 * A[2][1] + t12 * A[2][2] - t14 * A[2][1]);

178:   B[0][0] = (A[1][1] * A[2][2] - A[1][2] * A[2][1]) * t17;
179:   B[0][1] = -(A[0][1] * A[2][2] - A[0][2] * A[2][1]) * t17;
180:   B[0][2] = (A[0][1] * A[1][2] - A[0][2] * A[1][1]) * t17;
181:   B[1][0] = -(-A[2][0] * A[1][2] + A[1][0] * A[2][2]) * t17;
182:   B[1][1] = (-t6 + A[0][0] * A[2][2]) * t17;
183:   B[1][2] = -(-t10 + t14) * t17;
184:   B[2][0] = (-A[2][0] * A[1][1] + A[1][0] * A[2][1]) * t17;
185:   B[2][1] = -(-t4 + A[0][0] * A[2][1]) * t17;
186:   B[2][2] = (-t8 + t12) * t17;
187: }

189: static void ShapeFunctionQ13D_Evaluate_dx(PetscScalar GNi[][NODES_PER_EL],PetscScalar GNx[][NODES_PER_EL],PetscScalar coords[],PetscScalar *det_J)
190: {
191:   PetscScalar J00,J01,J02,J10,J11,J12,J20,J21,J22;
192:   PetscInt    n;
193:   PetscScalar iJ[3][3],JJ[3][3];

195:   J00 = J01 = J02 = 0.0;
196:   J10 = J11 = J12 = 0.0;
197:   J20 = J21 = J22 = 0.0;
198:   for (n=0; n<NODES_PER_EL; n++) {
199:     PetscScalar cx = coords[NSD*n + 0];
200:     PetscScalar cy = coords[NSD*n + 1];
201:     PetscScalar cz = coords[NSD*n + 2];

203:     /* J_ij = d(x_j) / d(xi_i) */ /* J_ij = \sum _I GNi[j][I} * x_i */
204:     J00 = J00 + GNi[0][n] * cx;   /* J_xx */
205:     J01 = J01 + GNi[0][n] * cy;   /* J_xy = dx/deta */
206:     J02 = J02 + GNi[0][n] * cz;   /* J_xz = dx/dzeta */

208:     J10 = J10 + GNi[1][n] * cx;   /* J_yx = dy/dxi */
209:     J11 = J11 + GNi[1][n] * cy;   /* J_yy */
210:     J12 = J12 + GNi[1][n] * cz;   /* J_yz */

212:     J20 = J20 + GNi[2][n] * cx;   /* J_zx */
213:     J21 = J21 + GNi[2][n] * cy;   /* J_zy */
214:     J22 = J22 + GNi[2][n] * cz;   /* J_zz */
215:   }

217:   JJ[0][0] = J00;      JJ[0][1] = J01;      JJ[0][2] = J02;
218:   JJ[1][0] = J10;      JJ[1][1] = J11;      JJ[1][2] = J12;
219:   JJ[2][0] = J20;      JJ[2][1] = J21;      JJ[2][2] = J22;

221:   matrix_inverse_3x3(JJ,iJ);

223:   *det_J = J00*J11*J22 - J00*J12*J21 - J10*J01*J22 + J10*J02*J21 + J20*J01*J12 - J20*J02*J11;

225:   for (n=0; n<NODES_PER_EL; n++) {
226:     GNx[0][n] = GNi[0][n]*iJ[0][0] + GNi[1][n]*iJ[0][1] + GNi[2][n]*iJ[0][2];
227:     GNx[1][n] = GNi[0][n]*iJ[1][0] + GNi[1][n]*iJ[1][1] + GNi[2][n]*iJ[1][2];
228:     GNx[2][n] = GNi[0][n]*iJ[2][0] + GNi[1][n]*iJ[2][1] + GNi[2][n]*iJ[2][2];
229:   }
230: }

232: static void ConstructGaussQuadrature3D(PetscInt *ngp,PetscScalar gp_xi[][NSD],PetscScalar gp_weight[])
233: {
234:   *ngp        = 8;
235:   gp_xi[0][0] = -0.57735026919; gp_xi[0][1] = -0.57735026919; gp_xi[0][2] = -0.57735026919;
236:   gp_xi[1][0] = -0.57735026919; gp_xi[1][1] =  0.57735026919; gp_xi[1][2] = -0.57735026919;
237:   gp_xi[2][0] =  0.57735026919; gp_xi[2][1] =  0.57735026919; gp_xi[2][2] = -0.57735026919;
238:   gp_xi[3][0] =  0.57735026919; gp_xi[3][1] = -0.57735026919; gp_xi[3][2] = -0.57735026919;

240:   gp_xi[4][0] = -0.57735026919; gp_xi[4][1] = -0.57735026919; gp_xi[4][2] =  0.57735026919;
241:   gp_xi[5][0] = -0.57735026919; gp_xi[5][1] =  0.57735026919; gp_xi[5][2] =  0.57735026919;
242:   gp_xi[6][0] =  0.57735026919; gp_xi[6][1] =  0.57735026919; gp_xi[6][2] =  0.57735026919;
243:   gp_xi[7][0] =  0.57735026919; gp_xi[7][1] = -0.57735026919; gp_xi[7][2] =  0.57735026919;

245:   gp_weight[0] = 1.0;
246:   gp_weight[1] = 1.0;
247:   gp_weight[2] = 1.0;
248:   gp_weight[3] = 1.0;

250:   gp_weight[4] = 1.0;
251:   gp_weight[5] = 1.0;
252:   gp_weight[6] = 1.0;
253:   gp_weight[7] = 1.0;
254: }

256: /* procs to the left claim the ghost node as their element */
259: static PetscErrorCode DMDAGetLocalElementSize(DM da,PetscInt *mxl,PetscInt *myl,PetscInt *mzl)
260: {
261:   PetscInt       m,n,p,M,N,P;
262:   PetscInt       sx,sy,sz;

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

269:   if (mxl) {
270:     *mxl = m;
271:     if ((sx+m) == M) *mxl = m-1;  /* last proc */
272:   }
273:   if (myl) {
274:     *myl = n;
275:     if ((sy+n) == N) *myl = n-1;  /* last proc */
276:   }
277:   if (mzl) {
278:     *mzl = p;
279:     if ((sz+p) == P) *mzl = p-1;  /* last proc */
280:   }
281:   return(0);
282: }

286: static PetscErrorCode DMDAGetElementCorners(DM da,PetscInt *sx,PetscInt *sy,PetscInt *sz,PetscInt *mx,PetscInt *my,PetscInt *mz)
287: {
288:   PetscInt       si,sj,sk;

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

294:   if (sx) {
295:     *sx = si;
296:     if (si != 0) *sx = si+1;
297:   }
298:   if (sy) {
299:     *sy = sj;
300:     if (sj != 0) *sy = sj+1;
301:   }
302:   if (sz) {
303:     *sz = sk;
304:     if (sk != 0) *sz = sk+1;
305:   }
306:   DMDAGetLocalElementSize(da,mx,my,mz);
307:   return(0);
308: }

310: /*
311:  i,j are the element indices
312:  The unknown is a vector quantity.
313:  The s[].c is used to indicate the degree of freedom.
314:  */
317: static PetscErrorCode DMDAGetElementEqnums3D_up(MatStencil s_u[],MatStencil s_p[],PetscInt i,PetscInt j,PetscInt k)
318: {
319:   PetscInt n;

322:   /* velocity */
323:   n = 0;
324:   /* node 0 */
325:   s_u[n].i = i; s_u[n].j = j; s_u[n].k = k; s_u[n].c = 0; n++; /* Vx0 */
326:   s_u[n].i = i; s_u[n].j = j; s_u[n].k = k; s_u[n].c = 1; n++; /* Vy0 */
327:   s_u[n].i = i; s_u[n].j = j; s_u[n].k = k; s_u[n].c = 2; n++; /* Vz0 */

329:   s_u[n].i = i; s_u[n].j = j+1; s_u[n].k = k; s_u[n].c = 0; n++;
330:   s_u[n].i = i; s_u[n].j = j+1; s_u[n].k = k; s_u[n].c = 1; n++;
331:   s_u[n].i = i; s_u[n].j = j+1; s_u[n].k = k; s_u[n].c = 2; n++;

333:   s_u[n].i = i+1; s_u[n].j = j+1; s_u[n].k = k; s_u[n].c = 0; n++;
334:   s_u[n].i = i+1; s_u[n].j = j+1; s_u[n].k = k; s_u[n].c = 1; n++;
335:   s_u[n].i = i+1; s_u[n].j = j+1; s_u[n].k = k; s_u[n].c = 2; n++;

337:   s_u[n].i = i+1; s_u[n].j = j; s_u[n].k = k; s_u[n].c = 0; n++;
338:   s_u[n].i = i+1; s_u[n].j = j; s_u[n].k = k; s_u[n].c = 1; n++;
339:   s_u[n].i = i+1; s_u[n].j = j; s_u[n].k = k; s_u[n].c = 2; n++;

341:   /* */
342:   s_u[n].i = i; s_u[n].j = j; s_u[n].k = k+1; s_u[n].c = 0; n++; /* Vx4 */
343:   s_u[n].i = i; s_u[n].j = j; s_u[n].k = k+1; s_u[n].c = 1; n++; /* Vy4 */
344:   s_u[n].i = i; s_u[n].j = j; s_u[n].k = k+1; s_u[n].c = 2; n++; /* Vz4 */

346:   s_u[n].i = i; s_u[n].j = j+1; s_u[n].k = k+1; s_u[n].c = 0; n++;
347:   s_u[n].i = i; s_u[n].j = j+1; s_u[n].k = k+1; s_u[n].c = 1; n++;
348:   s_u[n].i = i; s_u[n].j = j+1; s_u[n].k = k+1; s_u[n].c = 2; n++;

350:   s_u[n].i = i+1; s_u[n].j = j+1; s_u[n].k = k+1; s_u[n].c = 0; n++;
351:   s_u[n].i = i+1; s_u[n].j = j+1; s_u[n].k = k+1; s_u[n].c = 1; n++;
352:   s_u[n].i = i+1; s_u[n].j = j+1; s_u[n].k = k+1; s_u[n].c = 2; n++;

354:   s_u[n].i = i+1; s_u[n].j = j; s_u[n].k = k+1; s_u[n].c = 0; n++;
355:   s_u[n].i = i+1; s_u[n].j = j; s_u[n].k = k+1; s_u[n].c = 1; n++;
356:   s_u[n].i = i+1; s_u[n].j = j; s_u[n].k = k+1; s_u[n].c = 2; n++;

358:   /* pressure */
359:   n = 0;

361:   s_p[n].i = i;   s_p[n].j = j;   s_p[n].k = k; s_p[n].c = 3; n++; /* P0 */
362:   s_p[n].i = i;   s_p[n].j = j+1; s_p[n].k = k; s_p[n].c = 3; n++;
363:   s_p[n].i = i+1; s_p[n].j = j+1; s_p[n].k = k; s_p[n].c = 3; n++;
364:   s_p[n].i = i+1; s_p[n].j = j;   s_p[n].k = k; s_p[n].c = 3; n++;

366:   s_p[n].i = i;   s_p[n].j = j;   s_p[n].k = k+1; s_p[n].c = 3; n++; /* P0 */
367:   s_p[n].i = i;   s_p[n].j = j+1; s_p[n].k = k+1; s_p[n].c = 3; n++;
368:   s_p[n].i = i+1; s_p[n].j = j+1; s_p[n].k = k+1; s_p[n].c = 3; n++;
369:   s_p[n].i = i+1; s_p[n].j = j;   s_p[n].k = k+1; s_p[n].c = 3; n++;
370:   return(0);
371: }

375: static PetscErrorCode GetElementCoords3D(DMDACoor3d ***coords,PetscInt i,PetscInt j,PetscInt k,PetscScalar el_coord[])
376: {
378:   /* get coords for the element */
379:   el_coord[0] = coords[k][j][i].x;
380:   el_coord[1] = coords[k][j][i].y;
381:   el_coord[2] = coords[k][j][i].z;

383:   el_coord[3] = coords[k][j+1][i].x;
384:   el_coord[4] = coords[k][j+1][i].y;
385:   el_coord[5] = coords[k][j+1][i].z;

387:   el_coord[6] = coords[k][j+1][i+1].x;
388:   el_coord[7] = coords[k][j+1][i+1].y;
389:   el_coord[8] = coords[k][j+1][i+1].z;

391:   el_coord[9]  = coords[k][j][i+1].x;
392:   el_coord[10] = coords[k][j][i+1].y;
393:   el_coord[11] = coords[k][j][i+1].z;

395:   el_coord[12] = coords[k+1][j][i].x;
396:   el_coord[13] = coords[k+1][j][i].y;
397:   el_coord[14] = coords[k+1][j][i].z;

399:   el_coord[15] = coords[k+1][j+1][i].x;
400:   el_coord[16] = coords[k+1][j+1][i].y;
401:   el_coord[17] = coords[k+1][j+1][i].z;

403:   el_coord[18] = coords[k+1][j+1][i+1].x;
404:   el_coord[19] = coords[k+1][j+1][i+1].y;
405:   el_coord[20] = coords[k+1][j+1][i+1].z;

407:   el_coord[21] = coords[k+1][j][i+1].x;
408:   el_coord[22] = coords[k+1][j][i+1].y;
409:   el_coord[23] = coords[k+1][j][i+1].z;
410:   return(0);
411: }

415: static PetscErrorCode StokesDAGetNodalFields3D(StokesDOF ***field,PetscInt i,PetscInt j,PetscInt k,StokesDOF nodal_fields[])
416: {
418:   /* get the nodal fields for u */
419:   nodal_fields[0].u_dof = field[k][j][i].u_dof;
420:   nodal_fields[0].v_dof = field[k][j][i].v_dof;
421:   nodal_fields[0].w_dof = field[k][j][i].w_dof;

423:   nodal_fields[1].u_dof = field[k][j+1][i].u_dof;
424:   nodal_fields[1].v_dof = field[k][j+1][i].v_dof;
425:   nodal_fields[1].w_dof = field[k][j+1][i].w_dof;

427:   nodal_fields[2].u_dof = field[k][j+1][i+1].u_dof;
428:   nodal_fields[2].v_dof = field[k][j+1][i+1].v_dof;
429:   nodal_fields[2].w_dof = field[k][j+1][i+1].w_dof;

431:   nodal_fields[3].u_dof = field[k][j][i+1].u_dof;
432:   nodal_fields[3].v_dof = field[k][j][i+1].v_dof;
433:   nodal_fields[3].w_dof = field[k][j][i+1].w_dof;

435:   nodal_fields[4].u_dof = field[k+1][j][i].u_dof;
436:   nodal_fields[4].v_dof = field[k+1][j][i].v_dof;
437:   nodal_fields[4].w_dof = field[k+1][j][i].w_dof;

439:   nodal_fields[5].u_dof = field[k+1][j+1][i].u_dof;
440:   nodal_fields[5].v_dof = field[k+1][j+1][i].v_dof;
441:   nodal_fields[5].w_dof = field[k+1][j+1][i].w_dof;

443:   nodal_fields[6].u_dof = field[k+1][j+1][i+1].u_dof;
444:   nodal_fields[6].v_dof = field[k+1][j+1][i+1].v_dof;
445:   nodal_fields[6].w_dof = field[k+1][j+1][i+1].w_dof;

447:   nodal_fields[7].u_dof = field[k+1][j][i+1].u_dof;
448:   nodal_fields[7].v_dof = field[k+1][j][i+1].v_dof;
449:   nodal_fields[7].w_dof = field[k+1][j][i+1].w_dof;

451:   /* get the nodal fields for p */
452:   nodal_fields[0].p_dof = field[k][j][i].p_dof;
453:   nodal_fields[1].p_dof = field[k][j+1][i].p_dof;
454:   nodal_fields[2].p_dof = field[k][j+1][i+1].p_dof;
455:   nodal_fields[3].p_dof = field[k][j][i+1].p_dof;

457:   nodal_fields[4].p_dof = field[k+1][j][i].p_dof;
458:   nodal_fields[5].p_dof = field[k+1][j+1][i].p_dof;
459:   nodal_fields[6].p_dof = field[k+1][j+1][i+1].p_dof;
460:   nodal_fields[7].p_dof = field[k+1][j][i+1].p_dof;
461:   return(0);
462: }

464: 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)
465: {
466:   PetscInt              ij;
467:   PETSC_UNUSED PetscInt r,c,nr,nc;

469:   nr = w_NPE*w_dof;
470:   nc = u_NPE*u_dof;

472:   r = w_dof*wi+wd;
473:   c = u_dof*ui+ud;

475:   ij = r*nc+c;

477:   return ij;
478: }

482: static PetscErrorCode DMDASetValuesLocalStencil3D_ADD_VALUES(StokesDOF ***fields_F,MatStencil u_eqn[],MatStencil p_eqn[],PetscScalar Fe_u[],PetscScalar Fe_p[])
483: {
484:   PetscInt n,II,J,K;

487:   for (n = 0; n<NODES_PER_EL; n++) {
488:     II = u_eqn[NSD*n].i;
489:     J = u_eqn[NSD*n].j;
490:     K = u_eqn[NSD*n].k;

492:     fields_F[K][J][II].u_dof = fields_F[K][J][II].u_dof+Fe_u[NSD*n];

494:     II = u_eqn[NSD*n+1].i;
495:     J = u_eqn[NSD*n+1].j;
496:     K = u_eqn[NSD*n+1].k;

498:     fields_F[K][J][II].v_dof = fields_F[K][J][II].v_dof+Fe_u[NSD*n+1];

500:     II = u_eqn[NSD*n+2].i;
501:     J = u_eqn[NSD*n+2].j;
502:     K = u_eqn[NSD*n+2].k;
503:     fields_F[K][J][II].w_dof = fields_F[K][J][II].w_dof+Fe_u[NSD*n+2];

505:     II = p_eqn[n].i;
506:     J = p_eqn[n].j;
507:     K = p_eqn[n].k;

509:     fields_F[K][J][II].p_dof = fields_F[K][J][II].p_dof+Fe_p[n];

511:   }
512:   return(0);
513: }

515: static void FormStressOperatorQ13D(PetscScalar Ke[],PetscScalar coords[],PetscScalar eta[])
516: {
517:   PetscInt       ngp;
518:   PetscScalar    gp_xi[GAUSS_POINTS][NSD];
519:   PetscScalar    gp_weight[GAUSS_POINTS];
520:   PetscInt       p,i,j,k;
521:   PetscScalar    GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
522:   PetscScalar    J_p,tildeD[6];
523:   PetscScalar    B[6][U_DOFS*NODES_PER_EL];
524:   const PetscInt nvdof = U_DOFS*NODES_PER_EL;

526:   /* define quadrature rule */
527:   ConstructGaussQuadrature3D(&ngp,gp_xi,gp_weight);

529:   /* evaluate integral */
530:   for (p = 0; p < ngp; p++) {
531:     ShapeFunctionQ13D_Evaluate_dxi(gp_xi[p],GNi_p);
532:     ShapeFunctionQ13D_Evaluate_dx(GNi_p,GNx_p,coords,&J_p);

534:     for (i = 0; i < NODES_PER_EL; i++) {
535:       PetscScalar d_dx_i = GNx_p[0][i];
536:       PetscScalar d_dy_i = GNx_p[1][i];
537:       PetscScalar d_dz_i = GNx_p[2][i];

539:       B[0][3*i] = d_dx_i; B[0][3*i+1] = 0.0;     B[0][3*i+2] = 0.0;
540:       B[1][3*i] = 0.0;    B[1][3*i+1] = d_dy_i;  B[1][3*i+2] = 0.0;
541:       B[2][3*i] = 0.0;    B[2][3*i+1] = 0.0;     B[2][3*i+2] = d_dz_i;

543:       B[3][3*i] = d_dy_i; B[3][3*i+1] = d_dx_i;  B[3][3*i+2] = 0.0;   /* e_xy */
544:       B[4][3*i] = d_dz_i; B[4][3*i+1] = 0.0;     B[4][3*i+2] = d_dx_i; /* e_xz */
545:       B[5][3*i] = 0.0;    B[5][3*i+1] = d_dz_i;  B[5][3*i+2] = d_dy_i; /* e_yz */
546:     }


549:     tildeD[0] = 2.0*gp_weight[p]*J_p*eta[p];
550:     tildeD[1] = 2.0*gp_weight[p]*J_p*eta[p];
551:     tildeD[2] = 2.0*gp_weight[p]*J_p*eta[p];

553:     tildeD[3] =     gp_weight[p]*J_p*eta[p];
554:     tildeD[4] =     gp_weight[p]*J_p*eta[p];
555:     tildeD[5] =     gp_weight[p]*J_p*eta[p];

557:     /* form Bt tildeD B */
558:     /*
559:      Ke_ij = Bt_ik . D_kl . B_lj
560:      = B_ki . D_kl . B_lj
561:      = B_ki . D_kk . B_kj
562:      */
563:     for (i = 0; i < nvdof; i++) {
564:       for (j = i; j < nvdof; j++) {
565:         for (k = 0; k < 6; k++) {
566:           Ke[i*nvdof+j] += B[k][i]*tildeD[k]*B[k][j];
567:         }
568:       }
569:     }

571:   }
572:   /* fill lower triangular part */
573: #if defined(ASSEMBLE_LOWER_TRIANGULAR)
574:   for (i = 0; i < nvdof; i++) {
575:     for (j = i; j < nvdof; j++) {
576:       Ke[j*nvdof+i] = Ke[i*nvdof+j];
577:     }
578:   }
579: #endif
580: }

582: static void FormGradientOperatorQ13D(PetscScalar Ke[],PetscScalar coords[])
583: {
584:   PetscInt    ngp;
585:   PetscScalar gp_xi[GAUSS_POINTS][NSD];
586:   PetscScalar gp_weight[GAUSS_POINTS];
587:   PetscInt    p,i,j,di;
588:   PetscScalar Ni_p[NODES_PER_EL];
589:   PetscScalar GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
590:   PetscScalar J_p,fac;

592:   /* define quadrature rule */
593:   ConstructGaussQuadrature3D(&ngp,gp_xi,gp_weight);

595:   /* evaluate integral */
596:   for (p = 0; p < ngp; p++) {
597:     ShapeFunctionQ13D_Evaluate(gp_xi[p],Ni_p);
598:     ShapeFunctionQ13D_Evaluate_dxi(gp_xi[p],GNi_p);
599:     ShapeFunctionQ13D_Evaluate_dx(GNi_p,GNx_p,coords,&J_p);
600:     fac = gp_weight[p]*J_p;

602:     for (i = 0; i < NODES_PER_EL; i++) { /* u nodes */
603:       for (di = 0; di < NSD; di++) { /* u dofs */
604:         for (j = 0; j < NODES_PER_EL; j++) {  /* p nodes, p dofs = 1 (ie no loop) */
605:           PetscInt IJ;
606:           IJ = ASS_MAP_wIwDI_uJuDJ(i,di,NODES_PER_EL,3,j,0,NODES_PER_EL,1);

608:           Ke[IJ] = Ke[IJ]-GNx_p[di][i]*Ni_p[j]*fac;
609:         }
610:       }
611:     }
612:   }
613: }

615: static void FormDivergenceOperatorQ13D(PetscScalar De[],PetscScalar coords[])
616: {
617:   PetscScalar Ge[U_DOFS*NODES_PER_EL*P_DOFS*NODES_PER_EL];
618:   PetscInt    i,j;
619:   PetscInt    nr_g,nc_g;

621:   PetscMemzero(Ge,sizeof(PetscScalar)*U_DOFS*NODES_PER_EL*P_DOFS*NODES_PER_EL);
622:   FormGradientOperatorQ13D(Ge,coords);

624:   nr_g = U_DOFS*NODES_PER_EL;
625:   nc_g = P_DOFS*NODES_PER_EL;

627:   for (i = 0; i < nr_g; i++) {
628:     for (j = 0; j < nc_g; j++) {
629:       De[nr_g*j+i] = Ge[nc_g*i+j];
630:     }
631:   }
632: }

634: static void FormStabilisationOperatorQ13D(PetscScalar Ke[],PetscScalar coords[],PetscScalar eta[])
635: {
636:   PetscInt    ngp;
637:   PetscScalar gp_xi[GAUSS_POINTS][NSD];
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:   ConstructGaussQuadrature3D(&ngp,gp_xi,gp_weight);

647:   /* evaluate integral */
648:   for (p = 0; p < ngp; p++) {
649:     ShapeFunctionQ13D_Evaluate(gp_xi[p],Ni_p);
650:     ShapeFunctionQ13D_Evaluate_dxi(gp_xi[p],GNi_p);
651:     ShapeFunctionQ13D_Evaluate_dx(GNi_p,GNx_p,coords,&J_p);
652:     fac = gp_weight[p]*J_p;
653:     /*
654:      for (i = 0; i < NODES_PER_EL; i++) {
655:      for (j = i; j < NODES_PER_EL; j++) {
656:      Ke[NODES_PER_EL*i+j] -= fac*(Ni_p[i]*Ni_p[j]-0.015625);
657:      }
658:      }
659:      */

661:     for (i = 0; i < NODES_PER_EL; i++) {
662:       for (j = 0; j < NODES_PER_EL; j++) {
663:         Ke[NODES_PER_EL*i+j] -= fac*(Ni_p[i]*Ni_p[j]-0.015625);
664:       }
665:     }
666:   }

668:   /* scale */
669:   eta_avg = 0.0;
670:   for (p = 0; p < ngp; p++) eta_avg += eta[p];
671:   eta_avg = (1.0/((PetscScalar)ngp))*eta_avg;
672:   fac     = 1.0/eta_avg;

674:   /*
675:    for (i = 0; i < NODES_PER_EL; i++) {
676:    for (j = i; j < NODES_PER_EL; j++) {
677:    Ke[NODES_PER_EL*i+j] = fac*Ke[NODES_PER_EL*i+j];
678:    #if defined(ASSEMBLE_LOWER_TRIANGULAR)
679:    Ke[NODES_PER_EL*j+i] = Ke[NODES_PER_EL*i+j];
680:    #endif
681:    }
682:    }
683:    */

685:   for (i = 0; i < NODES_PER_EL; i++) {
686:     for (j = 0; j < NODES_PER_EL; j++) {
687:       Ke[NODES_PER_EL*i+j] = fac*Ke[NODES_PER_EL*i+j];
688:     }
689:   }
690: }

692: static void FormScaledMassMatrixOperatorQ13D(PetscScalar Ke[],PetscScalar coords[],PetscScalar eta[])
693: {
694:   PetscInt    ngp;
695:   PetscScalar gp_xi[GAUSS_POINTS][NSD];
696:   PetscScalar gp_weight[GAUSS_POINTS];
697:   PetscInt    p,i,j;
698:   PetscScalar Ni_p[NODES_PER_EL];
699:   PetscScalar GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
700:   PetscScalar J_p,fac,eta_avg;

702:   /* define quadrature rule */
703:   ConstructGaussQuadrature3D(&ngp,gp_xi,gp_weight);

705:   /* evaluate integral */
706:   for (p = 0; p < ngp; p++) {
707:     ShapeFunctionQ13D_Evaluate(gp_xi[p],Ni_p);
708:     ShapeFunctionQ13D_Evaluate_dxi(gp_xi[p],GNi_p);
709:     ShapeFunctionQ13D_Evaluate_dx(GNi_p,GNx_p,coords,&J_p);
710:     fac = gp_weight[p]*J_p;

712:     /*
713:      for (i = 0; i < NODES_PER_EL; i++) {
714:      for (j = i; j < NODES_PER_EL; j++) {
715:      Ke[NODES_PER_EL*i+j] = Ke[NODES_PER_EL*i+j]-fac*Ni_p[i]*Ni_p[j];
716:      }
717:      }
718:      */

720:     for (i = 0; i < NODES_PER_EL; i++) {
721:       for (j = 0; j < NODES_PER_EL; j++) {
722:         Ke[NODES_PER_EL*i+j] = Ke[NODES_PER_EL*i+j]-fac*Ni_p[i]*Ni_p[j];
723:       }
724:     }
725:   }

727:   /* scale */
728:   eta_avg = 0.0;
729:   for (p = 0; p < ngp; p++) eta_avg += eta[p];
730:   eta_avg = (1.0/((PetscScalar)ngp))*eta_avg;
731:   fac     = 1.0/eta_avg;
732:   /*
733:    for (i = 0; i < NODES_PER_EL; i++) {
734:    for (j = i; j < NODES_PER_EL; j++) {
735:    Ke[NODES_PER_EL*i+j] = fac*Ke[NODES_PER_EL*i+j];
736:    #if defined(ASSEMBLE_LOWER_TRIANGULAR)
737:    Ke[NODES_PER_EL*j+i] = Ke[NODES_PER_EL*i+j];
738:    #endif
739:    }
740:    }
741:    */

743:   for (i = 0; i < NODES_PER_EL; i++) {
744:     for (j = 0; j < NODES_PER_EL; j++) {
745:       Ke[NODES_PER_EL*i+j] = fac*Ke[NODES_PER_EL*i+j];
746:     }
747:   }
748: }

750: static void FormMomentumRhsQ13D(PetscScalar Fe[],PetscScalar coords[],PetscScalar fx[],PetscScalar fy[],PetscScalar fz[])
751: {
752:   PetscInt    ngp;
753:   PetscScalar gp_xi[GAUSS_POINTS][NSD];
754:   PetscScalar gp_weight[GAUSS_POINTS];
755:   PetscInt    p,i;
756:   PetscScalar Ni_p[NODES_PER_EL];
757:   PetscScalar GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
758:   PetscScalar J_p,fac;

760:   /* define quadrature rule */
761:   ConstructGaussQuadrature3D(&ngp,gp_xi,gp_weight);

763:   /* evaluate integral */
764:   for (p = 0; p < ngp; p++) {
765:     ShapeFunctionQ13D_Evaluate(gp_xi[p],Ni_p);
766:     ShapeFunctionQ13D_Evaluate_dxi(gp_xi[p],GNi_p);
767:     ShapeFunctionQ13D_Evaluate_dx(GNi_p,GNx_p,coords,&J_p);
768:     fac = gp_weight[p]*J_p;

770:     for (i = 0; i < NODES_PER_EL; i++) {
771:       Fe[NSD*i]   -= fac*Ni_p[i]*fx[p];
772:       Fe[NSD*i+1] -= fac*Ni_p[i]*fy[p];
773:       Fe[NSD*i+2] -= fac*Ni_p[i]*fz[p];
774:     }
775:   }
776: }

778: static void FormContinuityRhsQ13D(PetscScalar Fe[],PetscScalar coords[],PetscScalar hc[])
779: {
780:   PetscInt    ngp;
781:   PetscScalar gp_xi[GAUSS_POINTS][NSD];
782:   PetscScalar gp_weight[GAUSS_POINTS];
783:   PetscInt    p,i;
784:   PetscScalar Ni_p[NODES_PER_EL];
785:   PetscScalar GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
786:   PetscScalar J_p,fac;

788:   /* define quadrature rule */
789:   ConstructGaussQuadrature3D(&ngp,gp_xi,gp_weight);

791:   /* evaluate integral */
792:   for (p = 0; p < ngp; p++) {
793:     ShapeFunctionQ13D_Evaluate(gp_xi[p],Ni_p);
794:     ShapeFunctionQ13D_Evaluate_dxi(gp_xi[p],GNi_p);
795:     ShapeFunctionQ13D_Evaluate_dx(GNi_p,GNx_p,coords,&J_p);
796:     fac = gp_weight[p]*J_p;

798:     for (i = 0; i < NODES_PER_EL; i++) Fe[i] -= fac*Ni_p[i]*hc[p];
799:   }
800: }

802: #define _ZERO_ROWCOL_i(A,i) {                   \
803:     PetscInt    KK;                             \
804:     PetscScalar tmp = A[24*(i)+(i)];            \
805:     for (KK=0;KK<24;KK++) A[24*(i)+KK]=0.0;     \
806:     for (KK=0;KK<24;KK++) A[24*KK+(i)]=0.0;     \
807:     A[24*(i)+(i)] = tmp;}                       \

809: #define _ZERO_ROW_i(A,i) {                      \
810:     PetscInt KK;                                \
811:     for (KK=0;KK<8;KK++) A[8*(i)+KK]=0.0;}

813: #define _ZERO_COL_i(A,i) {                      \
814:     PetscInt KK;                                \
815:     for (KK=0;KK<8;KK++) A[24*KK+(i)]=0.0;}

819: static PetscErrorCode AssembleA_Stokes(Mat A,DM stokes_da,CellProperties cell_properties)
820: {
821:   DM                     cda;
822:   Vec                    coords;
823:   DMDACoor3d             ***_coords;
824:   MatStencil             u_eqn[NODES_PER_EL*U_DOFS];
825:   MatStencil             p_eqn[NODES_PER_EL*P_DOFS];
826:   PetscInt               sex,sey,sez,mx,my,mz;
827:   PetscInt               ei,ej,ek;
828:   PetscScalar            Ae[NODES_PER_EL*U_DOFS*NODES_PER_EL*U_DOFS];
829:   PetscScalar            Ge[NODES_PER_EL*U_DOFS*NODES_PER_EL*P_DOFS];
830:   PetscScalar            De[NODES_PER_EL*P_DOFS*NODES_PER_EL*U_DOFS];
831:   PetscScalar            Ce[NODES_PER_EL*P_DOFS*NODES_PER_EL*P_DOFS];
832:   PetscScalar            el_coords[NODES_PER_EL*NSD];
833:   GaussPointCoefficients *props;
834:   PetscScalar            *prop_eta;
835:   PetscInt               n,M,N,P;
836:   PetscErrorCode         ierr;

839:   DMDAGetInfo(stokes_da,0,&M,&N,&P,0,0,0, 0,0,0,0,0,0);
840:   /* setup for coords */
841:   DMGetCoordinateDM(stokes_da,&cda);
842:   DMGetCoordinatesLocal(stokes_da,&coords);
843:   DMDAVecGetArray(cda,coords,&_coords);

845:   DMDAGetElementCorners(stokes_da,&sex,&sey,&sez,&mx,&my,&mz);
846:   for (ek = sez; ek < sez+mz; ek++) {
847:     for (ej = sey; ej < sey+my; ej++) {
848:       for (ei = sex; ei < sex+mx; ei++) {
849:         /* get coords for the element */
850:         GetElementCoords3D(_coords,ei,ej,ek,el_coords);
851:         /* get cell properties */
852:         CellPropertiesGetCell(cell_properties,ei,ej,ek,&props);
853:         /* get coefficients for the element */
854:         prop_eta = props->eta;

856:         /* initialise element stiffness matrix */
857:         PetscMemzero(Ae,sizeof(PetscScalar)*NODES_PER_EL*U_DOFS*NODES_PER_EL*U_DOFS);
858:         PetscMemzero(Ge,sizeof(PetscScalar)*NODES_PER_EL*U_DOFS*NODES_PER_EL*P_DOFS);
859:         PetscMemzero(De,sizeof(PetscScalar)*NODES_PER_EL*P_DOFS*NODES_PER_EL*U_DOFS);
860:         PetscMemzero(Ce,sizeof(PetscScalar)*NODES_PER_EL*P_DOFS*NODES_PER_EL*P_DOFS);

862:         /* form element stiffness matrix */
863:         FormStressOperatorQ13D(Ae,el_coords,prop_eta);
864:         FormGradientOperatorQ13D(Ge,el_coords);
865:         /*#if defined(ASSEMBLE_LOWER_TRIANGULAR)*/
866:         FormDivergenceOperatorQ13D(De,el_coords);
867:         /*#endif*/
868:         FormStabilisationOperatorQ13D(Ce,el_coords,prop_eta);

870:         /* insert element matrix into global matrix */
871:         DMDAGetElementEqnums3D_up(u_eqn,p_eqn,ei,ej,ek);

873:         for (n=0; n<NODES_PER_EL; n++) {
874:           if ((u_eqn[3*n].i == 0) || (u_eqn[3*n].i == M-1)) {
875:             _ZERO_ROWCOL_i(Ae,3*n);
876:             _ZERO_ROW_i   (Ge,3*n);
877:             _ZERO_COL_i   (De,3*n);
878:           }

880:           if ((u_eqn[3*n+1].j == 0) || (u_eqn[3*n+1].j == N-1)) {
881:             _ZERO_ROWCOL_i(Ae,3*n+1);
882:             _ZERO_ROW_i   (Ge,3*n+1);
883:             _ZERO_COL_i   (De,3*n+1);
884:           }

886:           if ((u_eqn[3*n+2].k == 0) || (u_eqn[3*n+2].k == P-1)) {
887:             _ZERO_ROWCOL_i(Ae,3*n+2);
888:             _ZERO_ROW_i   (Ge,3*n+2);
889:             _ZERO_COL_i   (De,3*n+2);
890:           }
891:         }
892:         MatSetValuesStencil(A,NODES_PER_EL*U_DOFS,u_eqn,NODES_PER_EL*U_DOFS,u_eqn,Ae,ADD_VALUES);
893:         MatSetValuesStencil(A,NODES_PER_EL*U_DOFS,u_eqn,NODES_PER_EL*P_DOFS,p_eqn,Ge,ADD_VALUES);
894:         MatSetValuesStencil(A,NODES_PER_EL*P_DOFS,p_eqn,NODES_PER_EL*U_DOFS,u_eqn,De,ADD_VALUES);
895:         MatSetValuesStencil(A,NODES_PER_EL*P_DOFS,p_eqn,NODES_PER_EL*P_DOFS,p_eqn,Ce,ADD_VALUES);
896:       }
897:     }
898:   }
899:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
900:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

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

904:   return(0);
905: }

909: static PetscErrorCode AssembleA_PCStokes(Mat A,DM stokes_da,CellProperties cell_properties)
910: {
911:   DM                     cda;
912:   Vec                    coords;
913:   DMDACoor3d             ***_coords;
914:   MatStencil             u_eqn[NODES_PER_EL*U_DOFS];
915:   MatStencil             p_eqn[NODES_PER_EL*P_DOFS];
916:   PetscInt               sex,sey,sez,mx,my,mz;
917:   PetscInt               ei,ej,ek;
918:   PetscScalar            Ae[NODES_PER_EL*U_DOFS*NODES_PER_EL*U_DOFS];
919:   PetscScalar            Ge[NODES_PER_EL*U_DOFS*NODES_PER_EL*P_DOFS];
920:   PetscScalar            De[NODES_PER_EL*P_DOFS*NODES_PER_EL*U_DOFS];
921:   PetscScalar            Ce[NODES_PER_EL*P_DOFS*NODES_PER_EL*P_DOFS];
922:   PetscScalar            el_coords[NODES_PER_EL*NSD];
923:   GaussPointCoefficients *props;
924:   PetscScalar            *prop_eta;
925:   PetscInt               n,M,N,P;
926:   PetscErrorCode         ierr;

929:   DMDAGetInfo(stokes_da,0,&M,&N,&P,0,0,0, 0,0,0,0,0,0);
930:   /* setup for coords */
931:   DMGetCoordinateDM(stokes_da,&cda);
932:   DMGetCoordinatesLocal(stokes_da,&coords);
933:   DMDAVecGetArray(cda,coords,&_coords);

935:   DMDAGetElementCorners(stokes_da,&sex,&sey,&sez,&mx,&my,&mz);
936:   for (ek = sez; ek < sez+mz; ek++) {
937:     for (ej = sey; ej < sey+my; ej++) {
938:       for (ei = sex; ei < sex+mx; ei++) {
939:         /* get coords for the element */
940:         GetElementCoords3D(_coords,ei,ej,ek,el_coords);
941:         /* get cell properties */
942:         CellPropertiesGetCell(cell_properties,ei,ej,ek,&props);
943:         /* get coefficients for the element */
944:         prop_eta = props->eta;

946:         /* initialise element stiffness matrix */
947:         PetscMemzero(Ae,sizeof(PetscScalar)*NODES_PER_EL*U_DOFS*NODES_PER_EL*U_DOFS);
948:         PetscMemzero(Ge,sizeof(PetscScalar)*NODES_PER_EL*U_DOFS*NODES_PER_EL*P_DOFS);
949:         PetscMemzero(De,sizeof(PetscScalar)*NODES_PER_EL*P_DOFS*NODES_PER_EL*U_DOFS);
950:         PetscMemzero(Ce,sizeof(PetscScalar)*NODES_PER_EL*P_DOFS*NODES_PER_EL*P_DOFS);

952:         /* form element stiffness matrix */
953:         FormStressOperatorQ13D(Ae,el_coords,prop_eta);
954:         FormGradientOperatorQ13D(Ge,el_coords);
955:         /* FormDivergenceOperatorQ13D(De,el_coords); */
956:         FormScaledMassMatrixOperatorQ13D(Ce,el_coords,prop_eta);

958:         /* insert element matrix into global matrix */
959:         DMDAGetElementEqnums3D_up(u_eqn,p_eqn,ei,ej,ek);

961:         for (n=0; n<NODES_PER_EL; n++) {
962:           if ((u_eqn[3*n].i == 0) || (u_eqn[3*n].i == M-1)) {
963:             _ZERO_ROWCOL_i(Ae,3*n);
964:             _ZERO_ROW_i   (Ge,3*n);
965:             _ZERO_COL_i   (De,3*n);
966:           }

968:           if ((u_eqn[3*n+1].j == 0) || (u_eqn[3*n+1].j == N-1)) {
969:             _ZERO_ROWCOL_i(Ae,3*n+1);
970:             _ZERO_ROW_i   (Ge,3*n+1);
971:             _ZERO_COL_i   (De,3*n+1);
972:           }

974:           if ((u_eqn[3*n+2].k == 0) || (u_eqn[3*n+2].k == P-1)) {
975:             _ZERO_ROWCOL_i(Ae,3*n+2);
976:             _ZERO_ROW_i   (Ge,3*n+2);
977:             _ZERO_COL_i   (De,3*n+2);
978:           }
979:         }
980:         MatSetValuesStencil(A,NODES_PER_EL*U_DOFS,u_eqn,NODES_PER_EL*U_DOFS,u_eqn,Ae,ADD_VALUES);
981:         MatSetValuesStencil(A,NODES_PER_EL*U_DOFS,u_eqn,NODES_PER_EL*P_DOFS,p_eqn,Ge,ADD_VALUES);
982:         /*MatSetValuesStencil(A,NODES_PER_EL*P_DOFS,p_eqn,NODES_PER_EL*U_DOFS,u_eqn,De,ADD_VALUES);*/
983:         MatSetValuesStencil(A,NODES_PER_EL*P_DOFS,p_eqn,NODES_PER_EL*P_DOFS,p_eqn,Ce,ADD_VALUES);
984:       }
985:     }
986:   }
987:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
988:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

990:   DMDAVecRestoreArray(cda,coords,&_coords);
991:   return(0);
992: }

996: static PetscErrorCode AssembleF_Stokes(Vec F,DM stokes_da,CellProperties cell_properties)
997: {
998:   DM                     cda;
999:   Vec                    coords;
1000:   DMDACoor3d             ***_coords;
1001:   MatStencil             u_eqn[NODES_PER_EL*U_DOFS];
1002:   MatStencil             p_eqn[NODES_PER_EL*P_DOFS];
1003:   PetscInt               sex,sey,sez,mx,my,mz;
1004:   PetscInt               ei,ej,ek;
1005:   PetscScalar            Fe[NODES_PER_EL*U_DOFS];
1006:   PetscScalar            He[NODES_PER_EL*P_DOFS];
1007:   PetscScalar            el_coords[NODES_PER_EL*NSD];
1008:   GaussPointCoefficients *props;
1009:   PetscScalar            *prop_fx,*prop_fy,*prop_fz,*prop_hc;
1010:   Vec                    local_F;
1011:   StokesDOF              ***ff;
1012:   PetscInt               n,M,N,P;
1013:   PetscErrorCode         ierr;

1016:   DMDAGetInfo(stokes_da,0,&M,&N,&P,0,0,0, 0,0,0,0,0,0);
1017:   /* setup for coords */
1018:   DMGetCoordinateDM(stokes_da,&cda);
1019:   DMGetCoordinatesLocal(stokes_da,&coords);
1020:   DMDAVecGetArray(cda,coords,&_coords);

1022:   /* get acces to the vector */
1023:   DMGetLocalVector(stokes_da,&local_F);
1024:   VecZeroEntries(local_F);
1025:   DMDAVecGetArray(stokes_da,local_F,&ff);

1027:   DMDAGetElementCorners(stokes_da,&sex,&sey,&sez,&mx,&my,&mz);
1028:   for (ek = sez; ek < sez+mz; ek++) {
1029:     for (ej = sey; ej < sey+my; ej++) {
1030:       for (ei = sex; ei < sex+mx; ei++) {
1031:         /* get coords for the element */
1032:         GetElementCoords3D(_coords,ei,ej,ek,el_coords);
1033:         /* get cell properties */
1034:         CellPropertiesGetCell(cell_properties,ei,ej,ek,&props);
1035:         /* get coefficients for the element */
1036:         prop_fx = props->fx;
1037:         prop_fy = props->fy;
1038:         prop_fz = props->fz;
1039:         prop_hc = props->hc;

1041:         /* initialise element stiffness matrix */
1042:         PetscMemzero(Fe,sizeof(PetscScalar)*NODES_PER_EL*U_DOFS);
1043:         PetscMemzero(He,sizeof(PetscScalar)*NODES_PER_EL*P_DOFS);

1045:         /* form element stiffness matrix */
1046:         FormMomentumRhsQ13D(Fe,el_coords,prop_fx,prop_fy,prop_fz);
1047:         FormContinuityRhsQ13D(He,el_coords,prop_hc);

1049:         /* insert element matrix into global matrix */
1050:         DMDAGetElementEqnums3D_up(u_eqn,p_eqn,ei,ej,ek);

1052:         for (n=0; n<NODES_PER_EL; n++) {
1053:           if ((u_eqn[3*n].i == 0) || (u_eqn[3*n].i == M-1)) Fe[3*n] = 0.0;

1055:           if ((u_eqn[3*n+1].j == 0) || (u_eqn[3*n+1].j == N-1)) Fe[3*n+1] = 0.0;

1057:           if ((u_eqn[3*n+2].k == 0) || (u_eqn[3*n+2].k == P-1)) Fe[3*n+2] = 0.0;
1058:         }

1060:         DMDASetValuesLocalStencil3D_ADD_VALUES(ff,u_eqn,p_eqn,Fe,He);
1061:       }
1062:     }
1063:   }
1064:   DMDAVecRestoreArray(stokes_da,local_F,&ff);
1065:   DMLocalToGlobalBegin(stokes_da,local_F,ADD_VALUES,F);
1066:   DMLocalToGlobalEnd(stokes_da,local_F,ADD_VALUES,F);
1067:   DMRestoreLocalVector(stokes_da,&local_F);

1069:   DMDAVecRestoreArray(cda,coords,&_coords);
1070:   return(0);
1071: }

1073: static void evaluate_MS_FrankKamentski_constants(PetscReal *theta,PetscReal *MX,PetscReal *MY,PetscReal *MZ)
1074: {
1075:   *theta = 0.0;
1076:   *MX    = 2.0 * PETSC_PI;
1077:   *MY    = 2.0 * PETSC_PI;
1078:   *MZ    = 2.0 * PETSC_PI;
1079: }
1080: static void evaluate_MS_FrankKamentski(PetscReal pos[],PetscReal v[],PetscReal *p,PetscReal *eta,PetscReal Fm[],PetscReal *Fc)
1081: {
1082:   PetscReal x,y,z;
1083:   PetscReal theta,MX,MY,MZ;

1085:   evaluate_MS_FrankKamentski_constants(&theta,&MX,&MY,&MZ);
1086:   x = pos[0];
1087:   y = pos[1];
1088:   z = pos[2];
1089:   if (v) {
1090:     /*
1091:      v[0] = PetscPowRealInt(z,3)*PetscExpReal(y)*PetscSinReal(PETSC_PI*x);
1092:      v[1] = z*cos(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y);
1093:      v[2] = -(PetscPowRealInt(x,3) + PetscPowRealInt(y,3))*PetscSinReal(2.0*PETSC_PI*z);
1094:      */
1095:     /*
1096:      v[0] = PetscSinReal(PETSC_PI*x);
1097:      v[1] = PetscSinReal(PETSC_PI*y);
1098:      v[2] = PetscSinReal(PETSC_PI*z);
1099:      */
1100:     v[0] = PetscPowRealInt(z,3)*PetscExpReal(y)*PetscSinReal(PETSC_PI*x);
1101:     v[1] = z*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y);
1102:     v[2] = PetscPowRealInt(z,2)*(PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y)/2 - PETSC_PI*PetscCosReal(PETSC_PI*y)*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)/2) - PETSC_PI*PetscPowRealInt(z,4)*PetscCosReal(PETSC_PI*x)*PetscExpReal(y)/4;
1103:   }
1104:   if (p) *p = PetscPowRealInt(x,2) + PetscPowRealInt(y,2) + PetscPowRealInt(z,2);
1105:   if (eta) {
1106:     /**eta = PetscExpReal(-theta*(1.0 - y - 0.1*PetscCosReal(MX*x)*PetscCosReal(MZ*z)*PetscSinReal(MY*y)));*/
1107:     *eta = 1.0;
1108:   }
1109:   if (Fm) {
1110:     /*
1111:      Fm[0] = -2*x - 2.0*PetscPowRealInt(PETSC_PI,2)*PetscPowRealInt(z,3)*PetscExpReal(y)*PetscExpReal(-theta*(1 - y - 0.1*PetscCosReal(MX*x)*PetscCosReal(MZ*z)*PetscSinReal(MY*y)))*PetscSinReal(PETSC_PI*x) - 0.2*MZ*theta*(-1.5*PetscPowRealInt(x,2)*PetscSinReal(2.0*PETSC_PI*z) + 1.5*PetscPowRealInt(z,2)*PetscExpReal(y)*PetscSinReal(PETSC_PI*x))*PetscCosReal(MX*x)*PetscExpReal(-theta*(1 - y - 0.1*PetscCosReal(MX*x)*PetscCosReal(MZ*z)*PetscSinReal(MY*y)))*PetscSinReal(MY*y)*PetscSinReal(MZ*z) - 0.2*PETSC_PI*MX*theta*PetscPowRealInt(z,3)*PetscCosReal(PETSC_PI*x)*PetscCosReal(MZ*z)*PetscExpReal(y)*PetscExpReal(-theta*(1 - y - 0.1*PetscCosReal(MX*x)*PetscCosReal(MZ*z)*PetscSinReal(MY*y)))*PetscSinReal(MX*x)*PetscSinReal(MY*y) + 2.0*(3.0*z*PetscExpReal(y)*PetscSinReal(PETSC_PI*x) - 3.0*PETSC_PI*PetscPowRealInt(x,2)*PetscCosReal(2.0*PETSC_PI*z))*PetscExpReal(-theta*(1 - y - 0.1*PetscCosReal(MX*x)*PetscCosReal(MZ*z)*PetscSinReal(MY*y))) + 2.0*(0.5*PetscPowRealInt(z,3)*PetscExpReal(y)*PetscSinReal(PETSC_PI*x) + PETSC_PI*z*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y)*PetscSinReal(2.0*PETSC_PI*x) - 1.0*z*PetscPowRealInt(PETSC_PI,2)*PetscCosReal(PETSC_PI*y)*PetscExpReal(-y)*PetscSinReal(2.0*PETSC_PI*x))*PetscExpReal(-theta*(1 - y - 0.1*PetscCosReal(MX*x)*PetscCosReal(MZ*z)*PetscSinReal(MY*y))) + 2.0*theta*(1 + 0.1*MY*PetscCosReal(MX*x)*PetscCosReal(MY*y)*PetscCosReal(MZ*z))*(0.5*PetscPowRealInt(z,3)*PetscExpReal(y)*PetscSinReal(PETSC_PI*x) - 1.0*PETSC_PI*z*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y)*PetscSinReal(2.0*PETSC_PI*x))*PetscExpReal(-theta*(1 - y - 0.1*PetscCosReal(MX*x)*PetscCosReal(MZ*z)*PetscSinReal(MY*y))) ;
1112:      Fm[1] = -2*y - 0.2*MX*theta*(0.5*PetscPowRealInt(z,3)*PetscExpReal(y)*PetscSinReal(PETSC_PI*x) - 1.0*PETSC_PI*z*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y)*PetscSinReal(2.0*PETSC_PI*x))*PetscCosReal(MZ*z)*PetscExpReal(-theta*(1 - y - 0.1*PetscCosReal(MX*x)*PetscCosReal(MZ*z)*PetscSinReal(MY*y)))*PetscSinReal(MX*x)*PetscSinReal(MY*y) - 0.2*MZ*theta*(-1.5*PetscPowRealInt(y,2)*PetscSinReal(2.0*PETSC_PI*z) + 0.5*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y))*PetscCosReal(MX*x)*PetscExpReal(-theta*(1 - y - 0.1*PetscCosReal(MX*x)*PetscCosReal(MZ*z)*PetscSinReal(MY*y)))*PetscSinReal(MY*y)*PetscSinReal(MZ*z) + 2.0*(-2.0*z*PetscPowRealInt(PETSC_PI,2)*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y) + 0.5*PETSC_PI*PetscPowRealInt(z,3)*PetscCosReal(PETSC_PI*x)*PetscExpReal(y))*PetscExpReal(-theta*(1 - y - 0.1*PetscCosReal(MX*x)*PetscCosReal(MZ*z)*PetscSinReal(MY*y))) + 2.0*(z*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y) - z*PetscPowRealInt(PETSC_PI,2)*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y) - 2*PETSC_PI*z*PetscCosReal(PETSC_PI*y)*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y))*PetscExpReal(-theta*(1 - y - 0.1*PetscCosReal(MX*x)*PetscCosReal(MZ*z)*PetscSinReal(MY*y))) + 2.0*theta*(1 + 0.1*MY*PetscCosReal(MX*x)*PetscCosReal(MY*y)*PetscCosReal(MZ*z))*(-z*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y) + PETSC_PI*z*PetscCosReal(PETSC_PI*y)*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y))*PetscExpReal(-theta*(1 - y - 0.1*PetscCosReal(MX*x)*PetscCosReal(MZ*z)*PetscSinReal(MY*y))) - 6.0*PETSC_PI*PetscPowRealInt(y,2)*PetscCosReal(2.0*PETSC_PI*z)*PetscExpReal(-theta*(1 - y - 0.1*PetscCosReal(MX*x)*PetscCosReal(MZ*z)*PetscSinReal(MY*y)));
1113:      Fm[2] = -2*z + 8.0*PetscPowRealInt(PETSC_PI,2)*(PetscPowRealInt(x,3) + PetscPowRealInt(y,3))*PetscExpReal(-theta*(1 - y - 0.1*PetscCosReal(MX*x)*PetscCosReal(MZ*z)*PetscSinReal(MY*y)))*PetscSinReal(2.0*PETSC_PI*z) - 0.2*MX*theta*(-1.5*PetscPowRealInt(x,2)*PetscSinReal(2.0*PETSC_PI*z) + 1.5*PetscPowRealInt(z,2)*PetscExpReal(y)*PetscSinReal(PETSC_PI*x))*PetscCosReal(MZ*z)*PetscExpReal(-theta*(1 - y - 0.1*PetscCosReal(MX*x)*PetscCosReal(MZ*z)*PetscSinReal(MY*y)))*PetscSinReal(MX*x)*PetscSinReal(MY*y) + 0.4*PETSC_PI*MZ*theta*(PetscPowRealInt(x,3) + PetscPowRealInt(y,3))*PetscCosReal(MX*x)*PetscCosReal(2.0*PETSC_PI*z)*PetscExpReal(-theta*(1 - y - 0.1*PetscCosReal(MX*x)*PetscCosReal(MZ*z)*PetscSinReal(MY*y)))*PetscSinReal(MY*y)*PetscSinReal(MZ*z) + 2.0*(-3.0*x*PetscSinReal(2.0*PETSC_PI*z) + 1.5*PETSC_PI*PetscPowRealInt(z,2)*PetscCosReal(PETSC_PI*x)*PetscExpReal(y))*PetscExpReal(-theta*(1 - y - 0.1*PetscCosReal(MX*x)*PetscCosReal(MZ*z)*PetscSinReal(MY*y))) + 2.0*(-3.0*y*PetscSinReal(2.0*PETSC_PI*z) - 0.5*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y) + 0.5*PETSC_PI*PetscCosReal(PETSC_PI*y)*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y))*PetscExpReal(-theta*(1 - y - 0.1*PetscCosReal(MX*x)*PetscCosReal(MZ*z)*PetscSinReal(MY*y))) + 2.0*theta*(1 + 0.1*MY*PetscCosReal(MX*x)*PetscCosReal(MY*y)*PetscCosReal(MZ*z))*(-1.5*PetscPowRealInt(y,2)*PetscSinReal(2.0*PETSC_PI*z) + 0.5*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y))*PetscExpReal(-theta*(1 - y - 0.1*PetscCosReal(MX*x)*PetscCosReal(MZ*z)*PetscSinReal(MY*y)))  ;
1114:      */
1115:     /*
1116:      Fm[0]=-2*x - 2.0*PetscPowRealInt(PETSC_PI,2)*PetscSinReal(PETSC_PI*x);
1117:      Fm[1]=-2*y - 2.0*PetscPowRealInt(PETSC_PI,2)*PetscSinReal(PETSC_PI*y);
1118:      Fm[2]=-2*z - 2.0*PetscPowRealInt(PETSC_PI,2)*PetscSinReal(PETSC_PI*z);
1119:      */
1120:     /*
1121:      Fm[0] = -2*x + PetscPowRealInt(z,3)*PetscExpReal(y)*PetscSinReal(PETSC_PI*x) + 6.0*z*PetscExpReal(y)*PetscSinReal(PETSC_PI*x) - 6.0*PETSC_PI*PetscPowRealInt(x,2)*PetscCosReal(2.0*PETSC_PI*z) - 2.0*PetscPowRealInt(PETSC_PI,2)*PetscPowRealInt(z,3)*PetscExpReal(y)*PetscSinReal(PETSC_PI*x) + 2.0*PETSC_PI*z*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y)*PetscSinReal(2.0*PETSC_PI*x) - 2.0*z*PetscPowRealInt(PETSC_PI,2)*PetscCosReal(PETSC_PI*y)*PetscExpReal(-y)*PetscSinReal(2.0*PETSC_PI*x) ;
1122:      Fm[1] = -2*y - 6.0*PETSC_PI*PetscPowRealInt(y,2)*PetscCosReal(2.0*PETSC_PI*z) + 2.0*z*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y) - 6.0*z*PetscPowRealInt(PETSC_PI,2)*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y) + PETSC_PI*PetscPowRealInt(z,3)*PetscCosReal(PETSC_PI*x)*PetscExpReal(y) - 4.0*PETSC_PI*z*PetscCosReal(PETSC_PI*y)*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y);
1123:      Fm[2] = -2*z - 6.0*x*PetscSinReal(2.0*PETSC_PI*z) - 6.0*y*PetscSinReal(2.0*PETSC_PI*z) - PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y) + 8.0*PetscPowRealInt(PETSC_PI,2)*(PetscPowRealInt(x,3) + PetscPowRealInt(y,3))*PetscSinReal(2.0*PETSC_PI*z) + 3.0*PETSC_PI*PetscPowRealInt(z,2)*PetscCosReal(PETSC_PI*x)*PetscExpReal(y) + PETSC_PI*PetscCosReal(PETSC_PI*y)*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y) ;
1124:      */
1125:     Fm[0] = -2*x + 2*z*(PetscPowRealInt(PETSC_PI,2)*PetscCosReal(PETSC_PI*y)*PetscExpReal(-y)*PetscSinReal(2.0*PETSC_PI*x) - 1.0*PETSC_PI*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y)*PetscSinReal(2.0*PETSC_PI*x)) + PetscPowRealInt(z,3)*PetscExpReal(y)*PetscSinReal(PETSC_PI*x) + 6.0*z*PetscExpReal(y)*PetscSinReal(PETSC_PI*x) - 1.0*PetscPowRealInt(PETSC_PI,2)*PetscPowRealInt(z,3)*PetscExpReal(y)*PetscSinReal(PETSC_PI*x) + 2.0*PETSC_PI*z*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y)*PetscSinReal(2.0*PETSC_PI*x) - 2.0*z*PetscPowRealInt(PETSC_PI,2)*PetscCosReal(PETSC_PI*y)*PetscExpReal(-y)*PetscSinReal(2.0*PETSC_PI*x);
1126:     Fm[1] = -2*y + 2*z*(-PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y)/2 + PetscPowRealInt(PETSC_PI,2)*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y)/2 + PETSC_PI*PetscCosReal(PETSC_PI*y)*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)) + 2.0*z*PetscCosReal(2.0*PETSC_PI *x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y) - 6.0*z*PetscPowRealInt(PETSC_PI,2)*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y) - 4.0*PETSC_PI*z*PetscCosReal(PETSC_PI *y)*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y);
1127:     Fm[2] = -2*z + PetscPowRealInt(z,2)*(-2.0*PetscPowRealInt(PETSC_PI,2)*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y) + 2.0*PetscPowRealInt(PETSC_PI,3)*PetscCosReal(PETSC_PI*y)*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)) + PetscPowRealInt(z,2)*(PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y)/2 - 3*PetscPowRealInt(PETSC_PI,2)*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y)/2 + PetscPowRealInt(PETSC_PI,3)*PetscCosReal(PETSC_PI*y)*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)/2 - 3*PETSC_PI*PetscCosReal(PETSC_PI*y)*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)/2) + 1.0*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y) + 0.25*PetscPowRealInt(PETSC_PI,3)*PetscPowRealInt(z,4)*PetscCosReal(PETSC_PI*x)*PetscExpReal(y) - 0.25*PETSC_PI*PetscPowRealInt(z,4)*PetscCosReal(PETSC_PI*x)*PetscExpReal(y) - 3.0*PETSC_PI*PetscPowRealInt(z,2)*PetscCosReal(PETSC_PI*x)*PetscExpReal(y) - 1.0*PETSC_PI*PetscCosReal(PETSC_PI *y)*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y);
1128:   }
1129:   if (Fc) {
1130:     /**Fc = -2.0*PETSC_PI*(PetscPowRealInt(x,3) + PetscPowRealInt(y,3))*PetscCosReal(2.0*PETSC_PI*z) - z*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y) + PETSC_PI*PetscPowRealInt(z,3)*PetscCosReal(PETSC_PI*x)*PetscExpReal(y) + PETSC_PI*z*PetscCosReal(PETSC_PI*y)*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y) ;*/
1131:     /**Fc = PETSC_PI*PetscCosReal(PETSC_PI*x) + PETSC_PI*PetscCosReal(PETSC_PI*y) + PETSC_PI*PetscCosReal(PETSC_PI*z);*/
1132:     /**Fc = -2.0*PETSC_PI*(PetscPowRealInt(x,3) + PetscPowRealInt(y,3))*PetscCosReal(2.0*PETSC_PI*z) - z*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y) + PETSC_PI*PetscPowRealInt(z,3)*PetscCosReal(PETSC_PI*x)*PetscExpReal(y) + PETSC_PI*z*PetscCosReal(PETSC_PI*y)*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y);*/
1133:     *Fc = 0.0;
1134:   }
1135: }

1139: static PetscErrorCode DMDACreateManufacturedSolution(PetscInt mx,PetscInt my,PetscInt mz,DM *_da,Vec *_X)
1140: {
1141:   DM             da,cda;
1142:   Vec            X,local_X;
1143:   StokesDOF      ***_stokes;
1144:   Vec            coords;
1145:   DMDACoor3d     ***_coords;
1146:   PetscInt       si,sj,sk,ei,ej,ek,i,j,k;

1150:   DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,
1151:                       mx+1,my+1,mz+1,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,4,1,NULL,NULL,NULL,&da);
1152:   DMDASetFieldName(da,0,"anlytic_Vx");
1153:   DMDASetFieldName(da,1,"anlytic_Vy");
1154:   DMDASetFieldName(da,2,"anlytic_Vz");
1155:   DMDASetFieldName(da,3,"analytic_P");

1157:   DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,1.0);

1159:   DMGetCoordinatesLocal(da,&coords);
1160:   DMGetCoordinateDM(da,&cda);
1161:   DMDAVecGetArray(cda,coords,&_coords);

1163:   DMCreateGlobalVector(da,&X);
1164:   DMCreateLocalVector(da,&local_X);
1165:   DMDAVecGetArray(da,local_X,&_stokes);

1167:   DMDAGetGhostCorners(da,&si,&sj,&sk,&ei,&ej,&ek);
1168:   for (k = sk; k < sk+ek; k++) {
1169:     for (j = sj; j < sj+ej; j++) {
1170:       for (i = si; i < si+ei; i++) {
1171:         PetscReal pos[NSD],pressure,vel[NSD];

1173:         pos[0] = PetscRealPart(_coords[k][j][i].x);
1174:         pos[1] = PetscRealPart(_coords[k][j][i].y);
1175:         pos[2] = PetscRealPart(_coords[k][j][i].z);

1177:         evaluate_MS_FrankKamentski(pos,vel,&pressure,NULL,NULL,NULL);

1179:         _stokes[k][j][i].u_dof = vel[0];
1180:         _stokes[k][j][i].v_dof = vel[1];
1181:         _stokes[k][j][i].w_dof = vel[2];
1182:         _stokes[k][j][i].p_dof = pressure;
1183:       }
1184:     }
1185:   }
1186:   DMDAVecRestoreArray(da,local_X,&_stokes);
1187:   DMDAVecRestoreArray(cda,coords,&_coords);

1189:   DMLocalToGlobalBegin(da,local_X,INSERT_VALUES,X);
1190:   DMLocalToGlobalEnd(da,local_X,INSERT_VALUES,X);

1192:   VecDestroy(&local_X);

1194:   *_da = da;
1195:   *_X  = X;
1196:   return(0);
1197: }

1201: static PetscErrorCode DMDAIntegrateErrors3D(DM stokes_da,Vec X,Vec X_analytic)
1202: {
1203:   DM          cda;
1204:   Vec         coords,X_analytic_local,X_local;
1205:   DMDACoor3d  ***_coords;
1206:   PetscInt    sex,sey,sez,mx,my,mz;
1207:   PetscInt    ei,ej,ek;
1208:   PetscScalar el_coords[NODES_PER_EL*NSD];
1209:   StokesDOF   ***stokes_analytic,***stokes;
1210:   StokesDOF   stokes_analytic_e[NODES_PER_EL],stokes_e[NODES_PER_EL];

1212:   PetscScalar    GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
1213:   PetscScalar    Ni_p[NODES_PER_EL];
1214:   PetscInt       ngp;
1215:   PetscScalar    gp_xi[GAUSS_POINTS][NSD];
1216:   PetscScalar    gp_weight[GAUSS_POINTS];
1217:   PetscInt       p,i;
1218:   PetscScalar    J_p,fac;
1219:   PetscScalar    h,p_e_L2,u_e_L2,u_e_H1,p_L2,u_L2,u_H1,tp_L2,tu_L2,tu_H1;
1220:   PetscScalar    tint_p_ms,tint_p,int_p_ms,int_p;
1221:   PetscInt       M;
1222:   PetscReal      xymin[NSD],xymax[NSD];

1226:   /* define quadrature rule */
1227:   ConstructGaussQuadrature3D(&ngp,gp_xi,gp_weight);

1229:   /* setup for coords */
1230:   DMGetCoordinateDM(stokes_da,&cda);
1231:   DMGetCoordinatesLocal(stokes_da,&coords);
1232:   DMDAVecGetArray(cda,coords,&_coords);

1234:   /* setup for analytic */
1235:   DMCreateLocalVector(stokes_da,&X_analytic_local);
1236:   DMGlobalToLocalBegin(stokes_da,X_analytic,INSERT_VALUES,X_analytic_local);
1237:   DMGlobalToLocalEnd(stokes_da,X_analytic,INSERT_VALUES,X_analytic_local);
1238:   DMDAVecGetArray(stokes_da,X_analytic_local,&stokes_analytic);

1240:   /* setup for solution */
1241:   DMCreateLocalVector(stokes_da,&X_local);
1242:   DMGlobalToLocalBegin(stokes_da,X,INSERT_VALUES,X_local);
1243:   DMGlobalToLocalEnd(stokes_da,X,INSERT_VALUES,X_local);
1244:   DMDAVecGetArray(stokes_da,X_local,&stokes);

1246:   DMDAGetInfo(stokes_da,0,&M,0,0,0,0,0,0,0,0,0,0,0);
1247:   DMDAGetBoundingBox(stokes_da,xymin,xymax);

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

1251:   tp_L2     = tu_L2 = tu_H1 = 0.0;
1252:   tint_p_ms = tint_p = 0.0;

1254:   DMDAGetElementCorners(stokes_da,&sex,&sey,&sez,&mx,&my,&mz);

1256:   for (ek = sez; ek < sez+mz; ek++) {
1257:     for (ej = sey; ej < sey+my; ej++) {
1258:       for (ei = sex; ei < sex+mx; ei++) {
1259:         /* get coords for the element */
1260:         GetElementCoords3D(_coords,ei,ej,ek,el_coords);
1261:         StokesDAGetNodalFields3D(stokes,ei,ej,ek,stokes_e);
1262:         StokesDAGetNodalFields3D(stokes_analytic,ei,ej,ek,stokes_analytic_e);

1264:         /* evaluate integral */
1265:         for (p = 0; p < ngp; p++) {
1266:           ShapeFunctionQ13D_Evaluate(gp_xi[p],Ni_p);
1267:           ShapeFunctionQ13D_Evaluate_dxi(gp_xi[p],GNi_p);
1268:           ShapeFunctionQ13D_Evaluate_dx(GNi_p,GNx_p,el_coords,&J_p);
1269:           fac = gp_weight[p]*J_p;

1271:           for (i = 0; i < NODES_PER_EL; i++) {
1272:             tint_p_ms = tint_p_ms+fac*Ni_p[i]*stokes_analytic_e[i].p_dof;
1273:             tint_p    = tint_p   +fac*Ni_p[i]*stokes_e[i].p_dof;
1274:           }
1275:         }
1276:       }
1277:     }
1278:   }
1279:   MPI_Allreduce(&tint_p_ms,&int_p_ms,1,MPIU_SCALAR,MPIU_SUM,PETSC_COMM_WORLD);
1280:   MPI_Allreduce(&tint_p,&int_p,1,MPIU_SCALAR,MPIU_SUM,PETSC_COMM_WORLD);

1282:   PetscPrintf(PETSC_COMM_WORLD,"\\int P dv %1.4e (h)  %1.4e (ms)\n",PetscRealPart(int_p),PetscRealPart(int_p_ms));

1284:   /* remove mine and add manufacture one */
1285:   DMDAVecRestoreArray(stokes_da,X_analytic_local,&stokes_analytic);
1286:   DMDAVecRestoreArray(stokes_da,X_local,&stokes);

1288:   {
1289:     PetscInt    k,L,dof;
1290:     PetscScalar *fields;
1291:     DMDAGetInfo(stokes_da,0, 0,0,0, 0,0,0, &dof,0,0,0,0,0);

1293:     VecGetLocalSize(X_local,&L);
1294:     VecGetArray(X_local,&fields);
1295:     for (k=0; k<L/dof; k++) fields[dof*k+3] = fields[dof*k+3] - int_p + int_p_ms;
1296:     VecRestoreArray(X_local,&fields);

1298:     VecGetLocalSize(X,&L);
1299:     VecGetArray(X,&fields);
1300:     for (k=0; k<L/dof; k++) fields[dof*k+3] = fields[dof*k+3] - int_p + int_p_ms;
1301:     VecRestoreArray(X,&fields);
1302:   }

1304:   DMDAVecGetArray(stokes_da,X_local,&stokes);
1305:   DMDAVecGetArray(stokes_da,X_analytic_local,&stokes_analytic);

1307:   for (ek = sez; ek < sez+mz; ek++) {
1308:     for (ej = sey; ej < sey+my; ej++) {
1309:       for (ei = sex; ei < sex+mx; ei++) {
1310:         /* get coords for the element */
1311:         GetElementCoords3D(_coords,ei,ej,ek,el_coords);
1312:         StokesDAGetNodalFields3D(stokes,ei,ej,ek,stokes_e);
1313:         StokesDAGetNodalFields3D(stokes_analytic,ei,ej,ek,stokes_analytic_e);

1315:         /* evaluate integral */
1316:         p_e_L2 = 0.0;
1317:         u_e_L2 = 0.0;
1318:         u_e_H1 = 0.0;
1319:         for (p = 0; p < ngp; p++) {
1320:           ShapeFunctionQ13D_Evaluate(gp_xi[p],Ni_p);
1321:           ShapeFunctionQ13D_Evaluate_dxi(gp_xi[p],GNi_p);
1322:           ShapeFunctionQ13D_Evaluate_dx(GNi_p,GNx_p,el_coords,&J_p);
1323:           fac = gp_weight[p]*J_p;

1325:           for (i = 0; i < NODES_PER_EL; i++) {
1326:             PetscScalar u_error,v_error,w_error;

1328:             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);

1330:             u_error = stokes_e[i].u_dof-stokes_analytic_e[i].u_dof;
1331:             v_error = stokes_e[i].v_dof-stokes_analytic_e[i].v_dof;
1332:             w_error = stokes_e[i].w_dof-stokes_analytic_e[i].w_dof;
1333:             /*
1334:              if (p==0) {
1335:              printf("p=0: %d %d %d %1.4e,%1.4e,%1.4e \n", ei,ej,ek,u_error,v_error,w_error);
1336:              }
1337:              */
1338:             u_e_L2 += fac*Ni_p[i]*(u_error*u_error+v_error*v_error+w_error*w_error);

1340:             u_e_H1 = u_e_H1+fac*(GNx_p[0][i]*u_error*GNx_p[0][i]*u_error              /* du/dx */
1341:                                  +GNx_p[1][i]*u_error*GNx_p[1][i]*u_error
1342:                                  +GNx_p[2][i]*u_error*GNx_p[2][i]*u_error
1343:                                  +GNx_p[0][i]*v_error*GNx_p[0][i]*v_error               /* dv/dx */
1344:                                  +GNx_p[1][i]*v_error*GNx_p[1][i]*v_error
1345:                                  +GNx_p[2][i]*v_error*GNx_p[2][i]*v_error
1346:                                  +GNx_p[0][i]*w_error*GNx_p[0][i]*w_error               /* dw/dx */
1347:                                  +GNx_p[1][i]*w_error*GNx_p[1][i]*w_error
1348:                                  +GNx_p[2][i]*w_error*GNx_p[2][i]*w_error);
1349:           }
1350:         }

1352:         tp_L2 += p_e_L2;
1353:         tu_L2 += u_e_L2;
1354:         tu_H1 += u_e_H1;
1355:       }
1356:     }
1357:   }
1358:   MPI_Allreduce(&tp_L2,&p_L2,1,MPIU_SCALAR,MPIU_SUM,PETSC_COMM_WORLD);
1359:   MPI_Allreduce(&tu_L2,&u_L2,1,MPIU_SCALAR,MPIU_SUM,PETSC_COMM_WORLD);
1360:   MPI_Allreduce(&tu_H1,&u_H1,1,MPIU_SCALAR,MPIU_SUM,PETSC_COMM_WORLD);
1361:   p_L2 = PetscSqrtScalar(p_L2);
1362:   u_L2 = PetscSqrtScalar(u_L2);
1363:   u_H1 = PetscSqrtScalar(u_H1);

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


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

1370:   DMDAVecRestoreArray(stokes_da,X_analytic_local,&stokes_analytic);
1371:   VecDestroy(&X_analytic_local);
1372:   DMDAVecRestoreArray(stokes_da,X_local,&stokes);
1373:   VecDestroy(&X_local);
1374:   return(0);
1375: }

1379: PetscErrorCode DAView_3DVTK_StructuredGrid_appended(DM da,Vec FIELD,const char file_prefix[])
1380: {
1381:   char           vtk_filename[PETSC_MAX_PATH_LEN];
1382:   PetscMPIInt    rank;
1383:   MPI_Comm       comm;
1384:   FILE           *vtk_fp = NULL;
1385:   PetscInt       si,sj,sk,nx,ny,nz,i;
1386:   PetscInt       f,n_fields,N;
1387:   DM             cda;
1388:   Vec            coords;
1389:   Vec            l_FIELD;
1390:   PetscScalar    *_L_FIELD;
1391:   PetscInt       memory_offset;
1392:   PetscScalar    *buffer;


1397:   /* create file name */
1398:   PetscObjectGetComm((PetscObject)da,&comm);
1399:   MPI_Comm_rank(comm,&rank);
1400:   PetscSNPrintf(vtk_filename,sizeof(vtk_filename),"subdomain-%s-p%1.4d.vts",file_prefix,rank);

1402:   /* open file and write header */
1403:   vtk_fp = fopen(vtk_filename,"w");
1404:   if (!vtk_fp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SYS,"Cannot open file = %s \n",vtk_filename);

1406:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"<?xml version=\"1.0\"?>\n");

1408:   /* coords */
1409:   DMDAGetGhostCorners(da,&si,&sj,&sk,&nx,&ny,&nz);
1410:   N    = nx * ny * nz;

1412: #if defined(PETSC_WORDS_BIGENDIAN)
1413:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"BigEndian\">\n");
1414: #else
1415:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");
1416: #endif
1417:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"  <StructuredGrid WholeExtent=\"%D %D %D %D %D %D\">\n",si,si+nx-1,sj,sj+ny-1,sk,sk+nz-1);
1418:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"    <Piece Extent=\"%D %D %D %D %D %D\">\n",si,si+nx-1,sj,sj+ny-1,sk,sk+nz-1);

1420:   memory_offset = 0;

1422:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"      <CellData></CellData>\n");

1424:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"      <Points>\n");

1426:   /* copy coordinates */
1427:   DMGetCoordinateDM(da,&cda);
1428:   DMGetCoordinatesLocal(da,&coords);
1429:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"        <DataArray type=\"Float64\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%d\" />\n",memory_offset);
1430:   memory_offset = memory_offset + sizeof(PetscInt) + sizeof(PetscScalar)*N*3;

1432:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"      </Points>\n");

1434:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"      <PointData Scalars=\" ");
1435:   DMDAGetInfo(da,0,0,0,0,0,0,0,&n_fields,0,0,0,0,0);
1436:   for (f=0; f<n_fields; f++) {
1437:     const char *field_name;
1438:     DMDAGetFieldName(da,f,&field_name);
1439:     PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"%s ",field_name);
1440:   }
1441:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"\">\n");

1443:   for (f=0; f<n_fields; f++) {
1444:     const char *field_name;

1446:     DMDAGetFieldName(da,f,&field_name);
1447:     PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"        <DataArray type=\"Float64\" Name=\"%s\" format=\"appended\" offset=\"%d\"/>\n", field_name,memory_offset);
1448:     memory_offset = memory_offset + sizeof(PetscInt) + sizeof(PetscScalar)*N;
1449:   }

1451:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"      </PointData>\n");
1452:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"    </Piece>\n");
1453:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"  </StructuredGrid>\n");

1455:   PetscMalloc(sizeof(PetscScalar)*N,&buffer);
1456:   DMGetLocalVector(da,&l_FIELD);
1457:   DMGlobalToLocalBegin(da, FIELD,INSERT_VALUES,l_FIELD);
1458:   DMGlobalToLocalEnd(da,FIELD,INSERT_VALUES,l_FIELD);
1459:   VecGetArray(l_FIELD,&_L_FIELD);

1461:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"  <AppendedData encoding=\"raw\">\n");
1462:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"_");

1464:   /* write coordinates */
1465:   {
1466:     int         length = sizeof(PetscScalar)*N*3;
1467:     PetscScalar *allcoords;

1469:     fwrite(&length,sizeof(int),1,vtk_fp);
1470:     VecGetArray(coords,&allcoords);
1471:     fwrite(allcoords,sizeof(PetscScalar),3*N,vtk_fp);
1472:     VecRestoreArray(coords,&allcoords);
1473:   }
1474:   /* write fields */
1475:   for (f=0; f<n_fields; f++) {
1476:     int length = sizeof(PetscScalar)*N;
1477:     fwrite(&length,sizeof(int),1,vtk_fp);
1478:     /* load */
1479:     for (i=0; i<N; i++) buffer[i] = _L_FIELD[n_fields*i + f];

1481:     /* write */
1482:     fwrite(buffer,sizeof(PetscScalar),N,vtk_fp);
1483:   }
1484:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"\n  </AppendedData>\n");

1486:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"</VTKFile>\n");

1488:   PetscFree(buffer);
1489:   VecRestoreArray(l_FIELD,&_L_FIELD);
1490:   DMRestoreLocalVector(da,&l_FIELD);

1492:   if (vtk_fp) {
1493:     fclose(vtk_fp);
1494:     vtk_fp = NULL;
1495:   }

1497:   return(0);
1498: }

1502: PetscErrorCode DAViewVTK_write_PieceExtend(FILE *vtk_fp,PetscInt indent_level,DM da,const char local_file_prefix[])
1503: {
1504:   PetscMPIInt    size,rank;
1505:   MPI_Comm       comm;
1506:   const PetscInt *lx,*ly,*lz;
1507:   PetscInt       M,N,P,pM,pN,pP,sum,*olx,*oly,*olz;
1508:   PetscInt       *osx,*osy,*osz,*oex,*oey,*oez;
1509:   PetscInt       i,j,k,II,stencil;

1513:   /* create file name */
1514:   PetscObjectGetComm((PetscObject)da,&comm);
1515:   MPI_Comm_size(comm,&size);
1516:   MPI_Comm_rank(comm,&rank);

1518:   DMDAGetInfo(da,0,&M,&N,&P,&pM,&pN,&pP,0,&stencil,0,0,0,0);
1519:   DMDAGetOwnershipRanges(da,&lx,&ly,&lz);

1521:   /* generate start,end list */
1522:   PetscMalloc(sizeof(PetscInt)*(pM+1),&olx);
1523:   PetscMalloc(sizeof(PetscInt)*(pN+1),&oly);
1524:   PetscMalloc(sizeof(PetscInt)*(pP+1),&olz);
1525:   sum  = 0;
1526:   for (i=0; i<pM; i++) {
1527:     olx[i] = sum;
1528:     sum    = sum + lx[i];
1529:   }
1530:   olx[pM] = sum;
1531:   sum     = 0;
1532:   for (i=0; i<pN; i++) {
1533:     oly[i] = sum;
1534:     sum    = sum + ly[i];
1535:   }
1536:   oly[pN] = sum;
1537:   sum     = 0;
1538:   for (i=0; i<pP; i++) {
1539:     olz[i] = sum;
1540:     sum    = sum + lz[i];
1541:   }
1542:   olz[pP] = sum;

1544:   PetscMalloc(sizeof(PetscInt)*(pM),&osx);
1545:   PetscMalloc(sizeof(PetscInt)*(pN),&osy);
1546:   PetscMalloc(sizeof(PetscInt)*(pP),&osz);
1547:   PetscMalloc(sizeof(PetscInt)*(pM),&oex);
1548:   PetscMalloc(sizeof(PetscInt)*(pN),&oey);
1549:   PetscMalloc(sizeof(PetscInt)*(pP),&oez);
1550:   for (i=0; i<pM; i++) {
1551:     osx[i] = olx[i] - stencil;
1552:     oex[i] = olx[i] + lx[i] + stencil;
1553:     if (osx[i]<0) osx[i]=0;
1554:     if (oex[i]>M) oex[i]=M;
1555:   }

1557:   for (i=0; i<pN; i++) {
1558:     osy[i] = oly[i] - stencil;
1559:     oey[i] = oly[i] + ly[i] + stencil;
1560:     if (osy[i]<0)osy[i]=0;
1561:     if (oey[i]>M)oey[i]=N;
1562:   }
1563:   for (i=0; i<pP; i++) {
1564:     osz[i] = olz[i] - stencil;
1565:     oez[i] = olz[i] + lz[i] + stencil;
1566:     if (osz[i]<0) osz[i]=0;
1567:     if (oez[i]>P) oez[i]=P;
1568:   }

1570:   for (k=0; k<pP; k++) {
1571:     for (j=0; j<pN; j++) {
1572:       for (i=0; i<pM; i++) {
1573:         char     name[PETSC_MAX_PATH_LEN];
1574:         PetscInt procid = i + j*pM + k*pM*pN; /* convert proc(i,j,k) to pid */
1575:         PetscSNPrintf(name,sizeof(name),"subdomain-%s-p%1.4d.vts",local_file_prefix,procid);
1576:         for (II=0; II<indent_level; II++) PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"  ");

1578:         PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"<Piece Extent=\"%d %d %d %d %d %d\"      Source=\"%s\"/>\n",
1579:                      osx[i],oex[i]-1,
1580:                      osy[j],oey[j]-1,
1581:                      osz[k],oez[k]-1,name);
1582:       }
1583:     }
1584:   }
1585:   PetscFree(olx);
1586:   PetscFree(oly);
1587:   PetscFree(olz);
1588:   PetscFree(osx);
1589:   PetscFree(osy);
1590:   PetscFree(osz);
1591:   PetscFree(oex);
1592:   PetscFree(oey);
1593:   PetscFree(oez);
1594:   return(0);
1595: }

1599: PetscErrorCode DAView_3DVTK_PStructuredGrid(DM da,const char file_prefix[],const char local_file_prefix[])
1600: {
1601:   MPI_Comm       comm;
1602:   PetscMPIInt    size,rank;
1603:   char           vtk_filename[PETSC_MAX_PATH_LEN];
1604:   FILE           *vtk_fp = NULL;
1605:   PetscInt       M,N,P,si,sj,sk,nx,ny,nz;
1606:   PetscInt       i,dofs;

1610:   /* only master generates this file */
1611:   PetscObjectGetComm((PetscObject)da,&comm);
1612:   MPI_Comm_size(comm,&size);
1613:   MPI_Comm_rank(comm,&rank);

1615:   if (rank != 0) return(0);

1617:   /* create file name */
1618:   PetscSNPrintf(vtk_filename,sizeof(vtk_filename),"%s.pvts",file_prefix);
1619:   vtk_fp = fopen(vtk_filename,"w");
1620:   if (!vtk_fp) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SYS,"Cannot open file = %s \n",vtk_filename);

1622:   /* (VTK) generate pvts header */
1623:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"<?xml version=\"1.0\"?>\n");

1625: #if defined(PETSC_WORDS_BIGENDIAN)
1626:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"<VTKFile type=\"PStructuredGrid\" version=\"0.1\" byte_order=\"BigEndian\">\n");
1627: #else
1628:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"<VTKFile type=\"PStructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");
1629: #endif

1631:   /* define size of the nodal mesh based on the cell DM */
1632:   DMDAGetInfo(da,0,&M,&N,&P,0,0,0,&dofs,0,0,0,0,0);
1633:   DMDAGetGhostCorners(da,&si,&sj,&sk,&nx,&ny,&nz);
1634:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"  <PStructuredGrid GhostLevel=\"1\" WholeExtent=\"%d %d %d %d %d %d\">\n",0,M-1,0,N-1,0,P-1); /* note overlap = 1 for Q1 */

1636:   /* DUMP THE CELL REFERENCES */
1637:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"    <PCellData>\n");
1638:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"    </PCellData>\n");

1640:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"    <PPoints>\n");
1641:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"      <PDataArray type=\"Float64\" Name=\"Points\" NumberOfComponents=\"%d\"/>\n",NSD);
1642:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"    </PPoints>\n");

1644:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"    <PPointData>\n");
1645:   for (i=0; i<dofs; i++) {
1646:     const char *fieldname;
1647:     DMDAGetFieldName(da,i,&fieldname);
1648:     PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"      <PDataArray type=\"Float64\" Name=\"%s\" NumberOfComponents=\"1\"/>\n",fieldname);
1649:   }
1650:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"    </PPointData>\n");

1652:   /* write out the parallel information */
1653:   DAViewVTK_write_PieceExtend(vtk_fp,2,da,local_file_prefix);

1655:   /* close the file */
1656:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"  </PStructuredGrid>\n");
1657:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"</VTKFile>\n");

1659:   if (vtk_fp) {
1660:     fclose(vtk_fp);
1661:     vtk_fp = NULL;
1662:   }
1663:   return(0);
1664: }

1668: PetscErrorCode DAView3DPVTS(DM da, Vec x,const char NAME[])
1669: {
1670:   char           vts_filename[PETSC_MAX_PATH_LEN];
1671:   char           pvts_filename[PETSC_MAX_PATH_LEN];

1675:   PetscSNPrintf(vts_filename,sizeof(vts_filename),"%s-mesh",NAME);
1676:   DAView_3DVTK_StructuredGrid_appended(da,x,vts_filename);

1678:   PetscSNPrintf(pvts_filename,sizeof(pvts_filename),"%s-mesh",NAME);
1679:   DAView_3DVTK_PStructuredGrid(da,pvts_filename,vts_filename);
1680:   return(0);
1681: }

1685: PetscErrorCode KSPMonitorStokesBlocks(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy)
1686: {
1688:   PetscReal      norms[4];
1689:   Vec            Br,v,w;
1690:   Mat            A;

1693:   KSPGetOperators(ksp,&A,NULL);
1694:   MatGetVecs(A,&w,&v);

1696:   KSPBuildResidual(ksp,v,w,&Br);

1698:   VecStrideNorm(Br,0,NORM_2,&norms[0]);
1699:   VecStrideNorm(Br,1,NORM_2,&norms[1]);
1700:   VecStrideNorm(Br,2,NORM_2,&norms[2]);
1701:   VecStrideNorm(Br,3,NORM_2,&norms[3]);

1703:   VecDestroy(&v);
1704:   VecDestroy(&w);

1706:   PetscPrintf(PETSC_COMM_WORLD,"%3D KSP Component U,V,W,P residual norm [ %1.12e, %1.12e, %1.12e, %1.12e ]\n",n,norms[0],norms[1],norms[2],norms[3]);
1707:   return(0);
1708: }

1712: static PetscErrorCode PCMGSetupViaCoarsen(PC pc,DM da_fine)
1713: {
1714:   PetscInt       nlevels,k,PETSC_UNUSED finest;
1715:   DM             *da_list,*daclist;
1716:   Mat            R;

1720:   nlevels = 1;
1721:   PetscOptionsGetInt(NULL,"-levels",&nlevels,0);

1723:   PetscMalloc(sizeof(DM)*nlevels,&da_list);
1724:   for (k=0; k<nlevels; k++) da_list[k] = NULL;
1725:   PetscMalloc(sizeof(DM)*nlevels,&daclist);
1726:   for (k=0; k<nlevels; k++) daclist[k] = NULL;

1728:   /* finest grid is nlevels - 1 */
1729:   finest     = nlevels - 1;
1730:   daclist[0] = da_fine;
1731:   PetscObjectReference((PetscObject)da_fine);
1732:   DMCoarsenHierarchy(da_fine,nlevels-1,&daclist[1]);
1733:   for (k=0; k<nlevels; k++) {
1734:     da_list[k] = daclist[nlevels-1-k];
1735:     DMDASetUniformCoordinates(da_list[k],0.0,1.0,0.0,1.0,0.0,1.0);
1736:   }

1738:   PCMGSetLevels(pc,nlevels,NULL);
1739:   PCMGSetType(pc,PC_MG_MULTIPLICATIVE);
1740:   PCMGSetGalerkin(pc,PETSC_TRUE);

1742:   for (k=1; k<nlevels; k++) {
1743:     DMCreateInterpolation(da_list[k-1],da_list[k],&R,NULL);
1744:     PCMGSetInterpolation(pc,k,R);
1745:     MatDestroy(&R);
1746:   }

1748:   /* tidy up */
1749:   for (k=0; k<nlevels; k++) {
1750:     DMDestroy(&da_list[k]);
1751:   }
1752:   PetscFree(da_list);
1753:   PetscFree(daclist);
1754:   return(0);
1755: }

1759: static PetscErrorCode solve_stokes_3d_coupled(PetscInt mx,PetscInt my,PetscInt mz)
1760: {
1761:   DM             da_Stokes;
1762:   PetscInt       u_dof,p_dof,dof,stencil_width;
1763:   Mat            A,B;
1764:   PetscInt       mxl,myl,mzl;
1765:   DM             vel_cda;
1766:   Vec            vel_coords;
1767:   PetscInt       p;
1768:   Vec            f,X;
1769:   DMDACoor3d     ***_vel_coords;
1770:   PetscInt       its;
1771:   KSP            ksp_S;
1772:   PetscInt       model_definition = 0;
1773:   PetscInt       ei,ej,ek,sex,sey,sez,Imx,Jmy,Kmz;
1774:   CellProperties cell_properties;
1775:   PetscBool      write_output = PETSC_FALSE;

1779:   /* Generate the da for velocity and pressure */
1780:   /* Num nodes in each direction is mx+1, my+1, mz+1 */
1781:   u_dof         = U_DOFS; /* Vx, Vy - velocities */
1782:   p_dof         = P_DOFS; /* p - pressure */
1783:   dof           = u_dof+p_dof;
1784:   stencil_width = 1;
1785:   DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,
1786:                                mx+1,my+1,mz+1,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,dof,stencil_width,NULL,NULL,NULL,&da_Stokes);
1787:   DMDASetFieldName(da_Stokes,0,"Vx");
1788:   DMDASetFieldName(da_Stokes,1,"Vy");
1789:   DMDASetFieldName(da_Stokes,2,"Vz");
1790:   DMDASetFieldName(da_Stokes,3,"P");

1792:   /* unit box [0,1] x [0,1] x [0,1] */
1793:   DMDASetUniformCoordinates(da_Stokes,0.0,1.0,0.0,1.0,0.0,1.0);

1795:   /* local number of elements */
1796:   DMDAGetLocalElementSize(da_Stokes,&mxl,&myl,&mzl);

1798:   /* create quadrature point info for PDE coefficients */
1799:   CellPropertiesCreate(da_Stokes,&cell_properties);

1801:   /* interpolate the coordinates to quadrature points */
1802:   DMGetCoordinateDM(da_Stokes,&vel_cda);
1803:   DMGetCoordinatesLocal(da_Stokes,&vel_coords);
1804:   DMDAVecGetArray(vel_cda,vel_coords,&_vel_coords);
1805:   DMDAGetElementCorners(da_Stokes,&sex,&sey,&sez,&Imx,&Jmy,&Kmz);
1806:   for (ek = sez; ek < sez+Kmz; ek++) {
1807:     for (ej = sey; ej < sey+Jmy; ej++) {
1808:       for (ei = sex; ei < sex+Imx; ei++) {
1809:         /* get coords for the element */
1810:         PetscInt               ngp;
1811:         PetscScalar            gp_xi[GAUSS_POINTS][NSD],gp_weight[GAUSS_POINTS];
1812:         PetscScalar            el_coords[NSD*NODES_PER_EL];
1813:         GaussPointCoefficients *cell;

1815:         CellPropertiesGetCell(cell_properties,ei,ej,ek,&cell);
1816:         GetElementCoords3D(_vel_coords,ei,ej,ek,el_coords);
1817:         ConstructGaussQuadrature3D(&ngp,gp_xi,gp_weight);

1819:         for (p = 0; p < GAUSS_POINTS; p++) {
1820:           PetscScalar xi_p[NSD],Ni_p[NODES_PER_EL];
1821:           PetscScalar gp_x,gp_y,gp_z;
1822:           PetscInt    n;

1824:           xi_p[0] = gp_xi[p][0];
1825:           xi_p[1] = gp_xi[p][1];
1826:           xi_p[2] = gp_xi[p][2];
1827:           ShapeFunctionQ13D_Evaluate(xi_p,Ni_p);

1829:           gp_x = gp_y = gp_z = 0.0;
1830:           for (n = 0; n < NODES_PER_EL; n++) {
1831:             gp_x = gp_x+Ni_p[n]*el_coords[NSD*n];
1832:             gp_y = gp_y+Ni_p[n]*el_coords[NSD*n+1];
1833:             gp_z = gp_z+Ni_p[n]*el_coords[NSD*n+2];
1834:           }
1835:           cell->gp_coords[NSD*p]   = gp_x;
1836:           cell->gp_coords[NSD*p+1] = gp_y;
1837:           cell->gp_coords[NSD*p+2] = gp_z;
1838:         }
1839:       }
1840:     }
1841:   }

1843:   PetscOptionsGetInt(NULL,"-model",&model_definition,NULL);

1845:   switch (model_definition) {
1846:   case 0: /* isoviscous */
1847:     for (ek = sez; ek < sez+Kmz; ek++) {
1848:       for (ej = sey; ej < sey+Jmy; ej++) {
1849:         for (ei = sex; ei < sex+Imx; ei++) {
1850:           GaussPointCoefficients *cell;

1852:           CellPropertiesGetCell(cell_properties,ei,ej,ek,&cell);
1853:           for (p = 0; p < GAUSS_POINTS; p++) {
1854:             PetscReal coord_x = PetscRealPart(cell->gp_coords[NSD*p]);
1855:             PetscReal coord_y = PetscRealPart(cell->gp_coords[NSD*p+1]);
1856:             PetscReal coord_z = PetscRealPart(cell->gp_coords[NSD*p+2]);

1858:             cell->eta[p] = 1.0;

1860:             cell->fx[p] = 0.0*coord_x;
1861:             cell->fy[p] = -PetscSinReal(2.2*PETSC_PI*coord_y)*PetscCosReal(1.0*PETSC_PI*coord_x);
1862:             cell->fz[p] = 0.0*coord_z;
1863:             cell->hc[p] = 0.0;
1864:           }
1865:         }
1866:       }
1867:     }
1868:     break;

1870:   case 1: /* manufactured */
1871:     for (ek = sez; ek < sez+Kmz; ek++) {
1872:       for (ej = sey; ej < sey+Jmy; ej++) {
1873:         for (ei = sex; ei < sex+Imx; ei++) {
1874:           PetscReal              eta,Fm[NSD],Fc,pos[NSD];
1875:           GaussPointCoefficients *cell;

1877:           CellPropertiesGetCell(cell_properties,ei,ej,ek,&cell);
1878:           for (p = 0; p < GAUSS_POINTS; p++) {
1879:             PetscReal coord_x = PetscRealPart(cell->gp_coords[NSD*p]);
1880:             PetscReal coord_y = PetscRealPart(cell->gp_coords[NSD*p+1]);
1881:             PetscReal coord_z = PetscRealPart(cell->gp_coords[NSD*p+2]);

1883:             pos[0] = coord_x;
1884:             pos[1] = coord_y;
1885:             pos[2] = coord_z;

1887:             evaluate_MS_FrankKamentski(pos,NULL,NULL,&eta,Fm,&Fc);
1888:             cell->eta[p] = eta;

1890:             cell->fx[p] = Fm[0];
1891:             cell->fy[p] = Fm[1];
1892:             cell->fz[p] = Fm[2];
1893:             cell->hc[p] = Fc;
1894:           }
1895:         }
1896:       }
1897:     }
1898:     break;

1900:   case 2: /* solcx */
1901:     for (ek = sez; ek < sez+Kmz; ek++) {
1902:       for (ej = sey; ej < sey+Jmy; ej++) {
1903:         for (ei = sex; ei < sex+Imx; ei++) {
1904:           GaussPointCoefficients *cell;

1906:           CellPropertiesGetCell(cell_properties,ei,ej,ek,&cell);
1907:           for (p = 0; p < GAUSS_POINTS; p++) {
1908:             PetscReal coord_x = PetscRealPart(cell->gp_coords[NSD*p]);
1909:             PetscReal coord_y = PetscRealPart(cell->gp_coords[NSD*p+1]);
1910:             PetscReal coord_z = PetscRealPart(cell->gp_coords[NSD*p+2]);

1912:             cell->eta[p] = 1.0;

1914:             cell->fx[p] = 0.0;
1915:             cell->fy[p] = -PetscSinReal(3.0*PETSC_PI*coord_y)*PetscCosReal(1.0*PETSC_PI*coord_x);
1916:             cell->fz[p] = 0.0*coord_z;
1917:             cell->hc[p] = 0.0;
1918:           }
1919:         }
1920:       }
1921:     }
1922:     break;

1924:   case 3: /* sinker */
1925:     for (ek = sez; ek < sez+Kmz; ek++) {
1926:       for (ej = sey; ej < sey+Jmy; ej++) {
1927:         for (ei = sex; ei < sex+Imx; ei++) {
1928:           GaussPointCoefficients *cell;

1930:           CellPropertiesGetCell(cell_properties,ei,ej,ek,&cell);
1931:           for (p = 0; p < GAUSS_POINTS; p++) {
1932:             PetscReal xp = PetscRealPart(cell->gp_coords[NSD*p]);
1933:             PetscReal yp = PetscRealPart(cell->gp_coords[NSD*p+1]);
1934:             PetscReal zp = PetscRealPart(cell->gp_coords[NSD*p+2]);

1936:             cell->eta[p] = 1.0e-2;
1937:             cell->fx[p]  = 0.0;
1938:             cell->fy[p]  = 0.0;
1939:             cell->fz[p]  = 0.0;
1940:             cell->hc[p]  = 0.0;

1942:             if ((fabs(xp-0.5) < 0.2) &&
1943:                 (fabs(yp-0.5) < 0.2) &&
1944:                 (fabs(zp-0.5) < 0.2)) {
1945:               cell->eta[p] = 1.0;
1946:               cell->fz[p]  = 1.0;
1947:             }

1949:           }
1950:         }
1951:       }
1952:     }
1953:     break;

1955:   default:
1956:     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"No default model is supported. Choose either -model {0,1,2,3}");
1957:     break;
1958:   }

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

1962:   /* Generate a matrix with the correct non-zero pattern of type AIJ. This will work in parallel and serial */
1963:   DMSetMatType(da_Stokes,MATAIJ);
1964:   DMCreateMatrix(da_Stokes,&A);
1965:   DMCreateMatrix(da_Stokes,&B);
1966:   DMCreateGlobalVector(da_Stokes,&X);
1967:   DMCreateGlobalVector(da_Stokes,&f);

1969:   /* assemble A11 */
1970:   MatZeroEntries(A);
1971:   MatZeroEntries(B);
1972:   VecZeroEntries(f);

1974:   AssembleA_Stokes(A,da_Stokes,cell_properties);
1975:   AssembleA_PCStokes(B,da_Stokes,cell_properties);
1976:   /* build force vector */
1977:   AssembleF_Stokes(f,da_Stokes,cell_properties);

1979:   /* SOLVE */
1980:   KSPCreate(PETSC_COMM_WORLD,&ksp_S);
1981:   KSPSetOptionsPrefix(ksp_S,"stokes_"); /* stokes */
1982:   KSPSetOperators(ksp_S,A,B);
1983:   KSPSetFromOptions(ksp_S);

1985:   {
1986:     PC             pc;
1987:     const PetscInt ufields[] = {0,1,2},pfields[] = {3};
1988:     KSPGetPC(ksp_S,&pc);
1989:     PCFieldSplitSetBlockSize(pc,4);
1990:     PCFieldSplitSetFields(pc,"u",3,ufields,ufields);
1991:     PCFieldSplitSetFields(pc,"p",1,pfields,pfields);
1992:   }

1994:   {
1995:     PC        pc;
1996:     PetscBool same = PETSC_FALSE;
1997:     KSPGetPC(ksp_S,&pc);
1998:     PetscObjectTypeCompare((PetscObject)pc,PCMG,&same);
1999:     if (same) {
2000:       PCMGSetupViaCoarsen(pc,da_Stokes);
2001:     }
2002:   }

2004:   {
2005:     PetscBool stokes_monitor = PETSC_FALSE;
2006:     PetscOptionsGetBool(NULL,"-stokes_ksp_monitor_blocks",&stokes_monitor,0);
2007:     if (stokes_monitor) {
2008:       KSPMonitorSet(ksp_S,KSPMonitorStokesBlocks,NULL,NULL);
2009:     }
2010:   }
2011:   KSPSolve(ksp_S,f,X);

2013:   PetscOptionsGetBool(NULL,"-write_pvts",&write_output,NULL);
2014:   if (write_output) {DAView3DPVTS(da_Stokes,X,"up");}
2015:   {
2016:     PetscBool flg = PETSC_FALSE;
2017:     char      filename[PETSC_MAX_PATH_LEN];
2018:     PetscOptionsGetString(NULL,"-write_binary",filename,sizeof(filename),&flg);
2019:     if (flg) {
2020:       PetscViewer viewer;
2021:       /* PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename[0]?filename:"ex42-binaryoutput",FILE_MODE_WRITE,&viewer); */
2022:       PetscViewerVTKOpen(PETSC_COMM_WORLD,"ex42.vts",FILE_MODE_WRITE,&viewer);
2023:       VecView(X,viewer);
2024:       PetscViewerDestroy(&viewer);
2025:     }
2026:   }
2027:   KSPGetIterationNumber(ksp_S,&its);

2029:   /* verify */
2030:   if (model_definition == 1) {
2031:     DM  da_Stokes_analytic;
2032:     Vec X_analytic;

2034:     DMDACreateManufacturedSolution(mx,my,mz,&da_Stokes_analytic,&X_analytic);
2035:     if (write_output) {
2036:       DAView3DPVTS(da_Stokes_analytic,X_analytic,"ms");
2037:     }
2038:     DMDAIntegrateErrors3D(da_Stokes_analytic,X,X_analytic);
2039:     if (write_output) {
2040:       DAView3DPVTS(da_Stokes,X,"up2");
2041:     }
2042:     DMDestroy(&da_Stokes_analytic);
2043:     VecDestroy(&X_analytic);
2044:   }

2046:   KSPDestroy(&ksp_S);
2047:   VecDestroy(&X);
2048:   VecDestroy(&f);
2049:   MatDestroy(&A);
2050:   MatDestroy(&B);

2052:   CellPropertiesDestroy(&cell_properties);
2053:   DMDestroy(&da_Stokes);
2054:   return(0);
2055: }

2059: int main(int argc,char **args)
2060: {
2062:   PetscInt       mx,my,mz;

2064:   PetscInitialize(&argc,&args,(char*)0,help);

2066:   mx   = my = mz = 10;
2067:   PetscOptionsGetInt(NULL,"-mx",&mx,NULL);
2068:   my   = mx; mz = mx;
2069:   PetscOptionsGetInt(NULL,"-my",&my,NULL);
2070:   PetscOptionsGetInt(NULL,"-mz",&mz,NULL);

2072:   solve_stokes_3d_coupled(mx,my,mz);

2074:   PetscFinalize();
2075:   return 0;
2076: }