Actual source code: ex77.c

petsc-3.8.4 2018-03-24
Report Typos and Errors
  1: static char help[] = "Nonlinear elasticity problem in 3d with simplicial finite elements.\n\
  2: We solve a nonlinear elasticity problem, modelled as an incompressible Neo-Hookean solid, \n\
  3:  with pressure loading in a rectangular domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n";

  5: /*
  6: Nonlinear elasticity problem, which we discretize using the finite
  7: element method on an unstructured mesh. This uses both Dirichlet boundary conditions (fixed faces)
  8: and nonlinear Neumann boundary conditions (pressure loading).
  9: The Lagrangian density (modulo boundary conditions) for this problem is given by
 10: \begin{equation}
 11:   \frac{\mu}{2} (\mathrm{Tr}{C}-3) + J p + \frac{\kappa}{2} (J-1).
 12: \end{equation}

 14: Discretization:

 16: We use PetscFE to generate a tabulation of the finite element basis functions
 17: at quadrature points. We can currently generate an arbitrary order Lagrange
 18: element.

 20: Field Data:

 22:   DMPLEX data is organized by point, and the closure operation just stacks up the
 23: data from each sieve point in the closure. Thus, for a P_2-P_1 Stokes element, we
 24: have

 26:   cl{e} = {f e_0 e_1 e_2 v_0 v_1 v_2}
 27:   x     = [u_{e_0} v_{e_0} u_{e_1} v_{e_1} u_{e_2} v_{e_2} u_{v_0} v_{v_0} p_{v_0} u_{v_1} v_{v_1} p_{v_1} u_{v_2} v_{v_2} p_{v_2}]

 29: The problem here is that we would like to loop over each field separately for
 30: integration. Therefore, the closure visitor in DMPlexVecGetClosure() reorders
 31: the data so that each field is contiguous

 33:   x'    = [u_{e_0} v_{e_0} u_{e_1} v_{e_1} u_{e_2} v_{e_2} u_{v_0} v_{v_0} u_{v_1} v_{v_1} u_{v_2} v_{v_2} p_{v_0} p_{v_1} p_{v_2}]

 35: Likewise, DMPlexVecSetClosure() takes data partitioned by field, and correctly
 36: puts it into the Sieve ordering.

 38: */

 40:  #include <petscdmplex.h>
 41:  #include <petscsnes.h>
 42:  #include <petscds.h>

 44: typedef enum {RUN_FULL, RUN_TEST} RunType;

 46: typedef struct {
 47:   PetscInt      debug;             /* The debugging level */
 48:   RunType       runType;           /* Whether to run tests, or solve the full problem */
 49:   PetscLogEvent createMeshEvent;
 50:   PetscBool     showInitial, showSolution;
 51:   /* Domain and mesh definition */
 52:   PetscInt      dim;               /* The topological mesh dimension */
 53:   PetscBool     interpolate;       /* Generate intermediate mesh elements */
 54:   PetscBool     simplex;           /* Use simplices or tensor product cells */
 55:   PetscReal     refinementLimit;   /* The largest allowable cell volume */
 56:   PetscBool     testPartition;     /* Use a fixed partitioning for testing */
 57:   PetscReal     mu;                /* The shear modulus */
 58:   PetscReal     p_wall;            /* The wall pressure */
 59: } AppCtx;

 61: #if 0
 62: PETSC_STATIC_INLINE void Det2D(PetscReal *detJ, const PetscReal J[])
 63: {
 64:   *detJ = J[0]*J[3] - J[1]*J[2];
 65: }
 66: #endif

 68: PETSC_STATIC_INLINE void Det3D(PetscReal *detJ, const PetscScalar J[])
 69: {
 70:   *detJ = PetscRealPart(J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) +
 71:                         J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) +
 72:                         J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0]));
 73: }

 75: #if 0
 76: PETSC_STATIC_INLINE void Cof2D(PetscReal C[], const PetscReal A[])
 77: {
 78:   C[0] =  A[3];
 79:   C[1] = -A[2];
 80:   C[2] = -A[1];
 81:   C[3] =  A[0];
 82: }
 83: #endif

 85: PETSC_STATIC_INLINE void Cof3D(PetscReal C[], const PetscScalar A[])
 86: {
 87:   C[0*3+0] = PetscRealPart(A[1*3+1]*A[2*3+2] - A[1*3+2]*A[2*3+1]);
 88:   C[0*3+1] = PetscRealPart(A[1*3+2]*A[2*3+0] - A[1*3+0]*A[2*3+2]);
 89:   C[0*3+2] = PetscRealPart(A[1*3+0]*A[2*3+1] - A[1*3+1]*A[2*3+0]);
 90:   C[1*3+0] = PetscRealPart(A[0*3+2]*A[2*3+1] - A[0*3+1]*A[2*3+2]);
 91:   C[1*3+1] = PetscRealPart(A[0*3+0]*A[2*3+2] - A[0*3+2]*A[2*3+0]);
 92:   C[1*3+2] = PetscRealPart(A[0*3+1]*A[2*3+0] - A[0*3+0]*A[2*3+1]);
 93:   C[2*3+0] = PetscRealPart(A[0*3+1]*A[1*3+2] - A[0*3+2]*A[1*3+1]);
 94:   C[2*3+1] = PetscRealPart(A[0*3+2]*A[1*3+0] - A[0*3+0]*A[1*3+2]);
 95:   C[2*3+2] = PetscRealPart(A[0*3+0]*A[1*3+1] - A[0*3+1]*A[1*3+0]);
 96: }

 98: PetscErrorCode zero_scalar(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
 99: {
100:   u[0] = 0.0;
101:   return 0;
102: }

104: PetscErrorCode zero_vector(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
105: {
106:   const PetscInt Ncomp = dim;

108:   PetscInt       comp;
109:   for (comp = 0; comp < Ncomp; ++comp) u[comp] = 0.0;
110:   return 0;
111: }

113: PetscErrorCode coordinates(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
114: {
115:   const PetscInt Ncomp = dim;

117:   PetscInt       comp;
118:   for (comp = 0; comp < Ncomp; ++comp) u[comp] = x[comp];
119:   return 0;
120: }

122: PetscErrorCode elasticityMaterial(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
123: {
124:   AppCtx *user = (AppCtx *) ctx;
125:   u[0] = user->mu;
126:   return 0;
127: }

129: PetscErrorCode wallPressure(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
130: {
131:   AppCtx *user = (AppCtx *) ctx;
132:   u[0] = user->p_wall;
133:   return 0;
134: }

136: void f1_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
137:           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
138:           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
139:           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
140: {
141:   const PetscInt  Ncomp = dim;
142:   const PetscReal mu = PetscRealPart(a[0]), kappa = 3.0;
143:   PetscReal       cofu_x[9/*Ncomp*dim*/], detu_x, p = PetscRealPart(u[Ncomp]);
144:   PetscInt        comp, d;

146:   Cof3D(cofu_x, u_x);
147:   Det3D(&detu_x, u_x);
148:   p += kappa * (detu_x - 1.0);

150:   /* f1 is the first Piola-Kirchhoff tensor */
151:   for (comp = 0; comp < Ncomp; ++comp) {
152:     for (d = 0; d < dim; ++d) {
153:       f1[comp*dim+d] = mu * u_x[comp*dim+d] + p * cofu_x[comp*dim+d];
154:     }
155:   }
156: }

158: void g3_uu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
159:           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
160:           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
161:           PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
162: {
163:   const PetscInt  Ncomp = dim;
164:   const PetscReal mu = PetscRealPart(a[0]), kappa = 3.0;
165:   PetscReal       cofu_x[9/*Ncomp*dim*/], detu_x, pp, pm, p = PetscRealPart(u[Ncomp]);
166:   PetscInt        compI, compJ, d1, d2;

168:   Cof3D(cofu_x, u_x);
169:   Det3D(&detu_x, u_x);
170:   p += kappa * (detu_x - 1.0);
171:   pp = p/detu_x + kappa;
172:   pm = p/detu_x;

174:   /* g3 is the first elasticity tensor, i.e. A_i^I_j^J = S^{IJ} g_{ij} + C^{KILJ} F^k_K F^l_L g_{kj} g_{li} */
175:   for (compI = 0; compI < Ncomp; ++compI) {
176:     for (compJ = 0; compJ < Ncomp; ++compJ) {
177:       const PetscReal G = (compI == compJ) ? 1.0 : 0.0;

179:       for (d1 = 0; d1 < dim; ++d1) {
180:         for (d2 = 0; d2 < dim; ++d2) {
181:           const PetscReal g = (d1 == d2) ? 1.0 : 0.0;

183:           g3[((compI*Ncomp+compJ)*dim+d1)*dim+d2] = g * G * mu + pp * cofu_x[compI*dim+d1] * cofu_x[compJ*dim+d2] - pm * cofu_x[compI*dim+d2] * cofu_x[compJ*dim+d1];
184:         }
185:       }
186:     }
187:   }
188: }

190: void f0_bd_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
191:     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
192:     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
193:     PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
194: {
195:   const PetscInt    Ncomp = dim;
196:   const PetscScalar p = a[aOff[1]];
197:   PetscReal         cofu_x[9/*Ncomp*dim*/];
198:   PetscInt          comp, d;

200:   Cof3D(cofu_x, u_x);
201:   for (comp = 0; comp < Ncomp; ++comp) {
202:     for (d = 0, f0[comp] = 0.0; d < dim; ++d) f0[comp] += cofu_x[comp*dim+d] * n[d];
203:     f0[comp] *= p;
204:   }
205: }

207: void g1_bd_uu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
208:     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
209:     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
210:     PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
211: {
212:   const PetscInt Ncomp = dim;
213:   PetscScalar    p = a[aOff[1]];
214:   PetscReal      cofu_x[9/*Ncomp*dim*/], m[3/*Ncomp*/], detu_x;
215:   PetscInt       comp, compI, compJ, d;

217:   Cof3D(cofu_x, u_x);
218:   Det3D(&detu_x, u_x);
219:   p /= detu_x;

221:   for (comp = 0; comp < Ncomp; ++comp) for (d = 0, m[comp] = 0.0; d < dim; ++d) m[comp] += cofu_x[comp*dim+d] * n[d];
222:   for (compI = 0; compI < Ncomp; ++compI) {
223:     for (compJ = 0; compJ < Ncomp; ++compJ) {
224:       for (d = 0; d < dim; ++d) {
225:         g1[(compI*Ncomp+compJ)*dim+d] = p * (m[compI] * cofu_x[compJ*dim+d] - cofu_x[compI*dim+d] * m[compJ]);
226:       }
227:     }
228:   }
229: }

231: void f0_p_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
232:           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
233:           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
234:           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
235: {
236:   PetscReal detu_x;
237:   Det3D(&detu_x, u_x);
238:   f0[0] = detu_x - 1.0;
239: }

241: void g1_pu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
242:            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
243:            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
244:            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
245: {
246:   PetscReal cofu_x[9/*Ncomp*dim*/];
247:   PetscInt  compI, d;

249:   Cof3D(cofu_x, u_x);
250:   for (compI = 0; compI < dim; ++compI)
251:     for (d = 0; d < dim; ++d) g1[compI*dim+d] = cofu_x[compI*dim+d];
252: }

254: void g2_up_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux,
255:            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
256:            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
257:            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
258: {
259:   PetscReal cofu_x[9/*Ncomp*dim*/];
260:   PetscInt  compI, d;

262:   Cof3D(cofu_x, u_x);
263:   for (compI = 0; compI < dim; ++compI)
264:     for (d = 0; d < dim; ++d) g2[compI*dim+d] = cofu_x[compI*dim+d];
265: }

267: PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
268: {
269:   const char    *runTypes[2] = {"full", "test"};
270:   PetscInt       run;

274:   options->debug           = 0;
275:   options->runType         = RUN_FULL;
276:   options->dim             = 3;
277:   options->interpolate     = PETSC_FALSE;
278:   options->simplex         = PETSC_TRUE;
279:   options->refinementLimit = 0.0;
280:   options->mu              = 1.0;
281:   options->p_wall          = 0.4;
282:   options->testPartition   = PETSC_FALSE;
283:   options->showInitial     = PETSC_FALSE;
284:   options->showSolution    = PETSC_TRUE;

286:   PetscOptionsBegin(comm, "", "Nonlinear elasticity problem options", "DMPLEX");
287:   PetscOptionsInt("-debug", "The debugging level", "ex77.c", options->debug, &options->debug, NULL);
288:   run  = options->runType;
289:   PetscOptionsEList("-run_type", "The run type", "ex77.c", runTypes, 2, runTypes[options->runType], &run, NULL);

291:   options->runType = (RunType) run;

293:   PetscOptionsInt("-dim", "The topological mesh dimension", "ex77.c", options->dim, &options->dim, NULL);
294:   PetscOptionsBool("-interpolate", "Generate intermediate mesh elements", "ex77.c", options->interpolate, &options->interpolate, NULL);
295:   PetscOptionsBool("-simplex", "Use simplices or tensor product cells", "ex77.c", options->simplex, &options->simplex, NULL);
296:   PetscOptionsReal("-refinement_limit", "The largest allowable cell volume", "ex77.c", options->refinementLimit, &options->refinementLimit, NULL);
297:   PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex77.c", options->testPartition, &options->testPartition, NULL);
298:   PetscOptionsReal("-shear_modulus", "The shear modulus", "ex77.c", options->mu, &options->mu, NULL);
299:   PetscOptionsReal("-wall_pressure", "The wall pressure", "ex77.c", options->p_wall, &options->p_wall, NULL);

301:   PetscOptionsBool("-show_initial", "Output the initial guess for verification", "ex77.c", options->showInitial, &options->showInitial, NULL);
302:   PetscOptionsBool("-show_solution", "Output the solution for verification", "ex77.c", options->showSolution, &options->showSolution, NULL);
303:   PetscOptionsEnd();

305:   PetscLogEventRegister("CreateMesh", DM_CLASSID, &options->createMeshEvent);
306:   return(0);
307: }

309: PetscErrorCode DMVecViewLocal(DM dm, Vec v, PetscViewer viewer)
310: {
311:   Vec            lv;
312:   PetscInt       p;
313:   PetscMPIInt    rank, size;

317:   MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);
318:   MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size);
319:   DMGetLocalVector(dm, &lv);
320:   DMGlobalToLocalBegin(dm, v, INSERT_VALUES, lv);
321:   DMGlobalToLocalEnd(dm, v, INSERT_VALUES, lv);
322:   PetscPrintf(PETSC_COMM_WORLD, "Local function\n");
323:   for (p = 0; p < size; ++p) {
324:     if (p == rank) {VecView(lv, PETSC_VIEWER_STDOUT_SELF);}
325:     PetscBarrier((PetscObject) dm);
326:   }
327:   DMRestoreLocalVector(dm, &lv);
328:   return(0);
329: }

331: PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
332: {
333:   PetscInt       dim             = user->dim;
334:   PetscBool      interpolate     = user->interpolate;
335:   PetscReal      refinementLimit = user->refinementLimit;
336:   const PetscInt cells[3]        = {3, 3, 3};

340:   PetscLogEventBegin(user->createMeshEvent,0,0,0,0);
341:   if (user->simplex) {DMPlexCreateBoxMesh(comm, dim, dim == 2 ? 2 : 1, interpolate, dm);}
342:   else               {DMPlexCreateHexBoxMesh(comm, dim, cells, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, dm);}
343:   /* Label the faces (bit of a hack here, until it is properly implemented for simplices) */
344:   {
345:     DM              cdm;
346:     DMLabel         label;
347:     IS              is;
348:     PetscInt        d, dim = user->dim, b, f, Nf;
349:     const PetscInt *faces;
350:     PetscInt        csize;
351:     PetscScalar    *coords = NULL;
352:     PetscSection    cs;
353:     Vec             coordinates ;

355:     DMCreateLabel(*dm, "boundary");
356:     DMGetLabel(*dm, "boundary", &label);
357:     DMPlexMarkBoundaryFaces(*dm, label);
358:     DMGetStratumIS(*dm, "boundary", 1,  &is);
359:     if (is) {
360:       PetscReal faceCoord;
361:       PetscInt  v;

363:       ISGetLocalSize(is, &Nf);
364:       ISGetIndices(is, &faces);

366:       DMGetCoordinatesLocal(*dm, &coordinates);
367:       DMGetCoordinateDM(*dm, &cdm);
368:       DMGetDefaultSection(cdm, &cs);

370:       /* Check for each boundary face if any component of its centroid is either 0.0 or 1.0 */
371:       for (f = 0; f < Nf; ++f) {
372:         DMPlexVecGetClosure(cdm, cs, coordinates, faces[f], &csize, &coords);
373:         /* Calculate mean coordinate vector */
374:         for (d = 0; d < dim; ++d) {
375:           const PetscInt Nv = csize/dim;
376:           faceCoord = 0.0;
377:           for (v = 0; v < Nv; ++v) faceCoord += PetscRealPart(coords[v*dim+d]);
378:           faceCoord /= Nv;
379:           for (b = 0; b < 2; ++b) {
380:             if (PetscAbs(faceCoord - b*1.0) < PETSC_SMALL) {
381:               DMSetLabelValue(*dm, "Faces", faces[f], d*2+b+1);
382:             }
383:           }
384:         }
385:         DMPlexVecRestoreClosure(cdm, cs, coordinates, faces[f], &csize, &coords);
386:       }
387:       ISRestoreIndices(is, &faces);
388:     }
389:     ISDestroy(&is);
390:   }
391:   {
392:     DM refinedMesh     = NULL;
393:     DM distributedMesh = NULL;

395:     /* Refine mesh using a volume constraint */
396:     DMPlexSetRefinementLimit(*dm, refinementLimit);
397:     if (user->simplex) {DMRefine(*dm, comm, &refinedMesh);}
398:     if (refinedMesh) {
399:       DMDestroy(dm);
400:       *dm  = refinedMesh;
401:     }
402:     /* Setup test partitioning */
403:     if (user->testPartition) {
404:       PetscInt         triSizes_n2[2]       = {4, 4};
405:       PetscInt         triPoints_n2[8]      = {3, 5, 6, 7, 0, 1, 2, 4};
406:       PetscInt         triSizes_n3[3]       = {2, 3, 3};
407:       PetscInt         triPoints_n3[8]      = {3, 5, 1, 6, 7, 0, 2, 4};
408:       PetscInt         triSizes_n5[5]       = {1, 2, 2, 1, 2};
409:       PetscInt         triPoints_n5[8]      = {3, 5, 6, 4, 7, 0, 1, 2};
410:       PetscInt         triSizes_ref_n2[2]   = {8, 8};
411:       PetscInt         triPoints_ref_n2[16] = {1, 5, 6, 7, 10, 11, 14, 15, 0, 2, 3, 4, 8, 9, 12, 13};
412:       PetscInt         triSizes_ref_n3[3]   = {5, 6, 5};
413:       PetscInt         triPoints_ref_n3[16] = {1, 7, 10, 14, 15, 2, 6, 8, 11, 12, 13, 0, 3, 4, 5, 9};
414:       PetscInt         triSizes_ref_n5[5]   = {3, 4, 3, 3, 3};
415:       PetscInt         triPoints_ref_n5[16] = {1, 7, 10, 2, 11, 13, 14, 5, 6, 15, 0, 8, 9, 3, 4, 12};
416:       PetscInt         tetSizes_n2[2]       = {3, 3};
417:       PetscInt         tetPoints_n2[6]      = {1, 2, 3, 0, 4, 5};
418:       const PetscInt  *sizes = NULL;
419:       const PetscInt  *points = NULL;
420:       PetscPartitioner part;
421:       PetscInt         cEnd;
422:       PetscMPIInt      rank, size;

424:       MPI_Comm_rank(comm, &rank);
425:       MPI_Comm_size(comm, &size);
426:       DMPlexGetHeightStratum(*dm, 0, NULL, &cEnd);
427:       if (!rank) {
428:         if (dim == 2 && user->simplex && size == 2 && cEnd == 8) {
429:            sizes = triSizes_n2; points = triPoints_n2;
430:         } else if (dim == 2 && user->simplex && size == 3 && cEnd == 8) {
431:           sizes = triSizes_n3; points = triPoints_n3;
432:         } else if (dim == 2 && user->simplex && size == 5 && cEnd == 8) {
433:           sizes = triSizes_n5; points = triPoints_n5;
434:         } else if (dim == 2 && user->simplex && size == 2 && cEnd == 16) {
435:            sizes = triSizes_ref_n2; points = triPoints_ref_n2;
436:         } else if (dim == 2 && user->simplex && size == 3 && cEnd == 16) {
437:           sizes = triSizes_ref_n3; points = triPoints_ref_n3;
438:         } else if (dim == 2 && user->simplex && size == 5 && cEnd == 16) {
439:           sizes = triSizes_ref_n5; points = triPoints_ref_n5;
440:         } else if (dim == 3 && user->simplex && size == 2 && cEnd == 6) {
441:           sizes = tetSizes_n2; points = tetPoints_n2;
442:         } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "No stored partition matching run parameters");
443:       }
444:       DMPlexGetPartitioner(*dm, &part);
445:       PetscPartitionerSetType(part, PETSCPARTITIONERSHELL);
446:       PetscPartitionerShellSetPartition(part, size, sizes, points);
447:     }
448:     /* Distribute mesh over processes */
449:     DMPlexDistribute(*dm, 0, NULL, &distributedMesh);
450:     if (distributedMesh) {
451:       DMDestroy(dm);
452:       *dm  = distributedMesh;
453:     }
454:   }
455:   DMSetFromOptions(*dm);
456:   DMViewFromOptions(*dm, NULL, "-dm_view");
457: 
458:   PetscLogEventEnd(user->createMeshEvent,0,0,0,0);
459:   return(0);
460: }

462: PetscErrorCode SetupProblem(PetscDS prob, PetscInt dim, AppCtx *user)
463: {

467:   PetscDSSetResidual(prob, 0, NULL, f1_u_3d);
468:   PetscDSSetResidual(prob, 1, f0_p_3d, NULL);
469:   PetscDSSetJacobian(prob, 0, 0, NULL, NULL,  NULL,  g3_uu_3d);
470:   PetscDSSetJacobian(prob, 0, 1, NULL, NULL,  g2_up_3d, NULL);
471:   PetscDSSetJacobian(prob, 1, 0, NULL, g1_pu_3d, NULL,  NULL);

473:   PetscDSSetBdResidual(prob, 0, f0_bd_u_3d, NULL);
474:   PetscDSSetBdJacobian(prob, 0, 0, NULL, g1_bd_uu_3d, NULL, NULL);

476:   PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "fixed", "Faces", 0, 0, NULL, (void (*)(void)) coordinates, 0, NULL, user);
477:   PetscDSAddBoundary(prob, DM_BC_NATURAL, "pressure", "Faces", 0, 0, NULL, NULL, 0, NULL, user);
478:   return(0);
479: }

481: PetscErrorCode SetupMaterial(DM dm, DM dmAux, AppCtx *user)
482: {
483:   PetscErrorCode (*matFuncs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx) = {elasticityMaterial, wallPressure};
484:   Vec            A;
485:   void          *ctxs[2];

489:   ctxs[0] = user; ctxs[1] = user;
490:   DMCreateLocalVector(dmAux, &A);
491:   DMProjectFunctionLocal(dmAux, 0.0, matFuncs, ctxs, INSERT_ALL_VALUES, A);
492:   PetscObjectCompose((PetscObject) dm, "A", (PetscObject) A);
493:   VecDestroy(&A);
494:   return(0);
495: }

497: PetscErrorCode SetupNearNullSpace(DM dm, AppCtx *user)
498: {
499:   /* Set up the near null space (a.k.a. rigid body modes) that will be used by the multigrid preconditioner */
500:   DM             subdm;
501:   MatNullSpace   nearNullSpace;
502:   PetscInt       fields = 0;
503:   PetscObject    deformation;

507:   DMCreateSubDM(dm, 1, &fields, NULL, &subdm);
508:   DMPlexCreateRigidBody(subdm, &nearNullSpace);
509:   DMGetField(dm, 0, &deformation);
510:   PetscObjectCompose(deformation, "nearnullspace", (PetscObject) nearNullSpace);
511:   DMDestroy(&subdm);
512:   MatNullSpaceDestroy(&nearNullSpace);
513:   return(0);
514: }

516: PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
517: {
518:   DM              cdm   = dm, dmAux;
519:   const PetscInt  dim   = user->dim;
520:   const PetscBool simplex = user->simplex;
521:   PetscFE         fe[2], feAux[2];
522:   PetscQuadrature q, fq;
523:   PetscDS         prob, probAux;
524:   PetscErrorCode  ierr;

527:   /* Create finite element */
528:   PetscFECreateDefault(dm, dim, dim, simplex, "def_", PETSC_DEFAULT, &fe[0]);
529:   PetscObjectSetName((PetscObject) fe[0], "deformation");
530:   PetscFEGetQuadrature(fe[0], &q);
531:   PetscFEGetFaceQuadrature(fe[0], &fq);
532:   PetscFECreateDefault(dm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]);
533:   PetscFESetQuadrature(fe[1], q);
534:   PetscFESetFaceQuadrature(fe[1], fq);

536:   PetscObjectSetName((PetscObject) fe[1], "pressure");

538:   PetscFECreateDefault(dm, dim, 1, simplex, "elastMat_", PETSC_DEFAULT, &feAux[0]);
539:   PetscObjectSetName((PetscObject) feAux[0], "elasticityMaterial");
540:   PetscFESetQuadrature(feAux[0], q);
541:   PetscFESetFaceQuadrature(feAux[0], fq);
542:   /* It is not yet possible to define a field on a submesh (e.g. a boundary), so we will use a normal finite element */
543:   PetscFECreateDefault(dm, dim, 1, simplex, "wall_pres_", PETSC_DEFAULT, &feAux[1]);
544:   PetscObjectSetName((PetscObject) feAux[1], "wall_pressure");
545:   PetscFESetQuadrature(feAux[1], q);
546:   PetscFESetFaceQuadrature(feAux[1], fq);

548:   /* Set discretization and boundary conditions for each mesh */
549:   DMGetDS(dm, &prob);
550:   PetscDSSetDiscretization(prob, 0, (PetscObject) fe[0]);
551:   PetscDSSetDiscretization(prob, 1, (PetscObject) fe[1]);
552:   SetupProblem(prob, dim, user);
553:   PetscDSSetFromOptions(prob);
554:   PetscDSCreate(PetscObjectComm((PetscObject)dm),&probAux);
555:   PetscDSSetDiscretization(probAux, 0, (PetscObject) feAux[0]);
556:   PetscDSSetDiscretization(probAux, 1, (PetscObject) feAux[1]);
557:   PetscDSSetFromOptions(probAux);
558:   while (cdm) {
559:     DMSetDS(cdm, prob);
560:     DMClone(cdm, &dmAux);
561:     DMSetDS(dmAux, probAux);
562:     PetscObjectCompose((PetscObject) cdm, "dmAux", (PetscObject) dmAux);
563:     SetupMaterial(cdm, dmAux, user);
564:     DMDestroy(&dmAux);

566:     DMGetCoarseDM(cdm, &cdm);
567:   }
568:   PetscDSDestroy(&probAux);
569:   PetscFEDestroy(&fe[0]);
570:   PetscFEDestroy(&fe[1]);
571:   PetscFEDestroy(&feAux[0]);
572:   PetscFEDestroy(&feAux[1]);
573:   return(0);
574: }


577: int main(int argc, char **argv)
578: {
579:   SNES           snes;                 /* nonlinear solver */
580:   DM             dm;                   /* problem definition */
581:   Vec            u,r;                  /* solution, residual vectors */
582:   Mat            A,J;                  /* Jacobian matrix */
583:   AppCtx         user;                 /* user-defined work context */
584:   PetscInt       its;                  /* iterations for convergence */

587:   PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
588:   ProcessOptions(PETSC_COMM_WORLD, &user);
589:   SNESCreate(PETSC_COMM_WORLD, &snes);
590:   CreateMesh(PETSC_COMM_WORLD, &user, &dm);
591:   SNESSetDM(snes, dm);
592:   DMSetApplicationContext(dm, &user);

594:   SetupDiscretization(dm, &user);
595:   DMPlexCreateClosureIndex(dm, NULL);
596:   SetupNearNullSpace(dm, &user);

598:   DMCreateGlobalVector(dm, &u);
599:   VecDuplicate(u, &r);

601:   DMSetMatType(dm,MATAIJ);
602:   DMCreateMatrix(dm, &J);
603:   A = J;

605:   DMPlexSetSNESLocalFEM(dm,&user,&user,&user);
606:   SNESSetJacobian(snes, A, J, NULL, NULL);

608:   SNESSetFromOptions(snes);

610:   {
611:     PetscErrorCode (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void* ctx);
612:     initialGuess[0] = coordinates; initialGuess[1] = zero_scalar;
613:     DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u);
614:   }
615:   if (user.showInitial) {DMVecViewLocal(dm, u, PETSC_VIEWER_STDOUT_SELF);}

617:   if (user.runType == RUN_FULL) {
618:     if (user.debug) {
619:       PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n");
620:       VecView(u, PETSC_VIEWER_STDOUT_WORLD);
621:     }
622:     SNESSolve(snes, NULL, u);
623:     SNESGetIterationNumber(snes, &its);
624:     PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %D\n", its);
625:     if (user.showSolution) {
626:       PetscPrintf(PETSC_COMM_WORLD, "Solution\n");
627:       VecChop(u, 3.0e-9);
628:       VecView(u, PETSC_VIEWER_STDOUT_WORLD);
629:     }
630:   } else {
631:     PetscReal res = 0.0;

633:     /* Check initial guess */
634:     PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n");
635:     VecView(u, PETSC_VIEWER_STDOUT_WORLD);
636:     /* Check residual */
637:     SNESComputeFunction(snes, u, r);
638:     PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n");
639:     VecChop(r, 1.0e-10);
640:     VecView(r, PETSC_VIEWER_STDOUT_WORLD);
641:     VecNorm(r, NORM_2, &res);
642:     PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", res);
643:     /* Check Jacobian */
644:     {
645:       Vec b;

647:       SNESComputeJacobian(snes, u, A, A);
648:       VecDuplicate(u, &b);
649:       VecSet(r, 0.0);
650:       SNESComputeFunction(snes, r, b);
651:       MatMult(A, u, r);
652:       VecAXPY(r, 1.0, b);
653:       VecDestroy(&b);
654:       PetscPrintf(PETSC_COMM_WORLD, "Au - b = Au + F(0)\n");
655:       VecChop(r, 1.0e-10);
656:       VecView(r, PETSC_VIEWER_STDOUT_WORLD);
657:       VecNorm(r, NORM_2, &res);
658:       PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", res);
659:     }
660:   }
661:   VecViewFromOptions(u, NULL, "-sol_vec_view");

663:   if (A != J) {MatDestroy(&A);}
664:   MatDestroy(&J);
665:   VecDestroy(&u);
666:   VecDestroy(&r);
667:   SNESDestroy(&snes);
668:   DMDestroy(&dm);
669:   PetscFinalize();
670:   return ierr;
671: }

673: /*TEST
674:   build:
675:     requires: !complex

677:   test:
678:     suffix: 0
679:     requires: ctetgen
680:     args: -run_type full -dim 3 -dm_refine 3 -interpolate 1 -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4 -def_petscspace_order 2 -pres_petscspace_order 1 -elastMat_petscspace_order 0 -wall_pres_petscspace_order 0 -snes_rtol 1e-05 -ksp_type fgmres -ksp_rtol 1e-10 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper -fieldsplit_deformation_ksp_type preonly -fieldsplit_deformation_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi -snes_monitor_short -ksp_monitor_short -snes_converged_reason -ksp_converged_reason -show_solution 0
681:   test:
682:     suffix: 1
683:     requires: ctetgen superlu_dist
684:     nsize: 2
685:     args: -run_type full -dim 3 -dm_refine 0 -interpolate 1 -test_partition -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4 -def_petscspace_order 2 -pres_petscspace_order 1 -elastMat_petscspace_order 0 -wall_pres_petscspace_order 0 -snes_rtol 1e-05 -ksp_type fgmres -ksp_rtol 1e-10 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper -fieldsplit_deformation_ksp_type preonly -fieldsplit_deformation_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi -snes_monitor_short -ksp_monitor_short -snes_converged_reason -ksp_converged_reason -show_solution 0
686:   test:
687:     suffix: 2
688:     requires: ctetgen
689:     args: -run_type full -dim 3 -dm_refine 2 -interpolate 1 -bc_fixed 3,4,5,6 -bc_pressure 2 -wall_pressure 1.0 -def_petscspace_order 2 -pres_petscspace_order 1 -elastMat_petscspace_order 0 -wall_pres_petscspace_order 0 -snes_rtol 1e-05 -ksp_type fgmres -ksp_rtol 1e-10 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper -fieldsplit_deformation_ksp_type preonly -fieldsplit_deformation_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi -snes_monitor_short -ksp_monitor_short -snes_converged_reason -ksp_converged_reason -show_solution 0
690: TEST*/