Actual source code: ex42.c

  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), 4 (subdomain jumps) \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: /* elements */
 50: PetscErrorCode CellPropertiesCreate(DM da_stokes,CellProperties *C)
 51: {
 53:   CellProperties cells;
 54:   PetscInt       mx,my,mz,sex,sey,sez;

 57:   PetscNew(&cells);

 59:   DMDAGetElementsCorners(da_stokes,&sex,&sey,&sez);
 60:   DMDAGetElementsSizes(da_stokes,&mx,&my,&mz);
 61:   cells->mx     = mx;
 62:   cells->my     = my;
 63:   cells->mz     = mz;
 64:   cells->ncells = mx * my * mz;
 65:   cells->sex    = sex;
 66:   cells->sey    = sey;
 67:   cells->sez    = sez;

 69:   PetscMalloc1(mx*my*mz,&cells->gpc);

 71:   *C = cells;
 72:   return(0);
 73: }

 75: PetscErrorCode CellPropertiesDestroy(CellProperties *C)
 76: {
 78:   CellProperties cells;

 81:   if (!C) return(0);
 82:   cells = *C;
 83:   PetscFree(cells->gpc);
 84:   PetscFree(cells);
 85:   *C = NULL;
 86:   return(0);
 87: }

 89: PetscErrorCode CellPropertiesGetCell(CellProperties C,PetscInt II,PetscInt J,PetscInt K,GaussPointCoefficients **G)
 90: {
 92:   *G = &C->gpc[(II-C->sex) + (J-C->sey)*C->mx + (K-C->sez)*C->mx*C->my];
 93:   return(0);
 94: }

 96: /* FEM routines */
 97: /*
 98:  Element: Local basis function ordering
 99:  1-----2
100:  |     |
101:  |     |
102:  0-----3
103:  */
104: static void ShapeFunctionQ13D_Evaluate(PetscScalar _xi[],PetscScalar Ni[])
105: {
106:   PetscReal xi   = PetscRealPart(_xi[0]);
107:   PetscReal eta  = PetscRealPart(_xi[1]);
108:   PetscReal zeta = PetscRealPart(_xi[2]);

110:   Ni[0] = 0.125 * (1.0 - xi) * (1.0 - eta) * (1.0 - zeta);
111:   Ni[1] = 0.125 * (1.0 - xi) * (1.0 + eta) * (1.0 - zeta);
112:   Ni[2] = 0.125 * (1.0 + xi) * (1.0 + eta) * (1.0 - zeta);
113:   Ni[3] = 0.125 * (1.0 + xi) * (1.0 - eta) * (1.0 - zeta);

115:   Ni[4] = 0.125 * (1.0 - xi) * (1.0 - eta) * (1.0 + zeta);
116:   Ni[5] = 0.125 * (1.0 - xi) * (1.0 + eta) * (1.0 + zeta);
117:   Ni[6] = 0.125 * (1.0 + xi) * (1.0 + eta) * (1.0 + zeta);
118:   Ni[7] = 0.125 * (1.0 + xi) * (1.0 - eta) * (1.0 + zeta);
119: }

121: static void ShapeFunctionQ13D_Evaluate_dxi(PetscScalar _xi[],PetscScalar GNi[][NODES_PER_EL])
122: {
123:   PetscReal xi   = PetscRealPart(_xi[0]);
124:   PetscReal eta  = PetscRealPart(_xi[1]);
125:   PetscReal zeta = PetscRealPart(_xi[2]);
126:   /* xi deriv */
127:   GNi[0][0] = -0.125 * (1.0 - eta) * (1.0 - zeta);
128:   GNi[0][1] = -0.125 * (1.0 + eta) * (1.0 - zeta);
129:   GNi[0][2] =  0.125 * (1.0 + eta) * (1.0 - zeta);
130:   GNi[0][3] =  0.125 * (1.0 - eta) * (1.0 - zeta);

132:   GNi[0][4] = -0.125 * (1.0 - eta) * (1.0 + zeta);
133:   GNi[0][5] = -0.125 * (1.0 + eta) * (1.0 + zeta);
134:   GNi[0][6] =  0.125 * (1.0 + eta) * (1.0 + zeta);
135:   GNi[0][7] =  0.125 * (1.0 - eta) * (1.0 + zeta);
136:   /* eta deriv */
137:   GNi[1][0] = -0.125 * (1.0 - xi) * (1.0 - zeta);
138:   GNi[1][1] =  0.125 * (1.0 - xi) * (1.0 - zeta);
139:   GNi[1][2] =  0.125 * (1.0 + xi) * (1.0 - zeta);
140:   GNi[1][3] = -0.125 * (1.0 + xi) * (1.0 - zeta);

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

152:   GNi[2][4] = 0.125 * (1.0 - xi) * (1.0 - eta);
153:   GNi[2][5] = 0.125 * (1.0 - xi) * (1.0 + eta);
154:   GNi[2][6] = 0.125 * (1.0 + xi) * (1.0 + eta);
155:   GNi[2][7] = 0.125 * (1.0 + xi) * (1.0 - eta);
156: }

158: static void matrix_inverse_3x3(PetscScalar A[3][3],PetscScalar B[3][3])
159: {
160:   PetscScalar t4, t6, t8, t10, t12, t14, t17;

162:   t4  = A[2][0] * A[0][1];
163:   t6  = A[2][0] * A[0][2];
164:   t8  = A[1][0] * A[0][1];
165:   t10 = A[1][0] * A[0][2];
166:   t12 = A[0][0] * A[1][1];
167:   t14 = A[0][0] * A[1][2];
168:   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]);

170:   B[0][0] = (A[1][1] * A[2][2] - A[1][2] * A[2][1]) * t17;
171:   B[0][1] = -(A[0][1] * A[2][2] - A[0][2] * A[2][1]) * t17;
172:   B[0][2] = (A[0][1] * A[1][2] - A[0][2] * A[1][1]) * t17;
173:   B[1][0] = -(-A[2][0] * A[1][2] + A[1][0] * A[2][2]) * t17;
174:   B[1][1] = (-t6 + A[0][0] * A[2][2]) * t17;
175:   B[1][2] = -(-t10 + t14) * t17;
176:   B[2][0] = (-A[2][0] * A[1][1] + A[1][0] * A[2][1]) * t17;
177:   B[2][1] = -(-t4 + A[0][0] * A[2][1]) * t17;
178:   B[2][2] = (-t8 + t12) * t17;
179: }

181: static void ShapeFunctionQ13D_Evaluate_dx(PetscScalar GNi[][NODES_PER_EL],PetscScalar GNx[][NODES_PER_EL],PetscScalar coords[],PetscScalar *det_J)
182: {
183:   PetscScalar J00,J01,J02,J10,J11,J12,J20,J21,J22;
184:   PetscInt    n;
185:   PetscScalar iJ[3][3],JJ[3][3];

187:   J00 = J01 = J02 = 0.0;
188:   J10 = J11 = J12 = 0.0;
189:   J20 = J21 = J22 = 0.0;
190:   for (n=0; n<NODES_PER_EL; n++) {
191:     PetscScalar cx = coords[NSD*n + 0];
192:     PetscScalar cy = coords[NSD*n + 1];
193:     PetscScalar cz = coords[NSD*n + 2];

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

200:     J10 = J10 + GNi[1][n] * cx;   /* J_yx = dy/dxi */
201:     J11 = J11 + GNi[1][n] * cy;   /* J_yy */
202:     J12 = J12 + GNi[1][n] * cz;   /* J_yz */

204:     J20 = J20 + GNi[2][n] * cx;   /* J_zx */
205:     J21 = J21 + GNi[2][n] * cy;   /* J_zy */
206:     J22 = J22 + GNi[2][n] * cz;   /* J_zz */
207:   }

209:   JJ[0][0] = J00;      JJ[0][1] = J01;      JJ[0][2] = J02;
210:   JJ[1][0] = J10;      JJ[1][1] = J11;      JJ[1][2] = J12;
211:   JJ[2][0] = J20;      JJ[2][1] = J21;      JJ[2][2] = J22;

213:   matrix_inverse_3x3(JJ,iJ);

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

217:   for (n=0; n<NODES_PER_EL; n++) {
218:     GNx[0][n] = GNi[0][n]*iJ[0][0] + GNi[1][n]*iJ[0][1] + GNi[2][n]*iJ[0][2];
219:     GNx[1][n] = GNi[0][n]*iJ[1][0] + GNi[1][n]*iJ[1][1] + GNi[2][n]*iJ[1][2];
220:     GNx[2][n] = GNi[0][n]*iJ[2][0] + GNi[1][n]*iJ[2][1] + GNi[2][n]*iJ[2][2];
221:   }
222: }

224: static void ConstructGaussQuadrature3D(PetscInt *ngp,PetscScalar gp_xi[][NSD],PetscScalar gp_weight[])
225: {
226:   *ngp        = 8;
227:   gp_xi[0][0] = -0.57735026919; gp_xi[0][1] = -0.57735026919; gp_xi[0][2] = -0.57735026919;
228:   gp_xi[1][0] = -0.57735026919; gp_xi[1][1] =  0.57735026919; gp_xi[1][2] = -0.57735026919;
229:   gp_xi[2][0] =  0.57735026919; gp_xi[2][1] =  0.57735026919; gp_xi[2][2] = -0.57735026919;
230:   gp_xi[3][0] =  0.57735026919; gp_xi[3][1] = -0.57735026919; gp_xi[3][2] = -0.57735026919;

232:   gp_xi[4][0] = -0.57735026919; gp_xi[4][1] = -0.57735026919; gp_xi[4][2] =  0.57735026919;
233:   gp_xi[5][0] = -0.57735026919; gp_xi[5][1] =  0.57735026919; gp_xi[5][2] =  0.57735026919;
234:   gp_xi[6][0] =  0.57735026919; gp_xi[6][1] =  0.57735026919; gp_xi[6][2] =  0.57735026919;
235:   gp_xi[7][0] =  0.57735026919; gp_xi[7][1] = -0.57735026919; gp_xi[7][2] =  0.57735026919;

237:   gp_weight[0] = 1.0;
238:   gp_weight[1] = 1.0;
239:   gp_weight[2] = 1.0;
240:   gp_weight[3] = 1.0;

242:   gp_weight[4] = 1.0;
243:   gp_weight[5] = 1.0;
244:   gp_weight[6] = 1.0;
245:   gp_weight[7] = 1.0;
246: }

248: /*
249:  i,j are the element indices
250:  The unknown is a vector quantity.
251:  The s[].c is used to indicate the degree of freedom.
252:  */
253: static PetscErrorCode DMDAGetElementEqnums3D_up(MatStencil s_u[],MatStencil s_p[],PetscInt i,PetscInt j,PetscInt k)
254: {
255:   PetscInt n;

258:   /* velocity */
259:   n = 0;
260:   /* node 0 */
261:   s_u[n].i = i; s_u[n].j = j; s_u[n].k = k; s_u[n].c = 0; n++; /* Vx0 */
262:   s_u[n].i = i; s_u[n].j = j; s_u[n].k = k; s_u[n].c = 1; n++; /* Vy0 */
263:   s_u[n].i = i; s_u[n].j = j; s_u[n].k = k; s_u[n].c = 2; n++; /* Vz0 */

265:   s_u[n].i = i; s_u[n].j = j+1; s_u[n].k = k; s_u[n].c = 0; n++;
266:   s_u[n].i = i; s_u[n].j = j+1; s_u[n].k = k; s_u[n].c = 1; n++;
267:   s_u[n].i = i; s_u[n].j = j+1; s_u[n].k = k; s_u[n].c = 2; n++;

269:   s_u[n].i = i+1; s_u[n].j = j+1; s_u[n].k = k; s_u[n].c = 0; n++;
270:   s_u[n].i = i+1; s_u[n].j = j+1; s_u[n].k = k; s_u[n].c = 1; n++;
271:   s_u[n].i = i+1; s_u[n].j = j+1; s_u[n].k = k; s_u[n].c = 2; n++;

273:   s_u[n].i = i+1; s_u[n].j = j; s_u[n].k = k; s_u[n].c = 0; n++;
274:   s_u[n].i = i+1; s_u[n].j = j; s_u[n].k = k; s_u[n].c = 1; n++;
275:   s_u[n].i = i+1; s_u[n].j = j; s_u[n].k = k; s_u[n].c = 2; n++;

277:   /* */
278:   s_u[n].i = i; s_u[n].j = j; s_u[n].k = k+1; s_u[n].c = 0; n++; /* Vx4 */
279:   s_u[n].i = i; s_u[n].j = j; s_u[n].k = k+1; s_u[n].c = 1; n++; /* Vy4 */
280:   s_u[n].i = i; s_u[n].j = j; s_u[n].k = k+1; s_u[n].c = 2; n++; /* Vz4 */

282:   s_u[n].i = i; s_u[n].j = j+1; s_u[n].k = k+1; s_u[n].c = 0; n++;
283:   s_u[n].i = i; s_u[n].j = j+1; s_u[n].k = k+1; s_u[n].c = 1; n++;
284:   s_u[n].i = i; s_u[n].j = j+1; s_u[n].k = k+1; s_u[n].c = 2; n++;

286:   s_u[n].i = i+1; s_u[n].j = j+1; s_u[n].k = k+1; s_u[n].c = 0; n++;
287:   s_u[n].i = i+1; s_u[n].j = j+1; s_u[n].k = k+1; s_u[n].c = 1; n++;
288:   s_u[n].i = i+1; s_u[n].j = j+1; s_u[n].k = k+1; s_u[n].c = 2; n++;

290:   s_u[n].i = i+1; s_u[n].j = j; s_u[n].k = k+1; s_u[n].c = 0; n++;
291:   s_u[n].i = i+1; s_u[n].j = j; s_u[n].k = k+1; s_u[n].c = 1; n++;
292:   s_u[n].i = i+1; s_u[n].j = j; s_u[n].k = k+1; s_u[n].c = 2; n++;

294:   /* pressure */
295:   n = 0;

297:   s_p[n].i = i;   s_p[n].j = j;   s_p[n].k = k; s_p[n].c = 3; n++; /* P0 */
298:   s_p[n].i = i;   s_p[n].j = j+1; s_p[n].k = k; s_p[n].c = 3; n++;
299:   s_p[n].i = i+1; s_p[n].j = j+1; s_p[n].k = k; s_p[n].c = 3; n++;
300:   s_p[n].i = i+1; s_p[n].j = j;   s_p[n].k = k; s_p[n].c = 3; n++;

302:   s_p[n].i = i;   s_p[n].j = j;   s_p[n].k = k+1; s_p[n].c = 3; n++; /* P0 */
303:   s_p[n].i = i;   s_p[n].j = j+1; s_p[n].k = k+1; s_p[n].c = 3; n++;
304:   s_p[n].i = i+1; s_p[n].j = j+1; s_p[n].k = k+1; s_p[n].c = 3; n++;
305:   s_p[n].i = i+1; s_p[n].j = j;   s_p[n].k = k+1; s_p[n].c = 3; n++;
306:   return(0);
307: }

309: static PetscErrorCode GetElementCoords3D(DMDACoor3d ***coords,PetscInt i,PetscInt j,PetscInt k,PetscScalar el_coord[])
310: {
312:   /* get coords for the element */
313:   el_coord[0] = coords[k][j][i].x;
314:   el_coord[1] = coords[k][j][i].y;
315:   el_coord[2] = coords[k][j][i].z;

317:   el_coord[3] = coords[k][j+1][i].x;
318:   el_coord[4] = coords[k][j+1][i].y;
319:   el_coord[5] = coords[k][j+1][i].z;

321:   el_coord[6] = coords[k][j+1][i+1].x;
322:   el_coord[7] = coords[k][j+1][i+1].y;
323:   el_coord[8] = coords[k][j+1][i+1].z;

325:   el_coord[9]  = coords[k][j][i+1].x;
326:   el_coord[10] = coords[k][j][i+1].y;
327:   el_coord[11] = coords[k][j][i+1].z;

329:   el_coord[12] = coords[k+1][j][i].x;
330:   el_coord[13] = coords[k+1][j][i].y;
331:   el_coord[14] = coords[k+1][j][i].z;

333:   el_coord[15] = coords[k+1][j+1][i].x;
334:   el_coord[16] = coords[k+1][j+1][i].y;
335:   el_coord[17] = coords[k+1][j+1][i].z;

337:   el_coord[18] = coords[k+1][j+1][i+1].x;
338:   el_coord[19] = coords[k+1][j+1][i+1].y;
339:   el_coord[20] = coords[k+1][j+1][i+1].z;

341:   el_coord[21] = coords[k+1][j][i+1].x;
342:   el_coord[22] = coords[k+1][j][i+1].y;
343:   el_coord[23] = coords[k+1][j][i+1].z;
344:   return(0);
345: }

347: static PetscErrorCode StokesDAGetNodalFields3D(StokesDOF ***field,PetscInt i,PetscInt j,PetscInt k,StokesDOF nodal_fields[])
348: {
350:   /* get the nodal fields for u */
351:   nodal_fields[0].u_dof = field[k][j][i].u_dof;
352:   nodal_fields[0].v_dof = field[k][j][i].v_dof;
353:   nodal_fields[0].w_dof = field[k][j][i].w_dof;

355:   nodal_fields[1].u_dof = field[k][j+1][i].u_dof;
356:   nodal_fields[1].v_dof = field[k][j+1][i].v_dof;
357:   nodal_fields[1].w_dof = field[k][j+1][i].w_dof;

359:   nodal_fields[2].u_dof = field[k][j+1][i+1].u_dof;
360:   nodal_fields[2].v_dof = field[k][j+1][i+1].v_dof;
361:   nodal_fields[2].w_dof = field[k][j+1][i+1].w_dof;

363:   nodal_fields[3].u_dof = field[k][j][i+1].u_dof;
364:   nodal_fields[3].v_dof = field[k][j][i+1].v_dof;
365:   nodal_fields[3].w_dof = field[k][j][i+1].w_dof;

367:   nodal_fields[4].u_dof = field[k+1][j][i].u_dof;
368:   nodal_fields[4].v_dof = field[k+1][j][i].v_dof;
369:   nodal_fields[4].w_dof = field[k+1][j][i].w_dof;

371:   nodal_fields[5].u_dof = field[k+1][j+1][i].u_dof;
372:   nodal_fields[5].v_dof = field[k+1][j+1][i].v_dof;
373:   nodal_fields[5].w_dof = field[k+1][j+1][i].w_dof;

375:   nodal_fields[6].u_dof = field[k+1][j+1][i+1].u_dof;
376:   nodal_fields[6].v_dof = field[k+1][j+1][i+1].v_dof;
377:   nodal_fields[6].w_dof = field[k+1][j+1][i+1].w_dof;

379:   nodal_fields[7].u_dof = field[k+1][j][i+1].u_dof;
380:   nodal_fields[7].v_dof = field[k+1][j][i+1].v_dof;
381:   nodal_fields[7].w_dof = field[k+1][j][i+1].w_dof;

383:   /* get the nodal fields for p */
384:   nodal_fields[0].p_dof = field[k][j][i].p_dof;
385:   nodal_fields[1].p_dof = field[k][j+1][i].p_dof;
386:   nodal_fields[2].p_dof = field[k][j+1][i+1].p_dof;
387:   nodal_fields[3].p_dof = field[k][j][i+1].p_dof;

389:   nodal_fields[4].p_dof = field[k+1][j][i].p_dof;
390:   nodal_fields[5].p_dof = field[k+1][j+1][i].p_dof;
391:   nodal_fields[6].p_dof = field[k+1][j+1][i+1].p_dof;
392:   nodal_fields[7].p_dof = field[k+1][j][i+1].p_dof;
393:   return(0);
394: }

396: 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)
397: {
398:   PetscInt              ij;
399:   PETSC_UNUSED PetscInt r,c,nr,nc;

401:   nr = w_NPE*w_dof;
402:   nc = u_NPE*u_dof;

404:   r = w_dof*wi+wd;
405:   c = u_dof*ui+ud;

407:   ij = r*nc+c;

409:   return ij;
410: }

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

417:   for (n = 0; n<NODES_PER_EL; n++) {
418:     II = u_eqn[NSD*n].i;
419:     J = u_eqn[NSD*n].j;
420:     K = u_eqn[NSD*n].k;

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

424:     II = u_eqn[NSD*n+1].i;
425:     J = u_eqn[NSD*n+1].j;
426:     K = u_eqn[NSD*n+1].k;

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

430:     II = u_eqn[NSD*n+2].i;
431:     J = u_eqn[NSD*n+2].j;
432:     K = u_eqn[NSD*n+2].k;
433:     fields_F[K][J][II].w_dof = fields_F[K][J][II].w_dof+Fe_u[NSD*n+2];

435:     II = p_eqn[n].i;
436:     J = p_eqn[n].j;
437:     K = p_eqn[n].k;

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

441:   }
442:   return(0);
443: }

445: static void FormStressOperatorQ13D(PetscScalar Ke[],PetscScalar coords[],PetscScalar eta[])
446: {
447:   PetscInt       ngp;
448:   PetscScalar    gp_xi[GAUSS_POINTS][NSD];
449:   PetscScalar    gp_weight[GAUSS_POINTS];
450:   PetscInt       p,i,j,k;
451:   PetscScalar    GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
452:   PetscScalar    J_p,tildeD[6];
453:   PetscScalar    B[6][U_DOFS*NODES_PER_EL];
454:   const PetscInt nvdof = U_DOFS*NODES_PER_EL;

456:   /* define quadrature rule */
457:   ConstructGaussQuadrature3D(&ngp,gp_xi,gp_weight);

459:   /* evaluate integral */
460:   for (p = 0; p < ngp; p++) {
461:     ShapeFunctionQ13D_Evaluate_dxi(gp_xi[p],GNi_p);
462:     ShapeFunctionQ13D_Evaluate_dx(GNi_p,GNx_p,coords,&J_p);

464:     for (i = 0; i < NODES_PER_EL; i++) {
465:       PetscScalar d_dx_i = GNx_p[0][i];
466:       PetscScalar d_dy_i = GNx_p[1][i];
467:       PetscScalar d_dz_i = GNx_p[2][i];

469:       B[0][3*i] = d_dx_i; B[0][3*i+1] = 0.0;     B[0][3*i+2] = 0.0;
470:       B[1][3*i] = 0.0;    B[1][3*i+1] = d_dy_i;  B[1][3*i+2] = 0.0;
471:       B[2][3*i] = 0.0;    B[2][3*i+1] = 0.0;     B[2][3*i+2] = d_dz_i;

473:       B[3][3*i] = d_dy_i; B[3][3*i+1] = d_dx_i;  B[3][3*i+2] = 0.0;   /* e_xy */
474:       B[4][3*i] = d_dz_i; B[4][3*i+1] = 0.0;     B[4][3*i+2] = d_dx_i; /* e_xz */
475:       B[5][3*i] = 0.0;    B[5][3*i+1] = d_dz_i;  B[5][3*i+2] = d_dy_i; /* e_yz */
476:     }

478:     tildeD[0] = 2.0*gp_weight[p]*J_p*eta[p];
479:     tildeD[1] = 2.0*gp_weight[p]*J_p*eta[p];
480:     tildeD[2] = 2.0*gp_weight[p]*J_p*eta[p];

482:     tildeD[3] =     gp_weight[p]*J_p*eta[p];
483:     tildeD[4] =     gp_weight[p]*J_p*eta[p];
484:     tildeD[5] =     gp_weight[p]*J_p*eta[p];

486:     /* form Bt tildeD B */
487:     /*
488:      Ke_ij = Bt_ik . D_kl . B_lj
489:      = B_ki . D_kl . B_lj
490:      = B_ki . D_kk . B_kj
491:      */
492:     for (i = 0; i < nvdof; i++) {
493:       for (j = i; j < nvdof; j++) {
494:         for (k = 0; k < 6; k++) {
495:           Ke[i*nvdof+j] += B[k][i]*tildeD[k]*B[k][j];
496:         }
497:       }
498:     }

500:   }
501:   /* fill lower triangular part */
502: #if defined(ASSEMBLE_LOWER_TRIANGULAR)
503:   for (i = 0; i < nvdof; i++) {
504:     for (j = i; j < nvdof; j++) {
505:       Ke[j*nvdof+i] = Ke[i*nvdof+j];
506:     }
507:   }
508: #endif
509: }

511: static void FormGradientOperatorQ13D(PetscScalar Ke[],PetscScalar coords[])
512: {
513:   PetscInt    ngp;
514:   PetscScalar gp_xi[GAUSS_POINTS][NSD];
515:   PetscScalar gp_weight[GAUSS_POINTS];
516:   PetscInt    p,i,j,di;
517:   PetscScalar Ni_p[NODES_PER_EL];
518:   PetscScalar GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
519:   PetscScalar J_p,fac;

521:   /* define quadrature rule */
522:   ConstructGaussQuadrature3D(&ngp,gp_xi,gp_weight);

524:   /* evaluate integral */
525:   for (p = 0; p < ngp; p++) {
526:     ShapeFunctionQ13D_Evaluate(gp_xi[p],Ni_p);
527:     ShapeFunctionQ13D_Evaluate_dxi(gp_xi[p],GNi_p);
528:     ShapeFunctionQ13D_Evaluate_dx(GNi_p,GNx_p,coords,&J_p);
529:     fac = gp_weight[p]*J_p;

531:     for (i = 0; i < NODES_PER_EL; i++) { /* u nodes */
532:       for (di = 0; di < NSD; di++) { /* u dofs */
533:         for (j = 0; j < NODES_PER_EL; j++) {  /* p nodes, p dofs = 1 (ie no loop) */
534:           PetscInt IJ;
535:           IJ = ASS_MAP_wIwDI_uJuDJ(i,di,NODES_PER_EL,3,j,0,NODES_PER_EL,1);

537:           Ke[IJ] = Ke[IJ]-GNx_p[di][i]*Ni_p[j]*fac;
538:         }
539:       }
540:     }
541:   }
542: }

544: static void FormDivergenceOperatorQ13D(PetscScalar De[],PetscScalar coords[])
545: {
546:   PetscScalar Ge[U_DOFS*NODES_PER_EL*P_DOFS*NODES_PER_EL];
547:   PetscInt    i,j;
548:   PetscInt    nr_g,nc_g;

550:   PetscMemzero(Ge,sizeof(Ge));
551:   FormGradientOperatorQ13D(Ge,coords);

553:   nr_g = U_DOFS*NODES_PER_EL;
554:   nc_g = P_DOFS*NODES_PER_EL;

556:   for (i = 0; i < nr_g; i++) {
557:     for (j = 0; j < nc_g; j++) {
558:       De[nr_g*j+i] = Ge[nc_g*i+j];
559:     }
560:   }
561: }

563: static void FormStabilisationOperatorQ13D(PetscScalar Ke[],PetscScalar coords[],PetscScalar eta[])
564: {
565:   PetscInt    ngp;
566:   PetscScalar gp_xi[GAUSS_POINTS][NSD];
567:   PetscScalar gp_weight[GAUSS_POINTS];
568:   PetscInt    p,i,j;
569:   PetscScalar Ni_p[NODES_PER_EL];
570:   PetscScalar GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
571:   PetscScalar J_p,fac,eta_avg;

573:   /* define quadrature rule */
574:   ConstructGaussQuadrature3D(&ngp,gp_xi,gp_weight);

576:   /* evaluate integral */
577:   for (p = 0; p < ngp; p++) {
578:     ShapeFunctionQ13D_Evaluate(gp_xi[p],Ni_p);
579:     ShapeFunctionQ13D_Evaluate_dxi(gp_xi[p],GNi_p);
580:     ShapeFunctionQ13D_Evaluate_dx(GNi_p,GNx_p,coords,&J_p);
581:     fac = gp_weight[p]*J_p;
582:     /*
583:      for (i = 0; i < NODES_PER_EL; i++) {
584:      for (j = i; j < NODES_PER_EL; j++) {
585:      Ke[NODES_PER_EL*i+j] -= fac*(Ni_p[i]*Ni_p[j]-0.015625);
586:      }
587:      }
588:      */

590:     for (i = 0; i < NODES_PER_EL; i++) {
591:       for (j = 0; j < NODES_PER_EL; j++) {
592:         Ke[NODES_PER_EL*i+j] -= fac*(Ni_p[i]*Ni_p[j]-0.015625);
593:       }
594:     }
595:   }

597:   /* scale */
598:   eta_avg = 0.0;
599:   for (p = 0; p < ngp; p++) eta_avg += eta[p];
600:   eta_avg = (1.0/((PetscScalar)ngp))*eta_avg;
601:   fac     = 1.0/eta_avg;

603:   /*
604:    for (i = 0; i < NODES_PER_EL; i++) {
605:    for (j = i; j < NODES_PER_EL; j++) {
606:    Ke[NODES_PER_EL*i+j] = fac*Ke[NODES_PER_EL*i+j];
607:    #if defined(ASSEMBLE_LOWER_TRIANGULAR)
608:    Ke[NODES_PER_EL*j+i] = Ke[NODES_PER_EL*i+j];
609:    #endif
610:    }
611:    }
612:    */

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

621: static void FormScaledMassMatrixOperatorQ13D(PetscScalar Ke[],PetscScalar coords[],PetscScalar eta[])
622: {
623:   PetscInt    ngp;
624:   PetscScalar gp_xi[GAUSS_POINTS][NSD];
625:   PetscScalar gp_weight[GAUSS_POINTS];
626:   PetscInt    p,i,j;
627:   PetscScalar Ni_p[NODES_PER_EL];
628:   PetscScalar GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
629:   PetscScalar J_p,fac,eta_avg;

631:   /* define quadrature rule */
632:   ConstructGaussQuadrature3D(&ngp,gp_xi,gp_weight);

634:   /* evaluate integral */
635:   for (p = 0; p < ngp; p++) {
636:     ShapeFunctionQ13D_Evaluate(gp_xi[p],Ni_p);
637:     ShapeFunctionQ13D_Evaluate_dxi(gp_xi[p],GNi_p);
638:     ShapeFunctionQ13D_Evaluate_dx(GNi_p,GNx_p,coords,&J_p);
639:     fac = gp_weight[p]*J_p;

641:     /*
642:      for (i = 0; i < NODES_PER_EL; i++) {
643:      for (j = i; j < NODES_PER_EL; j++) {
644:      Ke[NODES_PER_EL*i+j] = Ke[NODES_PER_EL*i+j]-fac*Ni_p[i]*Ni_p[j];
645:      }
646:      }
647:      */

649:     for (i = 0; i < NODES_PER_EL; i++) {
650:       for (j = 0; j < NODES_PER_EL; j++) {
651:         Ke[NODES_PER_EL*i+j] = Ke[NODES_PER_EL*i+j]-fac*Ni_p[i]*Ni_p[j];
652:       }
653:     }
654:   }

656:   /* scale */
657:   eta_avg = 0.0;
658:   for (p = 0; p < ngp; p++) eta_avg += eta[p];
659:   eta_avg = (1.0/((PetscScalar)ngp))*eta_avg;
660:   fac     = 1.0/eta_avg;
661:   /*
662:    for (i = 0; i < NODES_PER_EL; i++) {
663:    for (j = i; j < NODES_PER_EL; j++) {
664:    Ke[NODES_PER_EL*i+j] = fac*Ke[NODES_PER_EL*i+j];
665:    #if defined(ASSEMBLE_LOWER_TRIANGULAR)
666:    Ke[NODES_PER_EL*j+i] = Ke[NODES_PER_EL*i+j];
667:    #endif
668:    }
669:    }
670:    */

672:   for (i = 0; i < NODES_PER_EL; i++) {
673:     for (j = 0; j < NODES_PER_EL; j++) {
674:       Ke[NODES_PER_EL*i+j] = fac*Ke[NODES_PER_EL*i+j];
675:     }
676:   }
677: }

679: static void FormMomentumRhsQ13D(PetscScalar Fe[],PetscScalar coords[],PetscScalar fx[],PetscScalar fy[],PetscScalar fz[])
680: {
681:   PetscInt    ngp;
682:   PetscScalar gp_xi[GAUSS_POINTS][NSD];
683:   PetscScalar gp_weight[GAUSS_POINTS];
684:   PetscInt    p,i;
685:   PetscScalar Ni_p[NODES_PER_EL];
686:   PetscScalar GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
687:   PetscScalar J_p,fac;

689:   /* define quadrature rule */
690:   ConstructGaussQuadrature3D(&ngp,gp_xi,gp_weight);

692:   /* evaluate integral */
693:   for (p = 0; p < ngp; p++) {
694:     ShapeFunctionQ13D_Evaluate(gp_xi[p],Ni_p);
695:     ShapeFunctionQ13D_Evaluate_dxi(gp_xi[p],GNi_p);
696:     ShapeFunctionQ13D_Evaluate_dx(GNi_p,GNx_p,coords,&J_p);
697:     fac = gp_weight[p]*J_p;

699:     for (i = 0; i < NODES_PER_EL; i++) {
700:       Fe[NSD*i]   -= fac*Ni_p[i]*fx[p];
701:       Fe[NSD*i+1] -= fac*Ni_p[i]*fy[p];
702:       Fe[NSD*i+2] -= fac*Ni_p[i]*fz[p];
703:     }
704:   }
705: }

707: static void FormContinuityRhsQ13D(PetscScalar Fe[],PetscScalar coords[],PetscScalar hc[])
708: {
709:   PetscInt    ngp;
710:   PetscScalar gp_xi[GAUSS_POINTS][NSD];
711:   PetscScalar gp_weight[GAUSS_POINTS];
712:   PetscInt    p,i;
713:   PetscScalar Ni_p[NODES_PER_EL];
714:   PetscScalar GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
715:   PetscScalar J_p,fac;

717:   /* define quadrature rule */
718:   ConstructGaussQuadrature3D(&ngp,gp_xi,gp_weight);

720:   /* evaluate integral */
721:   for (p = 0; p < ngp; p++) {
722:     ShapeFunctionQ13D_Evaluate(gp_xi[p],Ni_p);
723:     ShapeFunctionQ13D_Evaluate_dxi(gp_xi[p],GNi_p);
724:     ShapeFunctionQ13D_Evaluate_dx(GNi_p,GNx_p,coords,&J_p);
725:     fac = gp_weight[p]*J_p;

727:     for (i = 0; i < NODES_PER_EL; i++) Fe[i] -= fac*Ni_p[i]*hc[p];
728:   }
729: }

731: #define _ZERO_ROWCOL_i(A,i) {                   \
732:     PetscInt    KK;                             \
733:     PetscScalar tmp = A[24*(i)+(i)];            \
734:     for (KK=0;KK<24;KK++) A[24*(i)+KK]=0.0;     \
735:     for (KK=0;KK<24;KK++) A[24*KK+(i)]=0.0;     \
736:     A[24*(i)+(i)] = tmp;}                       \

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

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

746: static PetscErrorCode AssembleA_Stokes(Mat A,DM stokes_da,CellProperties cell_properties)
747: {
748:   DM                     cda;
749:   Vec                    coords;
750:   DMDACoor3d             ***_coords;
751:   MatStencil             u_eqn[NODES_PER_EL*U_DOFS];
752:   MatStencil             p_eqn[NODES_PER_EL*P_DOFS];
753:   PetscInt               sex,sey,sez,mx,my,mz;
754:   PetscInt               ei,ej,ek;
755:   PetscScalar            Ae[NODES_PER_EL*U_DOFS*NODES_PER_EL*U_DOFS];
756:   PetscScalar            Ge[NODES_PER_EL*U_DOFS*NODES_PER_EL*P_DOFS];
757:   PetscScalar            De[NODES_PER_EL*P_DOFS*NODES_PER_EL*U_DOFS];
758:   PetscScalar            Ce[NODES_PER_EL*P_DOFS*NODES_PER_EL*P_DOFS];
759:   PetscScalar            el_coords[NODES_PER_EL*NSD];
760:   GaussPointCoefficients *props;
761:   PetscScalar            *prop_eta;
762:   PetscInt               n,M,N,P;
763:   PetscErrorCode         ierr;

766:   DMDAGetInfo(stokes_da,0,&M,&N,&P,0,0,0, 0,0,0,0,0,0);
767:   /* setup for coords */
768:   DMGetCoordinateDM(stokes_da,&cda);
769:   DMGetCoordinatesLocal(stokes_da,&coords);
770:   DMDAVecGetArray(cda,coords,&_coords);
771:   DMDAGetElementsCorners(stokes_da,&sex,&sey,&sez);
772:   DMDAGetElementsSizes(stokes_da,&mx,&my,&mz);
773:   for (ek = sez; ek < sez+mz; ek++) {
774:     for (ej = sey; ej < sey+my; ej++) {
775:       for (ei = sex; ei < sex+mx; ei++) {
776:         /* get coords for the element */
777:         GetElementCoords3D(_coords,ei,ej,ek,el_coords);
778:         /* get cell properties */
779:         CellPropertiesGetCell(cell_properties,ei,ej,ek,&props);
780:         /* get coefficients for the element */
781:         prop_eta = props->eta;

783:         /* initialise element stiffness matrix */
784:         PetscMemzero(Ae,sizeof(Ae));
785:         PetscMemzero(Ge,sizeof(Ge));
786:         PetscMemzero(De,sizeof(De));
787:         PetscMemzero(Ce,sizeof(Ce));

789:         /* form element stiffness matrix */
790:         FormStressOperatorQ13D(Ae,el_coords,prop_eta);
791:         FormGradientOperatorQ13D(Ge,el_coords);
792:         /*#if defined(ASSEMBLE_LOWER_TRIANGULAR)*/
793:         FormDivergenceOperatorQ13D(De,el_coords);
794:         /*#endif*/
795:         FormStabilisationOperatorQ13D(Ce,el_coords,prop_eta);

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

800:         for (n=0; n<NODES_PER_EL; n++) {
801:           if ((u_eqn[3*n].i == 0) || (u_eqn[3*n].i == M-1)) {
802:             _ZERO_ROWCOL_i(Ae,3*n);
803:             _ZERO_ROW_i   (Ge,3*n);
804:             _ZERO_COL_i   (De,3*n);
805:           }

807:           if ((u_eqn[3*n+1].j == 0) || (u_eqn[3*n+1].j == N-1)) {
808:             _ZERO_ROWCOL_i(Ae,3*n+1);
809:             _ZERO_ROW_i   (Ge,3*n+1);
810:             _ZERO_COL_i   (De,3*n+1);
811:           }

813:           if ((u_eqn[3*n+2].k == 0) || (u_eqn[3*n+2].k == P-1)) {
814:             _ZERO_ROWCOL_i(Ae,3*n+2);
815:             _ZERO_ROW_i   (Ge,3*n+2);
816:             _ZERO_COL_i   (De,3*n+2);
817:           }
818:         }
819:         MatSetValuesStencil(A,NODES_PER_EL*U_DOFS,u_eqn,NODES_PER_EL*U_DOFS,u_eqn,Ae,ADD_VALUES);
820:         MatSetValuesStencil(A,NODES_PER_EL*U_DOFS,u_eqn,NODES_PER_EL*P_DOFS,p_eqn,Ge,ADD_VALUES);
821:         MatSetValuesStencil(A,NODES_PER_EL*P_DOFS,p_eqn,NODES_PER_EL*U_DOFS,u_eqn,De,ADD_VALUES);
822:         MatSetValuesStencil(A,NODES_PER_EL*P_DOFS,p_eqn,NODES_PER_EL*P_DOFS,p_eqn,Ce,ADD_VALUES);
823:       }
824:     }
825:   }
826:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
827:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

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

831:   return(0);
832: }

834: static PetscErrorCode AssembleA_PCStokes(Mat A,DM stokes_da,CellProperties cell_properties)
835: {
836:   DM                     cda;
837:   Vec                    coords;
838:   DMDACoor3d             ***_coords;
839:   MatStencil             u_eqn[NODES_PER_EL*U_DOFS];
840:   MatStencil             p_eqn[NODES_PER_EL*P_DOFS];
841:   PetscInt               sex,sey,sez,mx,my,mz;
842:   PetscInt               ei,ej,ek;
843:   PetscScalar            Ae[NODES_PER_EL*U_DOFS*NODES_PER_EL*U_DOFS];
844:   PetscScalar            Ge[NODES_PER_EL*U_DOFS*NODES_PER_EL*P_DOFS];
845:   PetscScalar            De[NODES_PER_EL*P_DOFS*NODES_PER_EL*U_DOFS];
846:   PetscScalar            Ce[NODES_PER_EL*P_DOFS*NODES_PER_EL*P_DOFS];
847:   PetscScalar            el_coords[NODES_PER_EL*NSD];
848:   GaussPointCoefficients *props;
849:   PetscScalar            *prop_eta;
850:   PetscInt               n,M,N,P;
851:   PetscErrorCode         ierr;

854:   DMDAGetInfo(stokes_da,0,&M,&N,&P,0,0,0, 0,0,0,0,0,0);
855:   /* setup for coords */
856:   DMGetCoordinateDM(stokes_da,&cda);
857:   DMGetCoordinatesLocal(stokes_da,&coords);
858:   DMDAVecGetArray(cda,coords,&_coords);

860:   DMDAGetElementsCorners(stokes_da,&sex,&sey,&sez);
861:   DMDAGetElementsSizes(stokes_da,&mx,&my,&mz);
862:   for (ek = sez; ek < sez+mz; ek++) {
863:     for (ej = sey; ej < sey+my; ej++) {
864:       for (ei = sex; ei < sex+mx; ei++) {
865:         /* get coords for the element */
866:         GetElementCoords3D(_coords,ei,ej,ek,el_coords);
867:         /* get cell properties */
868:         CellPropertiesGetCell(cell_properties,ei,ej,ek,&props);
869:         /* get coefficients for the element */
870:         prop_eta = props->eta;

872:         /* initialise element stiffness matrix */
873:         PetscMemzero(Ae,sizeof(Ae));
874:         PetscMemzero(Ge,sizeof(Ge));
875:         PetscMemzero(De,sizeof(De));
876:         PetscMemzero(Ce,sizeof(Ce));

878:         /* form element stiffness matrix */
879:         FormStressOperatorQ13D(Ae,el_coords,prop_eta);
880:         FormGradientOperatorQ13D(Ge,el_coords);
881:         /* FormDivergenceOperatorQ13D(De,el_coords); */
882:         FormScaledMassMatrixOperatorQ13D(Ce,el_coords,prop_eta);

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

887:         for (n=0; n<NODES_PER_EL; n++) {
888:           if ((u_eqn[3*n].i == 0) || (u_eqn[3*n].i == M-1)) {
889:             _ZERO_ROWCOL_i(Ae,3*n);
890:             _ZERO_ROW_i   (Ge,3*n);
891:             _ZERO_COL_i   (De,3*n);
892:           }

894:           if ((u_eqn[3*n+1].j == 0) || (u_eqn[3*n+1].j == N-1)) {
895:             _ZERO_ROWCOL_i(Ae,3*n+1);
896:             _ZERO_ROW_i   (Ge,3*n+1);
897:             _ZERO_COL_i   (De,3*n+1);
898:           }

900:           if ((u_eqn[3*n+2].k == 0) || (u_eqn[3*n+2].k == P-1)) {
901:             _ZERO_ROWCOL_i(Ae,3*n+2);
902:             _ZERO_ROW_i   (Ge,3*n+2);
903:             _ZERO_COL_i   (De,3*n+2);
904:           }
905:         }
906:         MatSetValuesStencil(A,NODES_PER_EL*U_DOFS,u_eqn,NODES_PER_EL*U_DOFS,u_eqn,Ae,ADD_VALUES);
907:         MatSetValuesStencil(A,NODES_PER_EL*U_DOFS,u_eqn,NODES_PER_EL*P_DOFS,p_eqn,Ge,ADD_VALUES);
908:         /*MatSetValuesStencil(A,NODES_PER_EL*P_DOFS,p_eqn,NODES_PER_EL*U_DOFS,u_eqn,De,ADD_VALUES);*/
909:         MatSetValuesStencil(A,NODES_PER_EL*P_DOFS,p_eqn,NODES_PER_EL*P_DOFS,p_eqn,Ce,ADD_VALUES);
910:       }
911:     }
912:   }
913:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
914:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

916:   DMDAVecRestoreArray(cda,coords,&_coords);
917:   return(0);
918: }

920: static PetscErrorCode AssembleF_Stokes(Vec F,DM stokes_da,CellProperties cell_properties)
921: {
922:   DM                     cda;
923:   Vec                    coords;
924:   DMDACoor3d             ***_coords;
925:   MatStencil             u_eqn[NODES_PER_EL*U_DOFS];
926:   MatStencil             p_eqn[NODES_PER_EL*P_DOFS];
927:   PetscInt               sex,sey,sez,mx,my,mz;
928:   PetscInt               ei,ej,ek;
929:   PetscScalar            Fe[NODES_PER_EL*U_DOFS];
930:   PetscScalar            He[NODES_PER_EL*P_DOFS];
931:   PetscScalar            el_coords[NODES_PER_EL*NSD];
932:   GaussPointCoefficients *props;
933:   PetscScalar            *prop_fx,*prop_fy,*prop_fz,*prop_hc;
934:   Vec                    local_F;
935:   StokesDOF              ***ff;
936:   PetscInt               n,M,N,P;
937:   PetscErrorCode         ierr;

940:   DMDAGetInfo(stokes_da,0,&M,&N,&P,0,0,0, 0,0,0,0,0,0);
941:   /* setup for coords */
942:   DMGetCoordinateDM(stokes_da,&cda);
943:   DMGetCoordinatesLocal(stokes_da,&coords);
944:   DMDAVecGetArray(cda,coords,&_coords);

946:   /* get access to the vector */
947:   DMGetLocalVector(stokes_da,&local_F);
948:   VecZeroEntries(local_F);
949:   DMDAVecGetArray(stokes_da,local_F,&ff);
950:   DMDAGetElementsCorners(stokes_da,&sex,&sey,&sez);
951:   DMDAGetElementsSizes(stokes_da,&mx,&my,&mz);
952:   for (ek = sez; ek < sez+mz; ek++) {
953:     for (ej = sey; ej < sey+my; ej++) {
954:       for (ei = sex; ei < sex+mx; ei++) {
955:         /* get coords for the element */
956:         GetElementCoords3D(_coords,ei,ej,ek,el_coords);
957:         /* get cell properties */
958:         CellPropertiesGetCell(cell_properties,ei,ej,ek,&props);
959:         /* get coefficients for the element */
960:         prop_fx = props->fx;
961:         prop_fy = props->fy;
962:         prop_fz = props->fz;
963:         prop_hc = props->hc;

965:         /* initialise element stiffness matrix */
966:         PetscMemzero(Fe,sizeof(Fe));
967:         PetscMemzero(He,sizeof(He));

969:         /* form element stiffness matrix */
970:         FormMomentumRhsQ13D(Fe,el_coords,prop_fx,prop_fy,prop_fz);
971:         FormContinuityRhsQ13D(He,el_coords,prop_hc);

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

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

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

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

984:         DMDASetValuesLocalStencil3D_ADD_VALUES(ff,u_eqn,p_eqn,Fe,He);
985:       }
986:     }
987:   }
988:   DMDAVecRestoreArray(stokes_da,local_F,&ff);
989:   DMLocalToGlobalBegin(stokes_da,local_F,ADD_VALUES,F);
990:   DMLocalToGlobalEnd(stokes_da,local_F,ADD_VALUES,F);
991:   DMRestoreLocalVector(stokes_da,&local_F);

993:   DMDAVecRestoreArray(cda,coords,&_coords);
994:   return(0);
995: }

997: static void evaluate_MS_FrankKamentski_constants(PetscReal *theta,PetscReal *MX,PetscReal *MY,PetscReal *MZ)
998: {
999:   *theta = 0.0;
1000:   *MX    = 2.0 * PETSC_PI;
1001:   *MY    = 2.0 * PETSC_PI;
1002:   *MZ    = 2.0 * PETSC_PI;
1003: }
1004: static void evaluate_MS_FrankKamentski(PetscReal pos[],PetscReal v[],PetscReal *p,PetscReal *eta,PetscReal Fm[],PetscReal *Fc)
1005: {
1006:   PetscReal x,y,z;
1007:   PetscReal theta,MX,MY,MZ;

1009:   evaluate_MS_FrankKamentski_constants(&theta,&MX,&MY,&MZ);
1010:   x = pos[0];
1011:   y = pos[1];
1012:   z = pos[2];
1013:   if (v) {
1014:     /*
1015:      v[0] = PetscPowRealInt(z,3)*PetscExpReal(y)*PetscSinReal(PETSC_PI*x);
1016:      v[1] = z*cos(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y);
1017:      v[2] = -(PetscPowRealInt(x,3) + PetscPowRealInt(y,3))*PetscSinReal(2.0*PETSC_PI*z);
1018:      */
1019:     /*
1020:      v[0] = PetscSinReal(PETSC_PI*x);
1021:      v[1] = PetscSinReal(PETSC_PI*y);
1022:      v[2] = PetscSinReal(PETSC_PI*z);
1023:      */
1024:     v[0] = PetscPowRealInt(z,3)*PetscExpReal(y)*PetscSinReal(PETSC_PI*x);
1025:     v[1] = z*PetscCosReal(2.0*PETSC_PI*x)*PetscExpReal(-y)*PetscSinReal(PETSC_PI*y);
1026:     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;
1027:   }
1028:   if (p) *p = PetscPowRealInt(x,2) + PetscPowRealInt(y,2) + PetscPowRealInt(z,2);
1029:   if (eta) {
1030:     /* eta = PetscExpReal(-theta*(1.0 - y - 0.1*PetscCosReal(MX*x)*PetscCosReal(MZ*z)*PetscSinReal(MY*y)));*/
1031:     *eta = 1.0;
1032:   }
1033:   if (Fm) {
1034:     /*
1035:      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))) ;
1036:      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)));
1037:      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)))  ;
1038:      */
1039:     /*
1040:      Fm[0]=-2*x - 2.0*PetscPowRealInt(PETSC_PI,2)*PetscSinReal(PETSC_PI*x);
1041:      Fm[1]=-2*y - 2.0*PetscPowRealInt(PETSC_PI,2)*PetscSinReal(PETSC_PI*y);
1042:      Fm[2]=-2*z - 2.0*PetscPowRealInt(PETSC_PI,2)*PetscSinReal(PETSC_PI*z);
1043:      */
1044:     /*
1045:      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) ;
1046:      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);
1047:      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) ;
1048:      */
1049:     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);
1050:     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);
1051:     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);
1052:   }
1053:   if (Fc) {
1054:     /* 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) ;*/
1055:     /* Fc = PETSC_PI*PetscCosReal(PETSC_PI*x) + PETSC_PI*PetscCosReal(PETSC_PI*y) + PETSC_PI*PetscCosReal(PETSC_PI*z);*/
1056:     /* 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);*/
1057:     *Fc = 0.0;
1058:   }
1059: }

1061: static PetscErrorCode DMDACreateManufacturedSolution(PetscInt mx,PetscInt my,PetscInt mz,DM *_da,Vec *_X)
1062: {
1063:   DM             da,cda;
1064:   Vec            X;
1065:   StokesDOF      ***_stokes;
1066:   Vec            coords;
1067:   DMDACoor3d     ***_coords;
1068:   PetscInt       si,sj,sk,ei,ej,ek,i,j,k;

1072:   DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,
1073:                       mx+1,my+1,mz+1,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,4,1,NULL,NULL,NULL,&da);
1074:   DMSetFromOptions(da);
1075:   DMSetUp(da);
1076:   DMDASetFieldName(da,0,"anlytic_Vx");
1077:   DMDASetFieldName(da,1,"anlytic_Vy");
1078:   DMDASetFieldName(da,2,"anlytic_Vz");
1079:   DMDASetFieldName(da,3,"analytic_P");

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

1083:   DMGetCoordinatesLocal(da,&coords);
1084:   DMGetCoordinateDM(da,&cda);
1085:   DMDAVecGetArray(cda,coords,&_coords);

1087:   DMCreateGlobalVector(da,&X);
1088:   DMDAVecGetArray(da,X,&_stokes);

1090:   DMDAGetCorners(da,&si,&sj,&sk,&ei,&ej,&ek);
1091:   for (k = sk; k < sk+ek; k++) {
1092:     for (j = sj; j < sj+ej; j++) {
1093:       for (i = si; i < si+ei; i++) {
1094:         PetscReal pos[NSD],pressure,vel[NSD];

1096:         pos[0] = PetscRealPart(_coords[k][j][i].x);
1097:         pos[1] = PetscRealPart(_coords[k][j][i].y);
1098:         pos[2] = PetscRealPart(_coords[k][j][i].z);

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

1102:         _stokes[k][j][i].u_dof = vel[0];
1103:         _stokes[k][j][i].v_dof = vel[1];
1104:         _stokes[k][j][i].w_dof = vel[2];
1105:         _stokes[k][j][i].p_dof = pressure;
1106:       }
1107:     }
1108:   }
1109:   DMDAVecRestoreArray(da,X,&_stokes);
1110:   DMDAVecRestoreArray(cda,coords,&_coords);

1112:   *_da = da;
1113:   *_X  = X;
1114:   return(0);
1115: }

1117: static PetscErrorCode DMDAIntegrateErrors3D(DM stokes_da,Vec X,Vec X_analytic)
1118: {
1119:   DM          cda;
1120:   Vec         coords,X_analytic_local,X_local;
1121:   DMDACoor3d  ***_coords;
1122:   PetscInt    sex,sey,sez,mx,my,mz;
1123:   PetscInt    ei,ej,ek;
1124:   PetscScalar el_coords[NODES_PER_EL*NSD];
1125:   StokesDOF   ***stokes_analytic,***stokes;
1126:   StokesDOF   stokes_analytic_e[NODES_PER_EL],stokes_e[NODES_PER_EL];

1128:   PetscScalar    GNi_p[NSD][NODES_PER_EL],GNx_p[NSD][NODES_PER_EL];
1129:   PetscScalar    Ni_p[NODES_PER_EL];
1130:   PetscInt       ngp;
1131:   PetscScalar    gp_xi[GAUSS_POINTS][NSD];
1132:   PetscScalar    gp_weight[GAUSS_POINTS];
1133:   PetscInt       p,i;
1134:   PetscScalar    J_p,fac;
1135:   PetscScalar    h,p_e_L2,u_e_L2,u_e_H1,p_L2,u_L2,u_H1,tp_L2,tu_L2,tu_H1;
1136:   PetscScalar    tint_p_ms,tint_p,int_p_ms,int_p;
1137:   PetscInt       M;
1138:   PetscReal      xymin[NSD],xymax[NSD];

1142:   /* define quadrature rule */
1143:   ConstructGaussQuadrature3D(&ngp,gp_xi,gp_weight);

1145:   /* setup for coords */
1146:   DMGetCoordinateDM(stokes_da,&cda);
1147:   DMGetCoordinatesLocal(stokes_da,&coords);
1148:   DMDAVecGetArray(cda,coords,&_coords);

1150:   /* setup for analytic */
1151:   DMCreateLocalVector(stokes_da,&X_analytic_local);
1152:   DMGlobalToLocalBegin(stokes_da,X_analytic,INSERT_VALUES,X_analytic_local);
1153:   DMGlobalToLocalEnd(stokes_da,X_analytic,INSERT_VALUES,X_analytic_local);
1154:   DMDAVecGetArray(stokes_da,X_analytic_local,&stokes_analytic);

1156:   /* setup for solution */
1157:   DMCreateLocalVector(stokes_da,&X_local);
1158:   DMGlobalToLocalBegin(stokes_da,X,INSERT_VALUES,X_local);
1159:   DMGlobalToLocalEnd(stokes_da,X,INSERT_VALUES,X_local);
1160:   DMDAVecGetArray(stokes_da,X_local,&stokes);

1162:   DMDAGetInfo(stokes_da,0,&M,0,0,0,0,0,0,0,0,0,0,0);
1163:   DMGetBoundingBox(stokes_da,xymin,xymax);

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

1167:   tp_L2     = tu_L2 = tu_H1 = 0.0;
1168:   tint_p_ms = tint_p = 0.0;

1170:   DMDAGetElementsCorners(stokes_da,&sex,&sey,&sez);
1171:   DMDAGetElementsSizes(stokes_da,&mx,&my,&mz);
1172:   for (ek = sez; ek < sez+mz; ek++) {
1173:     for (ej = sey; ej < sey+my; ej++) {
1174:       for (ei = sex; ei < sex+mx; ei++) {
1175:         /* get coords for the element */
1176:         GetElementCoords3D(_coords,ei,ej,ek,el_coords);
1177:         StokesDAGetNodalFields3D(stokes,ei,ej,ek,stokes_e);
1178:         StokesDAGetNodalFields3D(stokes_analytic,ei,ej,ek,stokes_analytic_e);

1180:         /* evaluate integral */
1181:         for (p = 0; p < ngp; p++) {
1182:           ShapeFunctionQ13D_Evaluate(gp_xi[p],Ni_p);
1183:           ShapeFunctionQ13D_Evaluate_dxi(gp_xi[p],GNi_p);
1184:           ShapeFunctionQ13D_Evaluate_dx(GNi_p,GNx_p,el_coords,&J_p);
1185:           fac = gp_weight[p]*J_p;

1187:           for (i = 0; i < NODES_PER_EL; i++) {
1188:             tint_p_ms = tint_p_ms+fac*Ni_p[i]*stokes_analytic_e[i].p_dof;
1189:             tint_p    = tint_p   +fac*Ni_p[i]*stokes_e[i].p_dof;
1190:           }
1191:         }
1192:       }
1193:     }
1194:   }
1195:   MPI_Allreduce(&tint_p_ms,&int_p_ms,1,MPIU_SCALAR,MPIU_SUM,PETSC_COMM_WORLD);
1196:   MPI_Allreduce(&tint_p,&int_p,1,MPIU_SCALAR,MPIU_SUM,PETSC_COMM_WORLD);

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

1200:   /* remove mine and add manufacture one */
1201:   DMDAVecRestoreArray(stokes_da,X_analytic_local,&stokes_analytic);
1202:   DMDAVecRestoreArray(stokes_da,X_local,&stokes);

1204:   {
1205:     PetscInt    k,L,dof;
1206:     PetscScalar *fields;
1207:     DMDAGetInfo(stokes_da,0, 0,0,0, 0,0,0, &dof,0,0,0,0,0);

1209:     VecGetLocalSize(X_local,&L);
1210:     VecGetArray(X_local,&fields);
1211:     for (k=0; k<L/dof; k++) fields[dof*k+3] = fields[dof*k+3] - int_p + int_p_ms;
1212:     VecRestoreArray(X_local,&fields);

1214:     VecGetLocalSize(X,&L);
1215:     VecGetArray(X,&fields);
1216:     for (k=0; k<L/dof; k++) fields[dof*k+3] = fields[dof*k+3] - int_p + int_p_ms;
1217:     VecRestoreArray(X,&fields);
1218:   }

1220:   DMDAVecGetArray(stokes_da,X_local,&stokes);
1221:   DMDAVecGetArray(stokes_da,X_analytic_local,&stokes_analytic);

1223:   for (ek = sez; ek < sez+mz; ek++) {
1224:     for (ej = sey; ej < sey+my; ej++) {
1225:       for (ei = sex; ei < sex+mx; ei++) {
1226:         /* get coords for the element */
1227:         GetElementCoords3D(_coords,ei,ej,ek,el_coords);
1228:         StokesDAGetNodalFields3D(stokes,ei,ej,ek,stokes_e);
1229:         StokesDAGetNodalFields3D(stokes_analytic,ei,ej,ek,stokes_analytic_e);

1231:         /* evaluate integral */
1232:         p_e_L2 = 0.0;
1233:         u_e_L2 = 0.0;
1234:         u_e_H1 = 0.0;
1235:         for (p = 0; p < ngp; p++) {
1236:           ShapeFunctionQ13D_Evaluate(gp_xi[p],Ni_p);
1237:           ShapeFunctionQ13D_Evaluate_dxi(gp_xi[p],GNi_p);
1238:           ShapeFunctionQ13D_Evaluate_dx(GNi_p,GNx_p,el_coords,&J_p);
1239:           fac = gp_weight[p]*J_p;

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

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

1246:             u_error = stokes_e[i].u_dof-stokes_analytic_e[i].u_dof;
1247:             v_error = stokes_e[i].v_dof-stokes_analytic_e[i].v_dof;
1248:             w_error = stokes_e[i].w_dof-stokes_analytic_e[i].w_dof;
1249:             u_e_L2 += fac*Ni_p[i]*(u_error*u_error+v_error*v_error+w_error*w_error);

1251:             u_e_H1 = u_e_H1+fac*(GNx_p[0][i]*u_error*GNx_p[0][i]*u_error              /* du/dx */
1252:                                  +GNx_p[1][i]*u_error*GNx_p[1][i]*u_error
1253:                                  +GNx_p[2][i]*u_error*GNx_p[2][i]*u_error
1254:                                  +GNx_p[0][i]*v_error*GNx_p[0][i]*v_error               /* dv/dx */
1255:                                  +GNx_p[1][i]*v_error*GNx_p[1][i]*v_error
1256:                                  +GNx_p[2][i]*v_error*GNx_p[2][i]*v_error
1257:                                  +GNx_p[0][i]*w_error*GNx_p[0][i]*w_error               /* dw/dx */
1258:                                  +GNx_p[1][i]*w_error*GNx_p[1][i]*w_error
1259:                                  +GNx_p[2][i]*w_error*GNx_p[2][i]*w_error);
1260:           }
1261:         }

1263:         tp_L2 += p_e_L2;
1264:         tu_L2 += u_e_L2;
1265:         tu_H1 += u_e_H1;
1266:       }
1267:     }
1268:   }
1269:   MPI_Allreduce(&tp_L2,&p_L2,1,MPIU_SCALAR,MPIU_SUM,PETSC_COMM_WORLD);
1270:   MPI_Allreduce(&tu_L2,&u_L2,1,MPIU_SCALAR,MPIU_SUM,PETSC_COMM_WORLD);
1271:   MPI_Allreduce(&tu_H1,&u_H1,1,MPIU_SCALAR,MPIU_SUM,PETSC_COMM_WORLD);
1272:   p_L2 = PetscSqrtScalar(p_L2);
1273:   u_L2 = PetscSqrtScalar(u_L2);
1274:   u_H1 = PetscSqrtScalar(u_H1);

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

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

1280:   DMDAVecRestoreArray(stokes_da,X_analytic_local,&stokes_analytic);
1281:   VecDestroy(&X_analytic_local);
1282:   DMDAVecRestoreArray(stokes_da,X_local,&stokes);
1283:   VecDestroy(&X_local);
1284:   return(0);
1285: }

1287: PetscErrorCode DAView_3DVTK_StructuredGrid_appended(DM da,Vec FIELD,const char file_prefix[])
1288: {
1289:   char           vtk_filename[PETSC_MAX_PATH_LEN];
1290:   PetscMPIInt    rank;
1291:   MPI_Comm       comm;
1292:   FILE           *vtk_fp = NULL;
1293:   const char     *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian";
1294:   PetscInt       si,sj,sk,nx,ny,nz,i;
1295:   PetscInt       f,n_fields,N;
1296:   DM             cda;
1297:   Vec            coords;
1298:   Vec            l_FIELD;
1299:   PetscScalar    *_L_FIELD;
1300:   PetscInt       memory_offset;
1301:   PetscScalar    *buffer;


1306:   /* create file name */
1307:   PetscObjectGetComm((PetscObject)da,&comm);
1308:   MPI_Comm_rank(comm,&rank);
1309:   PetscSNPrintf(vtk_filename,sizeof(vtk_filename),"subdomain-%s-p%1.4d.vts",file_prefix,rank);

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

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

1317:   /* coords */
1318:   DMDAGetGhostCorners(da,&si,&sj,&sk,&nx,&ny,&nz);
1319:   N    = nx * ny * nz;

1321:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"%s\">\n",byte_order);
1322:   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);
1323:   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);

1325:   memory_offset = 0;

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

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

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

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

1339:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"      <PointData Scalars=\" ");
1340:   DMDAGetInfo(da,0,0,0,0,0,0,0,&n_fields,0,0,0,0,0);
1341:   for (f=0; f<n_fields; f++) {
1342:     const char *field_name;
1343:     DMDAGetFieldName(da,f,&field_name);
1344:     PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"%s ",field_name);
1345:   }
1346:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"\">\n");

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

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

1356:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"      </PointData>\n");
1357:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"    </Piece>\n");
1358:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"  </StructuredGrid>\n");

1360:   PetscMalloc1(N,&buffer);
1361:   DMGetLocalVector(da,&l_FIELD);
1362:   DMGlobalToLocalBegin(da, FIELD,INSERT_VALUES,l_FIELD);
1363:   DMGlobalToLocalEnd(da,FIELD,INSERT_VALUES,l_FIELD);
1364:   VecGetArray(l_FIELD,&_L_FIELD);

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

1369:   /* write coordinates */
1370:   {
1371:     int         length = sizeof(PetscScalar)*N*3;
1372:     PetscScalar *allcoords;

1374:     fwrite(&length,sizeof(int),1,vtk_fp);
1375:     VecGetArray(coords,&allcoords);
1376:     fwrite(allcoords,sizeof(PetscScalar),3*N,vtk_fp);
1377:     VecRestoreArray(coords,&allcoords);
1378:   }
1379:   /* write fields */
1380:   for (f=0; f<n_fields; f++) {
1381:     int length = sizeof(PetscScalar)*N;
1382:     fwrite(&length,sizeof(int),1,vtk_fp);
1383:     /* load */
1384:     for (i=0; i<N; i++) buffer[i] = _L_FIELD[n_fields*i + f];

1386:     /* write */
1387:     fwrite(buffer,sizeof(PetscScalar),N,vtk_fp);
1388:   }
1389:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"\n  </AppendedData>\n");

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

1393:   PetscFree(buffer);
1394:   VecRestoreArray(l_FIELD,&_L_FIELD);
1395:   DMRestoreLocalVector(da,&l_FIELD);

1397:   if (vtk_fp) {
1398:     fclose(vtk_fp);
1399:     vtk_fp = NULL;
1400:   }

1402:   return(0);
1403: }

1405: PetscErrorCode DAViewVTK_write_PieceExtend(FILE *vtk_fp,PetscInt indent_level,DM da,const char local_file_prefix[])
1406: {
1407:   PetscMPIInt    size,rank;
1408:   MPI_Comm       comm;
1409:   const PetscInt *lx,*ly,*lz;
1410:   PetscInt       M,N,P,pM,pN,pP,sum,*olx,*oly,*olz;
1411:   PetscInt       *osx,*osy,*osz,*oex,*oey,*oez;
1412:   PetscInt       i,j,k,II,stencil;

1416:   /* create file name */
1417:   PetscObjectGetComm((PetscObject)da,&comm);
1418:   MPI_Comm_size(comm,&size);
1419:   MPI_Comm_rank(comm,&rank);

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

1424:   /* generate start,end list */
1425:   PetscMalloc1(pM+1,&olx);
1426:   PetscMalloc1(pN+1,&oly);
1427:   PetscMalloc1(pP+1,&olz);
1428:   sum  = 0;
1429:   for (i=0; i<pM; i++) {
1430:     olx[i] = sum;
1431:     sum    = sum + lx[i];
1432:   }
1433:   olx[pM] = sum;
1434:   sum     = 0;
1435:   for (i=0; i<pN; i++) {
1436:     oly[i] = sum;
1437:     sum    = sum + ly[i];
1438:   }
1439:   oly[pN] = sum;
1440:   sum     = 0;
1441:   for (i=0; i<pP; i++) {
1442:     olz[i] = sum;
1443:     sum    = sum + lz[i];
1444:   }
1445:   olz[pP] = sum;

1447:   PetscMalloc1(pM,&osx);
1448:   PetscMalloc1(pN,&osy);
1449:   PetscMalloc1(pP,&osz);
1450:   PetscMalloc1(pM,&oex);
1451:   PetscMalloc1(pN,&oey);
1452:   PetscMalloc1(pP,&oez);
1453:   for (i=0; i<pM; i++) {
1454:     osx[i] = olx[i] - stencil;
1455:     oex[i] = olx[i] + lx[i] + stencil;
1456:     if (osx[i]<0) osx[i]=0;
1457:     if (oex[i]>M) oex[i]=M;
1458:   }

1460:   for (i=0; i<pN; i++) {
1461:     osy[i] = oly[i] - stencil;
1462:     oey[i] = oly[i] + ly[i] + stencil;
1463:     if (osy[i]<0)osy[i]=0;
1464:     if (oey[i]>M)oey[i]=N;
1465:   }
1466:   for (i=0; i<pP; i++) {
1467:     osz[i] = olz[i] - stencil;
1468:     oez[i] = olz[i] + lz[i] + stencil;
1469:     if (osz[i]<0) osz[i]=0;
1470:     if (oez[i]>P) oez[i]=P;
1471:   }

1473:   for (k=0; k<pP; k++) {
1474:     for (j=0; j<pN; j++) {
1475:       for (i=0; i<pM; i++) {
1476:         char     name[PETSC_MAX_PATH_LEN];
1477:         PetscInt procid = i + j*pM + k*pM*pN; /* convert proc(i,j,k) to pid */
1478:         PetscSNPrintf(name,sizeof(name),"subdomain-%s-p%1.4d.vts",local_file_prefix,procid);
1479:         for (II=0; II<indent_level; II++) PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"  ");

1481:         PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"<Piece Extent=\"%d %d %d %d %d %d\"      Source=\"%s\"/>\n",
1482:                      osx[i],oex[i]-1,
1483:                      osy[j],oey[j]-1,
1484:                      osz[k],oez[k]-1,name);
1485:       }
1486:     }
1487:   }
1488:   PetscFree(olx);
1489:   PetscFree(oly);
1490:   PetscFree(olz);
1491:   PetscFree(osx);
1492:   PetscFree(osy);
1493:   PetscFree(osz);
1494:   PetscFree(oex);
1495:   PetscFree(oey);
1496:   PetscFree(oez);
1497:   return(0);
1498: }

1500: PetscErrorCode DAView_3DVTK_PStructuredGrid(DM da,const char file_prefix[],const char local_file_prefix[])
1501: {
1502:   MPI_Comm       comm;
1503:   PetscMPIInt    size,rank;
1504:   char           vtk_filename[PETSC_MAX_PATH_LEN];
1505:   FILE           *vtk_fp = NULL;
1506:   const char     *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian";
1507:   PetscInt       M,N,P,si,sj,sk,nx,ny,nz;
1508:   PetscInt       i,dofs;

1512:   /* only rank-0 generates this file */
1513:   PetscObjectGetComm((PetscObject)da,&comm);
1514:   MPI_Comm_size(comm,&size);
1515:   MPI_Comm_rank(comm,&rank);

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

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

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

1527:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"<VTKFile type=\"PStructuredGrid\" version=\"0.1\" byte_order=\"%s\">\n",byte_order);

1529:   /* define size of the nodal mesh based on the cell DM */
1530:   DMDAGetInfo(da,0,&M,&N,&P,0,0,0,&dofs,0,0,0,0,0);
1531:   DMDAGetGhostCorners(da,&si,&sj,&sk,&nx,&ny,&nz);
1532:   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 */

1534:   /* DUMP THE CELL REFERENCES */
1535:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"    <PCellData>\n");
1536:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"    </PCellData>\n");

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

1542:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"    <PPointData>\n");
1543:   for (i=0; i<dofs; i++) {
1544:     const char *fieldname;
1545:     DMDAGetFieldName(da,i,&fieldname);
1546:     PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"      <PDataArray type=\"Float64\" Name=\"%s\" NumberOfComponents=\"1\"/>\n",fieldname);
1547:   }
1548:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"    </PPointData>\n");

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

1553:   /* close the file */
1554:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"  </PStructuredGrid>\n");
1555:   PetscFPrintf(PETSC_COMM_SELF,vtk_fp,"</VTKFile>\n");

1557:   if (vtk_fp) {
1558:     fclose(vtk_fp);
1559:     vtk_fp = NULL;
1560:   }
1561:   return(0);
1562: }

1564: PetscErrorCode DAView3DPVTS(DM da, Vec x,const char NAME[])
1565: {
1566:   char           vts_filename[PETSC_MAX_PATH_LEN];
1567:   char           pvts_filename[PETSC_MAX_PATH_LEN];

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

1574:   PetscSNPrintf(pvts_filename,sizeof(pvts_filename),"%s-mesh",NAME);
1575:   DAView_3DVTK_PStructuredGrid(da,pvts_filename,vts_filename);
1576:   return(0);
1577: }

1579: PetscErrorCode KSPMonitorStokesBlocks(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy)
1580: {
1582:   PetscReal      norms[4];
1583:   Vec            Br,v,w;
1584:   Mat            A;

1587:   KSPGetOperators(ksp,&A,NULL);
1588:   MatCreateVecs(A,&w,&v);

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

1592:   VecStrideNorm(Br,0,NORM_2,&norms[0]);
1593:   VecStrideNorm(Br,1,NORM_2,&norms[1]);
1594:   VecStrideNorm(Br,2,NORM_2,&norms[2]);
1595:   VecStrideNorm(Br,3,NORM_2,&norms[3]);

1597:   VecDestroy(&v);
1598:   VecDestroy(&w);

1600:   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]);
1601:   return(0);
1602: }

1604: static PetscErrorCode PCMGSetupViaCoarsen(PC pc,DM da_fine)
1605: {
1606:   PetscInt              nlevels,k;
1607:   PETSC_UNUSED PetscInt finest;
1608:   DM                    *da_list,*daclist;
1609:   Mat                   R;
1610:   PetscErrorCode        ierr;

1613:   nlevels = 1;
1614:   PetscOptionsGetInt(NULL,NULL,"-levels",&nlevels,0);

1616:   PetscMalloc1(nlevels,&da_list);
1617:   for (k=0; k<nlevels; k++) da_list[k] = NULL;
1618:   PetscMalloc1(nlevels,&daclist);
1619:   for (k=0; k<nlevels; k++) daclist[k] = NULL;

1621:   /* finest grid is nlevels - 1 */
1622:   finest     = nlevels - 1;
1623:   daclist[0] = da_fine;
1624:   PetscObjectReference((PetscObject)da_fine);
1625:   DMCoarsenHierarchy(da_fine,nlevels-1,&daclist[1]);
1626:   for (k=0; k<nlevels; k++) {
1627:     da_list[k] = daclist[nlevels-1-k];
1628:     DMDASetUniformCoordinates(da_list[k],0.0,1.0,0.0,1.0,0.0,1.0);
1629:   }

1631:   PCMGSetLevels(pc,nlevels,NULL);
1632:   PCMGSetType(pc,PC_MG_MULTIPLICATIVE);
1633:   PCMGSetGalerkin(pc,PC_MG_GALERKIN_PMAT);

1635:   for (k=1; k<nlevels; k++) {
1636:     DMCreateInterpolation(da_list[k-1],da_list[k],&R,NULL);
1637:     PCMGSetInterpolation(pc,k,R);
1638:     MatDestroy(&R);
1639:   }

1641:   /* tidy up */
1642:   for (k=0; k<nlevels; k++) {
1643:     DMDestroy(&da_list[k]);
1644:   }
1645:   PetscFree(da_list);
1646:   PetscFree(daclist);
1647:   return(0);
1648: }

1650: static PetscErrorCode solve_stokes_3d_coupled(PetscInt mx,PetscInt my,PetscInt mz)
1651: {
1652:   DM             da_Stokes;
1653:   PetscInt       u_dof,p_dof,dof,stencil_width;
1654:   Mat            A,B;
1655:   DM             vel_cda;
1656:   Vec            vel_coords;
1657:   PetscInt       p;
1658:   Vec            f,X,X1;
1659:   DMDACoor3d     ***_vel_coords;
1660:   PetscInt       its;
1661:   KSP            ksp_S;
1662:   PetscInt       model_definition = 0;
1663:   PetscInt       ei,ej,ek,sex,sey,sez,Imx,Jmy,Kmz;
1664:   CellProperties cell_properties;
1665:   PetscBool      write_output = PETSC_FALSE,resolve= PETSC_FALSE;

1669:   PetscOptionsGetBool(NULL,NULL,"-resolve",&resolve,NULL);
1670:   /* Generate the da for velocity and pressure */
1671:   /* Num nodes in each direction is mx+1, my+1, mz+1 */
1672:   u_dof         = U_DOFS; /* Vx, Vy - velocities */
1673:   p_dof         = P_DOFS; /* p - pressure */
1674:   dof           = u_dof+p_dof;
1675:   stencil_width = 1;
1676:   DMDACreate3d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,
1677:                                mx+1,my+1,mz+1,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,dof,stencil_width,NULL,NULL,NULL,&da_Stokes);
1678:   DMSetMatType(da_Stokes,MATAIJ);
1679:   DMSetFromOptions(da_Stokes);
1680:   DMSetUp(da_Stokes);
1681:   DMDASetFieldName(da_Stokes,0,"Vx");
1682:   DMDASetFieldName(da_Stokes,1,"Vy");
1683:   DMDASetFieldName(da_Stokes,2,"Vz");
1684:   DMDASetFieldName(da_Stokes,3,"P");

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

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

1692:   /* interpolate the coordinates to quadrature points */
1693:   DMGetCoordinateDM(da_Stokes,&vel_cda);
1694:   DMGetCoordinatesLocal(da_Stokes,&vel_coords);
1695:   DMDAVecGetArray(vel_cda,vel_coords,&_vel_coords);
1696:   DMDAGetElementsCorners(da_Stokes,&sex,&sey,&sez);
1697:   DMDAGetElementsSizes(da_Stokes,&Imx,&Jmy,&Kmz);
1698:   for (ek = sez; ek < sez+Kmz; ek++) {
1699:     for (ej = sey; ej < sey+Jmy; ej++) {
1700:       for (ei = sex; ei < sex+Imx; ei++) {
1701:         /* get coords for the element */
1702:         PetscInt               ngp;
1703:         PetscScalar            gp_xi[GAUSS_POINTS][NSD],gp_weight[GAUSS_POINTS];
1704:         PetscScalar            el_coords[NSD*NODES_PER_EL];
1705:         GaussPointCoefficients *cell;

1707:         CellPropertiesGetCell(cell_properties,ei,ej,ek,&cell);
1708:         GetElementCoords3D(_vel_coords,ei,ej,ek,el_coords);
1709:         ConstructGaussQuadrature3D(&ngp,gp_xi,gp_weight);

1711:         for (p = 0; p < GAUSS_POINTS; p++) {
1712:           PetscScalar xi_p[NSD],Ni_p[NODES_PER_EL];
1713:           PetscScalar gp_x,gp_y,gp_z;
1714:           PetscInt    n;

1716:           xi_p[0] = gp_xi[p][0];
1717:           xi_p[1] = gp_xi[p][1];
1718:           xi_p[2] = gp_xi[p][2];
1719:           ShapeFunctionQ13D_Evaluate(xi_p,Ni_p);

1721:           gp_x = gp_y = gp_z = 0.0;
1722:           for (n = 0; n < NODES_PER_EL; n++) {
1723:             gp_x = gp_x+Ni_p[n]*el_coords[NSD*n];
1724:             gp_y = gp_y+Ni_p[n]*el_coords[NSD*n+1];
1725:             gp_z = gp_z+Ni_p[n]*el_coords[NSD*n+2];
1726:           }
1727:           cell->gp_coords[NSD*p]   = gp_x;
1728:           cell->gp_coords[NSD*p+1] = gp_y;
1729:           cell->gp_coords[NSD*p+2] = gp_z;
1730:         }
1731:       }
1732:     }
1733:   }

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

1737:   switch (model_definition) {
1738:   case 0: /* isoviscous */
1739:     for (ek = sez; ek < sez+Kmz; ek++) {
1740:       for (ej = sey; ej < sey+Jmy; ej++) {
1741:         for (ei = sex; ei < sex+Imx; ei++) {
1742:           GaussPointCoefficients *cell;

1744:           CellPropertiesGetCell(cell_properties,ei,ej,ek,&cell);
1745:           for (p = 0; p < GAUSS_POINTS; p++) {
1746:             PetscReal coord_x = PetscRealPart(cell->gp_coords[NSD*p]);
1747:             PetscReal coord_y = PetscRealPart(cell->gp_coords[NSD*p+1]);
1748:             PetscReal coord_z = PetscRealPart(cell->gp_coords[NSD*p+2]);

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

1752:             cell->fx[p] = 0.0*coord_x;
1753:             cell->fy[p] = -PetscSinReal(2.2*PETSC_PI*coord_y)*PetscCosReal(1.0*PETSC_PI*coord_x);
1754:             cell->fz[p] = 0.0*coord_z;
1755:             cell->hc[p] = 0.0;
1756:           }
1757:         }
1758:       }
1759:     }
1760:     break;

1762:   case 1: /* manufactured */
1763:     for (ek = sez; ek < sez+Kmz; ek++) {
1764:       for (ej = sey; ej < sey+Jmy; ej++) {
1765:         for (ei = sex; ei < sex+Imx; ei++) {
1766:           PetscReal              eta,Fm[NSD],Fc,pos[NSD];
1767:           GaussPointCoefficients *cell;

1769:           CellPropertiesGetCell(cell_properties,ei,ej,ek,&cell);
1770:           for (p = 0; p < GAUSS_POINTS; p++) {
1771:             PetscReal coord_x = PetscRealPart(cell->gp_coords[NSD*p]);
1772:             PetscReal coord_y = PetscRealPart(cell->gp_coords[NSD*p+1]);
1773:             PetscReal coord_z = PetscRealPart(cell->gp_coords[NSD*p+2]);

1775:             pos[0] = coord_x;
1776:             pos[1] = coord_y;
1777:             pos[2] = coord_z;

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

1782:             cell->fx[p] = Fm[0];
1783:             cell->fy[p] = Fm[1];
1784:             cell->fz[p] = Fm[2];
1785:             cell->hc[p] = Fc;
1786:           }
1787:         }
1788:       }
1789:     }
1790:     break;

1792:   case 2: /* solcx */
1793:     for (ek = sez; ek < sez+Kmz; ek++) {
1794:       for (ej = sey; ej < sey+Jmy; ej++) {
1795:         for (ei = sex; ei < sex+Imx; ei++) {
1796:           GaussPointCoefficients *cell;

1798:           CellPropertiesGetCell(cell_properties,ei,ej,ek,&cell);
1799:           for (p = 0; p < GAUSS_POINTS; p++) {
1800:             PetscReal coord_x = PetscRealPart(cell->gp_coords[NSD*p]);
1801:             PetscReal coord_y = PetscRealPart(cell->gp_coords[NSD*p+1]);
1802:             PetscReal coord_z = PetscRealPart(cell->gp_coords[NSD*p+2]);

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

1806:             cell->fx[p] = 0.0;
1807:             cell->fy[p] = -PetscSinReal(3.0*PETSC_PI*coord_y)*PetscCosReal(1.0*PETSC_PI*coord_x);
1808:             cell->fz[p] = 0.0*coord_z;
1809:             cell->hc[p] = 0.0;
1810:           }
1811:         }
1812:       }
1813:     }
1814:     break;

1816:   case 3: /* sinker */
1817:     for (ek = sez; ek < sez+Kmz; ek++) {
1818:       for (ej = sey; ej < sey+Jmy; ej++) {
1819:         for (ei = sex; ei < sex+Imx; ei++) {
1820:           GaussPointCoefficients *cell;

1822:           CellPropertiesGetCell(cell_properties,ei,ej,ek,&cell);
1823:           for (p = 0; p < GAUSS_POINTS; p++) {
1824:             PetscReal xp = PetscRealPart(cell->gp_coords[NSD*p]);
1825:             PetscReal yp = PetscRealPart(cell->gp_coords[NSD*p+1]);
1826:             PetscReal zp = PetscRealPart(cell->gp_coords[NSD*p+2]);

1828:             cell->eta[p] = 1.0e-2;
1829:             cell->fx[p]  = 0.0;
1830:             cell->fy[p]  = 0.0;
1831:             cell->fz[p]  = 0.0;
1832:             cell->hc[p]  = 0.0;

1834:             if ((PetscAbs(xp-0.5) < 0.2) && (PetscAbs(yp-0.5) < 0.2) && (PetscAbs(zp-0.5) < 0.2)) {
1835:               cell->eta[p] = 1.0;
1836:               cell->fz[p]  = 1.0;
1837:             }

1839:           }
1840:         }
1841:       }
1842:     }
1843:     break;

1845:   case 4: /* subdomain jumps */
1846:     for (ek = sez; ek < sez+Kmz; ek++) {
1847:       for (ej = sey; ej < sey+Jmy; ej++) {
1848:         for (ei = sex; ei < sex+Imx; ei++) {
1849:           PetscReal              opts_mag,opts_eta0;
1850:           PetscInt               px,py,pz;
1851:           PetscBool              jump;
1852:           PetscMPIInt            rr;
1853:           GaussPointCoefficients *cell;

1855:           opts_mag  = 1.0;
1856:           opts_eta0 = 1.e-2;

1858:           PetscOptionsGetReal(NULL,NULL,"-jump_eta0",&opts_eta0,NULL);
1859:           PetscOptionsGetReal(NULL,NULL,"-jump_magnitude",&opts_mag,NULL);
1860:           DMDAGetInfo(da_Stokes,NULL,NULL,NULL,NULL,&px,&py,&pz,NULL,NULL,NULL,NULL,NULL,NULL);
1861:           rr   = PetscGlobalRank%(px*py);
1862:           if (px%2) jump = (PetscBool)(rr%2);
1863:           else jump = (PetscBool)((rr/px)%2 ? rr%2 : !(rr%2));
1864:           rr = PetscGlobalRank/(px*py);
1865:           if (rr%2) jump = (PetscBool)!jump;
1866:           CellPropertiesGetCell(cell_properties,ei,ej,ek,&cell);
1867:           for (p = 0; p < GAUSS_POINTS; p++) {
1868:             PetscReal xp = PetscRealPart(cell->gp_coords[NSD*p]);
1869:             PetscReal yp = PetscRealPart(cell->gp_coords[NSD*p+1]);

1871:             cell->eta[p] = jump ? PetscPowReal(10.0,opts_mag) : opts_eta0;
1872:             cell->fx[p]  = 0.0;
1873:             cell->fy[p]  = -PetscSinReal(2.2*PETSC_PI*yp)*PetscCosReal(1.0*PETSC_PI*xp);
1874:             cell->fz[p]  = 0.0;
1875:             cell->hc[p]  = 0.0;
1876:           }
1877:         }
1878:       }
1879:     }
1880:     break;
1881:   default:
1882:     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"No default model is supported. Choose either -model {0,1,2,3}");
1883:   }

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

1887:   /* Generate a matrix with the correct non-zero pattern of type AIJ. This will work in parallel and serial */
1888:   DMCreateMatrix(da_Stokes,&A);
1889:   DMCreateMatrix(da_Stokes,&B);
1890:   DMCreateGlobalVector(da_Stokes,&X);
1891:   DMCreateGlobalVector(da_Stokes,&f);

1893:   /* assemble A11 */
1894:   MatZeroEntries(A);
1895:   MatZeroEntries(B);
1896:   VecZeroEntries(f);

1898:   AssembleA_Stokes(A,da_Stokes,cell_properties);
1899:   MatViewFromOptions(A,NULL,"-amat_view");
1900:   AssembleA_PCStokes(B,da_Stokes,cell_properties);
1901:   MatViewFromOptions(B,NULL,"-bmat_view");
1902:   /* build force vector */
1903:   AssembleF_Stokes(f,da_Stokes,cell_properties);

1905:   /* SOLVE */
1906:   KSPCreate(PETSC_COMM_WORLD,&ksp_S);
1907:   KSPSetOptionsPrefix(ksp_S,"stokes_"); /* stokes */ 
1908:   KSPSetOperators(ksp_S,A,B);
1909:   KSPSetFromOptions(ksp_S);

1911:   {
1912:     PC             pc;
1913:     const PetscInt ufields[] = {0,1,2},pfields[] = {3};
1914:     KSPGetPC(ksp_S,&pc);
1915:     PCFieldSplitSetBlockSize(pc,4);
1916:     PCFieldSplitSetFields(pc,"u",3,ufields,ufields);
1917:     PCFieldSplitSetFields(pc,"p",1,pfields,pfields);
1918:   }

1920:   {
1921:     PC        pc;
1922:     PetscBool same = PETSC_FALSE;
1923:     KSPGetPC(ksp_S,&pc);
1924:     PetscObjectTypeCompare((PetscObject)pc,PCMG,&same);
1925:     if (same) {
1926:       PCMGSetupViaCoarsen(pc,da_Stokes);
1927:     }
1928:   }

1930:   {
1931:     PC        pc;
1932:     PetscBool same = PETSC_FALSE;
1933:     KSPGetPC(ksp_S,&pc);
1934:     PetscObjectTypeCompare((PetscObject)pc,PCBDDC,&same);
1935:     if (same) {
1936:       KSPSetOperators(ksp_S,A,A);
1937:     }
1938:   }

1940:   {
1941:     PetscBool stokes_monitor = PETSC_FALSE;
1942:     PetscOptionsGetBool(NULL,NULL,"-stokes_ksp_monitor_blocks",&stokes_monitor,0);
1943:     if (stokes_monitor) {
1944:       KSPMonitorSet(ksp_S,KSPMonitorStokesBlocks,NULL,NULL);
1945:     }
1946:   }

1948:   if (resolve) {
1949:     /* Test changing matrix data structure and resolve */
1950:     VecDuplicate(X,&X1);
1951:   }

1953:   KSPSolve(ksp_S,f,X);
1954:   if (resolve) {
1955:     Mat C;
1956:     MatDuplicate(A,MAT_COPY_VALUES,&C);
1957:     KSPSetOperators(ksp_S,C,C);
1958:     KSPSolve(ksp_S,f,X1);
1959:     MatDestroy(&C);
1960:     VecDestroy(&X1);
1961:   }

1963:   PetscOptionsGetBool(NULL,NULL,"-write_pvts",&write_output,NULL);
1964:   if (write_output) {DAView3DPVTS(da_Stokes,X,"up");}
1965:   {
1966:     PetscBool flg = PETSC_FALSE;
1967:     char      filename[PETSC_MAX_PATH_LEN];
1968:     PetscOptionsGetString(NULL,NULL,"-write_binary",filename,sizeof(filename),&flg);
1969:     if (flg) {
1970:       PetscViewer viewer;
1971:       /* PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename[0]?filename:"ex42-binaryoutput",FILE_MODE_WRITE,&viewer); */
1972:       PetscViewerVTKOpen(PETSC_COMM_WORLD,"ex42.vts",FILE_MODE_WRITE,&viewer);
1973:       VecView(X,viewer);
1974:       PetscViewerDestroy(&viewer);
1975:     }
1976:   }
1977:   KSPGetIterationNumber(ksp_S,&its);

1979:   /* verify */
1980:   if (model_definition == 1) {
1981:     DM  da_Stokes_analytic;
1982:     Vec X_analytic;

1984:     DMDACreateManufacturedSolution(mx,my,mz,&da_Stokes_analytic,&X_analytic);
1985:     if (write_output) {
1986:       DAView3DPVTS(da_Stokes_analytic,X_analytic,"ms");
1987:     }
1988:     DMDAIntegrateErrors3D(da_Stokes_analytic,X,X_analytic);
1989:     if (write_output) {
1990:       DAView3DPVTS(da_Stokes,X,"up2");
1991:     }
1992:     DMDestroy(&da_Stokes_analytic);
1993:     VecDestroy(&X_analytic);
1994:   }

1996:   KSPDestroy(&ksp_S);
1997:   VecDestroy(&X);
1998:   VecDestroy(&f);
1999:   MatDestroy(&A);
2000:   MatDestroy(&B);

2002:   CellPropertiesDestroy(&cell_properties);
2003:   DMDestroy(&da_Stokes);
2004:   return(0);
2005: }

2007: int main(int argc,char **args)
2008: {
2010:   PetscInt       mx,my,mz;

2012:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
2013:   mx   = my = mz = 10;
2014:   PetscOptionsGetInt(NULL,NULL,"-mx",&mx,NULL);
2015:   my   = mx; mz = mx;
2016:   PetscOptionsGetInt(NULL,NULL,"-my",&my,NULL);
2017:   PetscOptionsGetInt(NULL,NULL,"-mz",&mz,NULL);
2018:   solve_stokes_3d_coupled(mx,my,mz);
2019:   PetscFinalize();
2020:   return ierr;
2021: }

2023: /*TEST

2025:    test:
2026:       args: -stokes_ksp_monitor_short -stokes_ksp_converged_reason -stokes_pc_type lu

2028:    test:
2029:       suffix: 2
2030:       nsize: 3
2031:       args: -stokes_ksp_monitor_short -stokes_ksp_converged_reason -stokes_pc_type redundant -stokes_redundant_pc_type lu

2033:    test:
2034:       suffix: bddc_stokes
2035:       nsize: 6
2036:       args: -mx 5 -my 4 -mz 3 -stokes_ksp_monitor_short -stokes_ksp_converged_reason -stokes_pc_type bddc -dm_mat_type is -stokes_pc_bddc_dirichlet_pc_type svd -stokes_pc_bddc_neumann_pc_type svd -stokes_pc_bddc_coarse_redundant_pc_type svd

2038:    test:
2039:       suffix: bddc_stokes_deluxe
2040:       nsize: 6
2041:       args: -mx 5 -my 4 -mz 3 -stokes_ksp_monitor_short -stokes_ksp_converged_reason -stokes_pc_type bddc -dm_mat_type is -stokes_pc_bddc_dirichlet_pc_type svd -stokes_pc_bddc_neumann_pc_type svd -stokes_pc_bddc_coarse_redundant_pc_type svd -stokes_pc_bddc_use_deluxe_scaling -stokes_sub_schurs_posdef 0 -stokes_sub_schurs_symmetric -stokes_sub_schurs_mat_solver_type petsc

2043:    test:
2044:       requires: !single
2045:       suffix: bddc_stokes_subdomainjump_deluxe
2046:       nsize: 9
2047:       args: -model 4 -jump_magnitude 4 -mx 6 -my 6 -mz 2 -stokes_ksp_monitor_short -stokes_ksp_converged_reason -stokes_pc_type bddc -dm_mat_type is -stokes_pc_bddc_use_deluxe_scaling -stokes_sub_schurs_posdef 0 -stokes_sub_schurs_symmetric -stokes_sub_schurs_mat_solver_type petsc -stokes_pc_bddc_schur_layers 1

2049:    test:
2050:       requires: !single
2051:       suffix: 3
2052:       args: -stokes_ksp_converged_reason -stokes_pc_type fieldsplit -resolve

2054:    test:
2055:       suffix: tut_1
2056:       nsize: 4
2057:       requires: !single
2058:       args: -stokes_ksp_monitor

2060:    test:
2061:       suffix: tut_2
2062:       nsize: 4
2063:       requires: !single
2064:       args: -stokes_ksp_monitor -stokes_pc_type fieldsplit -stokes_pc_fieldsplit_type schur

2066:    test:
2067:       suffix: tut_3
2068:       nsize: 4
2069:       requires: !single
2070:       args: -mx 20 -stokes_ksp_monitor -stokes_pc_type fieldsplit -stokes_pc_fieldsplit_type schur

2072: TEST*/