1c4762a1bSJed Brown static char help[] = "Linear elasticity in 2d and 3d with finite elements.\n\ 2c4762a1bSJed Brown We solve the elasticity problem in a rectangular\n\ 3c4762a1bSJed Brown domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\ 4c4762a1bSJed Brown This example supports automatic convergence estimation\n\ 5c4762a1bSJed Brown and eventually adaptivity.\n\n\n"; 6c4762a1bSJed Brown 7c4762a1bSJed Brown /* 8c4762a1bSJed Brown https://en.wikipedia.org/wiki/Linear_elasticity 9d12cd81dSMatthew G. Knepley 10d12cd81dSMatthew G. Knepley Converting elastic constants: 11d12cd81dSMatthew G. Knepley lambda = E nu / ((1 + nu) (1 - 2 nu)) 12d12cd81dSMatthew G. Knepley mu = E / (2 (1 + nu)) 13c4762a1bSJed Brown */ 14c4762a1bSJed Brown 15c4762a1bSJed Brown #include <petscdmplex.h> 16c4762a1bSJed Brown #include <petscsnes.h> 17c4762a1bSJed Brown #include <petscds.h> 18d12cd81dSMatthew G. Knepley #include <petscbag.h> 19c4762a1bSJed Brown #include <petscconvest.h> 20c4762a1bSJed Brown 219371c9d4SSatish Balay typedef enum { 229371c9d4SSatish Balay SOL_VLAP_QUADRATIC, 239371c9d4SSatish Balay SOL_ELAS_QUADRATIC, 249371c9d4SSatish Balay SOL_VLAP_TRIG, 259371c9d4SSatish Balay SOL_ELAS_TRIG, 269371c9d4SSatish Balay SOL_ELAS_AXIAL_DISP, 279371c9d4SSatish Balay SOL_ELAS_UNIFORM_STRAIN, 289371c9d4SSatish Balay SOL_ELAS_GE, 299371c9d4SSatish Balay SOL_MASS_QUADRATIC, 309371c9d4SSatish Balay NUM_SOLUTION_TYPES 319371c9d4SSatish Balay } SolutionType; 32d12cd81dSMatthew G. Knepley const char *solutionTypes[NUM_SOLUTION_TYPES + 1] = {"vlap_quad", "elas_quad", "vlap_trig", "elas_trig", "elas_axial_disp", "elas_uniform_strain", "elas_ge", "mass_quad", "unknown"}; 33482dd828SMatthew G. Knepley 349371c9d4SSatish Balay typedef enum { 359371c9d4SSatish Balay DEFORM_NONE, 369371c9d4SSatish Balay DEFORM_SHEAR, 379371c9d4SSatish Balay DEFORM_STEP, 389371c9d4SSatish Balay NUM_DEFORM_TYPES 399371c9d4SSatish Balay } DeformType; 40482dd828SMatthew G. Knepley const char *deformTypes[NUM_DEFORM_TYPES + 1] = {"none", "shear", "step", "unknown"}; 41c4762a1bSJed Brown 42c4762a1bSJed Brown typedef struct { 43d12cd81dSMatthew G. Knepley PetscScalar mu; /* shear modulus */ 44d12cd81dSMatthew G. Knepley PetscScalar lambda; /* Lame's first parameter */ 45d12cd81dSMatthew G. Knepley } Parameter; 46d12cd81dSMatthew G. Knepley 47d12cd81dSMatthew G. Knepley typedef struct { 48c4762a1bSJed Brown /* Domain and mesh definition */ 49c4762a1bSJed Brown char dmType[256]; /* DM type for the solve */ 50482dd828SMatthew G. Knepley DeformType deform; /* Domain deformation type */ 51c4762a1bSJed Brown /* Problem definition */ 52c4762a1bSJed Brown SolutionType solType; /* Type of exact solution */ 53d12cd81dSMatthew G. Knepley PetscBag bag; /* Problem parameters */ 54c4762a1bSJed Brown /* Solver definition */ 55c4762a1bSJed Brown PetscBool useNearNullspace; /* Use the rigid body modes as a near nullspace for AMG */ 56c4762a1bSJed Brown } AppCtx; 57c4762a1bSJed Brown 58d71ae5a4SJacob Faibussowitsch static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 59d71ae5a4SJacob Faibussowitsch { 60c4762a1bSJed Brown PetscInt d; 61c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[d] = 0.0; 62*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 63c4762a1bSJed Brown } 64c4762a1bSJed Brown 65d71ae5a4SJacob Faibussowitsch static PetscErrorCode ge_shift(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 66d71ae5a4SJacob Faibussowitsch { 67d12cd81dSMatthew G. Knepley PetscInt d; 68d12cd81dSMatthew G. Knepley u[0] = 0.1; 69d12cd81dSMatthew G. Knepley for (d = 1; d < dim; ++d) u[d] = 0.0; 70*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 71d12cd81dSMatthew G. Knepley } 72d12cd81dSMatthew G. Knepley 73d71ae5a4SJacob Faibussowitsch static PetscErrorCode quadratic_2d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 74d71ae5a4SJacob Faibussowitsch { 75c4762a1bSJed Brown u[0] = x[0] * x[0]; 76c4762a1bSJed Brown u[1] = x[1] * x[1] - 2.0 * x[0] * x[1]; 77*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 78c4762a1bSJed Brown } 79c4762a1bSJed Brown 80d71ae5a4SJacob Faibussowitsch static PetscErrorCode quadratic_3d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 81d71ae5a4SJacob Faibussowitsch { 82c4762a1bSJed Brown u[0] = x[0] * x[0]; 83c4762a1bSJed Brown u[1] = x[1] * x[1] - 2.0 * x[0] * x[1]; 84c4762a1bSJed Brown u[2] = x[2] * x[2] - 2.0 * x[1] * x[2]; 85*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 86c4762a1bSJed Brown } 87c4762a1bSJed Brown 88c4762a1bSJed Brown /* 89c4762a1bSJed Brown u = x^2 90c4762a1bSJed Brown v = y^2 - 2xy 91c4762a1bSJed Brown Delta <u,v> - f = <2, 2> - <2, 2> 92c4762a1bSJed Brown 93c4762a1bSJed Brown u = x^2 94c4762a1bSJed Brown v = y^2 - 2xy 95c4762a1bSJed Brown w = z^2 - 2yz 96c4762a1bSJed Brown Delta <u,v,w> - f = <2, 2, 2> - <2, 2, 2> 97c4762a1bSJed Brown */ 98d71ae5a4SJacob Faibussowitsch 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[]) 99d71ae5a4SJacob Faibussowitsch { 100c4762a1bSJed Brown PetscInt d; 101c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[d] += 2.0; 102c4762a1bSJed Brown } 103c4762a1bSJed Brown 104c4762a1bSJed Brown /* 105c4762a1bSJed Brown u = x^2 106c4762a1bSJed Brown v = y^2 - 2xy 107c4762a1bSJed Brown \varepsilon = / 2x -y \ 108c4762a1bSJed Brown \ -y 2y - 2x / 109c4762a1bSJed Brown Tr(\varepsilon) = div u = 2y 110c4762a1bSJed Brown div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij} 111c4762a1bSJed Brown = \lambda \partial_j (2y) + 2\mu < 2-1, 2 > 112c4762a1bSJed Brown = \lambda < 0, 2 > + \mu < 2, 4 > 113c4762a1bSJed Brown 114c4762a1bSJed Brown u = x^2 115c4762a1bSJed Brown v = y^2 - 2xy 116c4762a1bSJed Brown w = z^2 - 2yz 117c4762a1bSJed Brown \varepsilon = / 2x -y 0 \ 118c4762a1bSJed Brown | -y 2y - 2x -z | 119c4762a1bSJed Brown \ 0 -z 2z - 2y/ 120c4762a1bSJed Brown Tr(\varepsilon) = div u = 2z 121c4762a1bSJed Brown div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij} 122c4762a1bSJed Brown = \lambda \partial_j (2z) + 2\mu < 2-1, 2-1, 2 > 123c4762a1bSJed Brown = \lambda < 0, 0, 2 > + \mu < 2, 2, 4 > 124c4762a1bSJed Brown */ 125d71ae5a4SJacob Faibussowitsch 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[]) 126d71ae5a4SJacob Faibussowitsch { 127c4762a1bSJed Brown const PetscReal mu = 1.0; 128c4762a1bSJed Brown const PetscReal lambda = 1.0; 129c4762a1bSJed Brown PetscInt d; 130c4762a1bSJed Brown 131c4762a1bSJed Brown for (d = 0; d < dim - 1; ++d) f0[d] += 2.0 * mu; 132c4762a1bSJed Brown f0[dim - 1] += 2.0 * lambda + 4.0 * mu; 133c4762a1bSJed Brown } 134c4762a1bSJed Brown 135d71ae5a4SJacob Faibussowitsch 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[]) 136d71ae5a4SJacob Faibussowitsch { 137482dd828SMatthew G. Knepley if (dim == 2) { 138482dd828SMatthew G. Knepley f0[0] -= x[0] * x[0]; 139482dd828SMatthew G. Knepley f0[1] -= x[1] * x[1] - 2.0 * x[0] * x[1]; 140482dd828SMatthew G. Knepley } else { 141482dd828SMatthew G. Knepley f0[0] -= x[0] * x[0]; 142482dd828SMatthew G. Knepley f0[1] -= x[1] * x[1] - 2.0 * x[0] * x[1]; 143482dd828SMatthew G. Knepley f0[2] -= x[2] * x[2] - 2.0 * x[1] * x[2]; 144482dd828SMatthew G. Knepley } 145482dd828SMatthew G. Knepley } 146482dd828SMatthew G. Knepley 147d71ae5a4SJacob Faibussowitsch static PetscErrorCode trig_2d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 148d71ae5a4SJacob Faibussowitsch { 149c4762a1bSJed Brown u[0] = PetscSinReal(2.0 * PETSC_PI * x[0]); 150c4762a1bSJed Brown u[1] = PetscSinReal(2.0 * PETSC_PI * x[1]) - 2.0 * x[0] * x[1]; 151*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 152c4762a1bSJed Brown } 153c4762a1bSJed Brown 154d71ae5a4SJacob Faibussowitsch static PetscErrorCode trig_3d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 155d71ae5a4SJacob Faibussowitsch { 156c4762a1bSJed Brown u[0] = PetscSinReal(2.0 * PETSC_PI * x[0]); 157c4762a1bSJed Brown u[1] = PetscSinReal(2.0 * PETSC_PI * x[1]) - 2.0 * x[0] * x[1]; 158c4762a1bSJed Brown u[2] = PetscSinReal(2.0 * PETSC_PI * x[2]) - 2.0 * x[1] * x[2]; 159*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 160c4762a1bSJed Brown } 161c4762a1bSJed Brown 162c4762a1bSJed Brown /* 163c4762a1bSJed Brown u = sin(2 pi x) 164c4762a1bSJed Brown v = sin(2 pi y) - 2xy 165c4762a1bSJed 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)> 166c4762a1bSJed Brown 167c4762a1bSJed Brown u = sin(2 pi x) 168c4762a1bSJed Brown v = sin(2 pi y) - 2xy 169c4762a1bSJed Brown w = sin(2 pi z) - 2yz 170c4762a1bSJed 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)> 171c4762a1bSJed Brown */ 172d71ae5a4SJacob Faibussowitsch 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[]) 173d71ae5a4SJacob Faibussowitsch { 174c4762a1bSJed Brown PetscInt d; 175c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[d] += -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]); 176c4762a1bSJed Brown } 177c4762a1bSJed Brown 178c4762a1bSJed Brown /* 179c4762a1bSJed Brown u = sin(2 pi x) 180c4762a1bSJed Brown v = sin(2 pi y) - 2xy 181c4762a1bSJed Brown \varepsilon = / 2 pi cos(2 pi x) -y \ 182c4762a1bSJed Brown \ -y 2 pi cos(2 pi y) - 2x / 183c4762a1bSJed Brown Tr(\varepsilon) = div u = 2 pi (cos(2 pi x) + cos(2 pi y)) - 2 x 184c4762a1bSJed Brown div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij} 185c4762a1bSJed 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) > 186c4762a1bSJed 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) > 187c4762a1bSJed Brown 188c4762a1bSJed Brown u = sin(2 pi x) 189c4762a1bSJed Brown v = sin(2 pi y) - 2xy 190c4762a1bSJed Brown w = sin(2 pi z) - 2yz 191c4762a1bSJed Brown \varepsilon = / 2 pi cos(2 pi x) -y 0 \ 192c4762a1bSJed Brown | -y 2 pi cos(2 pi y) - 2x -z | 193c4762a1bSJed Brown \ 0 -z 2 pi cos(2 pi z) - 2y / 194c4762a1bSJed Brown Tr(\varepsilon) = div u = 2 pi (cos(2 pi x) + cos(2 pi y) + cos(2 pi z)) - 2 x - 2 y 195c4762a1bSJed Brown div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij} 196c4762a1bSJed 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) > 197c4762a1bSJed 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) > 198c4762a1bSJed Brown */ 199d71ae5a4SJacob Faibussowitsch 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[]) 200d71ae5a4SJacob Faibussowitsch { 201c4762a1bSJed Brown const PetscReal mu = 1.0; 202c4762a1bSJed Brown const PetscReal lambda = 1.0; 203c4762a1bSJed Brown const PetscReal fact = 4.0 * PetscSqr(PETSC_PI); 204c4762a1bSJed Brown PetscInt d; 205c4762a1bSJed Brown 206c4762a1bSJed 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); 207c4762a1bSJed Brown } 208c4762a1bSJed Brown 209d71ae5a4SJacob Faibussowitsch static PetscErrorCode axial_disp_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 210d71ae5a4SJacob Faibussowitsch { 211c4762a1bSJed Brown const PetscReal mu = 1.0; 212c4762a1bSJed Brown const PetscReal lambda = 1.0; 213c4762a1bSJed Brown const PetscReal N = 1.0; 214c4762a1bSJed Brown PetscInt d; 215c4762a1bSJed Brown 216c4762a1bSJed 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]; 217c4762a1bSJed Brown u[1] = -0.25 * lambda / mu / (lambda + mu) * N * x[1]; 218c4762a1bSJed Brown for (d = 2; d < dim; ++d) u[d] = 0.0; 219*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 220c4762a1bSJed Brown } 221c4762a1bSJed Brown 222c4762a1bSJed Brown /* 223c4762a1bSJed Brown We will pull/push on the right side of a block of linearly elastic material. The uniform traction conditions on the 224c4762a1bSJed Brown right side of the box will result in a uniform strain along x and y. The Neumann BC is given by 225c4762a1bSJed Brown 226c4762a1bSJed Brown n_i \sigma_{ij} = t_i 227c4762a1bSJed Brown 228c4762a1bSJed Brown u = (1/(2\mu) - 1) x 229c4762a1bSJed Brown v = -y 230c4762a1bSJed Brown f = 0 231c4762a1bSJed Brown t = <4\mu/\lambda (\lambda + \mu), 0> 232c4762a1bSJed Brown \varepsilon = / 1/(2\mu) - 1 0 \ 233c4762a1bSJed Brown \ 0 -1 / 234c4762a1bSJed Brown Tr(\varepsilon) = div u = 1/(2\mu) - 2 235c4762a1bSJed Brown div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij} 236c4762a1bSJed Brown = \lambda \partial_j (1/(2\mu) - 2) + 2\mu < 0, 0 > 237c4762a1bSJed Brown = \lambda < 0, 0 > + \mu < 0, 0 > = 0 238c4762a1bSJed Brown NBC = <1,0> . <4\mu/\lambda (\lambda + \mu), 0> = 4\mu/\lambda (\lambda + \mu) 239c4762a1bSJed Brown 240c4762a1bSJed Brown u = x - 1/2 241c4762a1bSJed Brown v = 0 242c4762a1bSJed Brown w = 0 243c4762a1bSJed Brown \varepsilon = / x 0 0 \ 244c4762a1bSJed Brown | 0 0 0 | 245c4762a1bSJed Brown \ 0 0 0 / 246c4762a1bSJed Brown Tr(\varepsilon) = div u = x 247c4762a1bSJed Brown div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij} 248c4762a1bSJed Brown = \lambda \partial_j x + 2\mu < 1, 0, 0 > 249c4762a1bSJed Brown = \lambda < 1, 0, 0 > + \mu < 2, 0, 0 > 250c4762a1bSJed Brown */ 251d71ae5a4SJacob Faibussowitsch 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[]) 252d71ae5a4SJacob Faibussowitsch { 253c4762a1bSJed Brown const PetscReal N = -1.0; 254c4762a1bSJed Brown 255c4762a1bSJed Brown f0[0] = N; 256c4762a1bSJed Brown } 257c4762a1bSJed Brown 258d71ae5a4SJacob Faibussowitsch static PetscErrorCode uniform_strain_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 259d71ae5a4SJacob Faibussowitsch { 260c4762a1bSJed Brown const PetscReal eps_xx = 0.1; 261c4762a1bSJed Brown const PetscReal eps_xy = 0.3; 262c4762a1bSJed Brown const PetscReal eps_yy = 0.25; 263c4762a1bSJed Brown PetscInt d; 264c4762a1bSJed Brown 265c4762a1bSJed Brown u[0] = eps_xx * x[0] + eps_xy * x[1]; 266c4762a1bSJed Brown u[1] = eps_xy * x[0] + eps_yy * x[1]; 267c4762a1bSJed Brown for (d = 2; d < dim; ++d) u[d] = 0.0; 268*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 269c4762a1bSJed Brown } 270c4762a1bSJed Brown 271d71ae5a4SJacob Faibussowitsch 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[]) 272d71ae5a4SJacob Faibussowitsch { 273482dd828SMatthew G. Knepley const PetscInt Nc = dim; 274482dd828SMatthew G. Knepley PetscInt c; 275482dd828SMatthew G. Knepley 276482dd828SMatthew G. Knepley for (c = 0; c < Nc; ++c) f0[c] = u[c]; 277482dd828SMatthew G. Knepley } 278482dd828SMatthew G. Knepley 279d71ae5a4SJacob Faibussowitsch 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[]) 280d71ae5a4SJacob Faibussowitsch { 281c4762a1bSJed Brown const PetscInt Nc = dim; 282c4762a1bSJed Brown PetscInt c, d; 283c4762a1bSJed Brown 2849371c9d4SSatish Balay for (c = 0; c < Nc; ++c) 2859371c9d4SSatish Balay for (d = 0; d < dim; ++d) f1[c * dim + d] += u_x[c * dim + d]; 286c4762a1bSJed Brown } 287c4762a1bSJed Brown 288d71ae5a4SJacob Faibussowitsch 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[]) 289d71ae5a4SJacob Faibussowitsch { 290c4762a1bSJed Brown const PetscInt Nc = dim; 291c4762a1bSJed Brown const PetscReal mu = 1.0; 292c4762a1bSJed Brown const PetscReal lambda = 1.0; 293c4762a1bSJed Brown PetscInt c, d; 294c4762a1bSJed Brown 295c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 296c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 297c4762a1bSJed Brown f1[c * dim + d] += mu * (u_x[c * dim + d] + u_x[d * dim + c]); 298c4762a1bSJed Brown f1[c * dim + c] += lambda * u_x[d * dim + d]; 299c4762a1bSJed Brown } 300c4762a1bSJed Brown } 301c4762a1bSJed Brown } 302c4762a1bSJed Brown 303d71ae5a4SJacob Faibussowitsch 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[]) 304d71ae5a4SJacob Faibussowitsch { 305482dd828SMatthew G. Knepley const PetscInt Nc = dim; 306482dd828SMatthew G. Knepley PetscInt c; 307482dd828SMatthew G. Knepley 308482dd828SMatthew G. Knepley for (c = 0; c < Nc; ++c) g0[c * Nc + c] = 1.0; 309482dd828SMatthew G. Knepley } 310482dd828SMatthew G. Knepley 311d71ae5a4SJacob Faibussowitsch 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[]) 312d71ae5a4SJacob Faibussowitsch { 313c4762a1bSJed Brown const PetscInt Nc = dim; 314c4762a1bSJed Brown PetscInt c, d; 315c4762a1bSJed Brown 316c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 317ad540459SPierre Jolivet for (d = 0; d < dim; ++d) g3[((c * Nc + c) * dim + d) * dim + d] = 1.0; 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 */ 327d71ae5a4SJacob Faibussowitsch 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[]) 328d71ae5a4SJacob Faibussowitsch { 329c4762a1bSJed Brown const PetscInt Nc = dim; 330c4762a1bSJed Brown const PetscReal mu = 1.0; 331c4762a1bSJed Brown const PetscReal lambda = 1.0; 332c4762a1bSJed Brown PetscInt c, d; 333c4762a1bSJed Brown 334c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 335c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 336c4762a1bSJed Brown g3[((c * Nc + c) * dim + d) * dim + d] += mu; 337c4762a1bSJed Brown g3[((c * Nc + d) * dim + d) * dim + c] += mu; 338c4762a1bSJed Brown g3[((c * Nc + d) * dim + c) * dim + d] += lambda; 339c4762a1bSJed Brown } 340c4762a1bSJed Brown } 341c4762a1bSJed Brown } 342c4762a1bSJed Brown 343d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 344d71ae5a4SJacob Faibussowitsch { 345482dd828SMatthew G. Knepley PetscInt sol = 0, def = 0; 346c4762a1bSJed Brown 347c4762a1bSJed Brown PetscFunctionBeginUser; 348482dd828SMatthew G. Knepley options->deform = DEFORM_NONE; 349c4762a1bSJed Brown options->solType = SOL_VLAP_QUADRATIC; 350c4762a1bSJed Brown options->useNearNullspace = PETSC_TRUE; 3519566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(options->dmType, DMPLEX, 256)); 352c4762a1bSJed Brown 353d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Linear Elasticity Problem Options", "DMPLEX"); 3549566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-deform_type", "Type of domain deformation", "ex17.c", deformTypes, NUM_DEFORM_TYPES, deformTypes[options->deform], &def, NULL)); 355482dd828SMatthew G. Knepley options->deform = (DeformType)def; 3569566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-sol_type", "Type of exact solution", "ex17.c", solutionTypes, NUM_SOLUTION_TYPES, solutionTypes[options->solType], &sol, NULL)); 357c4762a1bSJed Brown options->solType = (SolutionType)sol; 3589566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-near_nullspace", "Use the rigid body modes as an AMG near nullspace", "ex17.c", options->useNearNullspace, &options->useNearNullspace, NULL)); 3599566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-dm_type", "Convert DMPlex to another format", "ex17.c", DMList, options->dmType, options->dmType, 256, NULL)); 360d0609cedSBarry Smith PetscOptionsEnd(); 361*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 362c4762a1bSJed Brown } 363c4762a1bSJed Brown 364d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *ctx) 365d71ae5a4SJacob Faibussowitsch { 366d12cd81dSMatthew G. Knepley PetscBag bag; 367d12cd81dSMatthew G. Knepley Parameter *p; 368d12cd81dSMatthew G. Knepley 369d12cd81dSMatthew G. Knepley PetscFunctionBeginUser; 370d12cd81dSMatthew G. Knepley /* setup PETSc parameter bag */ 3719566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(ctx->bag, (void **)&p)); 3729566063dSJacob Faibussowitsch PetscCall(PetscBagSetName(ctx->bag, "par", "Elastic Parameters")); 373d12cd81dSMatthew G. Knepley bag = ctx->bag; 3749566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterScalar(bag, &p->mu, 1.0, "mu", "Shear Modulus, Pa")); 3759566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterScalar(bag, &p->lambda, 1.0, "lambda", "Lame's first parameter, Pa")); 3769566063dSJacob Faibussowitsch PetscCall(PetscBagSetFromOptions(bag)); 377d12cd81dSMatthew G. Knepley { 378d12cd81dSMatthew G. Knepley PetscViewer viewer; 379d12cd81dSMatthew G. Knepley PetscViewerFormat format; 380d12cd81dSMatthew G. Knepley PetscBool flg; 381d12cd81dSMatthew G. Knepley 3829566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(comm, NULL, NULL, "-param_view", &viewer, &format, &flg)); 383d12cd81dSMatthew G. Knepley if (flg) { 3849566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, format)); 3859566063dSJacob Faibussowitsch PetscCall(PetscBagView(bag, viewer)); 3869566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 3879566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 3889566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 389d12cd81dSMatthew G. Knepley } 390d12cd81dSMatthew G. Knepley } 391*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 392d12cd81dSMatthew G. Knepley } 393d12cd81dSMatthew G. Knepley 394d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexDistortGeometry(DM dm) 395d71ae5a4SJacob Faibussowitsch { 396482dd828SMatthew G. Knepley DM cdm; 397482dd828SMatthew G. Knepley DMLabel label; 398482dd828SMatthew G. Knepley Vec coordinates; 399482dd828SMatthew G. Knepley PetscScalar *coords; 400482dd828SMatthew G. Knepley PetscReal mid = 0.5; 401482dd828SMatthew G. Knepley PetscInt cdim, d, vStart, vEnd, v; 402482dd828SMatthew G. Knepley 403482dd828SMatthew G. Knepley PetscFunctionBeginUser; 4049566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 4059566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim)); 4069566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 4079566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 4089566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 4099566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(coordinates, &coords)); 410482dd828SMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 411482dd828SMatthew G. Knepley PetscScalar *pcoords, shift; 412482dd828SMatthew G. Knepley PetscInt val; 413482dd828SMatthew G. Knepley 4149566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(label, v, &val)); 415482dd828SMatthew G. Knepley if (val >= 0) continue; 4169566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRef(cdm, v, coords, &pcoords)); 417482dd828SMatthew G. Knepley shift = 0.2 * PetscAbsScalar(pcoords[0] - mid); 418482dd828SMatthew G. Knepley shift = PetscRealPart(pcoords[0]) > mid ? shift : -shift; 419482dd828SMatthew G. Knepley for (d = 1; d < cdim; ++d) pcoords[d] += shift; 420482dd828SMatthew G. Knepley } 4219566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(coordinates, &coords)); 422*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 423482dd828SMatthew G. Knepley } 424482dd828SMatthew G. Knepley 425d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 426d71ae5a4SJacob Faibussowitsch { 427c4762a1bSJed Brown PetscFunctionBeginUser; 4289566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 4299566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 4309566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 431482dd828SMatthew G. Knepley switch (user->deform) { 432d71ae5a4SJacob Faibussowitsch case DEFORM_NONE: 433d71ae5a4SJacob Faibussowitsch break; 434d71ae5a4SJacob Faibussowitsch case DEFORM_SHEAR: 435d71ae5a4SJacob Faibussowitsch PetscCall(DMPlexShearGeometry(*dm, DM_X, NULL)); 436d71ae5a4SJacob Faibussowitsch break; 437d71ae5a4SJacob Faibussowitsch case DEFORM_STEP: 438d71ae5a4SJacob Faibussowitsch PetscCall(DMPlexDistortGeometry(*dm)); 439d71ae5a4SJacob Faibussowitsch break; 440d71ae5a4SJacob Faibussowitsch default: 441d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid deformation type: %s (%d)", deformTypes[PetscMin(user->deform, NUM_DEFORM_TYPES)], user->deform); 442482dd828SMatthew G. Knepley } 4439566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*dm, user)); 4449566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 445*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 446c4762a1bSJed Brown } 447c4762a1bSJed Brown 448d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 449d71ae5a4SJacob Faibussowitsch { 450c4762a1bSJed Brown PetscErrorCode (*exact)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 451d12cd81dSMatthew G. Knepley Parameter *param; 452482dd828SMatthew G. Knepley PetscDS ds; 45345480ffeSMatthew G. Knepley PetscWeakForm wf; 45445480ffeSMatthew G. Knepley DMLabel label; 45545480ffeSMatthew G. Knepley PetscInt id, bd; 456c4762a1bSJed Brown PetscInt dim; 457c4762a1bSJed Brown 458c4762a1bSJed Brown PetscFunctionBeginUser; 4599566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 4609566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf)); 4619566063dSJacob Faibussowitsch PetscCall(PetscDSGetSpatialDimension(ds, &dim)); 4629566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **)¶m)); 463c4762a1bSJed Brown switch (user->solType) { 464482dd828SMatthew G. Knepley case SOL_MASS_QUADRATIC: 4659566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_mass_u, NULL)); 4669566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_mass_uu, NULL, NULL, NULL)); 4679566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, 0, 0, 1, f0_mass_quadratic_u, 0, NULL)); 468c4762a1bSJed Brown switch (dim) { 469d71ae5a4SJacob Faibussowitsch case 2: 470d71ae5a4SJacob Faibussowitsch exact = quadratic_2d_u; 471d71ae5a4SJacob Faibussowitsch break; 472d71ae5a4SJacob Faibussowitsch case 3: 473d71ae5a4SJacob Faibussowitsch exact = quadratic_3d_u; 474d71ae5a4SJacob Faibussowitsch break; 475d71ae5a4SJacob Faibussowitsch default: 476d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim); 477482dd828SMatthew G. Knepley } 478482dd828SMatthew G. Knepley break; 479482dd828SMatthew G. Knepley case SOL_VLAP_QUADRATIC: 4809566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_vlap_quadratic_u, f1_vlap_u)); 4819566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_vlap_uu)); 482482dd828SMatthew G. Knepley switch (dim) { 483d71ae5a4SJacob Faibussowitsch case 2: 484d71ae5a4SJacob Faibussowitsch exact = quadratic_2d_u; 485d71ae5a4SJacob Faibussowitsch break; 486d71ae5a4SJacob Faibussowitsch case 3: 487d71ae5a4SJacob Faibussowitsch exact = quadratic_3d_u; 488d71ae5a4SJacob Faibussowitsch break; 489d71ae5a4SJacob Faibussowitsch default: 490d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim); 491c4762a1bSJed Brown } 492c4762a1bSJed Brown break; 493c4762a1bSJed Brown case SOL_ELAS_QUADRATIC: 4949566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_elas_quadratic_u, f1_elas_u)); 4959566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu)); 496c4762a1bSJed Brown switch (dim) { 497d71ae5a4SJacob Faibussowitsch case 2: 498d71ae5a4SJacob Faibussowitsch exact = quadratic_2d_u; 499d71ae5a4SJacob Faibussowitsch break; 500d71ae5a4SJacob Faibussowitsch case 3: 501d71ae5a4SJacob Faibussowitsch exact = quadratic_3d_u; 502d71ae5a4SJacob Faibussowitsch break; 503d71ae5a4SJacob Faibussowitsch default: 504d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim); 505c4762a1bSJed Brown } 506c4762a1bSJed Brown break; 507c4762a1bSJed Brown case SOL_VLAP_TRIG: 5089566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_vlap_trig_u, f1_vlap_u)); 5099566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_vlap_uu)); 510c4762a1bSJed Brown switch (dim) { 511d71ae5a4SJacob Faibussowitsch case 2: 512d71ae5a4SJacob Faibussowitsch exact = trig_2d_u; 513d71ae5a4SJacob Faibussowitsch break; 514d71ae5a4SJacob Faibussowitsch case 3: 515d71ae5a4SJacob Faibussowitsch exact = trig_3d_u; 516d71ae5a4SJacob Faibussowitsch break; 517d71ae5a4SJacob Faibussowitsch default: 518d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim); 519c4762a1bSJed Brown } 520c4762a1bSJed Brown break; 521c4762a1bSJed Brown case SOL_ELAS_TRIG: 5229566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_elas_trig_u, f1_elas_u)); 5239566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu)); 524c4762a1bSJed Brown switch (dim) { 525d71ae5a4SJacob Faibussowitsch case 2: 526d71ae5a4SJacob Faibussowitsch exact = trig_2d_u; 527d71ae5a4SJacob Faibussowitsch break; 528d71ae5a4SJacob Faibussowitsch case 3: 529d71ae5a4SJacob Faibussowitsch exact = trig_3d_u; 530d71ae5a4SJacob Faibussowitsch break; 531d71ae5a4SJacob Faibussowitsch default: 532d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim); 533c4762a1bSJed Brown } 534c4762a1bSJed Brown break; 535c4762a1bSJed Brown case SOL_ELAS_AXIAL_DISP: 5369566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, NULL, f1_elas_u)); 5379566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu)); 53845480ffeSMatthew G. Knepley id = dim == 3 ? 5 : 2; 5399566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 5409566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "right", label, 1, &id, 0, 0, NULL, (void (*)(void))NULL, NULL, user, &bd)); 5419566063dSJacob Faibussowitsch PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 5429566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_elas_axial_disp_bd_u, 0, NULL)); 543c4762a1bSJed Brown exact = axial_disp_u; 544c4762a1bSJed Brown break; 545c4762a1bSJed Brown case SOL_ELAS_UNIFORM_STRAIN: 5469566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, NULL, f1_elas_u)); 5479566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu)); 548c4762a1bSJed Brown exact = uniform_strain_u; 549c4762a1bSJed Brown break; 550d12cd81dSMatthew G. Knepley case SOL_ELAS_GE: 5519566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, NULL, f1_elas_u)); 5529566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu)); 553d12cd81dSMatthew G. Knepley exact = zero; /* No exact solution available */ 554d12cd81dSMatthew G. Knepley break; 555d71ae5a4SJacob Faibussowitsch default: 556d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid solution type: %s (%d)", solutionTypes[PetscMin(user->solType, NUM_SOLUTION_TYPES)], user->solType); 557c4762a1bSJed Brown } 5589566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, exact, user)); 5599566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 560c4762a1bSJed Brown if (user->solType == SOL_ELAS_AXIAL_DISP) { 561c4762a1bSJed Brown PetscInt cmp; 562c4762a1bSJed Brown 563c4762a1bSJed Brown id = dim == 3 ? 6 : 4; 564c4762a1bSJed Brown cmp = 0; 5659566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "left", label, 1, &id, 0, 1, &cmp, (void (*)(void))zero, NULL, user, NULL)); 566c4762a1bSJed Brown cmp = dim == 3 ? 2 : 1; 567c4762a1bSJed Brown id = dim == 3 ? 1 : 1; 5689566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "bottom", label, 1, &id, 0, 1, &cmp, (void (*)(void))zero, NULL, user, NULL)); 569c4762a1bSJed Brown if (dim == 3) { 570c4762a1bSJed Brown cmp = 1; 571c4762a1bSJed Brown id = 3; 5729566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "front", label, 1, &id, 0, 1, &cmp, (void (*)(void))zero, NULL, user, NULL)); 573c4762a1bSJed Brown } 574d12cd81dSMatthew G. Knepley } else if (user->solType == SOL_ELAS_GE) { 575d12cd81dSMatthew G. Knepley PetscInt cmp; 576d12cd81dSMatthew G. Knepley 577d12cd81dSMatthew G. Knepley id = dim == 3 ? 6 : 4; 5789566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "left", label, 1, &id, 0, 0, NULL, (void (*)(void))zero, NULL, user, NULL)); 579d12cd81dSMatthew G. Knepley id = dim == 3 ? 5 : 2; 580d12cd81dSMatthew G. Knepley cmp = 0; 5819566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "right", label, 1, &id, 0, 1, &cmp, (void (*)(void))ge_shift, NULL, user, NULL)); 582c4762a1bSJed Brown } else { 583c4762a1bSJed Brown id = 1; 5849566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))exact, NULL, user, NULL)); 585c4762a1bSJed Brown } 586d12cd81dSMatthew G. Knepley /* Setup constants */ 587d12cd81dSMatthew G. Knepley { 588d12cd81dSMatthew G. Knepley PetscScalar constants[2]; 589d12cd81dSMatthew G. Knepley 590d12cd81dSMatthew G. Knepley constants[0] = param->mu; /* shear modulus, Pa */ 591d12cd81dSMatthew G. Knepley constants[1] = param->lambda; /* Lame's first parameter, Pa */ 5929566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(ds, 2, constants)); 593d12cd81dSMatthew G. Knepley } 594*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 595c4762a1bSJed Brown } 596c4762a1bSJed Brown 597d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateElasticityNullSpace(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullspace) 598d71ae5a4SJacob Faibussowitsch { 599c4762a1bSJed Brown PetscFunctionBegin; 6009566063dSJacob Faibussowitsch PetscCall(DMPlexCreateRigidBody(dm, origField, nullspace)); 601*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 602c4762a1bSJed Brown } 603c4762a1bSJed Brown 604d71ae5a4SJacob Faibussowitsch PetscErrorCode SetupFE(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), void *ctx) 605d71ae5a4SJacob Faibussowitsch { 606c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 607c4762a1bSJed Brown DM cdm = dm; 608c4762a1bSJed Brown PetscFE fe; 609c4762a1bSJed Brown char prefix[PETSC_MAX_PATH_LEN]; 61030602db0SMatthew G. Knepley DMPolytopeType ct; 61130602db0SMatthew G. Knepley PetscBool simplex; 61230602db0SMatthew G. Knepley PetscInt dim, cStart; 613c4762a1bSJed Brown 614c4762a1bSJed Brown PetscFunctionBegin; 615c4762a1bSJed Brown /* Create finite element */ 6169566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 6179566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 6189566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 61930602db0SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE; 6209566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name)); 6219566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, dim, simplex, name ? prefix : NULL, -1, &fe)); 6229566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, name)); 623c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 6249566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 6259566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 6269566063dSJacob Faibussowitsch PetscCall((*setup)(dm, user)); 627c4762a1bSJed Brown while (cdm) { 6289566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 6299566063dSJacob Faibussowitsch if (user->useNearNullspace) PetscCall(DMSetNearNullSpaceConstructor(cdm, 0, CreateElasticityNullSpace)); 630c4762a1bSJed Brown /* TODO: Check whether the boundary of coarse meshes is marked */ 6319566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 632c4762a1bSJed Brown } 6339566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 634*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 635c4762a1bSJed Brown } 636c4762a1bSJed Brown 637d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 638d71ae5a4SJacob Faibussowitsch { 639c4762a1bSJed Brown DM dm; /* Problem specification */ 640c4762a1bSJed Brown SNES snes; /* Nonlinear solver */ 641c4762a1bSJed Brown Vec u; /* Solutions */ 642c4762a1bSJed Brown AppCtx user; /* User-defined work context */ 643c4762a1bSJed Brown 644327415f7SBarry Smith PetscFunctionBeginUser; 6459566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 6469566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 6479566063dSJacob Faibussowitsch PetscCall(PetscBagCreate(PETSC_COMM_SELF, sizeof(Parameter), &user.bag)); 6489566063dSJacob Faibussowitsch PetscCall(SetupParameters(PETSC_COMM_WORLD, &user)); 649c4762a1bSJed Brown /* Primal system */ 6509566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 6519566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 6529566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 6539566063dSJacob Faibussowitsch PetscCall(SetupFE(dm, "displacement", SetupPrimalProblem, &user)); 6549566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 6559566063dSJacob Faibussowitsch PetscCall(VecSet(u, 0.0)); 6569566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "displacement")); 6579566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 6589566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 6599566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u)); 6609566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 6619566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &u)); 6629566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-displacement_view")); 663c4762a1bSJed Brown /* Cleanup */ 6649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 6659566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 6669566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 6679566063dSJacob Faibussowitsch PetscCall(PetscBagDestroy(&user.bag)); 6689566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 669b122ec5aSJacob Faibussowitsch return 0; 670c4762a1bSJed Brown } 671c4762a1bSJed Brown 672c4762a1bSJed Brown /*TEST 673c4762a1bSJed Brown 67430602db0SMatthew G. Knepley testset: 67530602db0SMatthew G. Knepley args: -dm_plex_box_faces 1,1,1 67630602db0SMatthew G. Knepley 677c4762a1bSJed Brown test: 678c4762a1bSJed Brown suffix: 2d_p1_quad_vlap 679c4762a1bSJed Brown requires: triangle 680c4762a1bSJed Brown args: -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate 681c4762a1bSJed Brown test: 682c4762a1bSJed Brown suffix: 2d_p2_quad_vlap 683c4762a1bSJed Brown requires: triangle 684c4762a1bSJed Brown args: -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001 685c4762a1bSJed Brown test: 686c4762a1bSJed Brown suffix: 2d_p3_quad_vlap 687c4762a1bSJed Brown requires: triangle 688c4762a1bSJed Brown args: -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001 689c4762a1bSJed Brown test: 690c4762a1bSJed Brown suffix: 2d_q1_quad_vlap 69130602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate 692c4762a1bSJed Brown test: 693c4762a1bSJed Brown suffix: 2d_q2_quad_vlap 69430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001 695c4762a1bSJed Brown test: 696c4762a1bSJed Brown suffix: 2d_q3_quad_vlap 697c4762a1bSJed Brown requires: !single 69830602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001 699c4762a1bSJed Brown test: 700c4762a1bSJed Brown suffix: 2d_p1_quad_elas 701c4762a1bSJed Brown requires: triangle 702c4762a1bSJed Brown args: -sol_type elas_quad -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate 703c4762a1bSJed Brown test: 704c4762a1bSJed Brown suffix: 2d_p2_quad_elas 705c4762a1bSJed Brown requires: triangle 706c4762a1bSJed Brown args: -sol_type elas_quad -displacement_petscspace_degree 2 -dmsnes_check .0001 707c4762a1bSJed Brown test: 708c4762a1bSJed Brown suffix: 2d_p3_quad_elas 709c4762a1bSJed Brown requires: triangle 710c4762a1bSJed Brown args: -sol_type elas_quad -displacement_petscspace_degree 3 -dmsnes_check .0001 711c4762a1bSJed Brown test: 712c4762a1bSJed Brown suffix: 2d_q1_quad_elas 71330602db0SMatthew 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 714c4762a1bSJed Brown test: 715c4762a1bSJed Brown suffix: 2d_q1_quad_elas_shear 716482dd828SMatthew 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 717c4762a1bSJed Brown test: 718c4762a1bSJed Brown suffix: 2d_q2_quad_elas 71930602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dmsnes_check .0001 720c4762a1bSJed Brown test: 721c4762a1bSJed Brown suffix: 2d_q2_quad_elas_shear 722482dd828SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_simplex 0 -deform_type shear -displacement_petscspace_degree 2 -dmsnes_check 723c4762a1bSJed Brown test: 724c4762a1bSJed Brown suffix: 2d_q3_quad_elas 72530602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dmsnes_check .0001 726c4762a1bSJed Brown test: 727c4762a1bSJed Brown suffix: 2d_q3_quad_elas_shear 728c4762a1bSJed Brown requires: !single 729482dd828SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_simplex 0 -deform_type shear -displacement_petscspace_degree 3 -dmsnes_check 730c4762a1bSJed Brown 731c4762a1bSJed Brown test: 732c4762a1bSJed Brown suffix: 3d_p1_quad_vlap 733c4762a1bSJed Brown requires: ctetgen 73430602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_refine 1 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate 735c4762a1bSJed Brown test: 736c4762a1bSJed Brown suffix: 3d_p2_quad_vlap 737c4762a1bSJed Brown requires: ctetgen 73830602db0SMatthew G. Knepley args: -dm_plex_dim 3 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001 739c4762a1bSJed Brown test: 740c4762a1bSJed Brown suffix: 3d_p3_quad_vlap 741c4762a1bSJed Brown requires: ctetgen 74230602db0SMatthew G. Knepley args: -dm_plex_dim 3 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001 743c4762a1bSJed Brown test: 744c4762a1bSJed Brown suffix: 3d_q1_quad_vlap 74530602db0SMatthew 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 746c4762a1bSJed Brown test: 747c4762a1bSJed Brown suffix: 3d_q2_quad_vlap 74830602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001 749c4762a1bSJed Brown test: 750c4762a1bSJed Brown suffix: 3d_q3_quad_vlap 75130602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001 752c4762a1bSJed Brown test: 753c4762a1bSJed Brown suffix: 3d_p1_quad_elas 754c4762a1bSJed Brown requires: ctetgen 75530602db0SMatthew 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 756c4762a1bSJed Brown test: 757c4762a1bSJed Brown suffix: 3d_p2_quad_elas 758c4762a1bSJed Brown requires: ctetgen 75930602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_dim 3 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001 760c4762a1bSJed Brown test: 761c4762a1bSJed Brown suffix: 3d_p3_quad_elas 762c4762a1bSJed Brown requires: ctetgen 76330602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_dim 3 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001 764c4762a1bSJed Brown test: 765c4762a1bSJed Brown suffix: 3d_q1_quad_elas 76630602db0SMatthew 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 767c4762a1bSJed Brown test: 768c4762a1bSJed Brown suffix: 3d_q2_quad_elas 76930602db0SMatthew 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 770c4762a1bSJed Brown test: 771c4762a1bSJed Brown suffix: 3d_q3_quad_elas 772c4762a1bSJed Brown requires: !single 77330602db0SMatthew 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 774c4762a1bSJed Brown 775c4762a1bSJed Brown test: 776c4762a1bSJed Brown suffix: 2d_p1_trig_vlap 777c4762a1bSJed Brown requires: triangle 778c4762a1bSJed Brown args: -sol_type vlap_trig -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 779c4762a1bSJed Brown test: 780c4762a1bSJed Brown suffix: 2d_p2_trig_vlap 781c4762a1bSJed Brown requires: triangle 782c4762a1bSJed Brown args: -sol_type vlap_trig -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 783c4762a1bSJed Brown test: 784c4762a1bSJed Brown suffix: 2d_p3_trig_vlap 785c4762a1bSJed Brown requires: triangle 786c4762a1bSJed Brown args: -sol_type vlap_trig -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 787c4762a1bSJed Brown test: 788c4762a1bSJed Brown suffix: 2d_q1_trig_vlap 78930602db0SMatthew 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 790c4762a1bSJed Brown test: 791c4762a1bSJed Brown suffix: 2d_q2_trig_vlap 79230602db0SMatthew 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 793c4762a1bSJed Brown test: 794c4762a1bSJed Brown suffix: 2d_q3_trig_vlap 79530602db0SMatthew 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 796c4762a1bSJed Brown test: 797c4762a1bSJed Brown suffix: 2d_p1_trig_elas 798c4762a1bSJed Brown requires: triangle 799c4762a1bSJed Brown args: -sol_type elas_trig -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 800c4762a1bSJed Brown test: 801c4762a1bSJed Brown suffix: 2d_p2_trig_elas 802c4762a1bSJed Brown requires: triangle 803c4762a1bSJed Brown args: -sol_type elas_trig -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 804c4762a1bSJed Brown test: 805c4762a1bSJed Brown suffix: 2d_p3_trig_elas 806c4762a1bSJed Brown requires: triangle 807c4762a1bSJed Brown args: -sol_type elas_trig -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 808c4762a1bSJed Brown test: 809c4762a1bSJed Brown suffix: 2d_q1_trig_elas 81030602db0SMatthew 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 811c4762a1bSJed Brown test: 812c4762a1bSJed Brown suffix: 2d_q1_trig_elas_shear 813482dd828SMatthew 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 814c4762a1bSJed Brown test: 815c4762a1bSJed Brown suffix: 2d_q2_trig_elas 81630602db0SMatthew 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 817c4762a1bSJed Brown test: 818c4762a1bSJed Brown suffix: 2d_q2_trig_elas_shear 819482dd828SMatthew 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 820c4762a1bSJed Brown test: 821c4762a1bSJed Brown suffix: 2d_q3_trig_elas 82230602db0SMatthew 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 823c4762a1bSJed Brown test: 824c4762a1bSJed Brown suffix: 2d_q3_trig_elas_shear 825482dd828SMatthew 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 826c4762a1bSJed Brown 827c4762a1bSJed Brown test: 828c4762a1bSJed Brown suffix: 3d_p1_trig_vlap 829c4762a1bSJed Brown requires: ctetgen 83030602db0SMatthew 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 831c4762a1bSJed Brown test: 832c4762a1bSJed Brown suffix: 3d_p2_trig_vlap 833c4762a1bSJed Brown requires: ctetgen 83430602db0SMatthew 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 835c4762a1bSJed Brown test: 836c4762a1bSJed Brown suffix: 3d_p3_trig_vlap 837c4762a1bSJed Brown requires: ctetgen 83830602db0SMatthew 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 839c4762a1bSJed Brown test: 840c4762a1bSJed Brown suffix: 3d_q1_trig_vlap 84130602db0SMatthew 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 842c4762a1bSJed Brown test: 843c4762a1bSJed Brown suffix: 3d_q2_trig_vlap 84430602db0SMatthew 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 845c4762a1bSJed Brown test: 846c4762a1bSJed Brown suffix: 3d_q3_trig_vlap 847c4762a1bSJed Brown requires: !__float128 84830602db0SMatthew 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 849c4762a1bSJed Brown test: 850c4762a1bSJed Brown suffix: 3d_p1_trig_elas 851c4762a1bSJed Brown requires: ctetgen 85230602db0SMatthew 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 853c4762a1bSJed Brown test: 854c4762a1bSJed Brown suffix: 3d_p2_trig_elas 855c4762a1bSJed Brown requires: ctetgen 85630602db0SMatthew 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 857c4762a1bSJed Brown test: 858c4762a1bSJed Brown suffix: 3d_p3_trig_elas 859c4762a1bSJed Brown requires: ctetgen 86030602db0SMatthew 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 861c4762a1bSJed Brown test: 862c4762a1bSJed Brown suffix: 3d_q1_trig_elas 86330602db0SMatthew 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 864c4762a1bSJed Brown test: 865c4762a1bSJed Brown suffix: 3d_q2_trig_elas 86630602db0SMatthew 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 867c4762a1bSJed Brown test: 868c4762a1bSJed Brown suffix: 3d_q3_trig_elas 869c4762a1bSJed Brown requires: !__float128 87030602db0SMatthew 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 871c4762a1bSJed Brown 872c4762a1bSJed Brown test: 873c4762a1bSJed Brown suffix: 2d_p1_axial_elas 874c4762a1bSJed Brown requires: triangle 875c4762a1bSJed Brown args: -sol_type elas_axial_disp -displacement_petscspace_degree 1 -dm_plex_separate_marker -dm_refine 2 -dmsnes_check .0001 -pc_type lu 876c4762a1bSJed Brown test: 877c4762a1bSJed Brown suffix: 2d_p2_axial_elas 878c4762a1bSJed Brown requires: triangle 879c4762a1bSJed Brown args: -sol_type elas_axial_disp -displacement_petscspace_degree 2 -dm_plex_separate_marker -dmsnes_check .0001 -pc_type lu 880c4762a1bSJed Brown test: 881c4762a1bSJed Brown suffix: 2d_p3_axial_elas 882c4762a1bSJed Brown requires: triangle 883c4762a1bSJed Brown args: -sol_type elas_axial_disp -displacement_petscspace_degree 3 -dm_plex_separate_marker -dmsnes_check .0001 -pc_type lu 884c4762a1bSJed Brown test: 885c4762a1bSJed Brown suffix: 2d_q1_axial_elas 88630602db0SMatthew 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 887c4762a1bSJed Brown test: 888c4762a1bSJed Brown suffix: 2d_q2_axial_elas 88930602db0SMatthew 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 890c4762a1bSJed Brown test: 891c4762a1bSJed Brown suffix: 2d_q3_axial_elas 89230602db0SMatthew 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 893c4762a1bSJed Brown 894c4762a1bSJed Brown test: 895c4762a1bSJed Brown suffix: 2d_p1_uniform_elas 896c4762a1bSJed Brown requires: triangle 897c4762a1bSJed Brown args: -sol_type elas_uniform_strain -displacement_petscspace_degree 1 -dm_refine 2 -dmsnes_check .0001 -pc_type lu 898c4762a1bSJed Brown test: 899c4762a1bSJed Brown suffix: 2d_p2_uniform_elas 900c4762a1bSJed Brown requires: triangle 901c4762a1bSJed Brown args: -sol_type elas_uniform_strain -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001 -pc_type lu 902c4762a1bSJed Brown test: 903c4762a1bSJed Brown suffix: 2d_p3_uniform_elas 904c4762a1bSJed Brown requires: triangle 905c4762a1bSJed Brown args: -sol_type elas_uniform_strain -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001 -pc_type lu 906c4762a1bSJed Brown test: 907c4762a1bSJed Brown suffix: 2d_q1_uniform_elas 90830602db0SMatthew 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 909c4762a1bSJed Brown test: 910c4762a1bSJed Brown suffix: 2d_q2_uniform_elas 911c4762a1bSJed Brown requires: !single 91230602db0SMatthew 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 913c4762a1bSJed Brown test: 914c4762a1bSJed Brown suffix: 2d_q3_uniform_elas 915c4762a1bSJed Brown requires: !single 91630602db0SMatthew 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 917482dd828SMatthew G. Knepley test: 918482dd828SMatthew G. Knepley suffix: 2d_p1_uniform_elas_step 919482dd828SMatthew G. Knepley requires: triangle 920482dd828SMatthew G. Knepley args: -sol_type elas_uniform_strain -deform_type step -displacement_petscspace_degree 1 -dm_refine 2 -dmsnes_check .0001 -pc_type lu 921482dd828SMatthew G. Knepley 922482dd828SMatthew G. Knepley testset: 923482dd828SMatthew 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 924482dd828SMatthew G. Knepley 925482dd828SMatthew G. Knepley test: 926482dd828SMatthew G. Knepley suffix: 2d_q1_uniform_elas_step 927482dd828SMatthew G. Knepley args: -sol_type elas_uniform_strain -dm_refine 2 928482dd828SMatthew G. Knepley test: 929482dd828SMatthew G. Knepley suffix: 2d_q1_quad_vlap_step 930482dd828SMatthew G. Knepley args: 931482dd828SMatthew G. Knepley test: 932482dd828SMatthew G. Knepley suffix: 2d_q2_quad_vlap_step 933482dd828SMatthew G. Knepley args: -displacement_petscspace_degree 2 934482dd828SMatthew G. Knepley test: 935482dd828SMatthew G. Knepley suffix: 2d_q1_quad_mass_step 936482dd828SMatthew G. Knepley args: -sol_type mass_quad 937c4762a1bSJed Brown 938d12cd81dSMatthew G. Knepley testset: 939908b9b43SStefano Zampini filter: grep -v "variant HERMITIAN" 940d12cd81dSMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower -5,-5,-0.25 -dm_plex_box_upper 5,5,0.25 \ 941d12cd81dSMatthew G. Knepley -dm_plex_box_faces 5,5,2 -dm_plex_separate_marker -dm_refine 0 -petscpartitioner_type simple \ 942908b9b43SStefano Zampini -sol_type elas_ge 943908b9b43SStefano Zampini 944908b9b43SStefano Zampini test: 945908b9b43SStefano Zampini suffix: ge_q1_0 946908b9b43SStefano Zampini args: -displacement_petscspace_degree 1 \ 947d12cd81dSMatthew G. Knepley -snes_max_it 2 -snes_rtol 1.e-10 \ 948d12cd81dSMatthew G. Knepley -ksp_type cg -ksp_rtol 1.e-10 -ksp_max_it 100 -ksp_norm_type unpreconditioned \ 949d12cd81dSMatthew G. Knepley -pc_type gamg -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 \ 950d12cd81dSMatthew G. Knepley -pc_gamg_coarse_eq_limit 10 -pc_gamg_reuse_interpolation true \ 951bae903cbSmarkadams4 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 \ 952d12cd81dSMatthew G. Knepley -mg_levels_ksp_max_it 2 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.1 -mg_levels_pc_type jacobi \ 953d12cd81dSMatthew G. Knepley -matptap_via scalable 954d12cd81dSMatthew G. Knepley test: 955908b9b43SStefano Zampini nsize: 5 956908b9b43SStefano Zampini suffix: ge_q1_gdsw 957908b9b43SStefano Zampini args: -snes_max_it 1 -ksp_type cg -ksp_norm_type natural -displacement_petscspace_degree 1 -snes_monitor_short -ksp_monitor_short -pc_type mg -pc_mg_adapt_interp_coarse_space gdsw -pc_mg_levels 2 -pc_mg_galerkin -mg_levels_pc_type bjacobi -mg_levels_esteig_ksp_type cg -mg_levels_sub_pc_type icc -mg_coarse_redundant_pc_type cholesky -ksp_view 958d12cd81dSMatthew G. Knepley 959c4762a1bSJed Brown TEST*/ 960