xref: /petsc/src/snes/tutorials/ex77.c (revision 6a5217c03994f2d95bb2e6dbd8bed42381aeb015)
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   PetscErrorCode ierr;
263 
264   PetscFunctionBeginUser;
265   options->runType      = RUN_FULL;
266   options->mu           = 1.0;
267   options->p_wall       = 0.4;
268 
269   ierr = PetscOptionsBegin(comm, "", "Nonlinear elasticity problem options", "DMPLEX");PetscCall(ierr);
270   run  = options->runType;
271   PetscCall(PetscOptionsEList("-run_type", "The run type", "ex77.c", runTypes, 2, runTypes[options->runType], &run, NULL));
272   options->runType = (RunType) run;
273   PetscCall(PetscOptionsReal("-shear_modulus", "The shear modulus", "ex77.c", options->mu, &options->mu, NULL));
274   PetscCall(PetscOptionsReal("-wall_pressure", "The wall pressure", "ex77.c", options->p_wall, &options->p_wall, NULL));
275   ierr = PetscOptionsEnd();PetscCall(ierr);
276   PetscFunctionReturn(0);
277 }
278 
279 PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
280 {
281   PetscFunctionBeginUser;
282   /* TODO The P1 coordinate space gives wrong results when compared to the affine version. Track this down */
283   if (0) {
284     PetscCall(DMPlexCreateBoxMesh(comm, 3, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm));
285   } else {
286     PetscCall(DMCreate(comm, dm));
287     PetscCall(DMSetType(*dm, DMPLEX));
288   }
289   PetscCall(DMSetFromOptions(*dm));
290   /* Label the faces (bit of a hack here, until it is properly implemented for simplices) */
291   PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view"));
292   {
293     DM              cdm;
294     DMLabel         label;
295     IS              is;
296     PetscInt        d, dim, b, f, Nf;
297     const PetscInt *faces;
298     PetscInt        csize;
299     PetscScalar    *coords = NULL;
300     PetscSection    cs;
301     Vec             coordinates ;
302 
303     PetscCall(DMGetDimension(*dm, &dim));
304     PetscCall(DMCreateLabel(*dm, "boundary"));
305     PetscCall(DMGetLabel(*dm, "boundary", &label));
306     PetscCall(DMPlexMarkBoundaryFaces(*dm, 1, label));
307     PetscCall(DMGetStratumIS(*dm, "boundary", 1,  &is));
308     if (is) {
309       PetscReal faceCoord;
310       PetscInt  v;
311 
312       PetscCall(ISGetLocalSize(is, &Nf));
313       PetscCall(ISGetIndices(is, &faces));
314 
315       PetscCall(DMGetCoordinatesLocal(*dm, &coordinates));
316       PetscCall(DMGetCoordinateDM(*dm, &cdm));
317       PetscCall(DMGetLocalSection(cdm, &cs));
318 
319       /* Check for each boundary face if any component of its centroid is either 0.0 or 1.0 */
320       for (f = 0; f < Nf; ++f) {
321         PetscCall(DMPlexVecGetClosure(cdm, cs, coordinates, faces[f], &csize, &coords));
322         /* Calculate mean coordinate vector */
323         for (d = 0; d < dim; ++d) {
324           const PetscInt Nv = csize/dim;
325           faceCoord = 0.0;
326           for (v = 0; v < Nv; ++v) faceCoord += PetscRealPart(coords[v*dim+d]);
327           faceCoord /= Nv;
328           for (b = 0; b < 2; ++b) {
329             if (PetscAbs(faceCoord - b*1.0) < PETSC_SMALL) {
330               PetscCall(DMSetLabelValue(*dm, "Faces", faces[f], d*2+b+1));
331             }
332           }
333         }
334         PetscCall(DMPlexVecRestoreClosure(cdm, cs, coordinates, faces[f], &csize, &coords));
335       }
336       PetscCall(ISRestoreIndices(is, &faces));
337     }
338     PetscCall(ISDestroy(&is));
339   }
340   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
341   PetscFunctionReturn(0);
342 }
343 
344 PetscErrorCode SetupProblem(DM dm, PetscInt dim, AppCtx *user)
345 {
346   PetscDS        ds;
347   PetscWeakForm  wf;
348   DMLabel        label;
349   PetscInt       bd;
350 
351   PetscFunctionBeginUser;
352   PetscCall(DMGetDS(dm, &ds));
353   PetscCall(PetscDSSetResidual(ds, 0, NULL, f1_u_3d));
354   PetscCall(PetscDSSetResidual(ds, 1, f0_p_3d, NULL));
355   PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL,  NULL,  g3_uu_3d));
356   PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL,  g2_up_3d, NULL));
357   PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu_3d, NULL,  NULL));
358 
359   PetscCall(DMGetLabel(dm, "Faces", &label));
360   PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "pressure", label, 0, NULL, 0, 0, NULL, NULL, NULL, user, &bd));
361   PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
362   PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, 0, f0_bd_u_3d, 0, NULL));
363   PetscCall(PetscWeakFormSetIndexBdJacobian(wf, label, 1, 0, 0, 0, 0, NULL, 0, g1_bd_uu_3d, 0, NULL, 0, NULL));
364 
365   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "fixed", label, 0, NULL, 0, 0, NULL, (void (*)(void)) coordinates, NULL, user, NULL));
366   PetscFunctionReturn(0);
367 }
368 
369 PetscErrorCode SetupMaterial(DM dm, DM dmAux, AppCtx *user)
370 {
371   PetscErrorCode (*matFuncs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx) = {elasticityMaterial, wallPressure};
372   Vec            A;
373   void          *ctxs[2];
374 
375   PetscFunctionBegin;
376   ctxs[0] = user; ctxs[1] = user;
377   PetscCall(DMCreateLocalVector(dmAux, &A));
378   PetscCall(DMProjectFunctionLocal(dmAux, 0.0, matFuncs, ctxs, INSERT_ALL_VALUES, A));
379   PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, A));
380   PetscCall(VecDestroy(&A));
381   PetscFunctionReturn(0);
382 }
383 
384 PetscErrorCode SetupNearNullSpace(DM dm, AppCtx *user)
385 {
386   /* Set up the near null space (a.k.a. rigid body modes) that will be used by the multigrid preconditioner */
387   DM             subdm;
388   MatNullSpace   nearNullSpace;
389   PetscInt       fields = 0;
390   PetscObject    deformation;
391 
392   PetscFunctionBeginUser;
393   PetscCall(DMCreateSubDM(dm, 1, &fields, NULL, &subdm));
394   PetscCall(DMPlexCreateRigidBody(subdm, 0, &nearNullSpace));
395   PetscCall(DMGetField(dm, 0, NULL, &deformation));
396   PetscCall(PetscObjectCompose(deformation, "nearnullspace", (PetscObject) nearNullSpace));
397   PetscCall(DMDestroy(&subdm));
398   PetscCall(MatNullSpaceDestroy(&nearNullSpace));
399   PetscFunctionReturn(0);
400 }
401 
402 static PetscErrorCode SetupAuxDM(DM dm, PetscInt NfAux, PetscFE feAux[], AppCtx *user)
403 {
404   DM             dmAux, coordDM;
405   PetscInt       f;
406 
407   PetscFunctionBegin;
408   /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */
409   PetscCall(DMGetCoordinateDM(dm, &coordDM));
410   if (!feAux) PetscFunctionReturn(0);
411   PetscCall(DMClone(dm, &dmAux));
412   PetscCall(DMSetCoordinateDM(dmAux, coordDM));
413   for (f = 0; f < NfAux; ++f) PetscCall(DMSetField(dmAux, f, NULL, (PetscObject) feAux[f]));
414   PetscCall(DMCreateDS(dmAux));
415   PetscCall(SetupMaterial(dm, dmAux, user));
416   PetscCall(DMDestroy(&dmAux));
417   PetscFunctionReturn(0);
418 }
419 
420 PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
421 {
422   DM              cdm = dm;
423   PetscFE         fe[2], feAux[2];
424   PetscBool       simplex;
425   PetscInt        dim;
426   MPI_Comm        comm;
427 
428   PetscFunctionBeginUser;
429   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
430   PetscCall(DMGetDimension(dm, &dim));
431   PetscCall(DMPlexIsSimplex(dm, &simplex));
432   /* Create finite element */
433   PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "def_", PETSC_DEFAULT, &fe[0]));
434   PetscCall(PetscObjectSetName((PetscObject) fe[0], "deformation"));
435   PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]));
436   PetscCall(PetscFECopyQuadrature(fe[0], fe[1]));
437 
438   PetscCall(PetscObjectSetName((PetscObject) fe[1], "pressure"));
439 
440   PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "elastMat_", PETSC_DEFAULT, &feAux[0]));
441   PetscCall(PetscObjectSetName((PetscObject) feAux[0], "elasticityMaterial"));
442   PetscCall(PetscFECopyQuadrature(fe[0], feAux[0]));
443   /* It is not yet possible to define a field on a submesh (e.g. a boundary), so we will use a normal finite element */
444   PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "wall_pres_", PETSC_DEFAULT, &feAux[1]));
445   PetscCall(PetscObjectSetName((PetscObject) feAux[1], "wall_pressure"));
446   PetscCall(PetscFECopyQuadrature(fe[0], feAux[1]));
447 
448   /* Set discretization and boundary conditions for each mesh */
449   PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe[0]));
450   PetscCall(DMSetField(dm, 1, NULL, (PetscObject) fe[1]));
451   PetscCall(DMCreateDS(dm));
452   PetscCall(SetupProblem(dm, dim, user));
453   while (cdm) {
454     PetscCall(SetupAuxDM(cdm, 2, feAux, user));
455     PetscCall(DMCopyDisc(dm, cdm));
456     PetscCall(DMGetCoarseDM(cdm, &cdm));
457   }
458   PetscCall(PetscFEDestroy(&fe[0]));
459   PetscCall(PetscFEDestroy(&fe[1]));
460   PetscCall(PetscFEDestroy(&feAux[0]));
461   PetscCall(PetscFEDestroy(&feAux[1]));
462   PetscFunctionReturn(0);
463 }
464 
465 int main(int argc, char **argv)
466 {
467   SNES           snes;                 /* nonlinear solver */
468   DM             dm;                   /* problem definition */
469   Vec            u,r;                  /* solution, residual vectors */
470   Mat            A,J;                  /* Jacobian matrix */
471   AppCtx         user;                 /* user-defined work context */
472   PetscInt       its;                  /* iterations for convergence */
473 
474   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
475   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
476   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
477   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
478   PetscCall(SNESSetDM(snes, dm));
479   PetscCall(DMSetApplicationContext(dm, &user));
480 
481   PetscCall(SetupDiscretization(dm, &user));
482   PetscCall(DMPlexCreateClosureIndex(dm, NULL));
483   PetscCall(SetupNearNullSpace(dm, &user));
484 
485   PetscCall(DMCreateGlobalVector(dm, &u));
486   PetscCall(VecDuplicate(u, &r));
487 
488   PetscCall(DMSetMatType(dm,MATAIJ));
489   PetscCall(DMCreateMatrix(dm, &J));
490   A = J;
491 
492   PetscCall(DMPlexSetSNESLocalFEM(dm,&user,&user,&user));
493   PetscCall(SNESSetJacobian(snes, A, J, NULL, NULL));
494 
495   PetscCall(SNESSetFromOptions(snes));
496 
497   {
498     PetscErrorCode (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void* ctx);
499     initialGuess[0] = coordinates; initialGuess[1] = zero_scalar;
500     PetscCall(DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u));
501   }
502 
503   if (user.runType == RUN_FULL) {
504     PetscCall(SNESSolve(snes, NULL, u));
505     PetscCall(SNESGetIterationNumber(snes, &its));
506     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %D\n", its));
507   } else {
508     PetscReal res = 0.0;
509 
510     /* Check initial guess */
511     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n"));
512     PetscCall(VecView(u, PETSC_VIEWER_STDOUT_WORLD));
513     /* Check residual */
514     PetscCall(SNESComputeFunction(snes, u, r));
515     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n"));
516     PetscCall(VecChop(r, 1.0e-10));
517     PetscCall(VecView(r, PETSC_VIEWER_STDOUT_WORLD));
518     PetscCall(VecNorm(r, NORM_2, &res));
519     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", (double)res));
520     /* Check Jacobian */
521     {
522       Vec b;
523 
524       PetscCall(SNESComputeJacobian(snes, u, A, A));
525       PetscCall(VecDuplicate(u, &b));
526       PetscCall(VecSet(r, 0.0));
527       PetscCall(SNESComputeFunction(snes, r, b));
528       PetscCall(MatMult(A, u, r));
529       PetscCall(VecAXPY(r, 1.0, b));
530       PetscCall(VecDestroy(&b));
531       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Au - b = Au + F(0)\n"));
532       PetscCall(VecChop(r, 1.0e-10));
533       PetscCall(VecView(r, PETSC_VIEWER_STDOUT_WORLD));
534       PetscCall(VecNorm(r, NORM_2, &res));
535       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", (double)res));
536     }
537   }
538   PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
539 
540   if (A != J) PetscCall(MatDestroy(&A));
541   PetscCall(MatDestroy(&J));
542   PetscCall(VecDestroy(&u));
543   PetscCall(VecDestroy(&r));
544   PetscCall(SNESDestroy(&snes));
545   PetscCall(DMDestroy(&dm));
546   PetscCall(PetscFinalize());
547   return 0;
548 }
549 
550 /*TEST
551 
552   build:
553     requires: !complex
554 
555   testset:
556     requires: ctetgen
557     args: -run_type full -dm_plex_dim 3 \
558           -def_petscspace_degree 2 -pres_petscspace_degree 1 -elastMat_petscspace_degree 0 -wall_pres_petscspace_degree 0 \
559           -snes_rtol 1e-05 -snes_monitor_short -snes_converged_reason \
560           -ksp_type fgmres -ksp_rtol 1e-10 -ksp_monitor_short -ksp_converged_reason \
561           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper \
562             -fieldsplit_deformation_ksp_type preonly -fieldsplit_deformation_pc_type lu \
563             -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
564 
565     test:
566       suffix: 0
567       requires: !single
568       args: -dm_refine 2 \
569             -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4
570 
571     test:
572       suffix: 1
573       requires: superlu_dist
574       nsize: 2
575       args: -dm_refine 0 -petscpartitioner_type simple \
576             -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4
577       timeoutfactor: 2
578 
579     test:
580       suffix: 4
581       requires: superlu_dist
582       nsize: 2
583       args: -dm_refine 0 -petscpartitioner_type simple \
584             -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4
585       output_file: output/ex77_1.out
586 
587     test:
588       suffix: 2
589       requires: !single
590       args: -dm_refine 2 \
591             -bc_fixed 3,4,5,6 -bc_pressure 2 -wall_pressure 1.0
592 
593     test:
594       suffix: 2_par
595       requires: superlu_dist !single
596       nsize: 4
597       args: -dm_refine 2 -petscpartitioner_type simple \
598             -bc_fixed 3,4,5,6 -bc_pressure 2 -wall_pressure 1.0
599       output_file: output/ex77_2.out
600 
601 TEST*/
602