xref: /petsc/src/snes/tutorials/ex77.c (revision 4302210d98c281cf4e8bc0fee406cf222e373a4c)
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";
4 
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}
13 
14 Discretization:
15 
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.
19 
20 Field Data:
21 
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
25 
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}]
28 
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
32 
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}]
34 
35 Likewise, DMPlexVecSetClosure() takes data partitioned by field, and correctly
36 puts it into the Sieve ordering.
37 
38 */
39 
40 #include <petscdmplex.h>
41 #include <petscsnes.h>
42 #include <petscds.h>
43 
44 typedef enum {RUN_FULL, RUN_TEST} RunType;
45 
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;
60 
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
67 
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 }
74 
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
84 
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 }
97 
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 }
103 
104 PetscErrorCode zero_vector(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
105 {
106   const PetscInt Ncomp = dim;
107 
108   PetscInt       comp;
109   for (comp = 0; comp < Ncomp; ++comp) u[comp] = 0.0;
110   return 0;
111 }
112 
113 PetscErrorCode coordinates(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
114 {
115   const PetscInt Ncomp = dim;
116 
117   PetscInt       comp;
118   for (comp = 0; comp < Ncomp; ++comp) u[comp] = x[comp];
119   return 0;
120 }
121 
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 }
128 
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 }
135 
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;
145 
146   Cof3D(cofu_x, u_x);
147   Det3D(&detu_x, u_x);
148   p += kappa * (detu_x - 1.0);
149 
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 }
157 
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;
167 
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;
173 
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;
178 
179       for (d1 = 0; d1 < dim; ++d1) {
180         for (d2 = 0; d2 < dim; ++d2) {
181           const PetscReal g = (d1 == d2) ? 1.0 : 0.0;
182 
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 }
189 
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;
199 
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 }
206 
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;
216 
217   Cof3D(cofu_x, u_x);
218   Det3D(&detu_x, u_x);
219   p /= detu_x;
220 
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 }
230 
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 }
240 
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;
248 
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 }
253 
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;
261 
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 }
266 
267 PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
268 {
269   const char    *runTypes[2] = {"full", "test"};
270   PetscInt       run;
271   PetscErrorCode ierr;
272 
273   PetscFunctionBeginUser;
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;
285 
286   ierr = PetscOptionsBegin(comm, "", "Nonlinear elasticity problem options", "DMPLEX");CHKERRQ(ierr);
287   ierr = PetscOptionsInt("-debug", "The debugging level", "ex77.c", options->debug, &options->debug, NULL);CHKERRQ(ierr);
288   run  = options->runType;
289   ierr = PetscOptionsEList("-run_type", "The run type", "ex77.c", runTypes, 2, runTypes[options->runType], &run, NULL);CHKERRQ(ierr);
290 
291   options->runType = (RunType) run;
292 
293   ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex77.c", options->dim, &options->dim, NULL);CHKERRQ(ierr);
294   ierr = PetscOptionsBool("-interpolate", "Generate intermediate mesh elements", "ex77.c", options->interpolate, &options->interpolate, NULL);CHKERRQ(ierr);
295   ierr = PetscOptionsBool("-simplex", "Use simplices or tensor product cells", "ex77.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr);
296   ierr = PetscOptionsReal("-refinement_limit", "The largest allowable cell volume", "ex77.c", options->refinementLimit, &options->refinementLimit, NULL);CHKERRQ(ierr);
297   ierr = PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex77.c", options->testPartition, &options->testPartition, NULL);CHKERRQ(ierr);
298   ierr = PetscOptionsReal("-shear_modulus", "The shear modulus", "ex77.c", options->mu, &options->mu, NULL);CHKERRQ(ierr);
299   ierr = PetscOptionsReal("-wall_pressure", "The wall pressure", "ex77.c", options->p_wall, &options->p_wall, NULL);CHKERRQ(ierr);
300 
301   ierr = PetscOptionsBool("-show_initial", "Output the initial guess for verification", "ex77.c", options->showInitial, &options->showInitial, NULL);CHKERRQ(ierr);
302   ierr = PetscOptionsBool("-show_solution", "Output the solution for verification", "ex77.c", options->showSolution, &options->showSolution, NULL);CHKERRQ(ierr);
303   ierr = PetscOptionsEnd();
304 
305   ierr = PetscLogEventRegister("CreateMesh", DM_CLASSID, &options->createMeshEvent);CHKERRQ(ierr);
306   PetscFunctionReturn(0);
307 }
308 
309 PetscErrorCode DMVecViewLocal(DM dm, Vec v, PetscViewer viewer)
310 {
311   Vec            lv;
312   PetscInt       p;
313   PetscMPIInt    rank, size;
314   PetscErrorCode ierr;
315 
316   PetscFunctionBeginUser;
317   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRMPI(ierr);
318   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size);CHKERRMPI(ierr);
319   ierr = DMGetLocalVector(dm, &lv);CHKERRQ(ierr);
320   ierr = DMGlobalToLocalBegin(dm, v, INSERT_VALUES, lv);CHKERRQ(ierr);
321   ierr = DMGlobalToLocalEnd(dm, v, INSERT_VALUES, lv);CHKERRQ(ierr);
322   ierr = PetscPrintf(PETSC_COMM_WORLD, "Local function\n");CHKERRQ(ierr);
323   for (p = 0; p < size; ++p) {
324     if (p == rank) {ierr = VecView(lv, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);}
325     ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr);
326   }
327   ierr = DMRestoreLocalVector(dm, &lv);CHKERRQ(ierr);
328   PetscFunctionReturn(0);
329 }
330 
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};
337   PetscErrorCode ierr;
338 
339   PetscFunctionBeginUser;
340   ierr = PetscLogEventBegin(user->createMeshEvent,0,0,0,0);CHKERRQ(ierr);
341   ierr = DMPlexCreateBoxMesh(comm, dim, user->simplex, user->simplex ? NULL : cells, NULL, NULL, NULL, interpolate, dm);CHKERRQ(ierr);
342   /* Label the faces (bit of a hack here, until it is properly implemented for simplices) */
343   {
344     DM              cdm;
345     DMLabel         label;
346     IS              is;
347     PetscInt        d, dim = user->dim, b, f, Nf;
348     const PetscInt *faces;
349     PetscInt        csize;
350     PetscScalar    *coords = NULL;
351     PetscSection    cs;
352     Vec             coordinates ;
353 
354     ierr = DMCreateLabel(*dm, "boundary");CHKERRQ(ierr);
355     ierr = DMGetLabel(*dm, "boundary", &label);CHKERRQ(ierr);
356     ierr = DMPlexMarkBoundaryFaces(*dm, 1, label);CHKERRQ(ierr);
357     ierr = DMGetStratumIS(*dm, "boundary", 1,  &is);CHKERRQ(ierr);
358     if (is) {
359       PetscReal faceCoord;
360       PetscInt  v;
361 
362       ierr = ISGetLocalSize(is, &Nf);CHKERRQ(ierr);
363       ierr = ISGetIndices(is, &faces);CHKERRQ(ierr);
364 
365       ierr = DMGetCoordinatesLocal(*dm, &coordinates);CHKERRQ(ierr);
366       ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr);
367       ierr = DMGetLocalSection(cdm, &cs);CHKERRQ(ierr);
368 
369       /* Check for each boundary face if any component of its centroid is either 0.0 or 1.0 */
370       for (f = 0; f < Nf; ++f) {
371         ierr = DMPlexVecGetClosure(cdm, cs, coordinates, faces[f], &csize, &coords);CHKERRQ(ierr);
372         /* Calculate mean coordinate vector */
373         for (d = 0; d < dim; ++d) {
374           const PetscInt Nv = csize/dim;
375           faceCoord = 0.0;
376           for (v = 0; v < Nv; ++v) faceCoord += PetscRealPart(coords[v*dim+d]);
377           faceCoord /= Nv;
378           for (b = 0; b < 2; ++b) {
379             if (PetscAbs(faceCoord - b*1.0) < PETSC_SMALL) {
380               ierr = DMSetLabelValue(*dm, "Faces", faces[f], d*2+b+1);CHKERRQ(ierr);
381             }
382           }
383         }
384         ierr = DMPlexVecRestoreClosure(cdm, cs, coordinates, faces[f], &csize, &coords);CHKERRQ(ierr);
385       }
386       ierr = ISRestoreIndices(is, &faces);CHKERRQ(ierr);
387     }
388     ierr = ISDestroy(&is);CHKERRQ(ierr);
389   }
390   {
391     DM refinedMesh     = NULL;
392     DM distributedMesh = NULL;
393     PetscPartitioner part;
394 
395     ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr);
396 
397     /* Refine mesh using a volume constraint */
398     ierr = DMPlexSetRefinementLimit(*dm, refinementLimit);CHKERRQ(ierr);
399     if (user->simplex) {ierr = DMRefine(*dm, comm, &refinedMesh);CHKERRQ(ierr);}
400     if (refinedMesh) {
401       ierr = DMDestroy(dm);CHKERRQ(ierr);
402       *dm  = refinedMesh;
403     }
404     /* Setup test partitioning */
405     if (user->testPartition) {
406       PetscInt         triSizes_n2[2]       = {4, 4};
407       PetscInt         triPoints_n2[8]      = {3, 5, 6, 7, 0, 1, 2, 4};
408       PetscInt         triSizes_n3[3]       = {2, 3, 3};
409       PetscInt         triPoints_n3[8]      = {3, 5, 1, 6, 7, 0, 2, 4};
410       PetscInt         triSizes_n5[5]       = {1, 2, 2, 1, 2};
411       PetscInt         triPoints_n5[8]      = {3, 5, 6, 4, 7, 0, 1, 2};
412       PetscInt         triSizes_ref_n2[2]   = {8, 8};
413       PetscInt         triPoints_ref_n2[16] = {1, 5, 6, 7, 10, 11, 14, 15, 0, 2, 3, 4, 8, 9, 12, 13};
414       PetscInt         triSizes_ref_n3[3]   = {5, 6, 5};
415       PetscInt         triPoints_ref_n3[16] = {1, 7, 10, 14, 15, 2, 6, 8, 11, 12, 13, 0, 3, 4, 5, 9};
416       PetscInt         triSizes_ref_n5[5]   = {3, 4, 3, 3, 3};
417       PetscInt         triPoints_ref_n5[16] = {1, 7, 10, 2, 11, 13, 14, 5, 6, 15, 0, 8, 9, 3, 4, 12};
418       PetscInt         tetSizes_n2[2]       = {3, 3};
419       PetscInt         tetPoints_n2[6]      = {1, 2, 3, 0, 4, 5};
420       const PetscInt  *sizes = NULL;
421       const PetscInt  *points = NULL;
422       PetscInt         cEnd;
423       PetscMPIInt      rank, size;
424 
425       ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
426       ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
427       ierr = DMPlexGetHeightStratum(*dm, 0, NULL, &cEnd);CHKERRQ(ierr);
428       if (!rank) {
429         if (dim == 2 && user->simplex && size == 2 && cEnd == 8) {
430            sizes = triSizes_n2; points = triPoints_n2;
431         } else if (dim == 2 && user->simplex && size == 3 && cEnd == 8) {
432           sizes = triSizes_n3; points = triPoints_n3;
433         } else if (dim == 2 && user->simplex && size == 5 && cEnd == 8) {
434           sizes = triSizes_n5; points = triPoints_n5;
435         } else if (dim == 2 && user->simplex && size == 2 && cEnd == 16) {
436            sizes = triSizes_ref_n2; points = triPoints_ref_n2;
437         } else if (dim == 2 && user->simplex && size == 3 && cEnd == 16) {
438           sizes = triSizes_ref_n3; points = triPoints_ref_n3;
439         } else if (dim == 2 && user->simplex && size == 5 && cEnd == 16) {
440           sizes = triSizes_ref_n5; points = triPoints_ref_n5;
441         } else if (dim == 3 && user->simplex && size == 2 && cEnd == 6) {
442           sizes = tetSizes_n2; points = tetPoints_n2;
443         } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "No stored partition matching run parameters");
444       }
445       ierr = PetscPartitionerSetType(part, PETSCPARTITIONERSHELL);CHKERRQ(ierr);
446       ierr = PetscPartitionerShellSetPartition(part, size, sizes, points);CHKERRQ(ierr);
447     } else {
448       ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
449     }
450     /* Distribute mesh over processes */
451     ierr = DMPlexDistribute(*dm, 0, NULL, &distributedMesh);CHKERRQ(ierr);
452     if (distributedMesh) {
453       ierr = DMDestroy(dm);CHKERRQ(ierr);
454       *dm  = distributedMesh;
455     }
456   }
457   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
458   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
459 
460   ierr = PetscLogEventEnd(user->createMeshEvent,0,0,0,0);CHKERRQ(ierr);
461   PetscFunctionReturn(0);
462 }
463 
464 PetscErrorCode SetupProblem(DM dm, PetscInt dim, AppCtx *user)
465 {
466   PetscDS        ds;
467   PetscWeakForm  wf;
468   DMLabel        label;
469   PetscInt       bd;
470   PetscErrorCode ierr;
471 
472   PetscFunctionBeginUser;
473   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
474   ierr = PetscDSSetResidual(ds, 0, NULL, f1_u_3d);CHKERRQ(ierr);
475   ierr = PetscDSSetResidual(ds, 1, f0_p_3d, NULL);CHKERRQ(ierr);
476   ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL,  NULL,  g3_uu_3d);CHKERRQ(ierr);
477   ierr = PetscDSSetJacobian(ds, 0, 1, NULL, NULL,  g2_up_3d, NULL);CHKERRQ(ierr);
478   ierr = PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu_3d, NULL,  NULL);CHKERRQ(ierr);
479 
480   ierr = DMGetLabel(dm, "Faces", &label);CHKERRQ(ierr);
481   ierr = DMAddBoundary(dm, DM_BC_NATURAL, "pressure", label, 0, NULL, 0, 0, NULL, NULL, NULL, user, &bd);CHKERRQ(ierr);
482   ierr = PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
483   ierr = PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, f0_bd_u_3d, 0, NULL);CHKERRQ(ierr);
484   ierr = PetscWeakFormSetIndexBdJacobian(wf, label, 1, 0, 0, 0, NULL, 0, g1_bd_uu_3d, 0, NULL, 0, NULL);CHKERRQ(ierr);
485 
486   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "fixed", label, 0, NULL, 0, 0, NULL, (void (*)(void)) coordinates, NULL, user, NULL);CHKERRQ(ierr);
487   PetscFunctionReturn(0);
488 }
489 
490 PetscErrorCode SetupMaterial(DM dm, DM dmAux, AppCtx *user)
491 {
492   PetscErrorCode (*matFuncs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx) = {elasticityMaterial, wallPressure};
493   Vec            A;
494   void          *ctxs[2];
495   PetscErrorCode ierr;
496 
497   PetscFunctionBegin;
498   ctxs[0] = user; ctxs[1] = user;
499   ierr = DMCreateLocalVector(dmAux, &A);CHKERRQ(ierr);
500   ierr = DMProjectFunctionLocal(dmAux, 0.0, matFuncs, ctxs, INSERT_ALL_VALUES, A);CHKERRQ(ierr);
501   ierr = DMSetAuxiliaryVec(dm, NULL, 0, A);CHKERRQ(ierr);
502   ierr = VecDestroy(&A);CHKERRQ(ierr);
503   PetscFunctionReturn(0);
504 }
505 
506 PetscErrorCode SetupNearNullSpace(DM dm, AppCtx *user)
507 {
508   /* Set up the near null space (a.k.a. rigid body modes) that will be used by the multigrid preconditioner */
509   DM             subdm;
510   MatNullSpace   nearNullSpace;
511   PetscInt       fields = 0;
512   PetscObject    deformation;
513   PetscErrorCode ierr;
514 
515   PetscFunctionBeginUser;
516   ierr = DMCreateSubDM(dm, 1, &fields, NULL, &subdm);CHKERRQ(ierr);
517   ierr = DMPlexCreateRigidBody(subdm, 0, &nearNullSpace);CHKERRQ(ierr);
518   ierr = DMGetField(dm, 0, NULL, &deformation);CHKERRQ(ierr);
519   ierr = PetscObjectCompose(deformation, "nearnullspace", (PetscObject) nearNullSpace);CHKERRQ(ierr);
520   ierr = DMDestroy(&subdm);CHKERRQ(ierr);
521   ierr = MatNullSpaceDestroy(&nearNullSpace);CHKERRQ(ierr);
522   PetscFunctionReturn(0);
523 }
524 
525 static PetscErrorCode SetupAuxDM(DM dm, PetscInt NfAux, PetscFE feAux[], AppCtx *user)
526 {
527   DM             dmAux, coordDM;
528   PetscInt       f;
529   PetscErrorCode ierr;
530 
531   PetscFunctionBegin;
532   /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */
533   ierr = DMGetCoordinateDM(dm, &coordDM);CHKERRQ(ierr);
534   if (!feAux) PetscFunctionReturn(0);
535   ierr = DMClone(dm, &dmAux);CHKERRQ(ierr);
536   ierr = DMSetCoordinateDM(dmAux, coordDM);CHKERRQ(ierr);
537   for (f = 0; f < NfAux; ++f) {ierr = DMSetField(dmAux, f, NULL, (PetscObject) feAux[f]);CHKERRQ(ierr);}
538   ierr = DMCreateDS(dmAux);CHKERRQ(ierr);
539   ierr = SetupMaterial(dm, dmAux, user);CHKERRQ(ierr);
540   ierr = DMDestroy(&dmAux);CHKERRQ(ierr);
541   PetscFunctionReturn(0);
542 }
543 
544 PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
545 {
546   DM              cdm   = dm;
547   const PetscInt  dim   = user->dim;
548   const PetscBool simplex = user->simplex;
549   PetscFE         fe[2], feAux[2];
550   MPI_Comm        comm;
551   PetscErrorCode  ierr;
552 
553   PetscFunctionBeginUser;
554   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
555   /* Create finite element */
556   ierr = PetscFECreateDefault(comm, dim, dim, simplex, "def_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr);
557   ierr = PetscObjectSetName((PetscObject) fe[0], "deformation");CHKERRQ(ierr);
558   ierr = PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr);
559   ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr);
560 
561   ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr);
562 
563   ierr = PetscFECreateDefault(comm, dim, 1, simplex, "elastMat_", PETSC_DEFAULT, &feAux[0]);CHKERRQ(ierr);
564   ierr = PetscObjectSetName((PetscObject) feAux[0], "elasticityMaterial");CHKERRQ(ierr);
565   ierr = PetscFECopyQuadrature(fe[0], feAux[0]);CHKERRQ(ierr);
566   /* It is not yet possible to define a field on a submesh (e.g. a boundary), so we will use a normal finite element */
567   ierr = PetscFECreateDefault(comm, dim, 1, simplex, "wall_pres_", PETSC_DEFAULT, &feAux[1]);CHKERRQ(ierr);
568   ierr = PetscObjectSetName((PetscObject) feAux[1], "wall_pressure");CHKERRQ(ierr);
569   ierr = PetscFECopyQuadrature(fe[0], feAux[1]);CHKERRQ(ierr);
570 
571   /* Set discretization and boundary conditions for each mesh */
572   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr);
573   ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr);
574   ierr = DMCreateDS(dm);CHKERRQ(ierr);
575   ierr = SetupProblem(dm, dim, user);CHKERRQ(ierr);
576   while (cdm) {
577     ierr = SetupAuxDM(cdm, 2, feAux, user);CHKERRQ(ierr);
578     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
579     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
580   }
581   ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr);
582   ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr);
583   ierr = PetscFEDestroy(&feAux[0]);CHKERRQ(ierr);
584   ierr = PetscFEDestroy(&feAux[1]);CHKERRQ(ierr);
585   PetscFunctionReturn(0);
586 }
587 
588 
589 int main(int argc, char **argv)
590 {
591   SNES           snes;                 /* nonlinear solver */
592   DM             dm;                   /* problem definition */
593   Vec            u,r;                  /* solution, residual vectors */
594   Mat            A,J;                  /* Jacobian matrix */
595   AppCtx         user;                 /* user-defined work context */
596   PetscInt       its;                  /* iterations for convergence */
597   PetscErrorCode ierr;
598 
599   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
600   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
601   ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
602   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
603   ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
604   ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr);
605 
606   ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr);
607   ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr);
608   ierr = SetupNearNullSpace(dm, &user);CHKERRQ(ierr);
609 
610   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
611   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
612 
613   ierr = DMSetMatType(dm,MATAIJ);CHKERRQ(ierr);
614   ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr);
615   A = J;
616 
617   ierr = DMPlexSetSNESLocalFEM(dm,&user,&user,&user);CHKERRQ(ierr);
618   ierr = SNESSetJacobian(snes, A, J, NULL, NULL);CHKERRQ(ierr);
619 
620   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
621 
622   {
623     PetscErrorCode (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void* ctx);
624     initialGuess[0] = coordinates; initialGuess[1] = zero_scalar;
625     ierr = DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u);CHKERRQ(ierr);
626   }
627   if (user.showInitial) {ierr = DMVecViewLocal(dm, u, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);}
628 
629   if (user.runType == RUN_FULL) {
630     if (user.debug) {
631       ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n");CHKERRQ(ierr);
632       ierr = VecView(u, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
633     }
634     ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
635     ierr = SNESGetIterationNumber(snes, &its);CHKERRQ(ierr);
636     ierr = PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %D\n", its);CHKERRQ(ierr);
637     if (user.showSolution) {
638       ierr = PetscPrintf(PETSC_COMM_WORLD, "Solution\n");CHKERRQ(ierr);
639       ierr = VecChop(u, 3.0e-9);CHKERRQ(ierr);
640       ierr = VecView(u, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
641     }
642   } else {
643     PetscReal res = 0.0;
644 
645     /* Check initial guess */
646     ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n");CHKERRQ(ierr);
647     ierr = VecView(u, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
648     /* Check residual */
649     ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr);
650     ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n");CHKERRQ(ierr);
651     ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
652     ierr = VecView(r, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
653     ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
654     ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", (double)res);CHKERRQ(ierr);
655     /* Check Jacobian */
656     {
657       Vec b;
658 
659       ierr = SNESComputeJacobian(snes, u, A, A);CHKERRQ(ierr);
660       ierr = VecDuplicate(u, &b);CHKERRQ(ierr);
661       ierr = VecSet(r, 0.0);CHKERRQ(ierr);
662       ierr = SNESComputeFunction(snes, r, b);CHKERRQ(ierr);
663       ierr = MatMult(A, u, r);CHKERRQ(ierr);
664       ierr = VecAXPY(r, 1.0, b);CHKERRQ(ierr);
665       ierr = VecDestroy(&b);CHKERRQ(ierr);
666       ierr = PetscPrintf(PETSC_COMM_WORLD, "Au - b = Au + F(0)\n");CHKERRQ(ierr);
667       ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr);
668       ierr = VecView(r, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
669       ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
670       ierr = PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", (double)res);CHKERRQ(ierr);
671     }
672   }
673   ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr);
674 
675   if (A != J) {ierr = MatDestroy(&A);CHKERRQ(ierr);}
676   ierr = MatDestroy(&J);CHKERRQ(ierr);
677   ierr = VecDestroy(&u);CHKERRQ(ierr);
678   ierr = VecDestroy(&r);CHKERRQ(ierr);
679   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
680   ierr = DMDestroy(&dm);CHKERRQ(ierr);
681   ierr = PetscFinalize();
682   return ierr;
683 }
684 
685 /*TEST
686 
687   build:
688     requires: !complex
689 
690   test:
691     suffix: 0
692     requires: ctetgen !single
693     args: -run_type full -dim 3 -dm_refine 2 -interpolate 1 -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4 -def_petscspace_degree 2 -pres_petscspace_degree 1 -elastMat_petscspace_degree 0 -wall_pres_petscspace_degree 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
694 
695   test:
696     suffix: 1
697     requires: ctetgen superlu_dist
698     nsize: 2
699     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_degree 2 -pres_petscspace_degree 1 -elastMat_petscspace_degree 0 -wall_pres_petscspace_degree 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
700     timeoutfactor: 2
701 
702   test:
703     suffix: 2
704     requires: ctetgen !single
705     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_degree 2 -pres_petscspace_degree 1 -elastMat_petscspace_degree 0 -wall_pres_petscspace_degree 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
706 
707   test:
708     requires: ctetgen superlu_dist
709     suffix: 4
710     nsize: 2
711     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_degree 2 -pres_petscspace_degree 1 -elastMat_petscspace_degree 0 -wall_pres_petscspace_degree 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
712     output_file: output/ex77_1.out
713 
714   test:
715     requires: ctetgen !single
716     suffix: 3
717     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_degree 2 -pres_petscspace_degree 1 -elastMat_petscspace_degree 0 -wall_pres_petscspace_degree 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
718     output_file: output/ex77_2.out
719 
720   #TODO: this example deadlocks for me when using ParMETIS
721   test:
722     requires: ctetgen superlu_dist !single
723     suffix: 2_par
724     nsize: 4
725     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_degree 2 -pres_petscspace_degree 1 -elastMat_petscspace_degree 0 -wall_pres_petscspace_degree 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 -petscpartitioner_type simple
726     output_file: output/ex77_2.out
727 
728 TEST*/
729