xref: /petsc/src/snes/tutorials/ex77.c (revision 21e3ffae2f3b73c0bd738cf6d0a809700fc04bb0)
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 {
45   RUN_FULL,
46   RUN_TEST
47 } RunType;
48 
49 typedef struct {
50   RunType   runType; /* Whether to run tests, or solve the full problem */
51   PetscReal mu;      /* The shear modulus */
52   PetscReal p_wall;  /* The wall pressure */
53 } AppCtx;
54 
55 #if 0
56 static inline void Det2D(PetscReal *detJ, const PetscReal J[])
57 {
58   *detJ = J[0]*J[3] - J[1]*J[2];
59 }
60 #endif
61 
62 static inline void Det3D(PetscReal *detJ, const PetscScalar J[])
63 {
64   *detJ = PetscRealPart(J[0 * 3 + 0] * (J[1 * 3 + 1] * J[2 * 3 + 2] - J[1 * 3 + 2] * J[2 * 3 + 1]) + J[0 * 3 + 1] * (J[1 * 3 + 2] * J[2 * 3 + 0] - J[1 * 3 + 0] * J[2 * 3 + 2]) + J[0 * 3 + 2] * (J[1 * 3 + 0] * J[2 * 3 + 1] - J[1 * 3 + 1] * J[2 * 3 + 0]));
65 }
66 
67 #if 0
68 static inline void Cof2D(PetscReal C[], const PetscReal A[])
69 {
70   C[0] =  A[3];
71   C[1] = -A[2];
72   C[2] = -A[1];
73   C[3] =  A[0];
74 }
75 #endif
76 
77 static inline void Cof3D(PetscReal C[], const PetscScalar A[])
78 {
79   C[0 * 3 + 0] = PetscRealPart(A[1 * 3 + 1] * A[2 * 3 + 2] - A[1 * 3 + 2] * A[2 * 3 + 1]);
80   C[0 * 3 + 1] = PetscRealPart(A[1 * 3 + 2] * A[2 * 3 + 0] - A[1 * 3 + 0] * A[2 * 3 + 2]);
81   C[0 * 3 + 2] = PetscRealPart(A[1 * 3 + 0] * A[2 * 3 + 1] - A[1 * 3 + 1] * A[2 * 3 + 0]);
82   C[1 * 3 + 0] = PetscRealPart(A[0 * 3 + 2] * A[2 * 3 + 1] - A[0 * 3 + 1] * A[2 * 3 + 2]);
83   C[1 * 3 + 1] = PetscRealPart(A[0 * 3 + 0] * A[2 * 3 + 2] - A[0 * 3 + 2] * A[2 * 3 + 0]);
84   C[1 * 3 + 2] = PetscRealPart(A[0 * 3 + 1] * A[2 * 3 + 0] - A[0 * 3 + 0] * A[2 * 3 + 1]);
85   C[2 * 3 + 0] = PetscRealPart(A[0 * 3 + 1] * A[1 * 3 + 2] - A[0 * 3 + 2] * A[1 * 3 + 1]);
86   C[2 * 3 + 1] = PetscRealPart(A[0 * 3 + 2] * A[1 * 3 + 0] - A[0 * 3 + 0] * A[1 * 3 + 2]);
87   C[2 * 3 + 2] = PetscRealPart(A[0 * 3 + 0] * A[1 * 3 + 1] - A[0 * 3 + 1] * A[1 * 3 + 0]);
88 }
89 
90 PetscErrorCode zero_scalar(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
91 {
92   u[0] = 0.0;
93   return PETSC_SUCCESS;
94 }
95 
96 PetscErrorCode zero_vector(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
97 {
98   const PetscInt Ncomp = dim;
99 
100   PetscInt comp;
101   for (comp = 0; comp < Ncomp; ++comp) u[comp] = 0.0;
102   return PETSC_SUCCESS;
103 }
104 
105 PetscErrorCode coordinates(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
106 {
107   const PetscInt Ncomp = dim;
108 
109   PetscInt comp;
110   for (comp = 0; comp < Ncomp; ++comp) u[comp] = x[comp];
111   return PETSC_SUCCESS;
112 }
113 
114 PetscErrorCode elasticityMaterial(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
115 {
116   AppCtx *user = (AppCtx *)ctx;
117   u[0]         = user->mu;
118   return PETSC_SUCCESS;
119 }
120 
121 PetscErrorCode wallPressure(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
122 {
123   AppCtx *user = (AppCtx *)ctx;
124   u[0]         = user->p_wall;
125   return PETSC_SUCCESS;
126 }
127 
128 void f1_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
129 {
130   const PetscInt  Ncomp = dim;
131   const PetscReal mu = PetscRealPart(a[0]), kappa = 3.0;
132   PetscReal       cofu_x[9 /* Ncomp*dim */], detu_x, p = PetscRealPart(u[Ncomp]);
133   PetscInt        comp, d;
134 
135   Cof3D(cofu_x, u_x);
136   Det3D(&detu_x, u_x);
137   p += kappa * (detu_x - 1.0);
138 
139   /* f1 is the first Piola-Kirchhoff tensor */
140   for (comp = 0; comp < Ncomp; ++comp) {
141     for (d = 0; d < dim; ++d) f1[comp * dim + d] = mu * u_x[comp * dim + d] + p * cofu_x[comp * dim + d];
142   }
143 }
144 
145 void g3_uu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
146 {
147   const PetscInt  Ncomp = dim;
148   const PetscReal mu = PetscRealPart(a[0]), kappa = 3.0;
149   PetscReal       cofu_x[9 /* Ncomp*dim */], detu_x, pp, pm, p = PetscRealPart(u[Ncomp]);
150   PetscInt        compI, compJ, d1, d2;
151 
152   Cof3D(cofu_x, u_x);
153   Det3D(&detu_x, u_x);
154   p += kappa * (detu_x - 1.0);
155   pp = p / detu_x + kappa;
156   pm = p / detu_x;
157 
158   /* 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} */
159   for (compI = 0; compI < Ncomp; ++compI) {
160     for (compJ = 0; compJ < Ncomp; ++compJ) {
161       const PetscReal G = (compI == compJ) ? 1.0 : 0.0;
162 
163       for (d1 = 0; d1 < dim; ++d1) {
164         for (d2 = 0; d2 < dim; ++d2) {
165           const PetscReal g = (d1 == d2) ? 1.0 : 0.0;
166 
167           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];
168         }
169       }
170     }
171   }
172 }
173 
174 void f0_bd_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
175 {
176   const PetscInt    Ncomp = dim;
177   const PetscScalar p     = a[aOff[1]];
178   PetscReal         cofu_x[9 /* Ncomp*dim */];
179   PetscInt          comp, d;
180 
181   Cof3D(cofu_x, u_x);
182   for (comp = 0; comp < Ncomp; ++comp) {
183     for (d = 0, f0[comp] = 0.0; d < dim; ++d) f0[comp] += cofu_x[comp * dim + d] * n[d];
184     f0[comp] *= p;
185   }
186 }
187 
188 void g1_bd_uu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
189 {
190   const PetscInt Ncomp = dim;
191   PetscScalar    p     = a[aOff[1]];
192   PetscReal      cofu_x[9 /* Ncomp*dim */], m[3 /* Ncomp */], detu_x;
193   PetscInt       comp, compI, compJ, d;
194 
195   Cof3D(cofu_x, u_x);
196   Det3D(&detu_x, u_x);
197   p /= detu_x;
198 
199   for (comp = 0; comp < Ncomp; ++comp)
200     for (d = 0, m[comp] = 0.0; d < dim; ++d) m[comp] += cofu_x[comp * dim + d] * n[d];
201   for (compI = 0; compI < Ncomp; ++compI) {
202     for (compJ = 0; compJ < Ncomp; ++compJ) {
203       for (d = 0; d < dim; ++d) g1[(compI * Ncomp + compJ) * dim + d] = p * (m[compI] * cofu_x[compJ * dim + d] - cofu_x[compI * dim + d] * m[compJ]);
204     }
205   }
206 }
207 
208 void f0_p_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
209 {
210   PetscReal detu_x;
211   Det3D(&detu_x, u_x);
212   f0[0] = detu_x - 1.0;
213 }
214 
215 void g1_pu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
216 {
217   PetscReal cofu_x[9 /* Ncomp*dim */];
218   PetscInt  compI, d;
219 
220   Cof3D(cofu_x, u_x);
221   for (compI = 0; compI < dim; ++compI)
222     for (d = 0; d < dim; ++d) g1[compI * dim + d] = cofu_x[compI * dim + d];
223 }
224 
225 void g2_up_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
226 {
227   PetscReal cofu_x[9 /* Ncomp*dim */];
228   PetscInt  compI, d;
229 
230   Cof3D(cofu_x, u_x);
231   for (compI = 0; compI < dim; ++compI)
232     for (d = 0; d < dim; ++d) g2[compI * dim + d] = cofu_x[compI * dim + d];
233 }
234 
235 PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
236 {
237   const char *runTypes[2] = {"full", "test"};
238   PetscInt    run;
239 
240   PetscFunctionBeginUser;
241   options->runType = RUN_FULL;
242   options->mu      = 1.0;
243   options->p_wall  = 0.4;
244   PetscOptionsBegin(comm, "", "Nonlinear elasticity problem options", "DMPLEX");
245   run = options->runType;
246   PetscCall(PetscOptionsEList("-run_type", "The run type", "ex77.c", runTypes, 2, runTypes[options->runType], &run, NULL));
247   options->runType = (RunType)run;
248   PetscCall(PetscOptionsReal("-shear_modulus", "The shear modulus", "ex77.c", options->mu, &options->mu, NULL));
249   PetscCall(PetscOptionsReal("-wall_pressure", "The wall pressure", "ex77.c", options->p_wall, &options->p_wall, NULL));
250   PetscOptionsEnd();
251   PetscFunctionReturn(PETSC_SUCCESS);
252 }
253 
254 PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
255 {
256   PetscFunctionBeginUser;
257   /* TODO The P1 coordinate space gives wrong results when compared to the affine version. Track this down */
258   if (0) {
259     PetscCall(DMPlexCreateBoxMesh(comm, 3, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm));
260   } else {
261     PetscCall(DMCreate(comm, dm));
262     PetscCall(DMSetType(*dm, DMPLEX));
263   }
264   PetscCall(DMSetFromOptions(*dm));
265   /* Label the faces (bit of a hack here, until it is properly implemented for simplices) */
266   PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view"));
267   {
268     DM              cdm;
269     DMLabel         label;
270     IS              is;
271     PetscInt        d, dim, b, f, Nf;
272     const PetscInt *faces;
273     PetscInt        csize;
274     PetscScalar    *coords = NULL;
275     PetscSection    cs;
276     Vec             coordinates;
277 
278     PetscCall(DMGetDimension(*dm, &dim));
279     PetscCall(DMCreateLabel(*dm, "boundary"));
280     PetscCall(DMGetLabel(*dm, "boundary", &label));
281     PetscCall(DMPlexMarkBoundaryFaces(*dm, 1, label));
282     PetscCall(DMGetStratumIS(*dm, "boundary", 1, &is));
283     if (is) {
284       PetscReal faceCoord;
285       PetscInt  v;
286 
287       PetscCall(ISGetLocalSize(is, &Nf));
288       PetscCall(ISGetIndices(is, &faces));
289 
290       PetscCall(DMGetCoordinatesLocal(*dm, &coordinates));
291       PetscCall(DMGetCoordinateDM(*dm, &cdm));
292       PetscCall(DMGetLocalSection(cdm, &cs));
293 
294       /* Check for each boundary face if any component of its centroid is either 0.0 or 1.0 */
295       for (f = 0; f < Nf; ++f) {
296         PetscCall(DMPlexVecGetClosure(cdm, cs, coordinates, faces[f], &csize, &coords));
297         /* Calculate mean coordinate vector */
298         for (d = 0; d < dim; ++d) {
299           const PetscInt Nv = csize / dim;
300           faceCoord         = 0.0;
301           for (v = 0; v < Nv; ++v) faceCoord += PetscRealPart(coords[v * dim + d]);
302           faceCoord /= Nv;
303           for (b = 0; b < 2; ++b) {
304             if (PetscAbs(faceCoord - b * 1.0) < PETSC_SMALL) PetscCall(DMSetLabelValue(*dm, "Faces", faces[f], d * 2 + b + 1));
305           }
306         }
307         PetscCall(DMPlexVecRestoreClosure(cdm, cs, coordinates, faces[f], &csize, &coords));
308       }
309       PetscCall(ISRestoreIndices(is, &faces));
310     }
311     PetscCall(ISDestroy(&is));
312   }
313   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
314   PetscFunctionReturn(PETSC_SUCCESS);
315 }
316 
317 PetscErrorCode SetupProblem(DM dm, PetscInt dim, AppCtx *user)
318 {
319   PetscDS       ds;
320   PetscWeakForm wf;
321   DMLabel       label;
322   PetscInt      bd;
323 
324   PetscFunctionBeginUser;
325   PetscCall(DMGetDS(dm, &ds));
326   PetscCall(PetscDSSetResidual(ds, 0, NULL, f1_u_3d));
327   PetscCall(PetscDSSetResidual(ds, 1, f0_p_3d, NULL));
328   PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu_3d));
329   PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up_3d, NULL));
330   PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu_3d, NULL, NULL));
331 
332   PetscCall(DMGetLabel(dm, "Faces", &label));
333   PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "pressure", label, 0, NULL, 0, 0, NULL, NULL, NULL, user, &bd));
334   PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
335   PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, 0, f0_bd_u_3d, 0, NULL));
336   PetscCall(PetscWeakFormSetIndexBdJacobian(wf, label, 1, 0, 0, 0, 0, NULL, 0, g1_bd_uu_3d, 0, NULL, 0, NULL));
337 
338   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "fixed", label, 0, NULL, 0, 0, NULL, (void (*)(void))coordinates, NULL, user, NULL));
339   PetscFunctionReturn(PETSC_SUCCESS);
340 }
341 
342 PetscErrorCode SetupMaterial(DM dm, DM dmAux, AppCtx *user)
343 {
344   PetscErrorCode (*matFuncs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx) = {elasticityMaterial, wallPressure};
345   Vec   A;
346   void *ctxs[2];
347 
348   PetscFunctionBegin;
349   ctxs[0] = user;
350   ctxs[1] = user;
351   PetscCall(DMCreateLocalVector(dmAux, &A));
352   PetscCall(DMProjectFunctionLocal(dmAux, 0.0, matFuncs, ctxs, INSERT_ALL_VALUES, A));
353   PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, A));
354   PetscCall(VecDestroy(&A));
355   PetscFunctionReturn(PETSC_SUCCESS);
356 }
357 
358 PetscErrorCode SetupNearNullSpace(DM dm, AppCtx *user)
359 {
360   /* Set up the near null space (a.k.a. rigid body modes) that will be used by the multigrid preconditioner */
361   DM           subdm;
362   MatNullSpace nearNullSpace;
363   PetscInt     fields = 0;
364   PetscObject  deformation;
365 
366   PetscFunctionBeginUser;
367   PetscCall(DMCreateSubDM(dm, 1, &fields, NULL, &subdm));
368   PetscCall(DMPlexCreateRigidBody(subdm, 0, &nearNullSpace));
369   PetscCall(DMGetField(dm, 0, NULL, &deformation));
370   PetscCall(PetscObjectCompose(deformation, "nearnullspace", (PetscObject)nearNullSpace));
371   PetscCall(DMDestroy(&subdm));
372   PetscCall(MatNullSpaceDestroy(&nearNullSpace));
373   PetscFunctionReturn(PETSC_SUCCESS);
374 }
375 
376 static PetscErrorCode SetupAuxDM(DM dm, PetscInt NfAux, PetscFE feAux[], AppCtx *user)
377 {
378   DM       dmAux, coordDM;
379   PetscInt f;
380 
381   PetscFunctionBegin;
382   /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */
383   PetscCall(DMGetCoordinateDM(dm, &coordDM));
384   if (!feAux) PetscFunctionReturn(PETSC_SUCCESS);
385   PetscCall(DMClone(dm, &dmAux));
386   PetscCall(DMSetCoordinateDM(dmAux, coordDM));
387   for (f = 0; f < NfAux; ++f) PetscCall(DMSetField(dmAux, f, NULL, (PetscObject)feAux[f]));
388   PetscCall(DMCreateDS(dmAux));
389   PetscCall(SetupMaterial(dm, dmAux, user));
390   PetscCall(DMDestroy(&dmAux));
391   PetscFunctionReturn(PETSC_SUCCESS);
392 }
393 
394 PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
395 {
396   DM        cdm = dm;
397   PetscFE   fe[2], feAux[2];
398   PetscBool simplex;
399   PetscInt  dim;
400   MPI_Comm  comm;
401 
402   PetscFunctionBeginUser;
403   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
404   PetscCall(DMGetDimension(dm, &dim));
405   PetscCall(DMPlexIsSimplex(dm, &simplex));
406   /* Create finite element */
407   PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "def_", PETSC_DEFAULT, &fe[0]));
408   PetscCall(PetscObjectSetName((PetscObject)fe[0], "deformation"));
409   PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]));
410   PetscCall(PetscFECopyQuadrature(fe[0], fe[1]));
411 
412   PetscCall(PetscObjectSetName((PetscObject)fe[1], "pressure"));
413 
414   PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "elastMat_", PETSC_DEFAULT, &feAux[0]));
415   PetscCall(PetscObjectSetName((PetscObject)feAux[0], "elasticityMaterial"));
416   PetscCall(PetscFECopyQuadrature(fe[0], feAux[0]));
417   /* It is not yet possible to define a field on a submesh (e.g. a boundary), so we will use a normal finite element */
418   PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "wall_pres_", PETSC_DEFAULT, &feAux[1]));
419   PetscCall(PetscObjectSetName((PetscObject)feAux[1], "wall_pressure"));
420   PetscCall(PetscFECopyQuadrature(fe[0], feAux[1]));
421 
422   /* Set discretization and boundary conditions for each mesh */
423   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe[0]));
424   PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe[1]));
425   PetscCall(DMCreateDS(dm));
426   PetscCall(SetupProblem(dm, dim, user));
427   while (cdm) {
428     PetscCall(SetupAuxDM(cdm, 2, feAux, user));
429     PetscCall(DMCopyDisc(dm, cdm));
430     PetscCall(DMGetCoarseDM(cdm, &cdm));
431   }
432   PetscCall(PetscFEDestroy(&fe[0]));
433   PetscCall(PetscFEDestroy(&fe[1]));
434   PetscCall(PetscFEDestroy(&feAux[0]));
435   PetscCall(PetscFEDestroy(&feAux[1]));
436   PetscFunctionReturn(PETSC_SUCCESS);
437 }
438 
439 int main(int argc, char **argv)
440 {
441   SNES     snes; /* nonlinear solver */
442   DM       dm;   /* problem definition */
443   Vec      u, r; /* solution, residual vectors */
444   Mat      A, J; /* Jacobian matrix */
445   AppCtx   user; /* user-defined work context */
446   PetscInt its;  /* iterations for convergence */
447 
448   PetscFunctionBeginUser;
449   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
450   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
451   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
452   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
453   PetscCall(SNESSetDM(snes, dm));
454   PetscCall(DMSetApplicationContext(dm, &user));
455 
456   PetscCall(SetupDiscretization(dm, &user));
457   PetscCall(DMPlexCreateClosureIndex(dm, NULL));
458   PetscCall(SetupNearNullSpace(dm, &user));
459 
460   PetscCall(DMCreateGlobalVector(dm, &u));
461   PetscCall(VecDuplicate(u, &r));
462 
463   PetscCall(DMSetMatType(dm, MATAIJ));
464   PetscCall(DMCreateMatrix(dm, &J));
465   A = J;
466 
467   PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user));
468   PetscCall(SNESSetJacobian(snes, A, J, NULL, NULL));
469 
470   PetscCall(SNESSetFromOptions(snes));
471 
472   {
473     PetscErrorCode (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
474     initialGuess[0] = coordinates;
475     initialGuess[1] = zero_scalar;
476     PetscCall(DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u));
477   }
478 
479   if (user.runType == RUN_FULL) {
480     PetscCall(SNESSolve(snes, NULL, u));
481     PetscCall(SNESGetIterationNumber(snes, &its));
482     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its));
483   } else {
484     PetscReal res = 0.0;
485 
486     /* Check initial guess */
487     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n"));
488     PetscCall(VecView(u, PETSC_VIEWER_STDOUT_WORLD));
489     /* Check residual */
490     PetscCall(SNESComputeFunction(snes, u, r));
491     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n"));
492     PetscCall(VecChop(r, 1.0e-10));
493     PetscCall(VecView(r, PETSC_VIEWER_STDOUT_WORLD));
494     PetscCall(VecNorm(r, NORM_2, &res));
495     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", (double)res));
496     /* Check Jacobian */
497     {
498       Vec b;
499 
500       PetscCall(SNESComputeJacobian(snes, u, A, A));
501       PetscCall(VecDuplicate(u, &b));
502       PetscCall(VecSet(r, 0.0));
503       PetscCall(SNESComputeFunction(snes, r, b));
504       PetscCall(MatMult(A, u, r));
505       PetscCall(VecAXPY(r, 1.0, b));
506       PetscCall(VecDestroy(&b));
507       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Au - b = Au + F(0)\n"));
508       PetscCall(VecChop(r, 1.0e-10));
509       PetscCall(VecView(r, PETSC_VIEWER_STDOUT_WORLD));
510       PetscCall(VecNorm(r, NORM_2, &res));
511       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", (double)res));
512     }
513   }
514   PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
515 
516   if (A != J) PetscCall(MatDestroy(&A));
517   PetscCall(MatDestroy(&J));
518   PetscCall(VecDestroy(&u));
519   PetscCall(VecDestroy(&r));
520   PetscCall(SNESDestroy(&snes));
521   PetscCall(DMDestroy(&dm));
522   PetscCall(PetscFinalize());
523   return 0;
524 }
525 
526 /*TEST
527 
528   build:
529     requires: !complex
530 
531   testset:
532     requires: ctetgen
533     args: -run_type full -dm_plex_dim 3 \
534           -def_petscspace_degree 2 -pres_petscspace_degree 1 -elastMat_petscspace_degree 0 -wall_pres_petscspace_degree 0 \
535           -snes_rtol 1e-05 -snes_monitor_short -snes_converged_reason \
536           -ksp_type fgmres -ksp_rtol 1e-10 -ksp_monitor_short -ksp_converged_reason \
537           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper \
538             -fieldsplit_deformation_ksp_type preonly -fieldsplit_deformation_pc_type lu \
539             -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
540 
541     test:
542       suffix: 0
543       requires: !single
544       args: -dm_refine 2 \
545             -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4
546 
547     test:
548       suffix: 1
549       requires: superlu_dist
550       nsize: 2
551       args: -dm_refine 0 -petscpartitioner_type simple \
552             -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4
553       timeoutfactor: 2
554 
555     test:
556       suffix: 4
557       requires: superlu_dist
558       nsize: 2
559       args: -dm_refine 0 -petscpartitioner_type simple \
560             -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4
561       output_file: output/ex77_1.out
562 
563     test:
564       suffix: 2
565       requires: !single
566       args: -dm_refine 2 \
567             -bc_fixed 3,4,5,6 -bc_pressure 2 -wall_pressure 1.0
568 
569     test:
570       suffix: 2_par
571       requires: superlu_dist !single
572       nsize: 4
573       args: -dm_refine 2 -petscpartitioner_type simple \
574             -bc_fixed 3,4,5,6 -bc_pressure 2 -wall_pressure 1.0
575       output_file: output/ex77_2.out
576 
577 TEST*/
578