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