xref: /petsc/src/snes/tutorials/ex17.c (revision 482dd828751c19448ab8eaa6131c34bae8c6237d)
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
9c4762a1bSJed Brown */
10c4762a1bSJed Brown 
11c4762a1bSJed Brown #include <petscdmplex.h>
12c4762a1bSJed Brown #include <petscsnes.h>
13c4762a1bSJed Brown #include <petscds.h>
14c4762a1bSJed Brown #include <petscconvest.h>
15c4762a1bSJed Brown 
16*482dd828SMatthew G. Knepley typedef enum {SOL_VLAP_QUADRATIC, SOL_ELAS_QUADRATIC, SOL_VLAP_TRIG, SOL_ELAS_TRIG, SOL_ELAS_AXIAL_DISP, SOL_ELAS_UNIFORM_STRAIN, SOL_MASS_QUADRATIC, NUM_SOLUTION_TYPES} SolutionType;
17*482dd828SMatthew G. Knepley const char *solutionTypes[NUM_SOLUTION_TYPES+1] = {"vlap_quad", "elas_quad", "vlap_trig", "elas_trig", "elas_axial_disp", "elas_uniform_strain", "mass_quad", "unknown"};
18*482dd828SMatthew G. Knepley 
19*482dd828SMatthew G. Knepley typedef enum {DEFORM_NONE, DEFORM_SHEAR, DEFORM_STEP, NUM_DEFORM_TYPES} DeformType;
20*482dd828SMatthew G. Knepley const char *deformTypes[NUM_DEFORM_TYPES+1] = {"none", "shear", "step", "unknown"};
21c4762a1bSJed Brown 
22c4762a1bSJed Brown typedef struct {
23c4762a1bSJed Brown   /* Domain and mesh definition */
24c4762a1bSJed Brown   char         dmType[256]; /* DM type for the solve */
25*482dd828SMatthew G. Knepley   DeformType   deform;      /* Domain deformation type */
26c4762a1bSJed Brown   /* Problem definition */
27c4762a1bSJed Brown   SolutionType solType;     /* Type of exact solution */
28c4762a1bSJed Brown   /* Solver definition */
29c4762a1bSJed Brown   PetscBool    useNearNullspace; /* Use the rigid body modes as a near nullspace for AMG */
30c4762a1bSJed Brown } AppCtx;
31c4762a1bSJed Brown 
32c4762a1bSJed Brown static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
33c4762a1bSJed Brown {
34c4762a1bSJed Brown   PetscInt d;
35c4762a1bSJed Brown   for (d = 0; d < dim; ++d) u[d] = 0.0;
36c4762a1bSJed Brown   return 0;
37c4762a1bSJed Brown }
38c4762a1bSJed Brown 
39c4762a1bSJed Brown static PetscErrorCode quadratic_2d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
40c4762a1bSJed Brown {
41c4762a1bSJed Brown   u[0] = x[0]*x[0];
42c4762a1bSJed Brown   u[1] = x[1]*x[1] - 2.0*x[0]*x[1];
43c4762a1bSJed Brown   return 0;
44c4762a1bSJed Brown }
45c4762a1bSJed Brown 
46c4762a1bSJed Brown static PetscErrorCode quadratic_3d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
47c4762a1bSJed Brown {
48c4762a1bSJed Brown   u[0] = x[0]*x[0];
49c4762a1bSJed Brown   u[1] = x[1]*x[1] - 2.0*x[0]*x[1];
50c4762a1bSJed Brown   u[2] = x[2]*x[2] - 2.0*x[1]*x[2];
51c4762a1bSJed Brown   return 0;
52c4762a1bSJed Brown }
53c4762a1bSJed Brown 
54c4762a1bSJed Brown /*
55c4762a1bSJed Brown   u = x^2
56c4762a1bSJed Brown   v = y^2 - 2xy
57c4762a1bSJed Brown   Delta <u,v> - f = <2, 2> - <2, 2>
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   u = x^2
60c4762a1bSJed Brown   v = y^2 - 2xy
61c4762a1bSJed Brown   w = z^2 - 2yz
62c4762a1bSJed Brown   Delta <u,v,w> - f = <2, 2, 2> - <2, 2, 2>
63c4762a1bSJed Brown */
64c4762a1bSJed Brown static void f0_vlap_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
65c4762a1bSJed Brown                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
66c4762a1bSJed Brown                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
67c4762a1bSJed Brown                                 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
68c4762a1bSJed Brown {
69c4762a1bSJed Brown   PetscInt d;
70c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f0[d] += 2.0;
71c4762a1bSJed Brown }
72c4762a1bSJed Brown 
73c4762a1bSJed Brown /*
74c4762a1bSJed Brown   u = x^2
75c4762a1bSJed Brown   v = y^2 - 2xy
76c4762a1bSJed Brown   \varepsilon = / 2x     -y    \
77c4762a1bSJed Brown                 \ -y   2y - 2x /
78c4762a1bSJed Brown   Tr(\varepsilon) = div u = 2y
79c4762a1bSJed Brown   div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij}
80c4762a1bSJed Brown     = \lambda \partial_j (2y) + 2\mu < 2-1, 2 >
81c4762a1bSJed Brown     = \lambda < 0, 2 > + \mu < 2, 4 >
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   u = x^2
84c4762a1bSJed Brown   v = y^2 - 2xy
85c4762a1bSJed Brown   w = z^2 - 2yz
86c4762a1bSJed Brown   \varepsilon = / 2x     -y       0   \
87c4762a1bSJed Brown                 | -y   2y - 2x   -z   |
88c4762a1bSJed Brown                 \  0     -z    2z - 2y/
89c4762a1bSJed Brown   Tr(\varepsilon) = div u = 2z
90c4762a1bSJed Brown   div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij}
91c4762a1bSJed Brown     = \lambda \partial_j (2z) + 2\mu < 2-1, 2-1, 2 >
92c4762a1bSJed Brown     = \lambda < 0, 0, 2 > + \mu < 2, 2, 4 >
93c4762a1bSJed Brown */
94c4762a1bSJed Brown static void f0_elas_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
95c4762a1bSJed Brown                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
96c4762a1bSJed Brown                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
97c4762a1bSJed Brown                                 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
98c4762a1bSJed Brown {
99c4762a1bSJed Brown   const PetscReal mu     = 1.0;
100c4762a1bSJed Brown   const PetscReal lambda = 1.0;
101c4762a1bSJed Brown   PetscInt        d;
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   for (d = 0; d < dim-1; ++d) f0[d] += 2.0*mu;
104c4762a1bSJed Brown   f0[dim-1] += 2.0*lambda + 4.0*mu;
105c4762a1bSJed Brown }
106c4762a1bSJed Brown 
107*482dd828SMatthew G. Knepley static void f0_mass_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
108*482dd828SMatthew G. Knepley                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
109*482dd828SMatthew G. Knepley                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
110*482dd828SMatthew G. Knepley                                 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
111*482dd828SMatthew G. Knepley {
112*482dd828SMatthew G. Knepley   if (dim == 2) {
113*482dd828SMatthew G. Knepley     f0[0] -= x[0]*x[0];
114*482dd828SMatthew G. Knepley     f0[1] -= x[1]*x[1] - 2.0*x[0]*x[1];
115*482dd828SMatthew G. Knepley   } else {
116*482dd828SMatthew G. Knepley     f0[0] -= x[0]*x[0];
117*482dd828SMatthew G. Knepley     f0[1] -= x[1]*x[1] - 2.0*x[0]*x[1];
118*482dd828SMatthew G. Knepley     f0[2] -= x[2]*x[2] - 2.0*x[1]*x[2];
119*482dd828SMatthew G. Knepley   }
120*482dd828SMatthew G. Knepley }
121*482dd828SMatthew G. Knepley 
122c4762a1bSJed Brown static PetscErrorCode trig_2d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
123c4762a1bSJed Brown {
124c4762a1bSJed Brown   u[0] = PetscSinReal(2.0*PETSC_PI*x[0]);
125c4762a1bSJed Brown   u[1] = PetscSinReal(2.0*PETSC_PI*x[1]) - 2.0*x[0]*x[1];
126c4762a1bSJed Brown   return 0;
127c4762a1bSJed Brown }
128c4762a1bSJed Brown 
129c4762a1bSJed Brown static PetscErrorCode trig_3d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
130c4762a1bSJed Brown {
131c4762a1bSJed Brown   u[0] = PetscSinReal(2.0*PETSC_PI*x[0]);
132c4762a1bSJed Brown   u[1] = PetscSinReal(2.0*PETSC_PI*x[1]) - 2.0*x[0]*x[1];
133c4762a1bSJed Brown   u[2] = PetscSinReal(2.0*PETSC_PI*x[2]) - 2.0*x[1]*x[2];
134c4762a1bSJed Brown   return 0;
135c4762a1bSJed Brown }
136c4762a1bSJed Brown 
137c4762a1bSJed Brown /*
138c4762a1bSJed Brown   u = sin(2 pi x)
139c4762a1bSJed Brown   v = sin(2 pi y) - 2xy
140c4762a1bSJed 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)>
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   u = sin(2 pi x)
143c4762a1bSJed Brown   v = sin(2 pi y) - 2xy
144c4762a1bSJed Brown   w = sin(2 pi z) - 2yz
145c4762a1bSJed 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)>
146c4762a1bSJed Brown */
147c4762a1bSJed Brown static void f0_vlap_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
148c4762a1bSJed Brown                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
149c4762a1bSJed Brown                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
150c4762a1bSJed Brown                            PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
151c4762a1bSJed Brown {
152c4762a1bSJed Brown   PetscInt d;
153c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f0[d] += -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]);
154c4762a1bSJed Brown }
155c4762a1bSJed Brown 
156c4762a1bSJed Brown /*
157c4762a1bSJed Brown   u = sin(2 pi x)
158c4762a1bSJed Brown   v = sin(2 pi y) - 2xy
159c4762a1bSJed Brown   \varepsilon = / 2 pi cos(2 pi x)             -y        \
160c4762a1bSJed Brown                 \      -y          2 pi cos(2 pi y) - 2x /
161c4762a1bSJed Brown   Tr(\varepsilon) = div u = 2 pi (cos(2 pi x) + cos(2 pi y)) - 2 x
162c4762a1bSJed Brown   div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij}
163c4762a1bSJed 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) >
164c4762a1bSJed 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) >
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   u = sin(2 pi x)
167c4762a1bSJed Brown   v = sin(2 pi y) - 2xy
168c4762a1bSJed Brown   w = sin(2 pi z) - 2yz
169c4762a1bSJed Brown   \varepsilon = / 2 pi cos(2 pi x)            -y                     0         \
170c4762a1bSJed Brown                 |         -y       2 pi cos(2 pi y) - 2x            -z         |
171c4762a1bSJed Brown                 \          0                  -z         2 pi cos(2 pi z) - 2y /
172c4762a1bSJed Brown   Tr(\varepsilon) = div u = 2 pi (cos(2 pi x) + cos(2 pi y) + cos(2 pi z)) - 2 x - 2 y
173c4762a1bSJed Brown   div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij}
174c4762a1bSJed 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) >
175c4762a1bSJed 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) >
176c4762a1bSJed Brown */
177c4762a1bSJed Brown static void f0_elas_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
178c4762a1bSJed Brown                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
179c4762a1bSJed Brown                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
180c4762a1bSJed Brown                            PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
181c4762a1bSJed Brown {
182c4762a1bSJed Brown   const PetscReal mu     = 1.0;
183c4762a1bSJed Brown   const PetscReal lambda = 1.0;
184c4762a1bSJed Brown   const PetscReal fact   = 4.0*PetscSqr(PETSC_PI);
185c4762a1bSJed Brown   PetscInt        d;
186c4762a1bSJed Brown 
187c4762a1bSJed 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);
188c4762a1bSJed Brown }
189c4762a1bSJed Brown 
190c4762a1bSJed Brown static PetscErrorCode axial_disp_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
191c4762a1bSJed Brown {
192c4762a1bSJed Brown   const PetscReal mu     = 1.0;
193c4762a1bSJed Brown   const PetscReal lambda = 1.0;
194c4762a1bSJed Brown   const PetscReal N      = 1.0;
195c4762a1bSJed Brown   PetscInt d;
196c4762a1bSJed Brown 
197c4762a1bSJed 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];
198c4762a1bSJed Brown   u[1] = -0.25*lambda/mu/(lambda+mu)*N*x[1];
199c4762a1bSJed Brown   for (d = 2; d < dim; ++d) u[d] = 0.0;
200c4762a1bSJed Brown   return 0;
201c4762a1bSJed Brown }
202c4762a1bSJed Brown 
203c4762a1bSJed Brown /*
204c4762a1bSJed Brown   We will pull/push on the right side of a block of linearly elastic material. The uniform traction conditions on the
205c4762a1bSJed Brown   right side of the box will result in a uniform strain along x and y. The Neumann BC is given by
206c4762a1bSJed Brown 
207c4762a1bSJed Brown      n_i \sigma_{ij} = t_i
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   u = (1/(2\mu) - 1) x
210c4762a1bSJed Brown   v = -y
211c4762a1bSJed Brown   f = 0
212c4762a1bSJed Brown   t = <4\mu/\lambda (\lambda + \mu), 0>
213c4762a1bSJed Brown   \varepsilon = / 1/(2\mu) - 1   0 \
214c4762a1bSJed Brown                 \ 0             -1 /
215c4762a1bSJed Brown   Tr(\varepsilon) = div u = 1/(2\mu) - 2
216c4762a1bSJed Brown   div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij}
217c4762a1bSJed Brown     = \lambda \partial_j (1/(2\mu) - 2) + 2\mu < 0, 0 >
218c4762a1bSJed Brown     = \lambda < 0, 0 > + \mu < 0, 0 > = 0
219c4762a1bSJed Brown   NBC =  <1,0> . <4\mu/\lambda (\lambda + \mu), 0> = 4\mu/\lambda (\lambda + \mu)
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   u = x - 1/2
222c4762a1bSJed Brown   v = 0
223c4762a1bSJed Brown   w = 0
224c4762a1bSJed Brown   \varepsilon = / x  0  0 \
225c4762a1bSJed Brown                 | 0  0  0 |
226c4762a1bSJed Brown                 \ 0  0  0 /
227c4762a1bSJed Brown   Tr(\varepsilon) = div u = x
228c4762a1bSJed Brown   div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij}
229c4762a1bSJed Brown     = \lambda \partial_j x + 2\mu < 1, 0, 0 >
230c4762a1bSJed Brown     = \lambda < 1, 0, 0 > + \mu < 2, 0, 0 >
231c4762a1bSJed Brown */
232c4762a1bSJed Brown static void f0_elas_axial_disp_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
233c4762a1bSJed Brown                                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
234c4762a1bSJed Brown                                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
235c4762a1bSJed Brown                                     PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
236c4762a1bSJed Brown {
237c4762a1bSJed Brown   const PetscReal N = -1.0;
238c4762a1bSJed Brown 
239c4762a1bSJed Brown   f0[0] = N;
240c4762a1bSJed Brown }
241c4762a1bSJed Brown 
242c4762a1bSJed Brown static PetscErrorCode uniform_strain_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
243c4762a1bSJed Brown {
244c4762a1bSJed Brown   const PetscReal eps_xx = 0.1;
245c4762a1bSJed Brown   const PetscReal eps_xy = 0.3;
246c4762a1bSJed Brown   const PetscReal eps_yy = 0.25;
247c4762a1bSJed Brown   PetscInt d;
248c4762a1bSJed Brown 
249c4762a1bSJed Brown   u[0] = eps_xx*x[0] + eps_xy*x[1];
250c4762a1bSJed Brown   u[1] = eps_xy*x[0] + eps_yy*x[1];
251c4762a1bSJed Brown   for (d = 2; d < dim; ++d) u[d] = 0.0;
252c4762a1bSJed Brown   return 0;
253c4762a1bSJed Brown }
254c4762a1bSJed Brown 
255*482dd828SMatthew G. Knepley static void f0_mass_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
256*482dd828SMatthew G. Knepley                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
257*482dd828SMatthew G. Knepley                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
258*482dd828SMatthew G. Knepley                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
259*482dd828SMatthew G. Knepley {
260*482dd828SMatthew G. Knepley   const PetscInt Nc = dim;
261*482dd828SMatthew G. Knepley   PetscInt       c;
262*482dd828SMatthew G. Knepley 
263*482dd828SMatthew G. Knepley   for (c = 0; c < Nc; ++c) f0[c] = u[c];
264*482dd828SMatthew G. Knepley }
265*482dd828SMatthew G. Knepley 
266c4762a1bSJed Brown static void f1_vlap_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
267c4762a1bSJed Brown                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
268c4762a1bSJed Brown                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
269c4762a1bSJed Brown                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
270c4762a1bSJed Brown {
271c4762a1bSJed Brown   const PetscInt Nc = dim;
272c4762a1bSJed Brown   PetscInt       c, d;
273c4762a1bSJed Brown 
274c4762a1bSJed Brown   for (c = 0; c < Nc; ++c) for (d = 0; d < dim; ++d) f1[c*dim+d] += u_x[c*dim+d];
275c4762a1bSJed Brown }
276c4762a1bSJed Brown 
277c4762a1bSJed Brown static void f1_elas_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
278c4762a1bSJed Brown                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
279c4762a1bSJed Brown                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
280c4762a1bSJed Brown                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
281c4762a1bSJed Brown {
282c4762a1bSJed Brown   const PetscInt  Nc     = dim;
283c4762a1bSJed Brown   const PetscReal mu     = 1.0;
284c4762a1bSJed Brown   const PetscReal lambda = 1.0;
285c4762a1bSJed Brown   PetscInt        c, d;
286c4762a1bSJed Brown 
287c4762a1bSJed Brown   for (c = 0; c < Nc; ++c) {
288c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
289c4762a1bSJed Brown       f1[c*dim+d] += mu*(u_x[c*dim+d] + u_x[d*dim+c]);
290c4762a1bSJed Brown       f1[c*dim+c] += lambda*u_x[d*dim+d];
291c4762a1bSJed Brown     }
292c4762a1bSJed Brown   }
293c4762a1bSJed Brown }
294c4762a1bSJed Brown 
295*482dd828SMatthew G. Knepley static void g0_mass_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
296*482dd828SMatthew G. Knepley                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
297*482dd828SMatthew G. Knepley                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
298*482dd828SMatthew G. Knepley                        PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
299*482dd828SMatthew G. Knepley {
300*482dd828SMatthew G. Knepley   const PetscInt Nc = dim;
301*482dd828SMatthew G. Knepley   PetscInt       c;
302*482dd828SMatthew G. Knepley 
303*482dd828SMatthew G. Knepley   for (c = 0; c < Nc; ++c) g0[c*Nc + c] = 1.0;
304*482dd828SMatthew G. Knepley }
305*482dd828SMatthew G. Knepley 
306c4762a1bSJed Brown static void g3_vlap_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
307c4762a1bSJed Brown                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
308c4762a1bSJed Brown                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
309c4762a1bSJed Brown                        PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
310c4762a1bSJed Brown {
311c4762a1bSJed Brown   const PetscInt Nc = dim;
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] = 1.0;
317c4762a1bSJed Brown     }
318c4762a1bSJed Brown   }
319c4762a1bSJed Brown }
320c4762a1bSJed Brown 
321c4762a1bSJed Brown /*
322c4762a1bSJed Brown   \partial_df \phi_fc g_{fc,gc,df,dg} \partial_dg \phi_gc
323c4762a1bSJed Brown 
324c4762a1bSJed Brown   \partial_df \phi_fc \lambda \delta_{fc,df} \sum_gc \partial_dg \phi_gc \delta_{gc,dg}
325c4762a1bSJed Brown   = \partial_fc \phi_fc \sum_gc \partial_gc \phi_gc
326c4762a1bSJed Brown */
327c4762a1bSJed Brown static void g3_elas_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
328c4762a1bSJed Brown                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
329c4762a1bSJed Brown                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
330c4762a1bSJed Brown                        PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
331c4762a1bSJed Brown {
332c4762a1bSJed Brown   const PetscInt  Nc     = dim;
333c4762a1bSJed Brown   const PetscReal mu     = 1.0;
334c4762a1bSJed Brown   const PetscReal lambda = 1.0;
335c4762a1bSJed Brown   PetscInt        c, d;
336c4762a1bSJed Brown 
337c4762a1bSJed Brown   for (c = 0; c < Nc; ++c) {
338c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
339c4762a1bSJed Brown       g3[((c*Nc + c)*dim + d)*dim + d] += mu;
340c4762a1bSJed Brown       g3[((c*Nc + d)*dim + d)*dim + c] += mu;
341c4762a1bSJed Brown       g3[((c*Nc + d)*dim + c)*dim + d] += lambda;
342c4762a1bSJed Brown     }
343c4762a1bSJed Brown   }
344c4762a1bSJed Brown }
345c4762a1bSJed Brown 
346c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
347c4762a1bSJed Brown {
348*482dd828SMatthew G. Knepley   PetscInt       sol = 0, def = 0;
349c4762a1bSJed Brown   PetscErrorCode ierr;
350c4762a1bSJed Brown 
351c4762a1bSJed Brown   PetscFunctionBeginUser;
352*482dd828SMatthew G. Knepley   options->deform           = DEFORM_NONE;
353c4762a1bSJed Brown   options->solType          = SOL_VLAP_QUADRATIC;
354c4762a1bSJed Brown   options->useNearNullspace = PETSC_TRUE;
355c4762a1bSJed Brown   ierr = PetscStrncpy(options->dmType, DMPLEX, 256);CHKERRQ(ierr);
356c4762a1bSJed Brown 
357c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Linear Elasticity Problem Options", "DMPLEX");CHKERRQ(ierr);
358*482dd828SMatthew G. Knepley   ierr = PetscOptionsEList("-deform_type", "Type of domain deformation", "ex17.c", deformTypes, NUM_DEFORM_TYPES, deformTypes[options->deform], &def, NULL);CHKERRQ(ierr);
359*482dd828SMatthew G. Knepley   options->deform = (DeformType) def;
360c4762a1bSJed Brown   ierr = PetscOptionsEList("-sol_type", "Type of exact solution", "ex17.c", solutionTypes, NUM_SOLUTION_TYPES, solutionTypes[options->solType], &sol, NULL);CHKERRQ(ierr);
361c4762a1bSJed Brown   options->solType = (SolutionType) sol;
362c4762a1bSJed Brown   ierr = PetscOptionsBool("-near_nullspace", "Use the rigid body modes as an AMG near nullspace", "ex17.c", options->useNearNullspace, &options->useNearNullspace, NULL);CHKERRQ(ierr);
363c4762a1bSJed Brown   ierr = PetscOptionsFList("-dm_type", "Convert DMPlex to another format", "ex17.c", DMList, options->dmType, options->dmType, 256, NULL);CHKERRQ(ierr);
364c4762a1bSJed Brown   ierr = PetscOptionsEnd();
365c4762a1bSJed Brown   PetscFunctionReturn(0);
366c4762a1bSJed Brown }
367c4762a1bSJed Brown 
368*482dd828SMatthew G. Knepley static PetscErrorCode DMPlexDistortGeometry(DM dm)
369*482dd828SMatthew G. Knepley {
370*482dd828SMatthew G. Knepley   DM             cdm;
371*482dd828SMatthew G. Knepley   DMLabel        label;
372*482dd828SMatthew G. Knepley   Vec            coordinates;
373*482dd828SMatthew G. Knepley   PetscScalar   *coords;
374*482dd828SMatthew G. Knepley   PetscReal      mid = 0.5;
375*482dd828SMatthew G. Knepley   PetscInt       cdim, d, vStart, vEnd, v;
376*482dd828SMatthew G. Knepley   PetscErrorCode ierr;
377*482dd828SMatthew G. Knepley 
378*482dd828SMatthew G. Knepley   PetscFunctionBeginUser;
379*482dd828SMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
380*482dd828SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr);
381*482dd828SMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
382*482dd828SMatthew G. Knepley   ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr);
383*482dd828SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
384*482dd828SMatthew G. Knepley   ierr = VecGetArrayWrite(coordinates, &coords);CHKERRQ(ierr);
385*482dd828SMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
386*482dd828SMatthew G. Knepley     PetscScalar *pcoords, shift;
387*482dd828SMatthew G. Knepley     PetscInt     val;
388*482dd828SMatthew G. Knepley 
389*482dd828SMatthew G. Knepley     ierr = DMLabelGetValue(label, v, &val);CHKERRQ(ierr);
390*482dd828SMatthew G. Knepley     if (val >= 0) continue;
391*482dd828SMatthew G. Knepley     ierr = DMPlexPointLocalRef(cdm, v, coords, &pcoords);CHKERRQ(ierr);
392*482dd828SMatthew G. Knepley     shift = 0.2 * PetscAbsScalar(pcoords[0] - mid);
393*482dd828SMatthew G. Knepley     shift = PetscRealPart(pcoords[0]) > mid ? shift : -shift;
394*482dd828SMatthew G. Knepley     for (d = 1; d < cdim; ++d) pcoords[d] += shift;
395*482dd828SMatthew G. Knepley   }
396*482dd828SMatthew G. Knepley   ierr = VecRestoreArrayWrite(coordinates, &coords);CHKERRQ(ierr);
397*482dd828SMatthew G. Knepley   PetscFunctionReturn(0);
398*482dd828SMatthew G. Knepley }
399*482dd828SMatthew G. Knepley 
400c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
401c4762a1bSJed Brown {
402c4762a1bSJed Brown   PetscErrorCode ierr;
403c4762a1bSJed Brown 
404c4762a1bSJed Brown   PetscFunctionBeginUser;
40530602db0SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
40630602db0SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
407c4762a1bSJed Brown   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
408*482dd828SMatthew G. Knepley   switch (user->deform) {
409*482dd828SMatthew G. Knepley     case DEFORM_NONE:  break;
410*482dd828SMatthew G. Knepley     case DEFORM_SHEAR: ierr = DMPlexShearGeometry(*dm, DM_X, NULL);CHKERRQ(ierr);break;
411*482dd828SMatthew G. Knepley     case DEFORM_STEP:  ierr = DMPlexDistortGeometry(*dm);CHKERRQ(ierr);break;
412*482dd828SMatthew G. Knepley     default: SETERRQ2(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid deformation type: %s (%D)", deformTypes[PetscMin(user->deform, NUM_DEFORM_TYPES)], user->deform);
413*482dd828SMatthew G. Knepley   }
41430602db0SMatthew G. Knepley   ierr = DMSetApplicationContext(*dm, user);CHKERRQ(ierr);
415c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
416c4762a1bSJed Brown   PetscFunctionReturn(0);
417c4762a1bSJed Brown }
418c4762a1bSJed Brown 
419c4762a1bSJed Brown static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
420c4762a1bSJed Brown {
421c4762a1bSJed Brown   PetscErrorCode (*exact)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
422*482dd828SMatthew G. Knepley   PetscDS        ds;
42345480ffeSMatthew G. Knepley   PetscWeakForm  wf;
42445480ffeSMatthew G. Knepley   DMLabel        label;
42545480ffeSMatthew G. Knepley   PetscInt       id, bd;
426c4762a1bSJed Brown   PetscInt       dim;
427c4762a1bSJed Brown   PetscErrorCode ierr;
428c4762a1bSJed Brown 
429c4762a1bSJed Brown   PetscFunctionBeginUser;
430*482dd828SMatthew G. Knepley   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
431*482dd828SMatthew G. Knepley   ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr);
432*482dd828SMatthew G. Knepley   ierr = PetscDSGetSpatialDimension(ds, &dim);CHKERRQ(ierr);
433c4762a1bSJed Brown   switch (user->solType) {
434*482dd828SMatthew G. Knepley   case SOL_MASS_QUADRATIC:
435*482dd828SMatthew G. Knepley     ierr = PetscDSSetResidual(ds, 0, f0_mass_u, NULL);CHKERRQ(ierr);
436*482dd828SMatthew G. Knepley     ierr = PetscDSSetJacobian(ds, 0, 0, g0_mass_uu, NULL, NULL, NULL);CHKERRQ(ierr);
437*482dd828SMatthew G. Knepley     ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, 0,  0, 1, f0_mass_quadratic_u, 0, NULL);CHKERRQ(ierr);
438c4762a1bSJed Brown     switch (dim) {
439c4762a1bSJed Brown     case 2: exact = quadratic_2d_u;break;
440c4762a1bSJed Brown     case 3: exact = quadratic_3d_u;break;
441*482dd828SMatthew G. Knepley     default: SETERRQ1(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %D", dim);
442*482dd828SMatthew G. Knepley     }
443*482dd828SMatthew G. Knepley     break;
444*482dd828SMatthew G. Knepley   case SOL_VLAP_QUADRATIC:
445*482dd828SMatthew G. Knepley     ierr = PetscDSSetResidual(ds, 0, f0_vlap_quadratic_u, f1_vlap_u);CHKERRQ(ierr);
446*482dd828SMatthew G. Knepley     ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_vlap_uu);CHKERRQ(ierr);
447*482dd828SMatthew G. Knepley     switch (dim) {
448*482dd828SMatthew G. Knepley     case 2: exact = quadratic_2d_u;break;
449*482dd828SMatthew G. Knepley     case 3: exact = quadratic_3d_u;break;
450*482dd828SMatthew G. Knepley     default: SETERRQ1(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %D", dim);
451c4762a1bSJed Brown     }
452c4762a1bSJed Brown     break;
453c4762a1bSJed Brown   case SOL_ELAS_QUADRATIC:
454*482dd828SMatthew G. Knepley     ierr = PetscDSSetResidual(ds, 0, f0_elas_quadratic_u, f1_elas_u);CHKERRQ(ierr);
455*482dd828SMatthew G. Knepley     ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu);CHKERRQ(ierr);
456c4762a1bSJed Brown     switch (dim) {
457c4762a1bSJed Brown     case 2: exact = quadratic_2d_u;break;
458c4762a1bSJed Brown     case 3: exact = quadratic_3d_u;break;
459*482dd828SMatthew G. Knepley     default: SETERRQ1(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %D", dim);
460c4762a1bSJed Brown     }
461c4762a1bSJed Brown     break;
462c4762a1bSJed Brown   case SOL_VLAP_TRIG:
463*482dd828SMatthew G. Knepley     ierr = PetscDSSetResidual(ds, 0, f0_vlap_trig_u, f1_vlap_u);CHKERRQ(ierr);
464*482dd828SMatthew G. Knepley     ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_vlap_uu);CHKERRQ(ierr);
465c4762a1bSJed Brown     switch (dim) {
466c4762a1bSJed Brown     case 2: exact = trig_2d_u;break;
467c4762a1bSJed Brown     case 3: exact = trig_3d_u;break;
468*482dd828SMatthew G. Knepley     default: SETERRQ1(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %D", dim);
469c4762a1bSJed Brown     }
470c4762a1bSJed Brown     break;
471c4762a1bSJed Brown   case SOL_ELAS_TRIG:
472*482dd828SMatthew G. Knepley     ierr = PetscDSSetResidual(ds, 0, f0_elas_trig_u, f1_elas_u);CHKERRQ(ierr);
473*482dd828SMatthew G. Knepley     ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu);CHKERRQ(ierr);
474c4762a1bSJed Brown     switch (dim) {
475c4762a1bSJed Brown     case 2: exact = trig_2d_u;break;
476c4762a1bSJed Brown     case 3: exact = trig_3d_u;break;
477*482dd828SMatthew G. Knepley     default: SETERRQ1(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %D", dim);
478c4762a1bSJed Brown     }
479c4762a1bSJed Brown     break;
480c4762a1bSJed Brown   case SOL_ELAS_AXIAL_DISP:
481*482dd828SMatthew G. Knepley     ierr = PetscDSSetResidual(ds, 0, NULL, f1_elas_u);CHKERRQ(ierr);
482*482dd828SMatthew G. Knepley     ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu);CHKERRQ(ierr);
48345480ffeSMatthew G. Knepley     id   = dim == 3 ? 5 : 2;
48445480ffeSMatthew G. Knepley     ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr);
48545480ffeSMatthew G. Knepley     ierr = DMAddBoundary(dm, DM_BC_NATURAL, "right", label, 1, &id, 0, 0, NULL, (void (*)(void)) NULL, NULL, user, &bd);CHKERRQ(ierr);
486*482dd828SMatthew G. Knepley     ierr = PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
48706ad1575SMatthew G. Knepley     ierr = PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_elas_axial_disp_bd_u, 0, NULL);CHKERRQ(ierr);
488c4762a1bSJed Brown     exact = axial_disp_u;
489c4762a1bSJed Brown     break;
490c4762a1bSJed Brown   case SOL_ELAS_UNIFORM_STRAIN:
491*482dd828SMatthew G. Knepley     ierr = PetscDSSetResidual(ds, 0, NULL, f1_elas_u);CHKERRQ(ierr);
492*482dd828SMatthew G. Knepley     ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu);CHKERRQ(ierr);
493c4762a1bSJed Brown     exact = uniform_strain_u;
494c4762a1bSJed Brown     break;
495*482dd828SMatthew G. Knepley   default: SETERRQ2(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Invalid solution type: %s (%D)", solutionTypes[PetscMin(user->solType, NUM_SOLUTION_TYPES)], user->solType);
496c4762a1bSJed Brown   }
497*482dd828SMatthew G. Knepley   ierr = PetscDSSetExactSolution(ds, 0, exact, user);CHKERRQ(ierr);
49845480ffeSMatthew G. Knepley   ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr);
499c4762a1bSJed Brown   if (user->solType == SOL_ELAS_AXIAL_DISP) {
500c4762a1bSJed Brown     PetscInt cmp;
501c4762a1bSJed Brown 
502c4762a1bSJed Brown     id   = dim == 3 ? 6 : 4;
503c4762a1bSJed Brown     cmp  = 0;
50445480ffeSMatthew G. Knepley     ierr = DMAddBoundary(dm,   DM_BC_ESSENTIAL, "left",   label, 1, &id, 0, 1, &cmp, (void (*)(void)) zero, NULL, user, NULL);CHKERRQ(ierr);
505c4762a1bSJed Brown     cmp  = dim == 3 ? 2 : 1;
506c4762a1bSJed Brown     id   = dim == 3 ? 1 : 1;
50745480ffeSMatthew G. Knepley     ierr = DMAddBoundary(dm,   DM_BC_ESSENTIAL, "bottom", label, 1, &id, 0, 1, &cmp, (void (*)(void)) zero, NULL, user, NULL);CHKERRQ(ierr);
508c4762a1bSJed Brown     if (dim == 3) {
509c4762a1bSJed Brown       cmp  = 1;
510c4762a1bSJed Brown       id   = 3;
51145480ffeSMatthew G. Knepley       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "front",  label, 1, &id, 0, 1, &cmp, (void (*)(void)) zero, NULL, user, NULL);CHKERRQ(ierr);
512c4762a1bSJed Brown     }
513c4762a1bSJed Brown   } else {
514c4762a1bSJed Brown     id = 1;
51545480ffeSMatthew G. Knepley     ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) exact, NULL, user, NULL);CHKERRQ(ierr);
516c4762a1bSJed Brown   }
517c4762a1bSJed Brown   PetscFunctionReturn(0);
518c4762a1bSJed Brown }
519c4762a1bSJed Brown 
5208cda7954SMatthew G. Knepley static PetscErrorCode CreateElasticityNullSpace(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullspace)
521c4762a1bSJed Brown {
522c4762a1bSJed Brown   PetscErrorCode ierr;
523c4762a1bSJed Brown 
524c4762a1bSJed Brown   PetscFunctionBegin;
5258cda7954SMatthew G. Knepley   ierr = DMPlexCreateRigidBody(dm, origField, nullspace);CHKERRQ(ierr);
526c4762a1bSJed Brown   PetscFunctionReturn(0);
527c4762a1bSJed Brown }
528c4762a1bSJed Brown 
52930602db0SMatthew G. Knepley PetscErrorCode SetupFE(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), void *ctx)
530c4762a1bSJed Brown {
531c4762a1bSJed Brown   AppCtx        *user = (AppCtx *) ctx;
532c4762a1bSJed Brown   DM             cdm  = dm;
533c4762a1bSJed Brown   PetscFE        fe;
534c4762a1bSJed Brown   char           prefix[PETSC_MAX_PATH_LEN];
53530602db0SMatthew G. Knepley   DMPolytopeType ct;
53630602db0SMatthew G. Knepley   PetscBool      simplex;
53730602db0SMatthew G. Knepley   PetscInt       dim, cStart;
538c4762a1bSJed Brown   PetscErrorCode ierr;
539c4762a1bSJed Brown 
540c4762a1bSJed Brown   PetscFunctionBegin;
541c4762a1bSJed Brown   /* Create finite element */
542c4762a1bSJed Brown   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
54330602db0SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr);
54430602db0SMatthew G. Knepley   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
54530602db0SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
546c4762a1bSJed Brown   ierr = PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name);CHKERRQ(ierr);
54730602db0SMatthew G. Knepley   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, dim, simplex, name ? prefix : NULL, -1, &fe);CHKERRQ(ierr);
548c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr);
549c4762a1bSJed Brown   /* Set discretization and boundary conditions for each mesh */
550c4762a1bSJed Brown   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr);
551c4762a1bSJed Brown   ierr = DMCreateDS(dm);CHKERRQ(ierr);
552c4762a1bSJed Brown   ierr = (*setup)(dm, user);CHKERRQ(ierr);
553c4762a1bSJed Brown   while (cdm) {
554c4762a1bSJed Brown     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
555c4762a1bSJed Brown     if (user->useNearNullspace) {ierr = DMSetNearNullSpaceConstructor(cdm, 0, CreateElasticityNullSpace);CHKERRQ(ierr);}
556c4762a1bSJed Brown     /* TODO: Check whether the boundary of coarse meshes is marked */
557c4762a1bSJed Brown     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
558c4762a1bSJed Brown   }
559c4762a1bSJed Brown   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
560c4762a1bSJed Brown   PetscFunctionReturn(0);
561c4762a1bSJed Brown }
562c4762a1bSJed Brown 
563c4762a1bSJed Brown int main(int argc, char **argv)
564c4762a1bSJed Brown {
565c4762a1bSJed Brown   DM             dm;   /* Problem specification */
566c4762a1bSJed Brown   SNES           snes; /* Nonlinear solver */
567c4762a1bSJed Brown   Vec            u;    /* Solutions */
568c4762a1bSJed Brown   AppCtx         user; /* User-defined work context */
569c4762a1bSJed Brown   PetscErrorCode ierr;
570c4762a1bSJed Brown 
571c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
572c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
573c4762a1bSJed Brown   /* Primal system */
574c4762a1bSJed Brown   ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
575c4762a1bSJed Brown   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
576c4762a1bSJed Brown   ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
57730602db0SMatthew G. Knepley   ierr = SetupFE(dm, "displacement", SetupPrimalProblem, &user);CHKERRQ(ierr);
578c4762a1bSJed Brown   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
579c4762a1bSJed Brown   ierr = VecSet(u, 0.0);CHKERRQ(ierr);
580c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) u, "displacement");CHKERRQ(ierr);
581c4762a1bSJed Brown   ierr = DMPlexSetSNESLocalFEM(dm, &user, &user, &user);CHKERRQ(ierr);
582c4762a1bSJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
583348a1646SMatthew G. Knepley   ierr = DMSNESCheckFromOptions(snes, u);CHKERRQ(ierr);
584c4762a1bSJed Brown   ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
585c4762a1bSJed Brown   ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr);
586c4762a1bSJed Brown   ierr = VecViewFromOptions(u, NULL, "-displacement_view");CHKERRQ(ierr);
587c4762a1bSJed Brown   /* Cleanup */
588c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
589c4762a1bSJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
590c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
591c4762a1bSJed Brown   ierr = PetscFinalize();
592c4762a1bSJed Brown   return ierr;
593c4762a1bSJed Brown }
594c4762a1bSJed Brown 
595c4762a1bSJed Brown /*TEST
596c4762a1bSJed Brown 
59730602db0SMatthew G. Knepley   testset:
59830602db0SMatthew G. Knepley     args: -dm_plex_box_faces 1,1,1
59930602db0SMatthew G. Knepley 
600c4762a1bSJed Brown     test:
601c4762a1bSJed Brown       suffix: 2d_p1_quad_vlap
602c4762a1bSJed Brown       requires: triangle
603c4762a1bSJed Brown       args: -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate
604c4762a1bSJed Brown     test:
605c4762a1bSJed Brown       suffix: 2d_p2_quad_vlap
606c4762a1bSJed Brown       requires: triangle
607c4762a1bSJed Brown       args: -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001
608c4762a1bSJed Brown     test:
609c4762a1bSJed Brown       suffix: 2d_p3_quad_vlap
610c4762a1bSJed Brown       requires: triangle
611c4762a1bSJed Brown       args: -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001
612c4762a1bSJed Brown     test:
613c4762a1bSJed Brown       suffix: 2d_q1_quad_vlap
61430602db0SMatthew G. Knepley       args: -dm_plex_simplex 0 -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate
615c4762a1bSJed Brown     test:
616c4762a1bSJed Brown       suffix: 2d_q2_quad_vlap
61730602db0SMatthew G. Knepley       args: -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001
618c4762a1bSJed Brown     test:
619c4762a1bSJed Brown       suffix: 2d_q3_quad_vlap
620c4762a1bSJed Brown       requires: !single
62130602db0SMatthew G. Knepley       args: -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001
622c4762a1bSJed Brown     test:
623c4762a1bSJed Brown       suffix: 2d_p1_quad_elas
624c4762a1bSJed Brown       requires: triangle
625c4762a1bSJed Brown       args: -sol_type elas_quad -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate
626c4762a1bSJed Brown     test:
627c4762a1bSJed Brown       suffix: 2d_p2_quad_elas
628c4762a1bSJed Brown       requires: triangle
629c4762a1bSJed Brown       args: -sol_type elas_quad -displacement_petscspace_degree 2 -dmsnes_check .0001
630c4762a1bSJed Brown     test:
631c4762a1bSJed Brown       suffix: 2d_p3_quad_elas
632c4762a1bSJed Brown       requires: triangle
633c4762a1bSJed Brown       args: -sol_type elas_quad -displacement_petscspace_degree 3 -dmsnes_check .0001
634c4762a1bSJed Brown     test:
635c4762a1bSJed Brown       suffix: 2d_q1_quad_elas
63630602db0SMatthew 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
637c4762a1bSJed Brown     test:
638c4762a1bSJed Brown       suffix: 2d_q1_quad_elas_shear
639*482dd828SMatthew 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
640c4762a1bSJed Brown     test:
641c4762a1bSJed Brown       suffix: 2d_q2_quad_elas
64230602db0SMatthew G. Knepley       args: -sol_type elas_quad -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dmsnes_check .0001
643c4762a1bSJed Brown     test:
644c4762a1bSJed Brown       suffix: 2d_q2_quad_elas_shear
645*482dd828SMatthew G. Knepley       args: -sol_type elas_quad -dm_plex_simplex 0 -deform_type shear -displacement_petscspace_degree 2 -dmsnes_check
646c4762a1bSJed Brown     test:
647c4762a1bSJed Brown       suffix: 2d_q3_quad_elas
64830602db0SMatthew G. Knepley       args: -sol_type elas_quad -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dmsnes_check .0001
649c4762a1bSJed Brown     test:
650c4762a1bSJed Brown       suffix: 2d_q3_quad_elas_shear
651c4762a1bSJed Brown       requires: !single
652*482dd828SMatthew G. Knepley       args: -sol_type elas_quad -dm_plex_simplex 0 -deform_type shear -displacement_petscspace_degree 3 -dmsnes_check
653c4762a1bSJed Brown 
654c4762a1bSJed Brown     test:
655c4762a1bSJed Brown       suffix: 3d_p1_quad_vlap
656c4762a1bSJed Brown       requires: ctetgen
65730602db0SMatthew G. Knepley       args: -dm_plex_dim 3 -dm_refine 1 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate
658c4762a1bSJed Brown     test:
659c4762a1bSJed Brown       suffix: 3d_p2_quad_vlap
660c4762a1bSJed Brown       requires: ctetgen
66130602db0SMatthew G. Knepley       args: -dm_plex_dim 3 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001
662c4762a1bSJed Brown     test:
663c4762a1bSJed Brown       suffix: 3d_p3_quad_vlap
664c4762a1bSJed Brown       requires: ctetgen
66530602db0SMatthew G. Knepley       args: -dm_plex_dim 3 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001
666c4762a1bSJed Brown     test:
667c4762a1bSJed Brown       suffix: 3d_q1_quad_vlap
66830602db0SMatthew 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
669c4762a1bSJed Brown     test:
670c4762a1bSJed Brown       suffix: 3d_q2_quad_vlap
67130602db0SMatthew G. Knepley       args: -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001
672c4762a1bSJed Brown     test:
673c4762a1bSJed Brown       suffix: 3d_q3_quad_vlap
67430602db0SMatthew G. Knepley       args: -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001
675c4762a1bSJed Brown     test:
676c4762a1bSJed Brown       suffix: 3d_p1_quad_elas
677c4762a1bSJed Brown       requires: ctetgen
67830602db0SMatthew 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
679c4762a1bSJed Brown     test:
680c4762a1bSJed Brown       suffix: 3d_p2_quad_elas
681c4762a1bSJed Brown       requires: ctetgen
68230602db0SMatthew G. Knepley       args: -sol_type elas_quad -dm_plex_dim 3 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001
683c4762a1bSJed Brown     test:
684c4762a1bSJed Brown       suffix: 3d_p3_quad_elas
685c4762a1bSJed Brown       requires: ctetgen
68630602db0SMatthew G. Knepley       args: -sol_type elas_quad -dm_plex_dim 3 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001
687c4762a1bSJed Brown     test:
688c4762a1bSJed Brown       suffix: 3d_q1_quad_elas
68930602db0SMatthew 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
690c4762a1bSJed Brown     test:
691c4762a1bSJed Brown       suffix: 3d_q2_quad_elas
69230602db0SMatthew 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
693c4762a1bSJed Brown     test:
694c4762a1bSJed Brown       suffix: 3d_q3_quad_elas
695c4762a1bSJed Brown       requires: !single
69630602db0SMatthew 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
697c4762a1bSJed Brown 
698c4762a1bSJed Brown     test:
699c4762a1bSJed Brown       suffix: 2d_p1_trig_vlap
700c4762a1bSJed Brown       requires: triangle
701c4762a1bSJed Brown       args: -sol_type vlap_trig -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
702c4762a1bSJed Brown     test:
703c4762a1bSJed Brown       suffix: 2d_p2_trig_vlap
704c4762a1bSJed Brown       requires: triangle
705c4762a1bSJed Brown       args: -sol_type vlap_trig -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
706c4762a1bSJed Brown     test:
707c4762a1bSJed Brown       suffix: 2d_p3_trig_vlap
708c4762a1bSJed Brown       requires: triangle
709c4762a1bSJed Brown       args: -sol_type vlap_trig -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
710c4762a1bSJed Brown     test:
711c4762a1bSJed Brown       suffix: 2d_q1_trig_vlap
71230602db0SMatthew 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
713c4762a1bSJed Brown     test:
714c4762a1bSJed Brown       suffix: 2d_q2_trig_vlap
71530602db0SMatthew 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
716c4762a1bSJed Brown     test:
717c4762a1bSJed Brown       suffix: 2d_q3_trig_vlap
71830602db0SMatthew 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
719c4762a1bSJed Brown     test:
720c4762a1bSJed Brown       suffix: 2d_p1_trig_elas
721c4762a1bSJed Brown       requires: triangle
722c4762a1bSJed Brown       args: -sol_type elas_trig -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
723c4762a1bSJed Brown     test:
724c4762a1bSJed Brown       suffix: 2d_p2_trig_elas
725c4762a1bSJed Brown       requires: triangle
726c4762a1bSJed Brown       args: -sol_type elas_trig -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
727c4762a1bSJed Brown     test:
728c4762a1bSJed Brown       suffix: 2d_p3_trig_elas
729c4762a1bSJed Brown       requires: triangle
730c4762a1bSJed Brown       args: -sol_type elas_trig -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
731c4762a1bSJed Brown     test:
732c4762a1bSJed Brown       suffix: 2d_q1_trig_elas
73330602db0SMatthew 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
734c4762a1bSJed Brown     test:
735c4762a1bSJed Brown       suffix: 2d_q1_trig_elas_shear
736*482dd828SMatthew 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
737c4762a1bSJed Brown     test:
738c4762a1bSJed Brown       suffix: 2d_q2_trig_elas
73930602db0SMatthew 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
740c4762a1bSJed Brown     test:
741c4762a1bSJed Brown       suffix: 2d_q2_trig_elas_shear
742*482dd828SMatthew 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
743c4762a1bSJed Brown     test:
744c4762a1bSJed Brown       suffix: 2d_q3_trig_elas
74530602db0SMatthew 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
746c4762a1bSJed Brown     test:
747c4762a1bSJed Brown       suffix: 2d_q3_trig_elas_shear
748*482dd828SMatthew 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
749c4762a1bSJed Brown 
750c4762a1bSJed Brown     test:
751c4762a1bSJed Brown       suffix: 3d_p1_trig_vlap
752c4762a1bSJed Brown       requires: ctetgen
75330602db0SMatthew 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
754c4762a1bSJed Brown     test:
755c4762a1bSJed Brown       suffix: 3d_p2_trig_vlap
756c4762a1bSJed Brown       requires: ctetgen
75730602db0SMatthew 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
758c4762a1bSJed Brown     test:
759c4762a1bSJed Brown       suffix: 3d_p3_trig_vlap
760c4762a1bSJed Brown       requires: ctetgen
76130602db0SMatthew 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
762c4762a1bSJed Brown     test:
763c4762a1bSJed Brown       suffix: 3d_q1_trig_vlap
76430602db0SMatthew 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
765c4762a1bSJed Brown     test:
766c4762a1bSJed Brown       suffix: 3d_q2_trig_vlap
76730602db0SMatthew 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
768c4762a1bSJed Brown     test:
769c4762a1bSJed Brown       suffix: 3d_q3_trig_vlap
770c4762a1bSJed Brown       requires: !__float128
77130602db0SMatthew 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
772c4762a1bSJed Brown     test:
773c4762a1bSJed Brown       suffix: 3d_p1_trig_elas
774c4762a1bSJed Brown       requires: ctetgen
77530602db0SMatthew 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
776c4762a1bSJed Brown     test:
777c4762a1bSJed Brown       suffix: 3d_p2_trig_elas
778c4762a1bSJed Brown       requires: ctetgen
77930602db0SMatthew 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
780c4762a1bSJed Brown     test:
781c4762a1bSJed Brown       suffix: 3d_p3_trig_elas
782c4762a1bSJed Brown       requires: ctetgen
78330602db0SMatthew 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
784c4762a1bSJed Brown     test:
785c4762a1bSJed Brown       suffix: 3d_q1_trig_elas
78630602db0SMatthew 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
787c4762a1bSJed Brown     test:
788c4762a1bSJed Brown       suffix: 3d_q2_trig_elas
78930602db0SMatthew 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
790c4762a1bSJed Brown     test:
791c4762a1bSJed Brown       suffix: 3d_q3_trig_elas
792c4762a1bSJed Brown       requires: !__float128
79330602db0SMatthew 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
794c4762a1bSJed Brown 
795c4762a1bSJed Brown     test:
796c4762a1bSJed Brown       suffix: 2d_p1_axial_elas
797c4762a1bSJed Brown       requires: triangle
798c4762a1bSJed Brown       args: -sol_type elas_axial_disp -displacement_petscspace_degree 1 -dm_plex_separate_marker -dm_refine 2 -dmsnes_check .0001 -pc_type lu
799c4762a1bSJed Brown     test:
800c4762a1bSJed Brown       suffix: 2d_p2_axial_elas
801c4762a1bSJed Brown       requires: triangle
802c4762a1bSJed Brown       args: -sol_type elas_axial_disp -displacement_petscspace_degree 2 -dm_plex_separate_marker -dmsnes_check .0001 -pc_type lu
803c4762a1bSJed Brown     test:
804c4762a1bSJed Brown       suffix: 2d_p3_axial_elas
805c4762a1bSJed Brown       requires: triangle
806c4762a1bSJed Brown       args: -sol_type elas_axial_disp -displacement_petscspace_degree 3 -dm_plex_separate_marker -dmsnes_check .0001 -pc_type lu
807c4762a1bSJed Brown     test:
808c4762a1bSJed Brown       suffix: 2d_q1_axial_elas
80930602db0SMatthew 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
810c4762a1bSJed Brown     test:
811c4762a1bSJed Brown       suffix: 2d_q2_axial_elas
81230602db0SMatthew 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
813c4762a1bSJed Brown     test:
814c4762a1bSJed Brown       suffix: 2d_q3_axial_elas
81530602db0SMatthew 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
816c4762a1bSJed Brown 
817c4762a1bSJed Brown     test:
818c4762a1bSJed Brown       suffix: 2d_p1_uniform_elas
819c4762a1bSJed Brown       requires: triangle
820c4762a1bSJed Brown       args: -sol_type elas_uniform_strain -displacement_petscspace_degree 1 -dm_refine 2 -dmsnes_check .0001 -pc_type lu
821c4762a1bSJed Brown     test:
822c4762a1bSJed Brown       suffix: 2d_p2_uniform_elas
823c4762a1bSJed Brown       requires: triangle
824c4762a1bSJed Brown       args: -sol_type elas_uniform_strain -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001 -pc_type lu
825c4762a1bSJed Brown     test:
826c4762a1bSJed Brown       suffix: 2d_p3_uniform_elas
827c4762a1bSJed Brown       requires: triangle
828c4762a1bSJed Brown       args: -sol_type elas_uniform_strain -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001 -pc_type lu
829c4762a1bSJed Brown     test:
830c4762a1bSJed Brown       suffix: 2d_q1_uniform_elas
83130602db0SMatthew 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
832c4762a1bSJed Brown     test:
833c4762a1bSJed Brown       suffix: 2d_q2_uniform_elas
834c4762a1bSJed Brown       requires: !single
83530602db0SMatthew 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
836c4762a1bSJed Brown     test:
837c4762a1bSJed Brown       suffix: 2d_q3_uniform_elas
838c4762a1bSJed Brown       requires: !single
83930602db0SMatthew 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
840*482dd828SMatthew G. Knepley     test:
841*482dd828SMatthew G. Knepley       suffix: 2d_p1_uniform_elas_step
842*482dd828SMatthew G. Knepley       requires: triangle
843*482dd828SMatthew G. Knepley       args: -sol_type elas_uniform_strain -deform_type step -displacement_petscspace_degree 1 -dm_refine 2 -dmsnes_check .0001 -pc_type lu
844*482dd828SMatthew G. Knepley 
845*482dd828SMatthew G. Knepley   testset:
846*482dd828SMatthew 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
847*482dd828SMatthew G. Knepley 
848*482dd828SMatthew G. Knepley     test:
849*482dd828SMatthew G. Knepley       suffix: 2d_q1_uniform_elas_step
850*482dd828SMatthew G. Knepley       args: -sol_type elas_uniform_strain -dm_refine 2
851*482dd828SMatthew G. Knepley     test:
852*482dd828SMatthew G. Knepley       suffix: 2d_q1_quad_vlap_step
853*482dd828SMatthew G. Knepley       args:
854*482dd828SMatthew G. Knepley     test:
855*482dd828SMatthew G. Knepley       suffix: 2d_q2_quad_vlap_step
856*482dd828SMatthew G. Knepley       args: -displacement_petscspace_degree 2
857*482dd828SMatthew G. Knepley     test:
858*482dd828SMatthew G. Knepley       suffix: 2d_q1_quad_mass_step
859*482dd828SMatthew G. Knepley       args: -sol_type mass_quad
860c4762a1bSJed Brown 
861c4762a1bSJed Brown TEST*/
862