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