Actual source code: ex17.c

  1: static char help[] = "Linear elasticity in 2d and 3d with finite elements.\n\
  2: We solve the elasticity problem in a rectangular\n\
  3: domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\
  4: This example supports automatic convergence estimation\n\
  5: and eventually adaptivity.\n\n\n";

  7: /*
  8:   https://en.wikipedia.org/wiki/Linear_elasticity

 10:   Converting elastic constants:
 11:     lambda = E nu / ((1 + nu) (1 - 2 nu))
 12:     mu     = E / (2 (1 + nu))
 13: */

 15: #include <petscdmplex.h>
 16: #include <petscsnes.h>
 17: #include <petscds.h>
 18: #include <petscbag.h>
 19: #include <petscconvest.h>

 21: typedef enum {
 22:   SOL_VLAP_QUADRATIC,
 23:   SOL_ELAS_QUADRATIC,
 24:   SOL_VLAP_TRIG,
 25:   SOL_ELAS_TRIG,
 26:   SOL_ELAS_AXIAL_DISP,
 27:   SOL_ELAS_UNIFORM_STRAIN,
 28:   SOL_ELAS_GE,
 29:   SOL_MASS_QUADRATIC,
 30:   NUM_SOLUTION_TYPES
 31: } SolutionType;
 32: 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"};

 34: typedef enum {
 35:   DEFORM_NONE,
 36:   DEFORM_SHEAR,
 37:   DEFORM_STEP,
 38:   NUM_DEFORM_TYPES
 39: } DeformType;
 40: const char *deformTypes[NUM_DEFORM_TYPES + 1] = {"none", "shear", "step", "unknown"};

 42: typedef struct {
 43:   PetscScalar mu;     /* shear modulus */
 44:   PetscScalar lambda; /* Lame's first parameter */
 45: } Parameter;

 47: typedef struct {
 48:   /* Domain and mesh definition */
 49:   char       dmType[256]; /* DM type for the solve */
 50:   DeformType deform;      /* Domain deformation type */
 51:   /* Problem definition */
 52:   SolutionType solType; /* Type of exact solution */
 53:   PetscBag     bag;     /* Problem parameters */
 54:   /* Solver definition */
 55:   PetscBool useNearNullspace; /* Use the rigid body modes as a near nullspace for AMG */
 56: } AppCtx;

 58: static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
 59: {
 60:   PetscInt d;
 61:   for (d = 0; d < dim; ++d) u[d] = 0.0;
 62:   return 0;
 63: }

 65: static PetscErrorCode ge_shift(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
 66: {
 67:   PetscInt d;
 68:   u[0] = 0.1;
 69:   for (d = 1; d < dim; ++d) u[d] = 0.0;
 70:   return 0;
 71: }

 73: static PetscErrorCode quadratic_2d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
 74: {
 75:   u[0] = x[0] * x[0];
 76:   u[1] = x[1] * x[1] - 2.0 * x[0] * x[1];
 77:   return 0;
 78: }

 80: static PetscErrorCode quadratic_3d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
 81: {
 82:   u[0] = x[0] * x[0];
 83:   u[1] = x[1] * x[1] - 2.0 * x[0] * x[1];
 84:   u[2] = x[2] * x[2] - 2.0 * x[1] * x[2];
 85:   return 0;
 86: }

 88: /*
 89:   u = x^2
 90:   v = y^2 - 2xy
 91:   Delta <u,v> - f = <2, 2> - <2, 2>

 93:   u = x^2
 94:   v = y^2 - 2xy
 95:   w = z^2 - 2yz
 96:   Delta <u,v,w> - f = <2, 2, 2> - <2, 2, 2>
 97: */
 98: 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[])
 99: {
100:   PetscInt d;
101:   for (d = 0; d < dim; ++d) f0[d] += 2.0;
102: }

104: /*
105:   u = x^2
106:   v = y^2 - 2xy
107:   \varepsilon = / 2x     -y    \
108:                 \ -y   2y - 2x /
109:   Tr(\varepsilon) = div u = 2y
110:   div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij}
111:     = \lambda \partial_j (2y) + 2\mu < 2-1, 2 >
112:     = \lambda < 0, 2 > + \mu < 2, 4 >

114:   u = x^2
115:   v = y^2 - 2xy
116:   w = z^2 - 2yz
117:   \varepsilon = / 2x     -y       0   \
118:                 | -y   2y - 2x   -z   |
119:                 \  0     -z    2z - 2y/
120:   Tr(\varepsilon) = div u = 2z
121:   div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij}
122:     = \lambda \partial_j (2z) + 2\mu < 2-1, 2-1, 2 >
123:     = \lambda < 0, 0, 2 > + \mu < 2, 2, 4 >
124: */
125: 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[])
126: {
127:   const PetscReal mu     = 1.0;
128:   const PetscReal lambda = 1.0;
129:   PetscInt        d;

131:   for (d = 0; d < dim - 1; ++d) f0[d] += 2.0 * mu;
132:   f0[dim - 1] += 2.0 * lambda + 4.0 * mu;
133: }

135: 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[])
136: {
137:   if (dim == 2) {
138:     f0[0] -= x[0] * x[0];
139:     f0[1] -= x[1] * x[1] - 2.0 * x[0] * x[1];
140:   } else {
141:     f0[0] -= x[0] * x[0];
142:     f0[1] -= x[1] * x[1] - 2.0 * x[0] * x[1];
143:     f0[2] -= x[2] * x[2] - 2.0 * x[1] * x[2];
144:   }
145: }

147: static PetscErrorCode trig_2d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
148: {
149:   u[0] = PetscSinReal(2.0 * PETSC_PI * x[0]);
150:   u[1] = PetscSinReal(2.0 * PETSC_PI * x[1]) - 2.0 * x[0] * x[1];
151:   return 0;
152: }

154: static PetscErrorCode trig_3d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
155: {
156:   u[0] = PetscSinReal(2.0 * PETSC_PI * x[0]);
157:   u[1] = PetscSinReal(2.0 * PETSC_PI * x[1]) - 2.0 * x[0] * x[1];
158:   u[2] = PetscSinReal(2.0 * PETSC_PI * x[2]) - 2.0 * x[1] * x[2];
159:   return 0;
160: }

162: /*
163:   u = sin(2 pi x)
164:   v = sin(2 pi y) - 2xy
165:   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)>

167:   u = sin(2 pi x)
168:   v = sin(2 pi y) - 2xy
169:   w = sin(2 pi z) - 2yz
170:   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)>
171: */
172: 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[])
173: {
174:   PetscInt d;
175:   for (d = 0; d < dim; ++d) f0[d] += -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]);
176: }

178: /*
179:   u = sin(2 pi x)
180:   v = sin(2 pi y) - 2xy
181:   \varepsilon = / 2 pi cos(2 pi x)             -y        \
182:                 \      -y          2 pi cos(2 pi y) - 2x /
183:   Tr(\varepsilon) = div u = 2 pi (cos(2 pi x) + cos(2 pi y)) - 2 x
184:   div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij}
185:     = \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) >
186:     = \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) >

188:   u = sin(2 pi x)
189:   v = sin(2 pi y) - 2xy
190:   w = sin(2 pi z) - 2yz
191:   \varepsilon = / 2 pi cos(2 pi x)            -y                     0         \
192:                 |         -y       2 pi cos(2 pi y) - 2x            -z         |
193:                 \          0                  -z         2 pi cos(2 pi z) - 2y /
194:   Tr(\varepsilon) = div u = 2 pi (cos(2 pi x) + cos(2 pi y) + cos(2 pi z)) - 2 x - 2 y
195:   div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij}
196:     = \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) >
197:     = \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) >
198: */
199: 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[])
200: {
201:   const PetscReal mu     = 1.0;
202:   const PetscReal lambda = 1.0;
203:   const PetscReal fact   = 4.0 * PetscSqr(PETSC_PI);
204:   PetscInt        d;

206:   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);
207: }

209: static PetscErrorCode axial_disp_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
210: {
211:   const PetscReal mu     = 1.0;
212:   const PetscReal lambda = 1.0;
213:   const PetscReal N      = 1.0;
214:   PetscInt        d;

216:   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];
217:   u[1] = -0.25 * lambda / mu / (lambda + mu) * N * x[1];
218:   for (d = 2; d < dim; ++d) u[d] = 0.0;
219:   return 0;
220: }

222: /*
223:   We will pull/push on the right side of a block of linearly elastic material. The uniform traction conditions on the
224:   right side of the box will result in a uniform strain along x and y. The Neumann BC is given by

226:      n_i \sigma_{ij} = t_i

228:   u = (1/(2\mu) - 1) x
229:   v = -y
230:   f = 0
231:   t = <4\mu/\lambda (\lambda + \mu), 0>
232:   \varepsilon = / 1/(2\mu) - 1   0 \
233:                 \ 0             -1 /
234:   Tr(\varepsilon) = div u = 1/(2\mu) - 2
235:   div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij}
236:     = \lambda \partial_j (1/(2\mu) - 2) + 2\mu < 0, 0 >
237:     = \lambda < 0, 0 > + \mu < 0, 0 > = 0
238:   NBC =  <1,0> . <4\mu/\lambda (\lambda + \mu), 0> = 4\mu/\lambda (\lambda + \mu)

240:   u = x - 1/2
241:   v = 0
242:   w = 0
243:   \varepsilon = / x  0  0 \
244:                 | 0  0  0 |
245:                 \ 0  0  0 /
246:   Tr(\varepsilon) = div u = x
247:   div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij}
248:     = \lambda \partial_j x + 2\mu < 1, 0, 0 >
249:     = \lambda < 1, 0, 0 > + \mu < 2, 0, 0 >
250: */
251: 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[])
252: {
253:   const PetscReal N = -1.0;

255:   f0[0] = N;
256: }

258: static PetscErrorCode uniform_strain_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
259: {
260:   const PetscReal eps_xx = 0.1;
261:   const PetscReal eps_xy = 0.3;
262:   const PetscReal eps_yy = 0.25;
263:   PetscInt        d;

265:   u[0] = eps_xx * x[0] + eps_xy * x[1];
266:   u[1] = eps_xy * x[0] + eps_yy * x[1];
267:   for (d = 2; d < dim; ++d) u[d] = 0.0;
268:   return 0;
269: }

271: 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[])
272: {
273:   const PetscInt Nc = dim;
274:   PetscInt       c;

276:   for (c = 0; c < Nc; ++c) f0[c] = u[c];
277: }

279: 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[])
280: {
281:   const PetscInt Nc = dim;
282:   PetscInt       c, d;

284:   for (c = 0; c < Nc; ++c)
285:     for (d = 0; d < dim; ++d) f1[c * dim + d] += u_x[c * dim + d];
286: }

288: 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[])
289: {
290:   const PetscInt  Nc     = dim;
291:   const PetscReal mu     = 1.0;
292:   const PetscReal lambda = 1.0;
293:   PetscInt        c, d;

295:   for (c = 0; c < Nc; ++c) {
296:     for (d = 0; d < dim; ++d) {
297:       f1[c * dim + d] += mu * (u_x[c * dim + d] + u_x[d * dim + c]);
298:       f1[c * dim + c] += lambda * u_x[d * dim + d];
299:     }
300:   }
301: }

303: 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[])
304: {
305:   const PetscInt Nc = dim;
306:   PetscInt       c;

308:   for (c = 0; c < Nc; ++c) g0[c * Nc + c] = 1.0;
309: }

311: 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[])
312: {
313:   const PetscInt Nc = dim;
314:   PetscInt       c, d;

316:   for (c = 0; c < Nc; ++c) {
317:     for (d = 0; d < dim; ++d) g3[((c * Nc + c) * dim + d) * dim + d] = 1.0;
318:   }
319: }

321: /*
322:   \partial_df \phi_fc g_{fc,gc,df,dg} \partial_dg \phi_gc

324:   \partial_df \phi_fc \lambda \delta_{fc,df} \sum_gc \partial_dg \phi_gc \delta_{gc,dg}
325:   = \partial_fc \phi_fc \sum_gc \partial_gc \phi_gc
326: */
327: 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[])
328: {
329:   const PetscInt  Nc     = dim;
330:   const PetscReal mu     = 1.0;
331:   const PetscReal lambda = 1.0;
332:   PetscInt        c, d;

334:   for (c = 0; c < Nc; ++c) {
335:     for (d = 0; d < dim; ++d) {
336:       g3[((c * Nc + c) * dim + d) * dim + d] += mu;
337:       g3[((c * Nc + d) * dim + d) * dim + c] += mu;
338:       g3[((c * Nc + d) * dim + c) * dim + d] += lambda;
339:     }
340:   }
341: }

343: static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
344: {
345:   PetscInt sol = 0, def = 0;

348:   options->deform           = DEFORM_NONE;
349:   options->solType          = SOL_VLAP_QUADRATIC;
350:   options->useNearNullspace = PETSC_TRUE;
351:   PetscStrncpy(options->dmType, DMPLEX, 256);

353:   PetscOptionsBegin(comm, "", "Linear Elasticity Problem Options", "DMPLEX");
354:   PetscOptionsEList("-deform_type", "Type of domain deformation", "ex17.c", deformTypes, NUM_DEFORM_TYPES, deformTypes[options->deform], &def, NULL);
355:   options->deform = (DeformType)def;
356:   PetscOptionsEList("-sol_type", "Type of exact solution", "ex17.c", solutionTypes, NUM_SOLUTION_TYPES, solutionTypes[options->solType], &sol, NULL);
357:   options->solType = (SolutionType)sol;
358:   PetscOptionsBool("-near_nullspace", "Use the rigid body modes as an AMG near nullspace", "ex17.c", options->useNearNullspace, &options->useNearNullspace, NULL);
359:   PetscOptionsFList("-dm_type", "Convert DMPlex to another format", "ex17.c", DMList, options->dmType, options->dmType, 256, NULL);
360:   PetscOptionsEnd();
361:   return 0;
362: }

364: static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *ctx)
365: {
366:   PetscBag   bag;
367:   Parameter *p;

370:   /* setup PETSc parameter bag */
371:   PetscBagGetData(ctx->bag, (void **)&p);
372:   PetscBagSetName(ctx->bag, "par", "Elastic Parameters");
373:   bag = ctx->bag;
374:   PetscBagRegisterScalar(bag, &p->mu, 1.0, "mu", "Shear Modulus, Pa");
375:   PetscBagRegisterScalar(bag, &p->lambda, 1.0, "lambda", "Lame's first parameter, Pa");
376:   PetscBagSetFromOptions(bag);
377:   {
378:     PetscViewer       viewer;
379:     PetscViewerFormat format;
380:     PetscBool         flg;

382:     PetscOptionsGetViewer(comm, NULL, NULL, "-param_view", &viewer, &format, &flg);
383:     if (flg) {
384:       PetscViewerPushFormat(viewer, format);
385:       PetscBagView(bag, viewer);
386:       PetscViewerFlush(viewer);
387:       PetscViewerPopFormat(viewer);
388:       PetscViewerDestroy(&viewer);
389:     }
390:   }
391:   return 0;
392: }

394: static PetscErrorCode DMPlexDistortGeometry(DM dm)
395: {
396:   DM           cdm;
397:   DMLabel      label;
398:   Vec          coordinates;
399:   PetscScalar *coords;
400:   PetscReal    mid = 0.5;
401:   PetscInt     cdim, d, vStart, vEnd, v;

404:   DMGetCoordinateDM(dm, &cdm);
405:   DMGetCoordinateDim(dm, &cdim);
406:   DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
407:   DMGetLabel(dm, "marker", &label);
408:   DMGetCoordinatesLocal(dm, &coordinates);
409:   VecGetArrayWrite(coordinates, &coords);
410:   for (v = vStart; v < vEnd; ++v) {
411:     PetscScalar *pcoords, shift;
412:     PetscInt     val;

414:     DMLabelGetValue(label, v, &val);
415:     if (val >= 0) continue;
416:     DMPlexPointLocalRef(cdm, v, coords, &pcoords);
417:     shift = 0.2 * PetscAbsScalar(pcoords[0] - mid);
418:     shift = PetscRealPart(pcoords[0]) > mid ? shift : -shift;
419:     for (d = 1; d < cdim; ++d) pcoords[d] += shift;
420:   }
421:   VecRestoreArrayWrite(coordinates, &coords);
422:   return 0;
423: }

425: static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
426: {
428:   DMCreate(comm, dm);
429:   DMSetType(*dm, DMPLEX);
430:   DMSetFromOptions(*dm);
431:   switch (user->deform) {
432:   case DEFORM_NONE:
433:     break;
434:   case DEFORM_SHEAR:
435:     DMPlexShearGeometry(*dm, DM_X, NULL);
436:     break;
437:   case DEFORM_STEP:
438:     DMPlexDistortGeometry(*dm);
439:     break;
440:   default:
441:     SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid deformation type: %s (%d)", deformTypes[PetscMin(user->deform, NUM_DEFORM_TYPES)], user->deform);
442:   }
443:   DMSetApplicationContext(*dm, user);
444:   DMViewFromOptions(*dm, NULL, "-dm_view");
445:   return 0;
446: }

448: static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
449: {
450:   PetscErrorCode (*exact)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
451:   Parameter    *param;
452:   PetscDS       ds;
453:   PetscWeakForm wf;
454:   DMLabel       label;
455:   PetscInt      id, bd;
456:   PetscInt      dim;

459:   DMGetDS(dm, &ds);
460:   PetscDSGetWeakForm(ds, &wf);
461:   PetscDSGetSpatialDimension(ds, &dim);
462:   PetscBagGetData(user->bag, (void **)&param);
463:   switch (user->solType) {
464:   case SOL_MASS_QUADRATIC:
465:     PetscDSSetResidual(ds, 0, f0_mass_u, NULL);
466:     PetscDSSetJacobian(ds, 0, 0, g0_mass_uu, NULL, NULL, NULL);
467:     PetscWeakFormSetIndexResidual(wf, NULL, 0, 0, 0, 1, f0_mass_quadratic_u, 0, NULL);
468:     switch (dim) {
469:     case 2:
470:       exact = quadratic_2d_u;
471:       break;
472:     case 3:
473:       exact = quadratic_3d_u;
474:       break;
475:     default:
476:       SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim);
477:     }
478:     break;
479:   case SOL_VLAP_QUADRATIC:
480:     PetscDSSetResidual(ds, 0, f0_vlap_quadratic_u, f1_vlap_u);
481:     PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_vlap_uu);
482:     switch (dim) {
483:     case 2:
484:       exact = quadratic_2d_u;
485:       break;
486:     case 3:
487:       exact = quadratic_3d_u;
488:       break;
489:     default:
490:       SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim);
491:     }
492:     break;
493:   case SOL_ELAS_QUADRATIC:
494:     PetscDSSetResidual(ds, 0, f0_elas_quadratic_u, f1_elas_u);
495:     PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu);
496:     switch (dim) {
497:     case 2:
498:       exact = quadratic_2d_u;
499:       break;
500:     case 3:
501:       exact = quadratic_3d_u;
502:       break;
503:     default:
504:       SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim);
505:     }
506:     break;
507:   case SOL_VLAP_TRIG:
508:     PetscDSSetResidual(ds, 0, f0_vlap_trig_u, f1_vlap_u);
509:     PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_vlap_uu);
510:     switch (dim) {
511:     case 2:
512:       exact = trig_2d_u;
513:       break;
514:     case 3:
515:       exact = trig_3d_u;
516:       break;
517:     default:
518:       SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim);
519:     }
520:     break;
521:   case SOL_ELAS_TRIG:
522:     PetscDSSetResidual(ds, 0, f0_elas_trig_u, f1_elas_u);
523:     PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu);
524:     switch (dim) {
525:     case 2:
526:       exact = trig_2d_u;
527:       break;
528:     case 3:
529:       exact = trig_3d_u;
530:       break;
531:     default:
532:       SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim);
533:     }
534:     break;
535:   case SOL_ELAS_AXIAL_DISP:
536:     PetscDSSetResidual(ds, 0, NULL, f1_elas_u);
537:     PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu);
538:     id = dim == 3 ? 5 : 2;
539:     DMGetLabel(dm, "marker", &label);
540:     DMAddBoundary(dm, DM_BC_NATURAL, "right", label, 1, &id, 0, 0, NULL, (void (*)(void))NULL, NULL, user, &bd);
541:     PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);
542:     PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_elas_axial_disp_bd_u, 0, NULL);
543:     exact = axial_disp_u;
544:     break;
545:   case SOL_ELAS_UNIFORM_STRAIN:
546:     PetscDSSetResidual(ds, 0, NULL, f1_elas_u);
547:     PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu);
548:     exact = uniform_strain_u;
549:     break;
550:   case SOL_ELAS_GE:
551:     PetscDSSetResidual(ds, 0, NULL, f1_elas_u);
552:     PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu);
553:     exact = zero; /* No exact solution available */
554:     break;
555:   default:
556:     SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid solution type: %s (%d)", solutionTypes[PetscMin(user->solType, NUM_SOLUTION_TYPES)], user->solType);
557:   }
558:   PetscDSSetExactSolution(ds, 0, exact, user);
559:   DMGetLabel(dm, "marker", &label);
560:   if (user->solType == SOL_ELAS_AXIAL_DISP) {
561:     PetscInt cmp;

563:     id  = dim == 3 ? 6 : 4;
564:     cmp = 0;
565:     DMAddBoundary(dm, DM_BC_ESSENTIAL, "left", label, 1, &id, 0, 1, &cmp, (void (*)(void))zero, NULL, user, NULL);
566:     cmp = dim == 3 ? 2 : 1;
567:     id  = dim == 3 ? 1 : 1;
568:     DMAddBoundary(dm, DM_BC_ESSENTIAL, "bottom", label, 1, &id, 0, 1, &cmp, (void (*)(void))zero, NULL, user, NULL);
569:     if (dim == 3) {
570:       cmp = 1;
571:       id  = 3;
572:       DMAddBoundary(dm, DM_BC_ESSENTIAL, "front", label, 1, &id, 0, 1, &cmp, (void (*)(void))zero, NULL, user, NULL);
573:     }
574:   } else if (user->solType == SOL_ELAS_GE) {
575:     PetscInt cmp;

577:     id = dim == 3 ? 6 : 4;
578:     DMAddBoundary(dm, DM_BC_ESSENTIAL, "left", label, 1, &id, 0, 0, NULL, (void (*)(void))zero, NULL, user, NULL);
579:     id  = dim == 3 ? 5 : 2;
580:     cmp = 0;
581:     DMAddBoundary(dm, DM_BC_ESSENTIAL, "right", label, 1, &id, 0, 1, &cmp, (void (*)(void))ge_shift, NULL, user, NULL);
582:   } else {
583:     id = 1;
584:     DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))exact, NULL, user, NULL);
585:   }
586:   /* Setup constants */
587:   {
588:     PetscScalar constants[2];

590:     constants[0] = param->mu;     /* shear modulus, Pa */
591:     constants[1] = param->lambda; /* Lame's first parameter, Pa */
592:     PetscDSSetConstants(ds, 2, constants);
593:   }
594:   return 0;
595: }

597: static PetscErrorCode CreateElasticityNullSpace(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullspace)
598: {
599:   DMPlexCreateRigidBody(dm, origField, nullspace);
600:   return 0;
601: }

603: PetscErrorCode SetupFE(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), void *ctx)
604: {
605:   AppCtx        *user = (AppCtx *)ctx;
606:   DM             cdm  = dm;
607:   PetscFE        fe;
608:   char           prefix[PETSC_MAX_PATH_LEN];
609:   DMPolytopeType ct;
610:   PetscBool      simplex;
611:   PetscInt       dim, cStart;

613:   /* Create finite element */
614:   DMGetDimension(dm, &dim);
615:   DMPlexGetHeightStratum(dm, 0, &cStart, NULL);
616:   DMPlexGetCellType(dm, cStart, &ct);
617:   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
618:   PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name);
619:   PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, dim, simplex, name ? prefix : NULL, -1, &fe);
620:   PetscObjectSetName((PetscObject)fe, name);
621:   /* Set discretization and boundary conditions for each mesh */
622:   DMSetField(dm, 0, NULL, (PetscObject)fe);
623:   DMCreateDS(dm);
624:   (*setup)(dm, user);
625:   while (cdm) {
626:     DMCopyDisc(dm, cdm);
627:     if (user->useNearNullspace) DMSetNearNullSpaceConstructor(cdm, 0, CreateElasticityNullSpace);
628:     /* TODO: Check whether the boundary of coarse meshes is marked */
629:     DMGetCoarseDM(cdm, &cdm);
630:   }
631:   PetscFEDestroy(&fe);
632:   return 0;
633: }

635: int main(int argc, char **argv)
636: {
637:   DM     dm;   /* Problem specification */
638:   SNES   snes; /* Nonlinear solver */
639:   Vec    u;    /* Solutions */
640:   AppCtx user; /* User-defined work context */

643:   PetscInitialize(&argc, &argv, NULL, help);
644:   ProcessOptions(PETSC_COMM_WORLD, &user);
645:   PetscBagCreate(PETSC_COMM_SELF, sizeof(Parameter), &user.bag);
646:   SetupParameters(PETSC_COMM_WORLD, &user);
647:   /* Primal system */
648:   SNESCreate(PETSC_COMM_WORLD, &snes);
649:   CreateMesh(PETSC_COMM_WORLD, &user, &dm);
650:   SNESSetDM(snes, dm);
651:   SetupFE(dm, "displacement", SetupPrimalProblem, &user);
652:   DMCreateGlobalVector(dm, &u);
653:   VecSet(u, 0.0);
654:   PetscObjectSetName((PetscObject)u, "displacement");
655:   DMPlexSetSNESLocalFEM(dm, &user, &user, &user);
656:   SNESSetFromOptions(snes);
657:   DMSNESCheckFromOptions(snes, u);
658:   SNESSolve(snes, NULL, u);
659:   SNESGetSolution(snes, &u);
660:   VecViewFromOptions(u, NULL, "-displacement_view");
661:   /* Cleanup */
662:   VecDestroy(&u);
663:   SNESDestroy(&snes);
664:   DMDestroy(&dm);
665:   PetscBagDestroy(&user.bag);
666:   PetscFinalize();
667:   return 0;
668: }

670: /*TEST

672:   testset:
673:     args: -dm_plex_box_faces 1,1,1

675:     test:
676:       suffix: 2d_p1_quad_vlap
677:       requires: triangle
678:       args: -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate
679:     test:
680:       suffix: 2d_p2_quad_vlap
681:       requires: triangle
682:       args: -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001
683:     test:
684:       suffix: 2d_p3_quad_vlap
685:       requires: triangle
686:       args: -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001
687:     test:
688:       suffix: 2d_q1_quad_vlap
689:       args: -dm_plex_simplex 0 -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate
690:     test:
691:       suffix: 2d_q2_quad_vlap
692:       args: -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001
693:     test:
694:       suffix: 2d_q3_quad_vlap
695:       requires: !single
696:       args: -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001
697:     test:
698:       suffix: 2d_p1_quad_elas
699:       requires: triangle
700:       args: -sol_type elas_quad -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate
701:     test:
702:       suffix: 2d_p2_quad_elas
703:       requires: triangle
704:       args: -sol_type elas_quad -displacement_petscspace_degree 2 -dmsnes_check .0001
705:     test:
706:       suffix: 2d_p3_quad_elas
707:       requires: triangle
708:       args: -sol_type elas_quad -displacement_petscspace_degree 3 -dmsnes_check .0001
709:     test:
710:       suffix: 2d_q1_quad_elas
711:       args: -sol_type elas_quad -dm_plex_simplex 0 -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
712:     test:
713:       suffix: 2d_q1_quad_elas_shear
714:       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
715:     test:
716:       suffix: 2d_q2_quad_elas
717:       args: -sol_type elas_quad -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dmsnes_check .0001
718:     test:
719:       suffix: 2d_q2_quad_elas_shear
720:       args: -sol_type elas_quad -dm_plex_simplex 0 -deform_type shear -displacement_petscspace_degree 2 -dmsnes_check
721:     test:
722:       suffix: 2d_q3_quad_elas
723:       args: -sol_type elas_quad -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dmsnes_check .0001
724:     test:
725:       suffix: 2d_q3_quad_elas_shear
726:       requires: !single
727:       args: -sol_type elas_quad -dm_plex_simplex 0 -deform_type shear -displacement_petscspace_degree 3 -dmsnes_check

729:     test:
730:       suffix: 3d_p1_quad_vlap
731:       requires: ctetgen
732:       args: -dm_plex_dim 3 -dm_refine 1 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate
733:     test:
734:       suffix: 3d_p2_quad_vlap
735:       requires: ctetgen
736:       args: -dm_plex_dim 3 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001
737:     test:
738:       suffix: 3d_p3_quad_vlap
739:       requires: ctetgen
740:       args: -dm_plex_dim 3 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001
741:     test:
742:       suffix: 3d_q1_quad_vlap
743:       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
744:     test:
745:       suffix: 3d_q2_quad_vlap
746:       args: -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001
747:     test:
748:       suffix: 3d_q3_quad_vlap
749:       args: -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001
750:     test:
751:       suffix: 3d_p1_quad_elas
752:       requires: ctetgen
753:       args: -sol_type elas_quad -dm_plex_dim 3 -dm_refine 1 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate
754:     test:
755:       suffix: 3d_p2_quad_elas
756:       requires: ctetgen
757:       args: -sol_type elas_quad -dm_plex_dim 3 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001
758:     test:
759:       suffix: 3d_p3_quad_elas
760:       requires: ctetgen
761:       args: -sol_type elas_quad -dm_plex_dim 3 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001
762:     test:
763:       suffix: 3d_q1_quad_elas
764:       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
765:     test:
766:       suffix: 3d_q2_quad_elas
767:       args: -sol_type elas_quad -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001
768:     test:
769:       suffix: 3d_q3_quad_elas
770:       requires: !single
771:       args: -sol_type elas_quad -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001

773:     test:
774:       suffix: 2d_p1_trig_vlap
775:       requires: triangle
776:       args: -sol_type vlap_trig -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
777:     test:
778:       suffix: 2d_p2_trig_vlap
779:       requires: triangle
780:       args: -sol_type vlap_trig -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
781:     test:
782:       suffix: 2d_p3_trig_vlap
783:       requires: triangle
784:       args: -sol_type vlap_trig -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
785:     test:
786:       suffix: 2d_q1_trig_vlap
787:       args: -sol_type vlap_trig -dm_plex_simplex 0 -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
788:     test:
789:       suffix: 2d_q2_trig_vlap
790:       args: -sol_type vlap_trig -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
791:     test:
792:       suffix: 2d_q3_trig_vlap
793:       args: -sol_type vlap_trig -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
794:     test:
795:       suffix: 2d_p1_trig_elas
796:       requires: triangle
797:       args: -sol_type elas_trig -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
798:     test:
799:       suffix: 2d_p2_trig_elas
800:       requires: triangle
801:       args: -sol_type elas_trig -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
802:     test:
803:       suffix: 2d_p3_trig_elas
804:       requires: triangle
805:       args: -sol_type elas_trig -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
806:     test:
807:       suffix: 2d_q1_trig_elas
808:       args: -sol_type elas_trig -dm_plex_simplex 0 -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
809:     test:
810:       suffix: 2d_q1_trig_elas_shear
811:       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
812:     test:
813:       suffix: 2d_q2_trig_elas
814:       args: -sol_type elas_trig -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
815:     test:
816:       suffix: 2d_q2_trig_elas_shear
817:       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
818:     test:
819:       suffix: 2d_q3_trig_elas
820:       args: -sol_type elas_trig -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate
821:     test:
822:       suffix: 2d_q3_trig_elas_shear
823:       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

825:     test:
826:       suffix: 3d_p1_trig_vlap
827:       requires: ctetgen
828:       args: -sol_type vlap_trig -dm_plex_dim 3 -dm_refine 1 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate
829:     test:
830:       suffix: 3d_p2_trig_vlap
831:       requires: ctetgen
832:       args: -sol_type vlap_trig -dm_plex_dim 3 -displacement_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate
833:     test:
834:       suffix: 3d_p3_trig_vlap
835:       requires: ctetgen
836:       args: -sol_type vlap_trig -dm_plex_dim 3 -displacement_petscspace_degree 3 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate
837:     test:
838:       suffix: 3d_q1_trig_vlap
839:       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
840:     test:
841:       suffix: 3d_q2_trig_vlap
842:       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
843:     test:
844:       suffix: 3d_q3_trig_vlap
845:       requires: !__float128
846:       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
847:     test:
848:       suffix: 3d_p1_trig_elas
849:       requires: ctetgen
850:       args: -sol_type elas_trig -dm_plex_dim 3 -dm_refine 1 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate
851:     test:
852:       suffix: 3d_p2_trig_elas
853:       requires: ctetgen
854:       args: -sol_type elas_trig -dm_plex_dim 3 -displacement_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate
855:     test:
856:       suffix: 3d_p3_trig_elas
857:       requires: ctetgen
858:       args: -sol_type elas_trig -dm_plex_dim 3 -displacement_petscspace_degree 3 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate
859:     test:
860:       suffix: 3d_q1_trig_elas
861:       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
862:     test:
863:       suffix: 3d_q2_trig_elas
864:       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
865:     test:
866:       suffix: 3d_q3_trig_elas
867:       requires: !__float128
868:       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

870:     test:
871:       suffix: 2d_p1_axial_elas
872:       requires: triangle
873:       args: -sol_type elas_axial_disp -displacement_petscspace_degree 1 -dm_plex_separate_marker -dm_refine 2 -dmsnes_check .0001 -pc_type lu
874:     test:
875:       suffix: 2d_p2_axial_elas
876:       requires: triangle
877:       args: -sol_type elas_axial_disp -displacement_petscspace_degree 2 -dm_plex_separate_marker -dmsnes_check .0001 -pc_type lu
878:     test:
879:       suffix: 2d_p3_axial_elas
880:       requires: triangle
881:       args: -sol_type elas_axial_disp -displacement_petscspace_degree 3 -dm_plex_separate_marker -dmsnes_check .0001 -pc_type lu
882:     test:
883:       suffix: 2d_q1_axial_elas
884:       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
885:     test:
886:       suffix: 2d_q2_axial_elas
887:       args: -sol_type elas_axial_disp -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_plex_separate_marker -dmsnes_check .0001 -pc_type   lu
888:     test:
889:       suffix: 2d_q3_axial_elas
890:       args: -sol_type elas_axial_disp -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_plex_separate_marker -dmsnes_check .0001 -pc_type   lu

892:     test:
893:       suffix: 2d_p1_uniform_elas
894:       requires: triangle
895:       args: -sol_type elas_uniform_strain -displacement_petscspace_degree 1 -dm_refine 2 -dmsnes_check .0001 -pc_type lu
896:     test:
897:       suffix: 2d_p2_uniform_elas
898:       requires: triangle
899:       args: -sol_type elas_uniform_strain -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001 -pc_type lu
900:     test:
901:       suffix: 2d_p3_uniform_elas
902:       requires: triangle
903:       args: -sol_type elas_uniform_strain -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001 -pc_type lu
904:     test:
905:       suffix: 2d_q1_uniform_elas
906:       args: -sol_type elas_uniform_strain -dm_plex_simplex 0 -displacement_petscspace_degree 1 -dm_refine 2 -dmsnes_check .0001 -pc_type lu
907:     test:
908:       suffix: 2d_q2_uniform_elas
909:       requires: !single
910:       args: -sol_type elas_uniform_strain -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001 -pc_type lu
911:     test:
912:       suffix: 2d_q3_uniform_elas
913:       requires: !single
914:       args: -sol_type elas_uniform_strain -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001 -pc_type lu
915:     test:
916:       suffix: 2d_p1_uniform_elas_step
917:       requires: triangle
918:       args: -sol_type elas_uniform_strain -deform_type step -displacement_petscspace_degree 1 -dm_refine 2 -dmsnes_check .0001 -pc_type lu

920:   testset:
921:     args: -dm_plex_simplex 0 -dm_plex_box_faces 3,3 -deform_type step -displacement_petscspace_degree 1 -dmsnes_check .0001 -pc_type lu

923:     test:
924:       suffix: 2d_q1_uniform_elas_step
925:       args: -sol_type elas_uniform_strain -dm_refine 2
926:     test:
927:       suffix: 2d_q1_quad_vlap_step
928:       args:
929:     test:
930:       suffix: 2d_q2_quad_vlap_step
931:       args: -displacement_petscspace_degree 2
932:     test:
933:       suffix: 2d_q1_quad_mass_step
934:       args: -sol_type mass_quad

936:   testset:
937:     filter: grep -v "variant HERMITIAN"
938:     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 \
939:           -dm_plex_box_faces 5,5,2 -dm_plex_separate_marker -dm_refine 0 -petscpartitioner_type simple \
940:           -sol_type elas_ge

942:     test:
943:       suffix: ge_q1_0
944:       args: -displacement_petscspace_degree 1 \
945:             -snes_max_it 2 -snes_rtol 1.e-10 \
946:             -ksp_type cg -ksp_rtol 1.e-10 -ksp_max_it 100 -ksp_norm_type unpreconditioned \
947:             -pc_type gamg -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 \
948:               -pc_gamg_coarse_eq_limit 10 -pc_gamg_reuse_interpolation true \
949:               -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 \
950:               -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 \
951:               -matptap_via scalable
952:     test:
953:       nsize: 5
954:       suffix: ge_q1_gdsw
955:       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

957: TEST*/