xref: /petsc/src/snes/tutorials/ex17.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown static char help[] = "Linear elasticity in 2d and 3d with finite elements.\n\
2c4762a1bSJed Brown We solve the elasticity problem in a rectangular\n\
3c4762a1bSJed Brown domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\
4c4762a1bSJed Brown This example supports automatic convergence estimation\n\
5c4762a1bSJed Brown and eventually adaptivity.\n\n\n";
6c4762a1bSJed Brown 
7c4762a1bSJed Brown /*
8c4762a1bSJed Brown   https://en.wikipedia.org/wiki/Linear_elasticity
9d12cd81dSMatthew G. Knepley 
10d12cd81dSMatthew G. Knepley   Converting elastic constants:
11d12cd81dSMatthew G. Knepley     lambda = E nu / ((1 + nu) (1 - 2 nu))
12d12cd81dSMatthew G. Knepley     mu     = E / (2 (1 + nu))
13c4762a1bSJed Brown */
14c4762a1bSJed Brown 
15c4762a1bSJed Brown #include <petscdmplex.h>
16c4762a1bSJed Brown #include <petscsnes.h>
17c4762a1bSJed Brown #include <petscds.h>
18d12cd81dSMatthew G. Knepley #include <petscbag.h>
19c4762a1bSJed Brown #include <petscconvest.h>
20c4762a1bSJed Brown 
21*9371c9d4SSatish Balay typedef enum {
22*9371c9d4SSatish Balay   SOL_VLAP_QUADRATIC,
23*9371c9d4SSatish Balay   SOL_ELAS_QUADRATIC,
24*9371c9d4SSatish Balay   SOL_VLAP_TRIG,
25*9371c9d4SSatish Balay   SOL_ELAS_TRIG,
26*9371c9d4SSatish Balay   SOL_ELAS_AXIAL_DISP,
27*9371c9d4SSatish Balay   SOL_ELAS_UNIFORM_STRAIN,
28*9371c9d4SSatish Balay   SOL_ELAS_GE,
29*9371c9d4SSatish Balay   SOL_MASS_QUADRATIC,
30*9371c9d4SSatish Balay   NUM_SOLUTION_TYPES
31*9371c9d4SSatish Balay } SolutionType;
32d12cd81dSMatthew G. Knepley const char *solutionTypes[NUM_SOLUTION_TYPES + 1] = {"vlap_quad", "elas_quad", "vlap_trig", "elas_trig", "elas_axial_disp", "elas_uniform_strain", "elas_ge", "mass_quad", "unknown"};
33482dd828SMatthew G. Knepley 
34*9371c9d4SSatish Balay typedef enum {
35*9371c9d4SSatish Balay   DEFORM_NONE,
36*9371c9d4SSatish Balay   DEFORM_SHEAR,
37*9371c9d4SSatish Balay   DEFORM_STEP,
38*9371c9d4SSatish Balay   NUM_DEFORM_TYPES
39*9371c9d4SSatish Balay } DeformType;
40482dd828SMatthew G. Knepley const char *deformTypes[NUM_DEFORM_TYPES + 1] = {"none", "shear", "step", "unknown"};
41c4762a1bSJed Brown 
42c4762a1bSJed Brown typedef struct {
43d12cd81dSMatthew G. Knepley   PetscScalar mu;     /* shear modulus */
44d12cd81dSMatthew G. Knepley   PetscScalar lambda; /* Lame's first parameter */
45d12cd81dSMatthew G. Knepley } Parameter;
46d12cd81dSMatthew G. Knepley 
47d12cd81dSMatthew G. Knepley typedef struct {
48c4762a1bSJed Brown   /* Domain and mesh definition */
49c4762a1bSJed Brown   char         dmType[256]; /* DM type for the solve */
50482dd828SMatthew G. Knepley   DeformType   deform;      /* Domain deformation type */
51c4762a1bSJed Brown   /* Problem definition */
52c4762a1bSJed Brown   SolutionType solType; /* Type of exact solution */
53d12cd81dSMatthew G. Knepley   PetscBag     bag;     /* Problem parameters */
54c4762a1bSJed Brown   /* Solver definition */
55c4762a1bSJed Brown   PetscBool    useNearNullspace; /* Use the rigid body modes as a near nullspace for AMG */
56c4762a1bSJed Brown } AppCtx;
57c4762a1bSJed Brown 
58*9371c9d4SSatish Balay static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
59c4762a1bSJed Brown   PetscInt d;
60c4762a1bSJed Brown   for (d = 0; d < dim; ++d) u[d] = 0.0;
61c4762a1bSJed Brown   return 0;
62c4762a1bSJed Brown }
63c4762a1bSJed Brown 
64*9371c9d4SSatish Balay static PetscErrorCode ge_shift(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
65d12cd81dSMatthew G. Knepley   PetscInt d;
66d12cd81dSMatthew G. Knepley   u[0] = 0.1;
67d12cd81dSMatthew G. Knepley   for (d = 1; d < dim; ++d) u[d] = 0.0;
68d12cd81dSMatthew G. Knepley   return 0;
69d12cd81dSMatthew G. Knepley }
70d12cd81dSMatthew G. Knepley 
71*9371c9d4SSatish Balay static PetscErrorCode quadratic_2d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
72c4762a1bSJed Brown   u[0] = x[0] * x[0];
73c4762a1bSJed Brown   u[1] = x[1] * x[1] - 2.0 * x[0] * x[1];
74c4762a1bSJed Brown   return 0;
75c4762a1bSJed Brown }
76c4762a1bSJed Brown 
77*9371c9d4SSatish Balay static PetscErrorCode quadratic_3d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
78c4762a1bSJed Brown   u[0] = x[0] * x[0];
79c4762a1bSJed Brown   u[1] = x[1] * x[1] - 2.0 * x[0] * x[1];
80c4762a1bSJed Brown   u[2] = x[2] * x[2] - 2.0 * x[1] * x[2];
81c4762a1bSJed Brown   return 0;
82c4762a1bSJed Brown }
83c4762a1bSJed Brown 
84c4762a1bSJed Brown /*
85c4762a1bSJed Brown   u = x^2
86c4762a1bSJed Brown   v = y^2 - 2xy
87c4762a1bSJed Brown   Delta <u,v> - f = <2, 2> - <2, 2>
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   u = x^2
90c4762a1bSJed Brown   v = y^2 - 2xy
91c4762a1bSJed Brown   w = z^2 - 2yz
92c4762a1bSJed Brown   Delta <u,v,w> - f = <2, 2, 2> - <2, 2, 2>
93c4762a1bSJed Brown */
94*9371c9d4SSatish Balay static void f0_vlap_quadratic_u(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[]) {
95c4762a1bSJed Brown   PetscInt d;
96c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f0[d] += 2.0;
97c4762a1bSJed Brown }
98c4762a1bSJed Brown 
99c4762a1bSJed Brown /*
100c4762a1bSJed Brown   u = x^2
101c4762a1bSJed Brown   v = y^2 - 2xy
102c4762a1bSJed Brown   \varepsilon = / 2x     -y    \
103c4762a1bSJed Brown                 \ -y   2y - 2x /
104c4762a1bSJed Brown   Tr(\varepsilon) = div u = 2y
105c4762a1bSJed Brown   div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij}
106c4762a1bSJed Brown     = \lambda \partial_j (2y) + 2\mu < 2-1, 2 >
107c4762a1bSJed Brown     = \lambda < 0, 2 > + \mu < 2, 4 >
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   u = x^2
110c4762a1bSJed Brown   v = y^2 - 2xy
111c4762a1bSJed Brown   w = z^2 - 2yz
112c4762a1bSJed Brown   \varepsilon = / 2x     -y       0   \
113c4762a1bSJed Brown                 | -y   2y - 2x   -z   |
114c4762a1bSJed Brown                 \  0     -z    2z - 2y/
115c4762a1bSJed Brown   Tr(\varepsilon) = div u = 2z
116c4762a1bSJed Brown   div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij}
117c4762a1bSJed Brown     = \lambda \partial_j (2z) + 2\mu < 2-1, 2-1, 2 >
118c4762a1bSJed Brown     = \lambda < 0, 0, 2 > + \mu < 2, 2, 4 >
119c4762a1bSJed Brown */
120*9371c9d4SSatish Balay static void f0_elas_quadratic_u(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[]) {
121c4762a1bSJed Brown   const PetscReal mu     = 1.0;
122c4762a1bSJed Brown   const PetscReal lambda = 1.0;
123c4762a1bSJed Brown   PetscInt        d;
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   for (d = 0; d < dim - 1; ++d) f0[d] += 2.0 * mu;
126c4762a1bSJed Brown   f0[dim - 1] += 2.0 * lambda + 4.0 * mu;
127c4762a1bSJed Brown }
128c4762a1bSJed Brown 
129*9371c9d4SSatish Balay static void f0_mass_quadratic_u(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[]) {
130482dd828SMatthew G. Knepley   if (dim == 2) {
131482dd828SMatthew G. Knepley     f0[0] -= x[0] * x[0];
132482dd828SMatthew G. Knepley     f0[1] -= x[1] * x[1] - 2.0 * x[0] * x[1];
133482dd828SMatthew G. Knepley   } else {
134482dd828SMatthew G. Knepley     f0[0] -= x[0] * x[0];
135482dd828SMatthew G. Knepley     f0[1] -= x[1] * x[1] - 2.0 * x[0] * x[1];
136482dd828SMatthew G. Knepley     f0[2] -= x[2] * x[2] - 2.0 * x[1] * x[2];
137482dd828SMatthew G. Knepley   }
138482dd828SMatthew G. Knepley }
139482dd828SMatthew G. Knepley 
140*9371c9d4SSatish Balay static PetscErrorCode trig_2d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
141c4762a1bSJed Brown   u[0] = PetscSinReal(2.0 * PETSC_PI * x[0]);
142c4762a1bSJed Brown   u[1] = PetscSinReal(2.0 * PETSC_PI * x[1]) - 2.0 * x[0] * x[1];
143c4762a1bSJed Brown   return 0;
144c4762a1bSJed Brown }
145c4762a1bSJed Brown 
146*9371c9d4SSatish Balay static PetscErrorCode trig_3d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
147c4762a1bSJed Brown   u[0] = PetscSinReal(2.0 * PETSC_PI * x[0]);
148c4762a1bSJed Brown   u[1] = PetscSinReal(2.0 * PETSC_PI * x[1]) - 2.0 * x[0] * x[1];
149c4762a1bSJed Brown   u[2] = PetscSinReal(2.0 * PETSC_PI * x[2]) - 2.0 * x[1] * x[2];
150c4762a1bSJed Brown   return 0;
151c4762a1bSJed Brown }
152c4762a1bSJed Brown 
153c4762a1bSJed Brown /*
154c4762a1bSJed Brown   u = sin(2 pi x)
155c4762a1bSJed Brown   v = sin(2 pi y) - 2xy
156c4762a1bSJed Brown   Delta <u,v> - f = <-4 pi^2 u, -4 pi^2 v> - <-4 pi^2 sin(2 pi x), -4 pi^2 sin(2 pi y)>
157c4762a1bSJed Brown 
158c4762a1bSJed Brown   u = sin(2 pi x)
159c4762a1bSJed Brown   v = sin(2 pi y) - 2xy
160c4762a1bSJed Brown   w = sin(2 pi z) - 2yz
161c4762a1bSJed Brown   Delta <u,v,2> - f = <-4 pi^2 u, -4 pi^2 v, -4 pi^2 w> - <-4 pi^2 sin(2 pi x), -4 pi^2 sin(2 pi y), -4 pi^2 sin(2 pi z)>
162c4762a1bSJed Brown */
163*9371c9d4SSatish Balay static void f0_vlap_trig_u(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[]) {
164c4762a1bSJed Brown   PetscInt d;
165c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f0[d] += -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]);
166c4762a1bSJed Brown }
167c4762a1bSJed Brown 
168c4762a1bSJed Brown /*
169c4762a1bSJed Brown   u = sin(2 pi x)
170c4762a1bSJed Brown   v = sin(2 pi y) - 2xy
171c4762a1bSJed Brown   \varepsilon = / 2 pi cos(2 pi x)             -y        \
172c4762a1bSJed Brown                 \      -y          2 pi cos(2 pi y) - 2x /
173c4762a1bSJed Brown   Tr(\varepsilon) = div u = 2 pi (cos(2 pi x) + cos(2 pi y)) - 2 x
174c4762a1bSJed Brown   div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij}
175c4762a1bSJed Brown     = \lambda \partial_j 2 pi (cos(2 pi x) + cos(2 pi y)) + 2\mu < -4 pi^2 sin(2 pi x) - 1, -4 pi^2 sin(2 pi y) >
176c4762a1bSJed Brown     = \lambda < -4 pi^2 sin(2 pi x) - 2, -4 pi^2 sin(2 pi y) > + \mu < -8 pi^2 sin(2 pi x) - 2, -8 pi^2 sin(2 pi y) >
177c4762a1bSJed Brown 
178c4762a1bSJed Brown   u = sin(2 pi x)
179c4762a1bSJed Brown   v = sin(2 pi y) - 2xy
180c4762a1bSJed Brown   w = sin(2 pi z) - 2yz
181c4762a1bSJed Brown   \varepsilon = / 2 pi cos(2 pi x)            -y                     0         \
182c4762a1bSJed Brown                 |         -y       2 pi cos(2 pi y) - 2x            -z         |
183c4762a1bSJed Brown                 \          0                  -z         2 pi cos(2 pi z) - 2y /
184c4762a1bSJed Brown   Tr(\varepsilon) = div u = 2 pi (cos(2 pi x) + cos(2 pi y) + cos(2 pi z)) - 2 x - 2 y
185c4762a1bSJed Brown   div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij}
186c4762a1bSJed Brown     = \lambda \partial_j (2 pi (cos(2 pi x) + cos(2 pi y) + cos(2 pi z)) - 2 x - 2 y) + 2\mu < -4 pi^2 sin(2 pi x) - 1, -4 pi^2 sin(2 pi y) - 1, -4 pi^2 sin(2 pi z) >
187c4762a1bSJed Brown     = \lambda < -4 pi^2 sin(2 pi x) - 2, -4 pi^2 sin(2 pi y) - 2, -4 pi^2 sin(2 pi z) > + 2\mu < -4 pi^2 sin(2 pi x) - 1, -4 pi^2 sin(2 pi y) - 1, -4 pi^2 sin(2 pi z) >
188c4762a1bSJed Brown */
189*9371c9d4SSatish Balay static void f0_elas_trig_u(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[]) {
190c4762a1bSJed Brown   const PetscReal mu     = 1.0;
191c4762a1bSJed Brown   const PetscReal lambda = 1.0;
192c4762a1bSJed Brown   const PetscReal fact   = 4.0 * PetscSqr(PETSC_PI);
193c4762a1bSJed Brown   PetscInt        d;
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f0[d] += -(2.0 * mu + lambda) * fact * PetscSinReal(2.0 * PETSC_PI * x[d]) - (d < dim - 1 ? 2.0 * (mu + lambda) : 0.0);
196c4762a1bSJed Brown }
197c4762a1bSJed Brown 
198*9371c9d4SSatish Balay static PetscErrorCode axial_disp_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
199c4762a1bSJed Brown   const PetscReal mu     = 1.0;
200c4762a1bSJed Brown   const PetscReal lambda = 1.0;
201c4762a1bSJed Brown   const PetscReal N      = 1.0;
202c4762a1bSJed Brown   PetscInt        d;
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   u[0] = (3. * lambda * lambda + 8. * lambda * mu + 4 * mu * mu) / (4 * mu * (3 * lambda * lambda + 5. * lambda * mu + 2 * mu * mu)) * N * x[0];
205c4762a1bSJed Brown   u[1] = -0.25 * lambda / mu / (lambda + mu) * N * x[1];
206c4762a1bSJed Brown   for (d = 2; d < dim; ++d) u[d] = 0.0;
207c4762a1bSJed Brown   return 0;
208c4762a1bSJed Brown }
209c4762a1bSJed Brown 
210c4762a1bSJed Brown /*
211c4762a1bSJed Brown   We will pull/push on the right side of a block of linearly elastic material. The uniform traction conditions on the
212c4762a1bSJed Brown   right side of the box will result in a uniform strain along x and y. The Neumann BC is given by
213c4762a1bSJed Brown 
214c4762a1bSJed Brown      n_i \sigma_{ij} = t_i
215c4762a1bSJed Brown 
216c4762a1bSJed Brown   u = (1/(2\mu) - 1) x
217c4762a1bSJed Brown   v = -y
218c4762a1bSJed Brown   f = 0
219c4762a1bSJed Brown   t = <4\mu/\lambda (\lambda + \mu), 0>
220c4762a1bSJed Brown   \varepsilon = / 1/(2\mu) - 1   0 \
221c4762a1bSJed Brown                 \ 0             -1 /
222c4762a1bSJed Brown   Tr(\varepsilon) = div u = 1/(2\mu) - 2
223c4762a1bSJed Brown   div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij}
224c4762a1bSJed Brown     = \lambda \partial_j (1/(2\mu) - 2) + 2\mu < 0, 0 >
225c4762a1bSJed Brown     = \lambda < 0, 0 > + \mu < 0, 0 > = 0
226c4762a1bSJed Brown   NBC =  <1,0> . <4\mu/\lambda (\lambda + \mu), 0> = 4\mu/\lambda (\lambda + \mu)
227c4762a1bSJed Brown 
228c4762a1bSJed Brown   u = x - 1/2
229c4762a1bSJed Brown   v = 0
230c4762a1bSJed Brown   w = 0
231c4762a1bSJed Brown   \varepsilon = / x  0  0 \
232c4762a1bSJed Brown                 | 0  0  0 |
233c4762a1bSJed Brown                 \ 0  0  0 /
234c4762a1bSJed Brown   Tr(\varepsilon) = div u = x
235c4762a1bSJed Brown   div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij}
236c4762a1bSJed Brown     = \lambda \partial_j x + 2\mu < 1, 0, 0 >
237c4762a1bSJed Brown     = \lambda < 1, 0, 0 > + \mu < 2, 0, 0 >
238c4762a1bSJed Brown */
239*9371c9d4SSatish Balay static void f0_elas_axial_disp_bd_u(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[]) {
240c4762a1bSJed Brown   const PetscReal N = -1.0;
241c4762a1bSJed Brown 
242c4762a1bSJed Brown   f0[0] = N;
243c4762a1bSJed Brown }
244c4762a1bSJed Brown 
245*9371c9d4SSatish Balay static PetscErrorCode uniform_strain_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
246c4762a1bSJed Brown   const PetscReal eps_xx = 0.1;
247c4762a1bSJed Brown   const PetscReal eps_xy = 0.3;
248c4762a1bSJed Brown   const PetscReal eps_yy = 0.25;
249c4762a1bSJed Brown   PetscInt        d;
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   u[0] = eps_xx * x[0] + eps_xy * x[1];
252c4762a1bSJed Brown   u[1] = eps_xy * x[0] + eps_yy * x[1];
253c4762a1bSJed Brown   for (d = 2; d < dim; ++d) u[d] = 0.0;
254c4762a1bSJed Brown   return 0;
255c4762a1bSJed Brown }
256c4762a1bSJed Brown 
257*9371c9d4SSatish Balay static void f0_mass_u(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[]) {
258482dd828SMatthew G. Knepley   const PetscInt Nc = dim;
259482dd828SMatthew G. Knepley   PetscInt       c;
260482dd828SMatthew G. Knepley 
261482dd828SMatthew G. Knepley   for (c = 0; c < Nc; ++c) f0[c] = u[c];
262482dd828SMatthew G. Knepley }
263482dd828SMatthew G. Knepley 
264*9371c9d4SSatish Balay static void f1_vlap_u(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[]) {
265c4762a1bSJed Brown   const PetscInt Nc = dim;
266c4762a1bSJed Brown   PetscInt       c, d;
267c4762a1bSJed Brown 
268*9371c9d4SSatish Balay   for (c = 0; c < Nc; ++c)
269*9371c9d4SSatish Balay     for (d = 0; d < dim; ++d) f1[c * dim + d] += u_x[c * dim + d];
270c4762a1bSJed Brown }
271c4762a1bSJed Brown 
272*9371c9d4SSatish Balay static void f1_elas_u(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[]) {
273c4762a1bSJed Brown   const PetscInt  Nc     = dim;
274c4762a1bSJed Brown   const PetscReal mu     = 1.0;
275c4762a1bSJed Brown   const PetscReal lambda = 1.0;
276c4762a1bSJed Brown   PetscInt        c, d;
277c4762a1bSJed Brown 
278c4762a1bSJed Brown   for (c = 0; c < Nc; ++c) {
279c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
280c4762a1bSJed Brown       f1[c * dim + d] += mu * (u_x[c * dim + d] + u_x[d * dim + c]);
281c4762a1bSJed Brown       f1[c * dim + c] += lambda * u_x[d * dim + d];
282c4762a1bSJed Brown     }
283c4762a1bSJed Brown   }
284c4762a1bSJed Brown }
285c4762a1bSJed Brown 
286*9371c9d4SSatish Balay static void g0_mass_uu(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 g0[]) {
287482dd828SMatthew G. Knepley   const PetscInt Nc = dim;
288482dd828SMatthew G. Knepley   PetscInt       c;
289482dd828SMatthew G. Knepley 
290482dd828SMatthew G. Knepley   for (c = 0; c < Nc; ++c) g0[c * Nc + c] = 1.0;
291482dd828SMatthew G. Knepley }
292482dd828SMatthew G. Knepley 
293*9371c9d4SSatish Balay static void g3_vlap_uu(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[]) {
294c4762a1bSJed Brown   const PetscInt Nc = dim;
295c4762a1bSJed Brown   PetscInt       c, d;
296c4762a1bSJed Brown 
297c4762a1bSJed Brown   for (c = 0; c < Nc; ++c) {
298*9371c9d4SSatish Balay     for (d = 0; d < dim; ++d) { g3[((c * Nc + c) * dim + d) * dim + d] = 1.0; }
299c4762a1bSJed Brown   }
300c4762a1bSJed Brown }
301c4762a1bSJed Brown 
302c4762a1bSJed Brown /*
303c4762a1bSJed Brown   \partial_df \phi_fc g_{fc,gc,df,dg} \partial_dg \phi_gc
304c4762a1bSJed Brown 
305c4762a1bSJed Brown   \partial_df \phi_fc \lambda \delta_{fc,df} \sum_gc \partial_dg \phi_gc \delta_{gc,dg}
306c4762a1bSJed Brown   = \partial_fc \phi_fc \sum_gc \partial_gc \phi_gc
307c4762a1bSJed Brown */
308*9371c9d4SSatish Balay static void g3_elas_uu(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[]) {
309c4762a1bSJed Brown   const PetscInt  Nc     = dim;
310c4762a1bSJed Brown   const PetscReal mu     = 1.0;
311c4762a1bSJed Brown   const PetscReal lambda = 1.0;
312c4762a1bSJed Brown   PetscInt        c, d;
313c4762a1bSJed Brown 
314c4762a1bSJed Brown   for (c = 0; c < Nc; ++c) {
315c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
316c4762a1bSJed Brown       g3[((c * Nc + c) * dim + d) * dim + d] += mu;
317c4762a1bSJed Brown       g3[((c * Nc + d) * dim + d) * dim + c] += mu;
318c4762a1bSJed Brown       g3[((c * Nc + d) * dim + c) * dim + d] += lambda;
319c4762a1bSJed Brown     }
320c4762a1bSJed Brown   }
321c4762a1bSJed Brown }
322c4762a1bSJed Brown 
323*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
324482dd828SMatthew G. Knepley   PetscInt sol = 0, def = 0;
325c4762a1bSJed Brown 
326c4762a1bSJed Brown   PetscFunctionBeginUser;
327482dd828SMatthew G. Knepley   options->deform           = DEFORM_NONE;
328c4762a1bSJed Brown   options->solType          = SOL_VLAP_QUADRATIC;
329c4762a1bSJed Brown   options->useNearNullspace = PETSC_TRUE;
3309566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(options->dmType, DMPLEX, 256));
331c4762a1bSJed Brown 
332d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Linear Elasticity Problem Options", "DMPLEX");
3339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-deform_type", "Type of domain deformation", "ex17.c", deformTypes, NUM_DEFORM_TYPES, deformTypes[options->deform], &def, NULL));
334482dd828SMatthew G. Knepley   options->deform = (DeformType)def;
3359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-sol_type", "Type of exact solution", "ex17.c", solutionTypes, NUM_SOLUTION_TYPES, solutionTypes[options->solType], &sol, NULL));
336c4762a1bSJed Brown   options->solType = (SolutionType)sol;
3379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-near_nullspace", "Use the rigid body modes as an AMG near nullspace", "ex17.c", options->useNearNullspace, &options->useNearNullspace, NULL));
3389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-dm_type", "Convert DMPlex to another format", "ex17.c", DMList, options->dmType, options->dmType, 256, NULL));
339d0609cedSBarry Smith   PetscOptionsEnd();
340c4762a1bSJed Brown   PetscFunctionReturn(0);
341c4762a1bSJed Brown }
342c4762a1bSJed Brown 
343*9371c9d4SSatish Balay static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *ctx) {
344d12cd81dSMatthew G. Knepley   PetscBag   bag;
345d12cd81dSMatthew G. Knepley   Parameter *p;
346d12cd81dSMatthew G. Knepley 
347d12cd81dSMatthew G. Knepley   PetscFunctionBeginUser;
348d12cd81dSMatthew G. Knepley   /* setup PETSc parameter bag */
3499566063dSJacob Faibussowitsch   PetscCall(PetscBagGetData(ctx->bag, (void **)&p));
3509566063dSJacob Faibussowitsch   PetscCall(PetscBagSetName(ctx->bag, "par", "Elastic Parameters"));
351d12cd81dSMatthew G. Knepley   bag = ctx->bag;
3529566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterScalar(bag, &p->mu, 1.0, "mu", "Shear Modulus, Pa"));
3539566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterScalar(bag, &p->lambda, 1.0, "lambda", "Lame's first parameter, Pa"));
3549566063dSJacob Faibussowitsch   PetscCall(PetscBagSetFromOptions(bag));
355d12cd81dSMatthew G. Knepley   {
356d12cd81dSMatthew G. Knepley     PetscViewer       viewer;
357d12cd81dSMatthew G. Knepley     PetscViewerFormat format;
358d12cd81dSMatthew G. Knepley     PetscBool         flg;
359d12cd81dSMatthew G. Knepley 
3609566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetViewer(comm, NULL, NULL, "-param_view", &viewer, &format, &flg));
361d12cd81dSMatthew G. Knepley     if (flg) {
3629566063dSJacob Faibussowitsch       PetscCall(PetscViewerPushFormat(viewer, format));
3639566063dSJacob Faibussowitsch       PetscCall(PetscBagView(bag, viewer));
3649566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
3659566063dSJacob Faibussowitsch       PetscCall(PetscViewerPopFormat(viewer));
3669566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&viewer));
367d12cd81dSMatthew G. Knepley     }
368d12cd81dSMatthew G. Knepley   }
369d12cd81dSMatthew G. Knepley   PetscFunctionReturn(0);
370d12cd81dSMatthew G. Knepley }
371d12cd81dSMatthew G. Knepley 
372*9371c9d4SSatish Balay static PetscErrorCode DMPlexDistortGeometry(DM dm) {
373482dd828SMatthew G. Knepley   DM           cdm;
374482dd828SMatthew G. Knepley   DMLabel      label;
375482dd828SMatthew G. Knepley   Vec          coordinates;
376482dd828SMatthew G. Knepley   PetscScalar *coords;
377482dd828SMatthew G. Knepley   PetscReal    mid = 0.5;
378482dd828SMatthew G. Knepley   PetscInt     cdim, d, vStart, vEnd, v;
379482dd828SMatthew G. Knepley 
380482dd828SMatthew G. Knepley   PetscFunctionBeginUser;
3819566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &cdm));
3829566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &cdim));
3839566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
3849566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
3859566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
3869566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(coordinates, &coords));
387482dd828SMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
388482dd828SMatthew G. Knepley     PetscScalar *pcoords, shift;
389482dd828SMatthew G. Knepley     PetscInt     val;
390482dd828SMatthew G. Knepley 
3919566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(label, v, &val));
392482dd828SMatthew G. Knepley     if (val >= 0) continue;
3939566063dSJacob Faibussowitsch     PetscCall(DMPlexPointLocalRef(cdm, v, coords, &pcoords));
394482dd828SMatthew G. Knepley     shift = 0.2 * PetscAbsScalar(pcoords[0] - mid);
395482dd828SMatthew G. Knepley     shift = PetscRealPart(pcoords[0]) > mid ? shift : -shift;
396482dd828SMatthew G. Knepley     for (d = 1; d < cdim; ++d) pcoords[d] += shift;
397482dd828SMatthew G. Knepley   }
3989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(coordinates, &coords));
399482dd828SMatthew G. Knepley   PetscFunctionReturn(0);
400482dd828SMatthew G. Knepley }
401482dd828SMatthew G. Knepley 
402*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) {
403c4762a1bSJed Brown   PetscFunctionBeginUser;
4049566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
4059566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
4069566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
407482dd828SMatthew G. Knepley   switch (user->deform) {
408482dd828SMatthew G. Knepley   case DEFORM_NONE: break;
4099566063dSJacob Faibussowitsch   case DEFORM_SHEAR: PetscCall(DMPlexShearGeometry(*dm, DM_X, NULL)); break;
4109566063dSJacob Faibussowitsch   case DEFORM_STEP: PetscCall(DMPlexDistortGeometry(*dm)); break;
41163a3b9bcSJacob Faibussowitsch   default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid deformation type: %s (%d)", deformTypes[PetscMin(user->deform, NUM_DEFORM_TYPES)], user->deform);
412482dd828SMatthew G. Knepley   }
4139566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(*dm, user));
4149566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
415c4762a1bSJed Brown   PetscFunctionReturn(0);
416c4762a1bSJed Brown }
417c4762a1bSJed Brown 
418*9371c9d4SSatish Balay static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) {
419c4762a1bSJed Brown   PetscErrorCode (*exact)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
420d12cd81dSMatthew G. Knepley   Parameter    *param;
421482dd828SMatthew G. Knepley   PetscDS       ds;
42245480ffeSMatthew G. Knepley   PetscWeakForm wf;
42345480ffeSMatthew G. Knepley   DMLabel       label;
42445480ffeSMatthew G. Knepley   PetscInt      id, bd;
425c4762a1bSJed Brown   PetscInt      dim;
426c4762a1bSJed Brown 
427c4762a1bSJed Brown   PetscFunctionBeginUser;
4289566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
4299566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
4309566063dSJacob Faibussowitsch   PetscCall(PetscDSGetSpatialDimension(ds, &dim));
4319566063dSJacob Faibussowitsch   PetscCall(PetscBagGetData(user->bag, (void **)&param));
432c4762a1bSJed Brown   switch (user->solType) {
433482dd828SMatthew G. Knepley   case SOL_MASS_QUADRATIC:
4349566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, f0_mass_u, NULL));
4359566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_mass_uu, NULL, NULL, NULL));
4369566063dSJacob Faibussowitsch     PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, 0, 0, 1, f0_mass_quadratic_u, 0, NULL));
437c4762a1bSJed Brown     switch (dim) {
438c4762a1bSJed Brown     case 2: exact = quadratic_2d_u; break;
439c4762a1bSJed Brown     case 3: exact = quadratic_3d_u; break;
44063a3b9bcSJacob Faibussowitsch     default: SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim);
441482dd828SMatthew G. Knepley     }
442482dd828SMatthew G. Knepley     break;
443482dd828SMatthew G. Knepley   case SOL_VLAP_QUADRATIC:
4449566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, f0_vlap_quadratic_u, f1_vlap_u));
4459566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_vlap_uu));
446482dd828SMatthew G. Knepley     switch (dim) {
447482dd828SMatthew G. Knepley     case 2: exact = quadratic_2d_u; break;
448482dd828SMatthew G. Knepley     case 3: exact = quadratic_3d_u; break;
44963a3b9bcSJacob Faibussowitsch     default: SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim);
450c4762a1bSJed Brown     }
451c4762a1bSJed Brown     break;
452c4762a1bSJed Brown   case SOL_ELAS_QUADRATIC:
4539566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, f0_elas_quadratic_u, f1_elas_u));
4549566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu));
455c4762a1bSJed Brown     switch (dim) {
456c4762a1bSJed Brown     case 2: exact = quadratic_2d_u; break;
457c4762a1bSJed Brown     case 3: exact = quadratic_3d_u; break;
45863a3b9bcSJacob Faibussowitsch     default: SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim);
459c4762a1bSJed Brown     }
460c4762a1bSJed Brown     break;
461c4762a1bSJed Brown   case SOL_VLAP_TRIG:
4629566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, f0_vlap_trig_u, f1_vlap_u));
4639566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_vlap_uu));
464c4762a1bSJed Brown     switch (dim) {
465c4762a1bSJed Brown     case 2: exact = trig_2d_u; break;
466c4762a1bSJed Brown     case 3: exact = trig_3d_u; break;
46763a3b9bcSJacob Faibussowitsch     default: SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim);
468c4762a1bSJed Brown     }
469c4762a1bSJed Brown     break;
470c4762a1bSJed Brown   case SOL_ELAS_TRIG:
4719566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, f0_elas_trig_u, f1_elas_u));
4729566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu));
473c4762a1bSJed Brown     switch (dim) {
474c4762a1bSJed Brown     case 2: exact = trig_2d_u; break;
475c4762a1bSJed Brown     case 3: exact = trig_3d_u; break;
47663a3b9bcSJacob Faibussowitsch     default: SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim);
477c4762a1bSJed Brown     }
478c4762a1bSJed Brown     break;
479c4762a1bSJed Brown   case SOL_ELAS_AXIAL_DISP:
4809566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, NULL, f1_elas_u));
4819566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu));
48245480ffeSMatthew G. Knepley     id = dim == 3 ? 5 : 2;
4839566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dm, "marker", &label));
4849566063dSJacob Faibussowitsch     PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "right", label, 1, &id, 0, 0, NULL, (void (*)(void))NULL, NULL, user, &bd));
4859566063dSJacob Faibussowitsch     PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
4869566063dSJacob Faibussowitsch     PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_elas_axial_disp_bd_u, 0, NULL));
487c4762a1bSJed Brown     exact = axial_disp_u;
488c4762a1bSJed Brown     break;
489c4762a1bSJed Brown   case SOL_ELAS_UNIFORM_STRAIN:
4909566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, NULL, f1_elas_u));
4919566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu));
492c4762a1bSJed Brown     exact = uniform_strain_u;
493c4762a1bSJed Brown     break;
494d12cd81dSMatthew G. Knepley   case SOL_ELAS_GE:
4959566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, NULL, f1_elas_u));
4969566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu));
497d12cd81dSMatthew G. Knepley     exact = zero; /* No exact solution available */
498d12cd81dSMatthew G. Knepley     break;
49963a3b9bcSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid solution type: %s (%d)", solutionTypes[PetscMin(user->solType, NUM_SOLUTION_TYPES)], user->solType);
500c4762a1bSJed Brown   }
5019566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(ds, 0, exact, user));
5029566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
503c4762a1bSJed Brown   if (user->solType == SOL_ELAS_AXIAL_DISP) {
504c4762a1bSJed Brown     PetscInt cmp;
505c4762a1bSJed Brown 
506c4762a1bSJed Brown     id  = dim == 3 ? 6 : 4;
507c4762a1bSJed Brown     cmp = 0;
5089566063dSJacob Faibussowitsch     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "left", label, 1, &id, 0, 1, &cmp, (void (*)(void))zero, NULL, user, NULL));
509c4762a1bSJed Brown     cmp = dim == 3 ? 2 : 1;
510c4762a1bSJed Brown     id  = dim == 3 ? 1 : 1;
5119566063dSJacob Faibussowitsch     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "bottom", label, 1, &id, 0, 1, &cmp, (void (*)(void))zero, NULL, user, NULL));
512c4762a1bSJed Brown     if (dim == 3) {
513c4762a1bSJed Brown       cmp = 1;
514c4762a1bSJed Brown       id  = 3;
5159566063dSJacob Faibussowitsch       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "front", label, 1, &id, 0, 1, &cmp, (void (*)(void))zero, NULL, user, NULL));
516c4762a1bSJed Brown     }
517d12cd81dSMatthew G. Knepley   } else if (user->solType == SOL_ELAS_GE) {
518d12cd81dSMatthew G. Knepley     PetscInt cmp;
519d12cd81dSMatthew G. Knepley 
520d12cd81dSMatthew G. Knepley     id = dim == 3 ? 6 : 4;
5219566063dSJacob Faibussowitsch     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "left", label, 1, &id, 0, 0, NULL, (void (*)(void))zero, NULL, user, NULL));
522d12cd81dSMatthew G. Knepley     id  = dim == 3 ? 5 : 2;
523d12cd81dSMatthew G. Knepley     cmp = 0;
5249566063dSJacob Faibussowitsch     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "right", label, 1, &id, 0, 1, &cmp, (void (*)(void))ge_shift, NULL, user, NULL));
525c4762a1bSJed Brown   } else {
526c4762a1bSJed Brown     id = 1;
5279566063dSJacob Faibussowitsch     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))exact, NULL, user, NULL));
528c4762a1bSJed Brown   }
529d12cd81dSMatthew G. Knepley   /* Setup constants */
530d12cd81dSMatthew G. Knepley   {
531d12cd81dSMatthew G. Knepley     PetscScalar constants[2];
532d12cd81dSMatthew G. Knepley 
533d12cd81dSMatthew G. Knepley     constants[0] = param->mu;     /* shear modulus, Pa */
534d12cd81dSMatthew G. Knepley     constants[1] = param->lambda; /* Lame's first parameter, Pa */
5359566063dSJacob Faibussowitsch     PetscCall(PetscDSSetConstants(ds, 2, constants));
536d12cd81dSMatthew G. Knepley   }
537c4762a1bSJed Brown   PetscFunctionReturn(0);
538c4762a1bSJed Brown }
539c4762a1bSJed Brown 
540*9371c9d4SSatish Balay static PetscErrorCode CreateElasticityNullSpace(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullspace) {
541c4762a1bSJed Brown   PetscFunctionBegin;
5429566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateRigidBody(dm, origField, nullspace));
543c4762a1bSJed Brown   PetscFunctionReturn(0);
544c4762a1bSJed Brown }
545c4762a1bSJed Brown 
546*9371c9d4SSatish Balay PetscErrorCode SetupFE(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), void *ctx) {
547c4762a1bSJed Brown   AppCtx        *user = (AppCtx *)ctx;
548c4762a1bSJed Brown   DM             cdm  = dm;
549c4762a1bSJed Brown   PetscFE        fe;
550c4762a1bSJed Brown   char           prefix[PETSC_MAX_PATH_LEN];
55130602db0SMatthew G. Knepley   DMPolytopeType ct;
55230602db0SMatthew G. Knepley   PetscBool      simplex;
55330602db0SMatthew G. Knepley   PetscInt       dim, cStart;
554c4762a1bSJed Brown 
555c4762a1bSJed Brown   PetscFunctionBegin;
556c4762a1bSJed Brown   /* Create finite element */
5579566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
5589566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
5599566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
56030602db0SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
5619566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
5629566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, dim, simplex, name ? prefix : NULL, -1, &fe));
5639566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe, name));
564c4762a1bSJed Brown   /* Set discretization and boundary conditions for each mesh */
5659566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
5669566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
5679566063dSJacob Faibussowitsch   PetscCall((*setup)(dm, user));
568c4762a1bSJed Brown   while (cdm) {
5699566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm, cdm));
5709566063dSJacob Faibussowitsch     if (user->useNearNullspace) PetscCall(DMSetNearNullSpaceConstructor(cdm, 0, CreateElasticityNullSpace));
571c4762a1bSJed Brown     /* TODO: Check whether the boundary of coarse meshes is marked */
5729566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm, &cdm));
573c4762a1bSJed Brown   }
5749566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
575c4762a1bSJed Brown   PetscFunctionReturn(0);
576c4762a1bSJed Brown }
577c4762a1bSJed Brown 
578*9371c9d4SSatish Balay int main(int argc, char **argv) {
579c4762a1bSJed Brown   DM     dm;   /* Problem specification */
580c4762a1bSJed Brown   SNES   snes; /* Nonlinear solver */
581c4762a1bSJed Brown   Vec    u;    /* Solutions */
582c4762a1bSJed Brown   AppCtx user; /* User-defined work context */
583c4762a1bSJed Brown 
584327415f7SBarry Smith   PetscFunctionBeginUser;
5859566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
5869566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
5879566063dSJacob Faibussowitsch   PetscCall(PetscBagCreate(PETSC_COMM_SELF, sizeof(Parameter), &user.bag));
5889566063dSJacob Faibussowitsch   PetscCall(SetupParameters(PETSC_COMM_WORLD, &user));
589c4762a1bSJed Brown   /* Primal system */
5909566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
5919566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
5929566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, dm));
5939566063dSJacob Faibussowitsch   PetscCall(SetupFE(dm, "displacement", SetupPrimalProblem, &user));
5949566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &u));
5959566063dSJacob Faibussowitsch   PetscCall(VecSet(u, 0.0));
5969566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)u, "displacement"));
5979566063dSJacob Faibussowitsch   PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user));
5989566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
5999566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckFromOptions(snes, u));
6009566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, u));
6019566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes, &u));
6029566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-displacement_view"));
603c4762a1bSJed Brown   /* Cleanup */
6049566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
6059566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
6069566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
6079566063dSJacob Faibussowitsch   PetscCall(PetscBagDestroy(&user.bag));
6089566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
609b122ec5aSJacob Faibussowitsch   return 0;
610c4762a1bSJed Brown }
611c4762a1bSJed Brown 
612c4762a1bSJed Brown /*TEST
613c4762a1bSJed Brown 
61430602db0SMatthew G. Knepley   testset:
61530602db0SMatthew G. Knepley     args: -dm_plex_box_faces 1,1,1
61630602db0SMatthew G. Knepley 
617c4762a1bSJed Brown     test:
618c4762a1bSJed Brown       suffix: 2d_p1_quad_vlap
619c4762a1bSJed Brown       requires: triangle
620c4762a1bSJed Brown       args: -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate
621c4762a1bSJed Brown     test:
622c4762a1bSJed Brown       suffix: 2d_p2_quad_vlap
623c4762a1bSJed Brown       requires: triangle
624c4762a1bSJed Brown       args: -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001
625c4762a1bSJed Brown     test:
626c4762a1bSJed Brown       suffix: 2d_p3_quad_vlap
627c4762a1bSJed Brown       requires: triangle
628c4762a1bSJed Brown       args: -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001
629c4762a1bSJed Brown     test:
630c4762a1bSJed Brown       suffix: 2d_q1_quad_vlap
63130602db0SMatthew G. Knepley       args: -dm_plex_simplex 0 -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate
632c4762a1bSJed Brown     test:
633c4762a1bSJed Brown       suffix: 2d_q2_quad_vlap
63430602db0SMatthew G. Knepley       args: -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001
635c4762a1bSJed Brown     test:
636c4762a1bSJed Brown       suffix: 2d_q3_quad_vlap
637c4762a1bSJed Brown       requires: !single
63830602db0SMatthew G. Knepley       args: -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001
639c4762a1bSJed Brown     test:
640c4762a1bSJed Brown       suffix: 2d_p1_quad_elas
641c4762a1bSJed Brown       requires: triangle
642c4762a1bSJed Brown       args: -sol_type elas_quad -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate
643c4762a1bSJed Brown     test:
644c4762a1bSJed Brown       suffix: 2d_p2_quad_elas
645c4762a1bSJed Brown       requires: triangle
646c4762a1bSJed Brown       args: -sol_type elas_quad -displacement_petscspace_degree 2 -dmsnes_check .0001
647c4762a1bSJed Brown     test:
648c4762a1bSJed Brown       suffix: 2d_p3_quad_elas
649c4762a1bSJed Brown       requires: triangle
650c4762a1bSJed Brown       args: -sol_type elas_quad -displacement_petscspace_degree 3 -dmsnes_check .0001
651c4762a1bSJed Brown     test:
652c4762a1bSJed Brown       suffix: 2d_q1_quad_elas
65330602db0SMatthew G. Knepley       args: -sol_type elas_quad -dm_plex_simplex 0 -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
654c4762a1bSJed Brown     test:
655c4762a1bSJed Brown       suffix: 2d_q1_quad_elas_shear
656482dd828SMatthew G. Knepley       args: -sol_type elas_quad -dm_plex_simplex 0 -deform_type shear -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
657c4762a1bSJed Brown     test:
658c4762a1bSJed Brown       suffix: 2d_q2_quad_elas
65930602db0SMatthew G. Knepley       args: -sol_type elas_quad -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dmsnes_check .0001
660c4762a1bSJed Brown     test:
661c4762a1bSJed Brown       suffix: 2d_q2_quad_elas_shear
662482dd828SMatthew G. Knepley       args: -sol_type elas_quad -dm_plex_simplex 0 -deform_type shear -displacement_petscspace_degree 2 -dmsnes_check
663c4762a1bSJed Brown     test:
664c4762a1bSJed Brown       suffix: 2d_q3_quad_elas
66530602db0SMatthew G. Knepley       args: -sol_type elas_quad -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dmsnes_check .0001
666c4762a1bSJed Brown     test:
667c4762a1bSJed Brown       suffix: 2d_q3_quad_elas_shear
668c4762a1bSJed Brown       requires: !single
669482dd828SMatthew G. Knepley       args: -sol_type elas_quad -dm_plex_simplex 0 -deform_type shear -displacement_petscspace_degree 3 -dmsnes_check
670c4762a1bSJed Brown 
671c4762a1bSJed Brown     test:
672c4762a1bSJed Brown       suffix: 3d_p1_quad_vlap
673c4762a1bSJed Brown       requires: ctetgen
67430602db0SMatthew G. Knepley       args: -dm_plex_dim 3 -dm_refine 1 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate
675c4762a1bSJed Brown     test:
676c4762a1bSJed Brown       suffix: 3d_p2_quad_vlap
677c4762a1bSJed Brown       requires: ctetgen
67830602db0SMatthew G. Knepley       args: -dm_plex_dim 3 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001
679c4762a1bSJed Brown     test:
680c4762a1bSJed Brown       suffix: 3d_p3_quad_vlap
681c4762a1bSJed Brown       requires: ctetgen
68230602db0SMatthew G. Knepley       args: -dm_plex_dim 3 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001
683c4762a1bSJed Brown     test:
684c4762a1bSJed Brown       suffix: 3d_q1_quad_vlap
68530602db0SMatthew G. Knepley       args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_plex_simplex 0 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate
686c4762a1bSJed Brown     test:
687c4762a1bSJed Brown       suffix: 3d_q2_quad_vlap
68830602db0SMatthew G. Knepley       args: -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001
689c4762a1bSJed Brown     test:
690c4762a1bSJed Brown       suffix: 3d_q3_quad_vlap
69130602db0SMatthew G. Knepley       args: -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001
692c4762a1bSJed Brown     test:
693c4762a1bSJed Brown       suffix: 3d_p1_quad_elas
694c4762a1bSJed Brown       requires: ctetgen
69530602db0SMatthew G. Knepley       args: -sol_type elas_quad -dm_plex_dim 3 -dm_refine 1 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate
696c4762a1bSJed Brown     test:
697c4762a1bSJed Brown       suffix: 3d_p2_quad_elas
698c4762a1bSJed Brown       requires: ctetgen
69930602db0SMatthew G. Knepley       args: -sol_type elas_quad -dm_plex_dim 3 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001
700c4762a1bSJed Brown     test:
701c4762a1bSJed Brown       suffix: 3d_p3_quad_elas
702c4762a1bSJed Brown       requires: ctetgen
70330602db0SMatthew G. Knepley       args: -sol_type elas_quad -dm_plex_dim 3 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001
704c4762a1bSJed Brown     test:
705c4762a1bSJed Brown       suffix: 3d_q1_quad_elas
70630602db0SMatthew G. Knepley       args: -sol_type elas_quad -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_plex_simplex 0 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate
707c4762a1bSJed Brown     test:
708c4762a1bSJed Brown       suffix: 3d_q2_quad_elas
70930602db0SMatthew G. Knepley       args: -sol_type elas_quad -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001
710c4762a1bSJed Brown     test:
711c4762a1bSJed Brown       suffix: 3d_q3_quad_elas
712c4762a1bSJed Brown       requires: !single
71330602db0SMatthew G. Knepley       args: -sol_type elas_quad -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001
714c4762a1bSJed Brown 
715c4762a1bSJed Brown     test:
716c4762a1bSJed Brown       suffix: 2d_p1_trig_vlap
717c4762a1bSJed Brown       requires: triangle
718c4762a1bSJed Brown       args: -sol_type vlap_trig -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
719c4762a1bSJed Brown     test:
720c4762a1bSJed Brown       suffix: 2d_p2_trig_vlap
721c4762a1bSJed Brown       requires: triangle
722c4762a1bSJed Brown       args: -sol_type vlap_trig -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
723c4762a1bSJed Brown     test:
724c4762a1bSJed Brown       suffix: 2d_p3_trig_vlap
725c4762a1bSJed Brown       requires: triangle
726c4762a1bSJed Brown       args: -sol_type vlap_trig -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
727c4762a1bSJed Brown     test:
728c4762a1bSJed Brown       suffix: 2d_q1_trig_vlap
72930602db0SMatthew G. Knepley       args: -sol_type vlap_trig -dm_plex_simplex 0 -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
730c4762a1bSJed Brown     test:
731c4762a1bSJed Brown       suffix: 2d_q2_trig_vlap
73230602db0SMatthew G. Knepley       args: -sol_type vlap_trig -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
733c4762a1bSJed Brown     test:
734c4762a1bSJed Brown       suffix: 2d_q3_trig_vlap
73530602db0SMatthew G. Knepley       args: -sol_type vlap_trig -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
736c4762a1bSJed Brown     test:
737c4762a1bSJed Brown       suffix: 2d_p1_trig_elas
738c4762a1bSJed Brown       requires: triangle
739c4762a1bSJed Brown       args: -sol_type elas_trig -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
740c4762a1bSJed Brown     test:
741c4762a1bSJed Brown       suffix: 2d_p2_trig_elas
742c4762a1bSJed Brown       requires: triangle
743c4762a1bSJed Brown       args: -sol_type elas_trig -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
744c4762a1bSJed Brown     test:
745c4762a1bSJed Brown       suffix: 2d_p3_trig_elas
746c4762a1bSJed Brown       requires: triangle
747c4762a1bSJed Brown       args: -sol_type elas_trig -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
748c4762a1bSJed Brown     test:
749c4762a1bSJed Brown       suffix: 2d_q1_trig_elas
75030602db0SMatthew G. Knepley       args: -sol_type elas_trig -dm_plex_simplex 0 -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
751c4762a1bSJed Brown     test:
752c4762a1bSJed Brown       suffix: 2d_q1_trig_elas_shear
753482dd828SMatthew G. Knepley       args: -sol_type elas_trig -dm_plex_simplex 0 -deform_type shear -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
754c4762a1bSJed Brown     test:
755c4762a1bSJed Brown       suffix: 2d_q2_trig_elas
75630602db0SMatthew G. Knepley       args: -sol_type elas_trig -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
757c4762a1bSJed Brown     test:
758c4762a1bSJed Brown       suffix: 2d_q2_trig_elas_shear
759482dd828SMatthew G. Knepley       args: -sol_type elas_trig -dm_plex_simplex 0 -deform_type shear -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
760c4762a1bSJed Brown     test:
761c4762a1bSJed Brown       suffix: 2d_q3_trig_elas
76230602db0SMatthew G. Knepley       args: -sol_type elas_trig -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
763c4762a1bSJed Brown     test:
764c4762a1bSJed Brown       suffix: 2d_q3_trig_elas_shear
765482dd828SMatthew G. Knepley       args: -sol_type elas_trig -dm_plex_simplex 0 -deform_type shear -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
766c4762a1bSJed Brown 
767c4762a1bSJed Brown     test:
768c4762a1bSJed Brown       suffix: 3d_p1_trig_vlap
769c4762a1bSJed Brown       requires: ctetgen
77030602db0SMatthew G. Knepley       args: -sol_type vlap_trig -dm_plex_dim 3 -dm_refine 1 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate
771c4762a1bSJed Brown     test:
772c4762a1bSJed Brown       suffix: 3d_p2_trig_vlap
773c4762a1bSJed Brown       requires: ctetgen
77430602db0SMatthew G. Knepley       args: -sol_type vlap_trig -dm_plex_dim 3 -displacement_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate
775c4762a1bSJed Brown     test:
776c4762a1bSJed Brown       suffix: 3d_p3_trig_vlap
777c4762a1bSJed Brown       requires: ctetgen
77830602db0SMatthew G. Knepley       args: -sol_type vlap_trig -dm_plex_dim 3 -displacement_petscspace_degree 3 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate
779c4762a1bSJed Brown     test:
780c4762a1bSJed Brown       suffix: 3d_q1_trig_vlap
78130602db0SMatthew G. Knepley       args: -sol_type vlap_trig -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_plex_simplex 0 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate
782c4762a1bSJed Brown     test:
783c4762a1bSJed Brown       suffix: 3d_q2_trig_vlap
78430602db0SMatthew G. Knepley       args: -sol_type vlap_trig -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate
785c4762a1bSJed Brown     test:
786c4762a1bSJed Brown       suffix: 3d_q3_trig_vlap
787c4762a1bSJed Brown       requires: !__float128
78830602db0SMatthew G. Knepley       args: -sol_type vlap_trig -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate
789c4762a1bSJed Brown     test:
790c4762a1bSJed Brown       suffix: 3d_p1_trig_elas
791c4762a1bSJed Brown       requires: ctetgen
79230602db0SMatthew G. Knepley       args: -sol_type elas_trig -dm_plex_dim 3 -dm_refine 1 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate
793c4762a1bSJed Brown     test:
794c4762a1bSJed Brown       suffix: 3d_p2_trig_elas
795c4762a1bSJed Brown       requires: ctetgen
79630602db0SMatthew G. Knepley       args: -sol_type elas_trig -dm_plex_dim 3 -displacement_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate
797c4762a1bSJed Brown     test:
798c4762a1bSJed Brown       suffix: 3d_p3_trig_elas
799c4762a1bSJed Brown       requires: ctetgen
80030602db0SMatthew G. Knepley       args: -sol_type elas_trig -dm_plex_dim 3 -displacement_petscspace_degree 3 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate
801c4762a1bSJed Brown     test:
802c4762a1bSJed Brown       suffix: 3d_q1_trig_elas
80330602db0SMatthew G. Knepley       args: -sol_type elas_trig -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_plex_simplex 0 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate
804c4762a1bSJed Brown     test:
805c4762a1bSJed Brown       suffix: 3d_q2_trig_elas
80630602db0SMatthew G. Knepley       args: -sol_type elas_trig -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate
807c4762a1bSJed Brown     test:
808c4762a1bSJed Brown       suffix: 3d_q3_trig_elas
809c4762a1bSJed Brown       requires: !__float128
81030602db0SMatthew G. Knepley       args: -sol_type elas_trig -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate
811c4762a1bSJed Brown 
812c4762a1bSJed Brown     test:
813c4762a1bSJed Brown       suffix: 2d_p1_axial_elas
814c4762a1bSJed Brown       requires: triangle
815c4762a1bSJed Brown       args: -sol_type elas_axial_disp -displacement_petscspace_degree 1 -dm_plex_separate_marker -dm_refine 2 -dmsnes_check .0001 -pc_type lu
816c4762a1bSJed Brown     test:
817c4762a1bSJed Brown       suffix: 2d_p2_axial_elas
818c4762a1bSJed Brown       requires: triangle
819c4762a1bSJed Brown       args: -sol_type elas_axial_disp -displacement_petscspace_degree 2 -dm_plex_separate_marker -dmsnes_check .0001 -pc_type lu
820c4762a1bSJed Brown     test:
821c4762a1bSJed Brown       suffix: 2d_p3_axial_elas
822c4762a1bSJed Brown       requires: triangle
823c4762a1bSJed Brown       args: -sol_type elas_axial_disp -displacement_petscspace_degree 3 -dm_plex_separate_marker -dmsnes_check .0001 -pc_type lu
824c4762a1bSJed Brown     test:
825c4762a1bSJed Brown       suffix: 2d_q1_axial_elas
82630602db0SMatthew G. Knepley       args: -sol_type elas_axial_disp -dm_plex_simplex 0 -displacement_petscspace_degree 1 -dm_plex_separate_marker -dm_refine 1 -dmsnes_check   .0001 -pc_type lu
827c4762a1bSJed Brown     test:
828c4762a1bSJed Brown       suffix: 2d_q2_axial_elas
82930602db0SMatthew G. Knepley       args: -sol_type elas_axial_disp -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_plex_separate_marker -dmsnes_check .0001 -pc_type   lu
830c4762a1bSJed Brown     test:
831c4762a1bSJed Brown       suffix: 2d_q3_axial_elas
83230602db0SMatthew G. Knepley       args: -sol_type elas_axial_disp -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_plex_separate_marker -dmsnes_check .0001 -pc_type   lu
833c4762a1bSJed Brown 
834c4762a1bSJed Brown     test:
835c4762a1bSJed Brown       suffix: 2d_p1_uniform_elas
836c4762a1bSJed Brown       requires: triangle
837c4762a1bSJed Brown       args: -sol_type elas_uniform_strain -displacement_petscspace_degree 1 -dm_refine 2 -dmsnes_check .0001 -pc_type lu
838c4762a1bSJed Brown     test:
839c4762a1bSJed Brown       suffix: 2d_p2_uniform_elas
840c4762a1bSJed Brown       requires: triangle
841c4762a1bSJed Brown       args: -sol_type elas_uniform_strain -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001 -pc_type lu
842c4762a1bSJed Brown     test:
843c4762a1bSJed Brown       suffix: 2d_p3_uniform_elas
844c4762a1bSJed Brown       requires: triangle
845c4762a1bSJed Brown       args: -sol_type elas_uniform_strain -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001 -pc_type lu
846c4762a1bSJed Brown     test:
847c4762a1bSJed Brown       suffix: 2d_q1_uniform_elas
84830602db0SMatthew G. Knepley       args: -sol_type elas_uniform_strain -dm_plex_simplex 0 -displacement_petscspace_degree 1 -dm_refine 2 -dmsnes_check .0001 -pc_type lu
849c4762a1bSJed Brown     test:
850c4762a1bSJed Brown       suffix: 2d_q2_uniform_elas
851c4762a1bSJed Brown       requires: !single
85230602db0SMatthew G. Knepley       args: -sol_type elas_uniform_strain -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001 -pc_type lu
853c4762a1bSJed Brown     test:
854c4762a1bSJed Brown       suffix: 2d_q3_uniform_elas
855c4762a1bSJed Brown       requires: !single
85630602db0SMatthew G. Knepley       args: -sol_type elas_uniform_strain -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001 -pc_type lu
857482dd828SMatthew G. Knepley     test:
858482dd828SMatthew G. Knepley       suffix: 2d_p1_uniform_elas_step
859482dd828SMatthew G. Knepley       requires: triangle
860482dd828SMatthew G. Knepley       args: -sol_type elas_uniform_strain -deform_type step -displacement_petscspace_degree 1 -dm_refine 2 -dmsnes_check .0001 -pc_type lu
861482dd828SMatthew G. Knepley 
862482dd828SMatthew G. Knepley   testset:
863482dd828SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 3,3 -deform_type step -displacement_petscspace_degree 1 -dmsnes_check .0001 -pc_type lu
864482dd828SMatthew G. Knepley 
865482dd828SMatthew G. Knepley     test:
866482dd828SMatthew G. Knepley       suffix: 2d_q1_uniform_elas_step
867482dd828SMatthew G. Knepley       args: -sol_type elas_uniform_strain -dm_refine 2
868482dd828SMatthew G. Knepley     test:
869482dd828SMatthew G. Knepley       suffix: 2d_q1_quad_vlap_step
870482dd828SMatthew G. Knepley       args:
871482dd828SMatthew G. Knepley     test:
872482dd828SMatthew G. Knepley       suffix: 2d_q2_quad_vlap_step
873482dd828SMatthew G. Knepley       args: -displacement_petscspace_degree 2
874482dd828SMatthew G. Knepley     test:
875482dd828SMatthew G. Knepley       suffix: 2d_q1_quad_mass_step
876482dd828SMatthew G. Knepley       args: -sol_type mass_quad
877c4762a1bSJed Brown 
878d12cd81dSMatthew G. Knepley   testset:
879908b9b43SStefano Zampini     filter: grep -v "variant HERMITIAN"
880d12cd81dSMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower -5,-5,-0.25 -dm_plex_box_upper 5,5,0.25 \
881d12cd81dSMatthew G. Knepley           -dm_plex_box_faces 5,5,2 -dm_plex_separate_marker -dm_refine 0 -petscpartitioner_type simple \
882908b9b43SStefano Zampini           -sol_type elas_ge
883908b9b43SStefano Zampini 
884908b9b43SStefano Zampini     test:
885908b9b43SStefano Zampini       suffix: ge_q1_0
886908b9b43SStefano Zampini       args: -displacement_petscspace_degree 1 \
887d12cd81dSMatthew G. Knepley             -snes_max_it 2 -snes_rtol 1.e-10 \
888d12cd81dSMatthew G. Knepley             -ksp_type cg -ksp_rtol 1.e-10 -ksp_max_it 100 -ksp_norm_type unpreconditioned \
889d12cd81dSMatthew G. Knepley             -pc_type gamg -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 \
890d12cd81dSMatthew G. Knepley               -pc_gamg_coarse_eq_limit 10 -pc_gamg_reuse_interpolation true \
891bae903cbSmarkadams4               -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 \
892d12cd81dSMatthew G. Knepley               -mg_levels_ksp_max_it 2 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.1 -mg_levels_pc_type jacobi \
893d12cd81dSMatthew G. Knepley               -matptap_via scalable
894d12cd81dSMatthew G. Knepley     test:
895908b9b43SStefano Zampini       nsize: 5
896908b9b43SStefano Zampini       suffix: ge_q1_gdsw
897908b9b43SStefano Zampini       args: -snes_max_it 1 -ksp_type cg -ksp_norm_type natural -displacement_petscspace_degree 1 -snes_monitor_short -ksp_monitor_short -pc_type mg -pc_mg_adapt_interp_coarse_space gdsw -pc_mg_levels 2 -pc_mg_galerkin -mg_levels_pc_type bjacobi -mg_levels_esteig_ksp_type cg -mg_levels_sub_pc_type icc -mg_coarse_redundant_pc_type cholesky -ksp_view
898d12cd81dSMatthew G. Knepley 
899c4762a1bSJed Brown TEST*/
900