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 21d12cd81dSMatthew G. Knepley typedef enum {SOL_VLAP_QUADRATIC, SOL_ELAS_QUADRATIC, SOL_VLAP_TRIG, SOL_ELAS_TRIG, SOL_ELAS_AXIAL_DISP, SOL_ELAS_UNIFORM_STRAIN, SOL_ELAS_GE, SOL_MASS_QUADRATIC, NUM_SOLUTION_TYPES} SolutionType; 22d12cd81dSMatthew 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"}; 23482dd828SMatthew G. Knepley 24482dd828SMatthew G. Knepley typedef enum {DEFORM_NONE, DEFORM_SHEAR, DEFORM_STEP, NUM_DEFORM_TYPES} DeformType; 25482dd828SMatthew G. Knepley const char *deformTypes[NUM_DEFORM_TYPES+1] = {"none", "shear", "step", "unknown"}; 26c4762a1bSJed Brown 27c4762a1bSJed Brown typedef struct { 28d12cd81dSMatthew G. Knepley PetscScalar mu; /* shear modulus */ 29d12cd81dSMatthew G. Knepley PetscScalar lambda; /* Lame's first parameter */ 30d12cd81dSMatthew G. Knepley } Parameter; 31d12cd81dSMatthew G. Knepley 32d12cd81dSMatthew G. Knepley typedef struct { 33c4762a1bSJed Brown /* Domain and mesh definition */ 34c4762a1bSJed Brown char dmType[256]; /* DM type for the solve */ 35482dd828SMatthew G. Knepley DeformType deform; /* Domain deformation type */ 36c4762a1bSJed Brown /* Problem definition */ 37c4762a1bSJed Brown SolutionType solType; /* Type of exact solution */ 38d12cd81dSMatthew G. Knepley PetscBag bag; /* Problem parameters */ 39c4762a1bSJed Brown /* Solver definition */ 40c4762a1bSJed Brown PetscBool useNearNullspace; /* Use the rigid body modes as a near nullspace for AMG */ 41c4762a1bSJed Brown } AppCtx; 42c4762a1bSJed Brown 43c4762a1bSJed Brown static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 44c4762a1bSJed Brown { 45c4762a1bSJed Brown PetscInt d; 46c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[d] = 0.0; 47c4762a1bSJed Brown return 0; 48c4762a1bSJed Brown } 49c4762a1bSJed Brown 50d12cd81dSMatthew G. Knepley static PetscErrorCode ge_shift(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 51d12cd81dSMatthew G. Knepley { 52d12cd81dSMatthew G. Knepley PetscInt d; 53d12cd81dSMatthew G. Knepley u[0] = 0.1; 54d12cd81dSMatthew G. Knepley for (d = 1; d < dim; ++d) u[d] = 0.0; 55d12cd81dSMatthew G. Knepley return 0; 56d12cd81dSMatthew G. Knepley } 57d12cd81dSMatthew G. Knepley 58c4762a1bSJed Brown static PetscErrorCode quadratic_2d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 59c4762a1bSJed Brown { 60c4762a1bSJed Brown u[0] = x[0]*x[0]; 61c4762a1bSJed Brown u[1] = x[1]*x[1] - 2.0*x[0]*x[1]; 62c4762a1bSJed Brown return 0; 63c4762a1bSJed Brown } 64c4762a1bSJed Brown 65c4762a1bSJed Brown static PetscErrorCode quadratic_3d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 66c4762a1bSJed Brown { 67c4762a1bSJed Brown u[0] = x[0]*x[0]; 68c4762a1bSJed Brown u[1] = x[1]*x[1] - 2.0*x[0]*x[1]; 69c4762a1bSJed Brown u[2] = x[2]*x[2] - 2.0*x[1]*x[2]; 70c4762a1bSJed Brown return 0; 71c4762a1bSJed Brown } 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* 74c4762a1bSJed Brown u = x^2 75c4762a1bSJed Brown v = y^2 - 2xy 76c4762a1bSJed Brown Delta <u,v> - f = <2, 2> - <2, 2> 77c4762a1bSJed Brown 78c4762a1bSJed Brown u = x^2 79c4762a1bSJed Brown v = y^2 - 2xy 80c4762a1bSJed Brown w = z^2 - 2yz 81c4762a1bSJed Brown Delta <u,v,w> - f = <2, 2, 2> - <2, 2, 2> 82c4762a1bSJed Brown */ 83c4762a1bSJed Brown static void f0_vlap_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 84c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 85c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 86c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 87c4762a1bSJed Brown { 88c4762a1bSJed Brown PetscInt d; 89c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[d] += 2.0; 90c4762a1bSJed Brown } 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* 93c4762a1bSJed Brown u = x^2 94c4762a1bSJed Brown v = y^2 - 2xy 95c4762a1bSJed Brown \varepsilon = / 2x -y \ 96c4762a1bSJed Brown \ -y 2y - 2x / 97c4762a1bSJed Brown Tr(\varepsilon) = div u = 2y 98c4762a1bSJed Brown div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij} 99c4762a1bSJed Brown = \lambda \partial_j (2y) + 2\mu < 2-1, 2 > 100c4762a1bSJed Brown = \lambda < 0, 2 > + \mu < 2, 4 > 101c4762a1bSJed Brown 102c4762a1bSJed Brown u = x^2 103c4762a1bSJed Brown v = y^2 - 2xy 104c4762a1bSJed Brown w = z^2 - 2yz 105c4762a1bSJed Brown \varepsilon = / 2x -y 0 \ 106c4762a1bSJed Brown | -y 2y - 2x -z | 107c4762a1bSJed Brown \ 0 -z 2z - 2y/ 108c4762a1bSJed Brown Tr(\varepsilon) = div u = 2z 109c4762a1bSJed Brown div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij} 110c4762a1bSJed Brown = \lambda \partial_j (2z) + 2\mu < 2-1, 2-1, 2 > 111c4762a1bSJed Brown = \lambda < 0, 0, 2 > + \mu < 2, 2, 4 > 112c4762a1bSJed Brown */ 113c4762a1bSJed Brown static void f0_elas_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 114c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 115c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 116c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 117c4762a1bSJed Brown { 118c4762a1bSJed Brown const PetscReal mu = 1.0; 119c4762a1bSJed Brown const PetscReal lambda = 1.0; 120c4762a1bSJed Brown PetscInt d; 121c4762a1bSJed Brown 122c4762a1bSJed Brown for (d = 0; d < dim-1; ++d) f0[d] += 2.0*mu; 123c4762a1bSJed Brown f0[dim-1] += 2.0*lambda + 4.0*mu; 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126482dd828SMatthew G. Knepley static void f0_mass_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 127482dd828SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 128482dd828SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 129482dd828SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 130482dd828SMatthew G. Knepley { 131482dd828SMatthew G. Knepley if (dim == 2) { 132482dd828SMatthew G. Knepley f0[0] -= x[0]*x[0]; 133482dd828SMatthew G. Knepley f0[1] -= x[1]*x[1] - 2.0*x[0]*x[1]; 134482dd828SMatthew G. Knepley } else { 135482dd828SMatthew G. Knepley f0[0] -= x[0]*x[0]; 136482dd828SMatthew G. Knepley f0[1] -= x[1]*x[1] - 2.0*x[0]*x[1]; 137482dd828SMatthew G. Knepley f0[2] -= x[2]*x[2] - 2.0*x[1]*x[2]; 138482dd828SMatthew G. Knepley } 139482dd828SMatthew G. Knepley } 140482dd828SMatthew G. Knepley 141c4762a1bSJed Brown static PetscErrorCode trig_2d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 142c4762a1bSJed Brown { 143c4762a1bSJed Brown u[0] = PetscSinReal(2.0*PETSC_PI*x[0]); 144c4762a1bSJed Brown u[1] = PetscSinReal(2.0*PETSC_PI*x[1]) - 2.0*x[0]*x[1]; 145c4762a1bSJed Brown return 0; 146c4762a1bSJed Brown } 147c4762a1bSJed Brown 148c4762a1bSJed Brown static PetscErrorCode trig_3d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 149c4762a1bSJed Brown { 150c4762a1bSJed Brown u[0] = PetscSinReal(2.0*PETSC_PI*x[0]); 151c4762a1bSJed Brown u[1] = PetscSinReal(2.0*PETSC_PI*x[1]) - 2.0*x[0]*x[1]; 152c4762a1bSJed Brown u[2] = PetscSinReal(2.0*PETSC_PI*x[2]) - 2.0*x[1]*x[2]; 153c4762a1bSJed Brown return 0; 154c4762a1bSJed Brown } 155c4762a1bSJed Brown 156c4762a1bSJed Brown /* 157c4762a1bSJed Brown u = sin(2 pi x) 158c4762a1bSJed Brown v = sin(2 pi y) - 2xy 159c4762a1bSJed 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)> 160c4762a1bSJed Brown 161c4762a1bSJed Brown u = sin(2 pi x) 162c4762a1bSJed Brown v = sin(2 pi y) - 2xy 163c4762a1bSJed Brown w = sin(2 pi z) - 2yz 164c4762a1bSJed 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)> 165c4762a1bSJed Brown */ 166c4762a1bSJed Brown static void f0_vlap_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 167c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 168c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 169c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 170c4762a1bSJed Brown { 171c4762a1bSJed Brown PetscInt d; 172c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[d] += -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]); 173c4762a1bSJed Brown } 174c4762a1bSJed Brown 175c4762a1bSJed Brown /* 176c4762a1bSJed Brown u = sin(2 pi x) 177c4762a1bSJed Brown v = sin(2 pi y) - 2xy 178c4762a1bSJed Brown \varepsilon = / 2 pi cos(2 pi x) -y \ 179c4762a1bSJed Brown \ -y 2 pi cos(2 pi y) - 2x / 180c4762a1bSJed Brown Tr(\varepsilon) = div u = 2 pi (cos(2 pi x) + cos(2 pi y)) - 2 x 181c4762a1bSJed Brown div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij} 182c4762a1bSJed 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) > 183c4762a1bSJed 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) > 184c4762a1bSJed Brown 185c4762a1bSJed Brown u = sin(2 pi x) 186c4762a1bSJed Brown v = sin(2 pi y) - 2xy 187c4762a1bSJed Brown w = sin(2 pi z) - 2yz 188c4762a1bSJed Brown \varepsilon = / 2 pi cos(2 pi x) -y 0 \ 189c4762a1bSJed Brown | -y 2 pi cos(2 pi y) - 2x -z | 190c4762a1bSJed Brown \ 0 -z 2 pi cos(2 pi z) - 2y / 191c4762a1bSJed Brown Tr(\varepsilon) = div u = 2 pi (cos(2 pi x) + cos(2 pi y) + cos(2 pi z)) - 2 x - 2 y 192c4762a1bSJed Brown div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij} 193c4762a1bSJed 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) > 194c4762a1bSJed 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) > 195c4762a1bSJed Brown */ 196c4762a1bSJed Brown static void f0_elas_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 197c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 198c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 199c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 200c4762a1bSJed Brown { 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 209c4762a1bSJed Brown static PetscErrorCode axial_disp_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 210c4762a1bSJed Brown { 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; 219c4762a1bSJed Brown return 0; 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 */ 251c4762a1bSJed Brown static void f0_elas_axial_disp_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 252c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 253c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 254c4762a1bSJed Brown PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 255c4762a1bSJed Brown { 256c4762a1bSJed Brown const PetscReal N = -1.0; 257c4762a1bSJed Brown 258c4762a1bSJed Brown f0[0] = N; 259c4762a1bSJed Brown } 260c4762a1bSJed Brown 261c4762a1bSJed Brown static PetscErrorCode uniform_strain_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 262c4762a1bSJed Brown { 263c4762a1bSJed Brown const PetscReal eps_xx = 0.1; 264c4762a1bSJed Brown const PetscReal eps_xy = 0.3; 265c4762a1bSJed Brown const PetscReal eps_yy = 0.25; 266c4762a1bSJed Brown PetscInt d; 267c4762a1bSJed Brown 268c4762a1bSJed Brown u[0] = eps_xx*x[0] + eps_xy*x[1]; 269c4762a1bSJed Brown u[1] = eps_xy*x[0] + eps_yy*x[1]; 270c4762a1bSJed Brown for (d = 2; d < dim; ++d) u[d] = 0.0; 271c4762a1bSJed Brown return 0; 272c4762a1bSJed Brown } 273c4762a1bSJed Brown 274482dd828SMatthew G. Knepley static void f0_mass_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 275482dd828SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 276482dd828SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 277482dd828SMatthew G. Knepley PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 278482dd828SMatthew G. Knepley { 279482dd828SMatthew G. Knepley const PetscInt Nc = dim; 280482dd828SMatthew G. Knepley PetscInt c; 281482dd828SMatthew G. Knepley 282482dd828SMatthew G. Knepley for (c = 0; c < Nc; ++c) f0[c] = u[c]; 283482dd828SMatthew G. Knepley } 284482dd828SMatthew G. Knepley 285c4762a1bSJed Brown static void f1_vlap_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 286c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 287c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 288c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 289c4762a1bSJed Brown { 290c4762a1bSJed Brown const PetscInt Nc = dim; 291c4762a1bSJed Brown PetscInt c, d; 292c4762a1bSJed Brown 293c4762a1bSJed Brown for (c = 0; c < Nc; ++c) for (d = 0; d < dim; ++d) f1[c*dim+d] += u_x[c*dim+d]; 294c4762a1bSJed Brown } 295c4762a1bSJed Brown 296c4762a1bSJed Brown static void f1_elas_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 297c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 298c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 299c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 300c4762a1bSJed Brown { 301c4762a1bSJed Brown const PetscInt Nc = dim; 302c4762a1bSJed Brown const PetscReal mu = 1.0; 303c4762a1bSJed Brown const PetscReal lambda = 1.0; 304c4762a1bSJed Brown PetscInt c, d; 305c4762a1bSJed Brown 306c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 307c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 308c4762a1bSJed Brown f1[c*dim+d] += mu*(u_x[c*dim+d] + u_x[d*dim+c]); 309c4762a1bSJed Brown f1[c*dim+c] += lambda*u_x[d*dim+d]; 310c4762a1bSJed Brown } 311c4762a1bSJed Brown } 312c4762a1bSJed Brown } 313c4762a1bSJed Brown 314482dd828SMatthew G. Knepley static void g0_mass_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 315482dd828SMatthew G. Knepley const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 316482dd828SMatthew G. Knepley const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 317482dd828SMatthew G. Knepley PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 318482dd828SMatthew G. Knepley { 319482dd828SMatthew G. Knepley const PetscInt Nc = dim; 320482dd828SMatthew G. Knepley PetscInt c; 321482dd828SMatthew G. Knepley 322482dd828SMatthew G. Knepley for (c = 0; c < Nc; ++c) g0[c*Nc + c] = 1.0; 323482dd828SMatthew G. Knepley } 324482dd828SMatthew G. Knepley 325c4762a1bSJed Brown static void g3_vlap_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 326c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 327c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 328c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 329c4762a1bSJed Brown { 330c4762a1bSJed Brown const PetscInt Nc = dim; 331c4762a1bSJed Brown PetscInt c, d; 332c4762a1bSJed Brown 333c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 334c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 335c4762a1bSJed Brown g3[((c*Nc + c)*dim + d)*dim + d] = 1.0; 336c4762a1bSJed Brown } 337c4762a1bSJed Brown } 338c4762a1bSJed Brown } 339c4762a1bSJed Brown 340c4762a1bSJed Brown /* 341c4762a1bSJed Brown \partial_df \phi_fc g_{fc,gc,df,dg} \partial_dg \phi_gc 342c4762a1bSJed Brown 343c4762a1bSJed Brown \partial_df \phi_fc \lambda \delta_{fc,df} \sum_gc \partial_dg \phi_gc \delta_{gc,dg} 344c4762a1bSJed Brown = \partial_fc \phi_fc \sum_gc \partial_gc \phi_gc 345c4762a1bSJed Brown */ 346c4762a1bSJed Brown static void g3_elas_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 347c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 348c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 349c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 350c4762a1bSJed Brown { 351c4762a1bSJed Brown const PetscInt Nc = dim; 352c4762a1bSJed Brown const PetscReal mu = 1.0; 353c4762a1bSJed Brown const PetscReal lambda = 1.0; 354c4762a1bSJed Brown PetscInt c, d; 355c4762a1bSJed Brown 356c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 357c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 358c4762a1bSJed Brown g3[((c*Nc + c)*dim + d)*dim + d] += mu; 359c4762a1bSJed Brown g3[((c*Nc + d)*dim + d)*dim + c] += mu; 360c4762a1bSJed Brown g3[((c*Nc + d)*dim + c)*dim + d] += lambda; 361c4762a1bSJed Brown } 362c4762a1bSJed Brown } 363c4762a1bSJed Brown } 364c4762a1bSJed Brown 365c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 366c4762a1bSJed Brown { 367482dd828SMatthew G. Knepley PetscInt sol = 0, def = 0; 368c4762a1bSJed Brown PetscErrorCode ierr; 369c4762a1bSJed Brown 370c4762a1bSJed Brown PetscFunctionBeginUser; 371482dd828SMatthew G. Knepley options->deform = DEFORM_NONE; 372c4762a1bSJed Brown options->solType = SOL_VLAP_QUADRATIC; 373c4762a1bSJed Brown options->useNearNullspace = PETSC_TRUE; 374*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncpy(options->dmType, DMPLEX, 256)); 375c4762a1bSJed Brown 376c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Linear Elasticity Problem Options", "DMPLEX");CHKERRQ(ierr); 377*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEList("-deform_type", "Type of domain deformation", "ex17.c", deformTypes, NUM_DEFORM_TYPES, deformTypes[options->deform], &def, NULL)); 378482dd828SMatthew G. Knepley options->deform = (DeformType) def; 379*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEList("-sol_type", "Type of exact solution", "ex17.c", solutionTypes, NUM_SOLUTION_TYPES, solutionTypes[options->solType], &sol, NULL)); 380c4762a1bSJed Brown options->solType = (SolutionType) sol; 381*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-near_nullspace", "Use the rigid body modes as an AMG near nullspace", "ex17.c", options->useNearNullspace, &options->useNearNullspace, NULL)); 382*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsFList("-dm_type", "Convert DMPlex to another format", "ex17.c", DMList, options->dmType, options->dmType, 256, NULL)); 3831e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 384c4762a1bSJed Brown PetscFunctionReturn(0); 385c4762a1bSJed Brown } 386c4762a1bSJed Brown 387d12cd81dSMatthew G. Knepley static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *ctx) 388d12cd81dSMatthew G. Knepley { 389d12cd81dSMatthew G. Knepley PetscBag bag; 390d12cd81dSMatthew G. Knepley Parameter *p; 391d12cd81dSMatthew G. Knepley 392d12cd81dSMatthew G. Knepley PetscFunctionBeginUser; 393d12cd81dSMatthew G. Knepley /* setup PETSc parameter bag */ 394*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagGetData(ctx->bag,(void**)&p)); 395*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagSetName(ctx->bag,"par","Elastic Parameters")); 396d12cd81dSMatthew G. Knepley bag = ctx->bag; 397*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagRegisterScalar(bag, &p->mu, 1.0, "mu", "Shear Modulus, Pa")); 398*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagRegisterScalar(bag, &p->lambda, 1.0, "lambda", "Lame's first parameter, Pa")); 399*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagSetFromOptions(bag)); 400d12cd81dSMatthew G. Knepley { 401d12cd81dSMatthew G. Knepley PetscViewer viewer; 402d12cd81dSMatthew G. Knepley PetscViewerFormat format; 403d12cd81dSMatthew G. Knepley PetscBool flg; 404d12cd81dSMatthew G. Knepley 405*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetViewer(comm, NULL, NULL, "-param_view", &viewer, &format, &flg)); 406d12cd81dSMatthew G. Knepley if (flg) { 407*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(viewer, format)); 408*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagView(bag, viewer)); 409*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFlush(viewer)); 410*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPopFormat(viewer)); 411*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 412d12cd81dSMatthew G. Knepley } 413d12cd81dSMatthew G. Knepley } 414d12cd81dSMatthew G. Knepley PetscFunctionReturn(0); 415d12cd81dSMatthew G. Knepley } 416d12cd81dSMatthew G. Knepley 417482dd828SMatthew G. Knepley static PetscErrorCode DMPlexDistortGeometry(DM dm) 418482dd828SMatthew G. Knepley { 419482dd828SMatthew G. Knepley DM cdm; 420482dd828SMatthew G. Knepley DMLabel label; 421482dd828SMatthew G. Knepley Vec coordinates; 422482dd828SMatthew G. Knepley PetscScalar *coords; 423482dd828SMatthew G. Knepley PetscReal mid = 0.5; 424482dd828SMatthew G. Knepley PetscInt cdim, d, vStart, vEnd, v; 425482dd828SMatthew G. Knepley 426482dd828SMatthew G. Knepley PetscFunctionBeginUser; 427*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(dm, &cdm)); 428*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDim(dm, &cdim)); 429*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 430*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "marker", &label)); 431*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocal(dm, &coordinates)); 432*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayWrite(coordinates, &coords)); 433482dd828SMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 434482dd828SMatthew G. Knepley PetscScalar *pcoords, shift; 435482dd828SMatthew G. Knepley PetscInt val; 436482dd828SMatthew G. Knepley 437*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(label, v, &val)); 438482dd828SMatthew G. Knepley if (val >= 0) continue; 439*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexPointLocalRef(cdm, v, coords, &pcoords)); 440482dd828SMatthew G. Knepley shift = 0.2 * PetscAbsScalar(pcoords[0] - mid); 441482dd828SMatthew G. Knepley shift = PetscRealPart(pcoords[0]) > mid ? shift : -shift; 442482dd828SMatthew G. Knepley for (d = 1; d < cdim; ++d) pcoords[d] += shift; 443482dd828SMatthew G. Knepley } 444*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayWrite(coordinates, &coords)); 445482dd828SMatthew G. Knepley PetscFunctionReturn(0); 446482dd828SMatthew G. Knepley } 447482dd828SMatthew G. Knepley 448c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 449c4762a1bSJed Brown { 450c4762a1bSJed Brown PetscFunctionBeginUser; 451*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 452*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 453*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 454482dd828SMatthew G. Knepley switch (user->deform) { 455482dd828SMatthew G. Knepley case DEFORM_NONE: break; 456*5f80ce2aSJacob Faibussowitsch case DEFORM_SHEAR: CHKERRQ(DMPlexShearGeometry(*dm, DM_X, NULL));break; 457*5f80ce2aSJacob Faibussowitsch case DEFORM_STEP: CHKERRQ(DMPlexDistortGeometry(*dm));break; 45898921bdaSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid deformation type: %s (%D)", deformTypes[PetscMin(user->deform, NUM_DEFORM_TYPES)], user->deform); 459482dd828SMatthew G. Knepley } 460*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetApplicationContext(*dm, user)); 461*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 462c4762a1bSJed Brown PetscFunctionReturn(0); 463c4762a1bSJed Brown } 464c4762a1bSJed Brown 465c4762a1bSJed Brown static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 466c4762a1bSJed Brown { 467c4762a1bSJed Brown PetscErrorCode (*exact)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 468d12cd81dSMatthew G. Knepley Parameter *param; 469482dd828SMatthew G. Knepley PetscDS ds; 47045480ffeSMatthew G. Knepley PetscWeakForm wf; 47145480ffeSMatthew G. Knepley DMLabel label; 47245480ffeSMatthew G. Knepley PetscInt id, bd; 473c4762a1bSJed Brown PetscInt dim; 474c4762a1bSJed Brown 475c4762a1bSJed Brown PetscFunctionBeginUser; 476*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm, &ds)); 477*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetWeakForm(ds, &wf)); 478*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetSpatialDimension(ds, &dim)); 479*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagGetData(user->bag, (void **) ¶m)); 480c4762a1bSJed Brown switch (user->solType) { 481482dd828SMatthew G. Knepley case SOL_MASS_QUADRATIC: 482*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, f0_mass_u, NULL)); 483*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, g0_mass_uu, NULL, NULL, NULL)); 484*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscWeakFormSetIndexResidual(wf, NULL, 0, 0, 0, 1, f0_mass_quadratic_u, 0, NULL)); 485c4762a1bSJed Brown switch (dim) { 486c4762a1bSJed Brown case 2: exact = quadratic_2d_u;break; 487c4762a1bSJed Brown case 3: exact = quadratic_3d_u;break; 48898921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %D", dim); 489482dd828SMatthew G. Knepley } 490482dd828SMatthew G. Knepley break; 491482dd828SMatthew G. Knepley case SOL_VLAP_QUADRATIC: 492*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, f0_vlap_quadratic_u, f1_vlap_u)); 493*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_vlap_uu)); 494482dd828SMatthew G. Knepley switch (dim) { 495482dd828SMatthew G. Knepley case 2: exact = quadratic_2d_u;break; 496482dd828SMatthew G. Knepley case 3: exact = quadratic_3d_u;break; 49798921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %D", dim); 498c4762a1bSJed Brown } 499c4762a1bSJed Brown break; 500c4762a1bSJed Brown case SOL_ELAS_QUADRATIC: 501*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, f0_elas_quadratic_u, f1_elas_u)); 502*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu)); 503c4762a1bSJed Brown switch (dim) { 504c4762a1bSJed Brown case 2: exact = quadratic_2d_u;break; 505c4762a1bSJed Brown case 3: exact = quadratic_3d_u;break; 50698921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %D", dim); 507c4762a1bSJed Brown } 508c4762a1bSJed Brown break; 509c4762a1bSJed Brown case SOL_VLAP_TRIG: 510*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, f0_vlap_trig_u, f1_vlap_u)); 511*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_vlap_uu)); 512c4762a1bSJed Brown switch (dim) { 513c4762a1bSJed Brown case 2: exact = trig_2d_u;break; 514c4762a1bSJed Brown case 3: exact = trig_3d_u;break; 51598921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %D", dim); 516c4762a1bSJed Brown } 517c4762a1bSJed Brown break; 518c4762a1bSJed Brown case SOL_ELAS_TRIG: 519*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, f0_elas_trig_u, f1_elas_u)); 520*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu)); 521c4762a1bSJed Brown switch (dim) { 522c4762a1bSJed Brown case 2: exact = trig_2d_u;break; 523c4762a1bSJed Brown case 3: exact = trig_3d_u;break; 52498921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %D", dim); 525c4762a1bSJed Brown } 526c4762a1bSJed Brown break; 527c4762a1bSJed Brown case SOL_ELAS_AXIAL_DISP: 528*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, NULL, f1_elas_u)); 529*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu)); 53045480ffeSMatthew G. Knepley id = dim == 3 ? 5 : 2; 531*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "marker", &label)); 532*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddBoundary(dm, DM_BC_NATURAL, "right", label, 1, &id, 0, 0, NULL, (void (*)(void)) NULL, NULL, user, &bd)); 533*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 534*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_elas_axial_disp_bd_u, 0, NULL)); 535c4762a1bSJed Brown exact = axial_disp_u; 536c4762a1bSJed Brown break; 537c4762a1bSJed Brown case SOL_ELAS_UNIFORM_STRAIN: 538*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, NULL, f1_elas_u)); 539*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu)); 540c4762a1bSJed Brown exact = uniform_strain_u; 541c4762a1bSJed Brown break; 542d12cd81dSMatthew G. Knepley case SOL_ELAS_GE: 543*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetResidual(ds, 0, NULL, f1_elas_u)); 544*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu)); 545d12cd81dSMatthew G. Knepley exact = zero; /* No exact solution available */ 546d12cd81dSMatthew G. Knepley break; 54798921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Invalid solution type: %s (%D)", solutionTypes[PetscMin(user->solType, NUM_SOLUTION_TYPES)], user->solType); 548c4762a1bSJed Brown } 549*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetExactSolution(ds, 0, exact, user)); 550*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "marker", &label)); 551c4762a1bSJed Brown if (user->solType == SOL_ELAS_AXIAL_DISP) { 552c4762a1bSJed Brown PetscInt cmp; 553c4762a1bSJed Brown 554c4762a1bSJed Brown id = dim == 3 ? 6 : 4; 555c4762a1bSJed Brown cmp = 0; 556*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "left", label, 1, &id, 0, 1, &cmp, (void (*)(void)) zero, NULL, user, NULL)); 557c4762a1bSJed Brown cmp = dim == 3 ? 2 : 1; 558c4762a1bSJed Brown id = dim == 3 ? 1 : 1; 559*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "bottom", label, 1, &id, 0, 1, &cmp, (void (*)(void)) zero, NULL, user, NULL)); 560c4762a1bSJed Brown if (dim == 3) { 561c4762a1bSJed Brown cmp = 1; 562c4762a1bSJed Brown id = 3; 563*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "front", label, 1, &id, 0, 1, &cmp, (void (*)(void)) zero, NULL, user, NULL)); 564c4762a1bSJed Brown } 565d12cd81dSMatthew G. Knepley } else if (user->solType == SOL_ELAS_GE) { 566d12cd81dSMatthew G. Knepley PetscInt cmp; 567d12cd81dSMatthew G. Knepley 568d12cd81dSMatthew G. Knepley id = dim == 3 ? 6 : 4; 569*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "left", label, 1, &id, 0, 0, NULL, (void (*)(void)) zero, NULL, user, NULL)); 570d12cd81dSMatthew G. Knepley id = dim == 3 ? 5 : 2; 571d12cd81dSMatthew G. Knepley cmp = 0; 572*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "right", label, 1, &id, 0, 1, &cmp, (void (*)(void)) ge_shift, NULL, user, NULL)); 573c4762a1bSJed Brown } else { 574c4762a1bSJed Brown id = 1; 575*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) exact, NULL, user, NULL)); 576c4762a1bSJed Brown } 577d12cd81dSMatthew G. Knepley /* Setup constants */ 578d12cd81dSMatthew G. Knepley { 579d12cd81dSMatthew G. Knepley PetscScalar constants[2]; 580d12cd81dSMatthew G. Knepley 581d12cd81dSMatthew G. Knepley constants[0] = param->mu; /* shear modulus, Pa */ 582d12cd81dSMatthew G. Knepley constants[1] = param->lambda; /* Lame's first parameter, Pa */ 583*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetConstants(ds, 2, constants)); 584d12cd81dSMatthew G. Knepley } 585c4762a1bSJed Brown PetscFunctionReturn(0); 586c4762a1bSJed Brown } 587c4762a1bSJed Brown 5888cda7954SMatthew G. Knepley static PetscErrorCode CreateElasticityNullSpace(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullspace) 589c4762a1bSJed Brown { 590c4762a1bSJed Brown PetscFunctionBegin; 591*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateRigidBody(dm, origField, nullspace)); 592c4762a1bSJed Brown PetscFunctionReturn(0); 593c4762a1bSJed Brown } 594c4762a1bSJed Brown 59530602db0SMatthew G. Knepley PetscErrorCode SetupFE(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), void *ctx) 596c4762a1bSJed Brown { 597c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 598c4762a1bSJed Brown DM cdm = dm; 599c4762a1bSJed Brown PetscFE fe; 600c4762a1bSJed Brown char prefix[PETSC_MAX_PATH_LEN]; 60130602db0SMatthew G. Knepley DMPolytopeType ct; 60230602db0SMatthew G. Knepley PetscBool simplex; 60330602db0SMatthew G. Knepley PetscInt dim, cStart; 604c4762a1bSJed Brown 605c4762a1bSJed Brown PetscFunctionBegin; 606c4762a1bSJed Brown /* Create finite element */ 607*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 608*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 609*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCellType(dm, cStart, &ct)); 61030602db0SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 611*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name)); 612*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, dim, simplex, name ? prefix : NULL, -1, &fe)); 613*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) fe, name)); 614c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 615*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe)); 616*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(dm)); 617*5f80ce2aSJacob Faibussowitsch CHKERRQ((*setup)(dm, user)); 618c4762a1bSJed Brown while (cdm) { 619*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCopyDisc(dm, cdm)); 620*5f80ce2aSJacob Faibussowitsch if (user->useNearNullspace) CHKERRQ(DMSetNearNullSpaceConstructor(cdm, 0, CreateElasticityNullSpace)); 621c4762a1bSJed Brown /* TODO: Check whether the boundary of coarse meshes is marked */ 622*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoarseDM(cdm, &cdm)); 623c4762a1bSJed Brown } 624*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe)); 625c4762a1bSJed Brown PetscFunctionReturn(0); 626c4762a1bSJed Brown } 627c4762a1bSJed Brown 628c4762a1bSJed Brown int main(int argc, char **argv) 629c4762a1bSJed Brown { 630c4762a1bSJed Brown DM dm; /* Problem specification */ 631c4762a1bSJed Brown SNES snes; /* Nonlinear solver */ 632c4762a1bSJed Brown Vec u; /* Solutions */ 633c4762a1bSJed Brown AppCtx user; /* User-defined work context */ 634c4762a1bSJed Brown PetscErrorCode ierr; 635c4762a1bSJed Brown 636c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 637*5f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user)); 638*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagCreate(PETSC_COMM_SELF, sizeof(Parameter), &user.bag)); 639*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetupParameters(PETSC_COMM_WORLD, &user)); 640c4762a1bSJed Brown /* Primal system */ 641*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESCreate(PETSC_COMM_WORLD, &snes)); 642*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 643*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetDM(snes, dm)); 644*5f80ce2aSJacob Faibussowitsch CHKERRQ(SetupFE(dm, "displacement", SetupPrimalProblem, &user)); 645*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(dm, &u)); 646*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(u, 0.0)); 647*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) u, "displacement")); 648*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 649*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFromOptions(snes)); 650*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSNESCheckFromOptions(snes, u)); 651*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSolve(snes, NULL, u)); 652*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetSolution(snes, &u)); 653*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(u, NULL, "-displacement_view")); 654c4762a1bSJed Brown /* Cleanup */ 655*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 656*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESDestroy(&snes)); 657*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 658*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBagDestroy(&user.bag)); 659c4762a1bSJed Brown ierr = PetscFinalize(); 660c4762a1bSJed Brown return ierr; 661c4762a1bSJed Brown } 662c4762a1bSJed Brown 663c4762a1bSJed Brown /*TEST 664c4762a1bSJed Brown 66530602db0SMatthew G. Knepley testset: 66630602db0SMatthew G. Knepley args: -dm_plex_box_faces 1,1,1 66730602db0SMatthew G. Knepley 668c4762a1bSJed Brown test: 669c4762a1bSJed Brown suffix: 2d_p1_quad_vlap 670c4762a1bSJed Brown requires: triangle 671c4762a1bSJed Brown args: -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate 672c4762a1bSJed Brown test: 673c4762a1bSJed Brown suffix: 2d_p2_quad_vlap 674c4762a1bSJed Brown requires: triangle 675c4762a1bSJed Brown args: -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001 676c4762a1bSJed Brown test: 677c4762a1bSJed Brown suffix: 2d_p3_quad_vlap 678c4762a1bSJed Brown requires: triangle 679c4762a1bSJed Brown args: -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001 680c4762a1bSJed Brown test: 681c4762a1bSJed Brown suffix: 2d_q1_quad_vlap 68230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate 683c4762a1bSJed Brown test: 684c4762a1bSJed Brown suffix: 2d_q2_quad_vlap 68530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001 686c4762a1bSJed Brown test: 687c4762a1bSJed Brown suffix: 2d_q3_quad_vlap 688c4762a1bSJed Brown requires: !single 68930602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001 690c4762a1bSJed Brown test: 691c4762a1bSJed Brown suffix: 2d_p1_quad_elas 692c4762a1bSJed Brown requires: triangle 693c4762a1bSJed Brown args: -sol_type elas_quad -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate 694c4762a1bSJed Brown test: 695c4762a1bSJed Brown suffix: 2d_p2_quad_elas 696c4762a1bSJed Brown requires: triangle 697c4762a1bSJed Brown args: -sol_type elas_quad -displacement_petscspace_degree 2 -dmsnes_check .0001 698c4762a1bSJed Brown test: 699c4762a1bSJed Brown suffix: 2d_p3_quad_elas 700c4762a1bSJed Brown requires: triangle 701c4762a1bSJed Brown args: -sol_type elas_quad -displacement_petscspace_degree 3 -dmsnes_check .0001 702c4762a1bSJed Brown test: 703c4762a1bSJed Brown suffix: 2d_q1_quad_elas 70430602db0SMatthew 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 705c4762a1bSJed Brown test: 706c4762a1bSJed Brown suffix: 2d_q1_quad_elas_shear 707482dd828SMatthew 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 708c4762a1bSJed Brown test: 709c4762a1bSJed Brown suffix: 2d_q2_quad_elas 71030602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dmsnes_check .0001 711c4762a1bSJed Brown test: 712c4762a1bSJed Brown suffix: 2d_q2_quad_elas_shear 713482dd828SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_simplex 0 -deform_type shear -displacement_petscspace_degree 2 -dmsnes_check 714c4762a1bSJed Brown test: 715c4762a1bSJed Brown suffix: 2d_q3_quad_elas 71630602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dmsnes_check .0001 717c4762a1bSJed Brown test: 718c4762a1bSJed Brown suffix: 2d_q3_quad_elas_shear 719c4762a1bSJed Brown requires: !single 720482dd828SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_simplex 0 -deform_type shear -displacement_petscspace_degree 3 -dmsnes_check 721c4762a1bSJed Brown 722c4762a1bSJed Brown test: 723c4762a1bSJed Brown suffix: 3d_p1_quad_vlap 724c4762a1bSJed Brown requires: ctetgen 72530602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_refine 1 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate 726c4762a1bSJed Brown test: 727c4762a1bSJed Brown suffix: 3d_p2_quad_vlap 728c4762a1bSJed Brown requires: ctetgen 72930602db0SMatthew G. Knepley args: -dm_plex_dim 3 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001 730c4762a1bSJed Brown test: 731c4762a1bSJed Brown suffix: 3d_p3_quad_vlap 732c4762a1bSJed Brown requires: ctetgen 73330602db0SMatthew G. Knepley args: -dm_plex_dim 3 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001 734c4762a1bSJed Brown test: 735c4762a1bSJed Brown suffix: 3d_q1_quad_vlap 73630602db0SMatthew 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 737c4762a1bSJed Brown test: 738c4762a1bSJed Brown suffix: 3d_q2_quad_vlap 73930602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001 740c4762a1bSJed Brown test: 741c4762a1bSJed Brown suffix: 3d_q3_quad_vlap 74230602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001 743c4762a1bSJed Brown test: 744c4762a1bSJed Brown suffix: 3d_p1_quad_elas 745c4762a1bSJed Brown requires: ctetgen 74630602db0SMatthew 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 747c4762a1bSJed Brown test: 748c4762a1bSJed Brown suffix: 3d_p2_quad_elas 749c4762a1bSJed Brown requires: ctetgen 75030602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_dim 3 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001 751c4762a1bSJed Brown test: 752c4762a1bSJed Brown suffix: 3d_p3_quad_elas 753c4762a1bSJed Brown requires: ctetgen 75430602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_dim 3 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001 755c4762a1bSJed Brown test: 756c4762a1bSJed Brown suffix: 3d_q1_quad_elas 75730602db0SMatthew 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 758c4762a1bSJed Brown test: 759c4762a1bSJed Brown suffix: 3d_q2_quad_elas 76030602db0SMatthew 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 761c4762a1bSJed Brown test: 762c4762a1bSJed Brown suffix: 3d_q3_quad_elas 763c4762a1bSJed Brown requires: !single 76430602db0SMatthew 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 765c4762a1bSJed Brown 766c4762a1bSJed Brown test: 767c4762a1bSJed Brown suffix: 2d_p1_trig_vlap 768c4762a1bSJed Brown requires: triangle 769c4762a1bSJed Brown args: -sol_type vlap_trig -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 770c4762a1bSJed Brown test: 771c4762a1bSJed Brown suffix: 2d_p2_trig_vlap 772c4762a1bSJed Brown requires: triangle 773c4762a1bSJed Brown args: -sol_type vlap_trig -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 774c4762a1bSJed Brown test: 775c4762a1bSJed Brown suffix: 2d_p3_trig_vlap 776c4762a1bSJed Brown requires: triangle 777c4762a1bSJed Brown args: -sol_type vlap_trig -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 778c4762a1bSJed Brown test: 779c4762a1bSJed Brown suffix: 2d_q1_trig_vlap 78030602db0SMatthew 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 781c4762a1bSJed Brown test: 782c4762a1bSJed Brown suffix: 2d_q2_trig_vlap 78330602db0SMatthew 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 784c4762a1bSJed Brown test: 785c4762a1bSJed Brown suffix: 2d_q3_trig_vlap 78630602db0SMatthew 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 787c4762a1bSJed Brown test: 788c4762a1bSJed Brown suffix: 2d_p1_trig_elas 789c4762a1bSJed Brown requires: triangle 790c4762a1bSJed Brown args: -sol_type elas_trig -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 791c4762a1bSJed Brown test: 792c4762a1bSJed Brown suffix: 2d_p2_trig_elas 793c4762a1bSJed Brown requires: triangle 794c4762a1bSJed Brown args: -sol_type elas_trig -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 795c4762a1bSJed Brown test: 796c4762a1bSJed Brown suffix: 2d_p3_trig_elas 797c4762a1bSJed Brown requires: triangle 798c4762a1bSJed Brown args: -sol_type elas_trig -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 799c4762a1bSJed Brown test: 800c4762a1bSJed Brown suffix: 2d_q1_trig_elas 80130602db0SMatthew 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 802c4762a1bSJed Brown test: 803c4762a1bSJed Brown suffix: 2d_q1_trig_elas_shear 804482dd828SMatthew 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 805c4762a1bSJed Brown test: 806c4762a1bSJed Brown suffix: 2d_q2_trig_elas 80730602db0SMatthew 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 808c4762a1bSJed Brown test: 809c4762a1bSJed Brown suffix: 2d_q2_trig_elas_shear 810482dd828SMatthew 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 811c4762a1bSJed Brown test: 812c4762a1bSJed Brown suffix: 2d_q3_trig_elas 81330602db0SMatthew 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 814c4762a1bSJed Brown test: 815c4762a1bSJed Brown suffix: 2d_q3_trig_elas_shear 816482dd828SMatthew 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 817c4762a1bSJed Brown 818c4762a1bSJed Brown test: 819c4762a1bSJed Brown suffix: 3d_p1_trig_vlap 820c4762a1bSJed Brown requires: ctetgen 82130602db0SMatthew 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 822c4762a1bSJed Brown test: 823c4762a1bSJed Brown suffix: 3d_p2_trig_vlap 824c4762a1bSJed Brown requires: ctetgen 82530602db0SMatthew 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 826c4762a1bSJed Brown test: 827c4762a1bSJed Brown suffix: 3d_p3_trig_vlap 828c4762a1bSJed Brown requires: ctetgen 82930602db0SMatthew 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 830c4762a1bSJed Brown test: 831c4762a1bSJed Brown suffix: 3d_q1_trig_vlap 83230602db0SMatthew 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 833c4762a1bSJed Brown test: 834c4762a1bSJed Brown suffix: 3d_q2_trig_vlap 83530602db0SMatthew 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 836c4762a1bSJed Brown test: 837c4762a1bSJed Brown suffix: 3d_q3_trig_vlap 838c4762a1bSJed Brown requires: !__float128 83930602db0SMatthew 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 840c4762a1bSJed Brown test: 841c4762a1bSJed Brown suffix: 3d_p1_trig_elas 842c4762a1bSJed Brown requires: ctetgen 84330602db0SMatthew 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 844c4762a1bSJed Brown test: 845c4762a1bSJed Brown suffix: 3d_p2_trig_elas 846c4762a1bSJed Brown requires: ctetgen 84730602db0SMatthew 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 848c4762a1bSJed Brown test: 849c4762a1bSJed Brown suffix: 3d_p3_trig_elas 850c4762a1bSJed Brown requires: ctetgen 85130602db0SMatthew 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 852c4762a1bSJed Brown test: 853c4762a1bSJed Brown suffix: 3d_q1_trig_elas 85430602db0SMatthew 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 855c4762a1bSJed Brown test: 856c4762a1bSJed Brown suffix: 3d_q2_trig_elas 85730602db0SMatthew 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 858c4762a1bSJed Brown test: 859c4762a1bSJed Brown suffix: 3d_q3_trig_elas 860c4762a1bSJed Brown requires: !__float128 86130602db0SMatthew 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 862c4762a1bSJed Brown 863c4762a1bSJed Brown test: 864c4762a1bSJed Brown suffix: 2d_p1_axial_elas 865c4762a1bSJed Brown requires: triangle 866c4762a1bSJed Brown args: -sol_type elas_axial_disp -displacement_petscspace_degree 1 -dm_plex_separate_marker -dm_refine 2 -dmsnes_check .0001 -pc_type lu 867c4762a1bSJed Brown test: 868c4762a1bSJed Brown suffix: 2d_p2_axial_elas 869c4762a1bSJed Brown requires: triangle 870c4762a1bSJed Brown args: -sol_type elas_axial_disp -displacement_petscspace_degree 2 -dm_plex_separate_marker -dmsnes_check .0001 -pc_type lu 871c4762a1bSJed Brown test: 872c4762a1bSJed Brown suffix: 2d_p3_axial_elas 873c4762a1bSJed Brown requires: triangle 874c4762a1bSJed Brown args: -sol_type elas_axial_disp -displacement_petscspace_degree 3 -dm_plex_separate_marker -dmsnes_check .0001 -pc_type lu 875c4762a1bSJed Brown test: 876c4762a1bSJed Brown suffix: 2d_q1_axial_elas 87730602db0SMatthew 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 878c4762a1bSJed Brown test: 879c4762a1bSJed Brown suffix: 2d_q2_axial_elas 88030602db0SMatthew 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 881c4762a1bSJed Brown test: 882c4762a1bSJed Brown suffix: 2d_q3_axial_elas 88330602db0SMatthew 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 884c4762a1bSJed Brown 885c4762a1bSJed Brown test: 886c4762a1bSJed Brown suffix: 2d_p1_uniform_elas 887c4762a1bSJed Brown requires: triangle 888c4762a1bSJed Brown args: -sol_type elas_uniform_strain -displacement_petscspace_degree 1 -dm_refine 2 -dmsnes_check .0001 -pc_type lu 889c4762a1bSJed Brown test: 890c4762a1bSJed Brown suffix: 2d_p2_uniform_elas 891c4762a1bSJed Brown requires: triangle 892c4762a1bSJed Brown args: -sol_type elas_uniform_strain -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001 -pc_type lu 893c4762a1bSJed Brown test: 894c4762a1bSJed Brown suffix: 2d_p3_uniform_elas 895c4762a1bSJed Brown requires: triangle 896c4762a1bSJed Brown args: -sol_type elas_uniform_strain -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001 -pc_type lu 897c4762a1bSJed Brown test: 898c4762a1bSJed Brown suffix: 2d_q1_uniform_elas 89930602db0SMatthew 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 900c4762a1bSJed Brown test: 901c4762a1bSJed Brown suffix: 2d_q2_uniform_elas 902c4762a1bSJed Brown requires: !single 90330602db0SMatthew 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 904c4762a1bSJed Brown test: 905c4762a1bSJed Brown suffix: 2d_q3_uniform_elas 906c4762a1bSJed Brown requires: !single 90730602db0SMatthew 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 908482dd828SMatthew G. Knepley test: 909482dd828SMatthew G. Knepley suffix: 2d_p1_uniform_elas_step 910482dd828SMatthew G. Knepley requires: triangle 911482dd828SMatthew G. Knepley args: -sol_type elas_uniform_strain -deform_type step -displacement_petscspace_degree 1 -dm_refine 2 -dmsnes_check .0001 -pc_type lu 912482dd828SMatthew G. Knepley 913482dd828SMatthew G. Knepley testset: 914482dd828SMatthew 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 915482dd828SMatthew G. Knepley 916482dd828SMatthew G. Knepley test: 917482dd828SMatthew G. Knepley suffix: 2d_q1_uniform_elas_step 918482dd828SMatthew G. Knepley args: -sol_type elas_uniform_strain -dm_refine 2 919482dd828SMatthew G. Knepley test: 920482dd828SMatthew G. Knepley suffix: 2d_q1_quad_vlap_step 921482dd828SMatthew G. Knepley args: 922482dd828SMatthew G. Knepley test: 923482dd828SMatthew G. Knepley suffix: 2d_q2_quad_vlap_step 924482dd828SMatthew G. Knepley args: -displacement_petscspace_degree 2 925482dd828SMatthew G. Knepley test: 926482dd828SMatthew G. Knepley suffix: 2d_q1_quad_mass_step 927482dd828SMatthew G. Knepley args: -sol_type mass_quad 928c4762a1bSJed Brown 929d12cd81dSMatthew G. Knepley testset: 930d12cd81dSMatthew 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 \ 931d12cd81dSMatthew G. Knepley -dm_plex_box_faces 5,5,2 -dm_plex_separate_marker -dm_refine 0 -petscpartitioner_type simple \ 932d12cd81dSMatthew G. Knepley -sol_type elas_ge \ 933d12cd81dSMatthew G. Knepley -snes_max_it 2 -snes_rtol 1.e-10 \ 934d12cd81dSMatthew G. Knepley -ksp_type cg -ksp_rtol 1.e-10 -ksp_max_it 100 -ksp_norm_type unpreconditioned \ 935d12cd81dSMatthew G. Knepley -pc_type gamg -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 \ 936d12cd81dSMatthew G. Knepley -pc_gamg_coarse_eq_limit 10 -pc_gamg_reuse_interpolation true \ 937d12cd81dSMatthew G. Knepley -pc_gamg_square_graph 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 \ 938d12cd81dSMatthew 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 \ 939d12cd81dSMatthew G. Knepley -matptap_via scalable 940d12cd81dSMatthew G. Knepley test: 941d12cd81dSMatthew G. Knepley suffix: ge_q1_0 942d12cd81dSMatthew G. Knepley args: -displacement_petscspace_degree 1 943d12cd81dSMatthew G. Knepley 944c4762a1bSJed Brown TEST*/ 945