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 369c4762a1bSJed Brown PetscFunctionBeginUser; 370482dd828SMatthew G. Knepley options->deform = DEFORM_NONE; 371c4762a1bSJed Brown options->solType = SOL_VLAP_QUADRATIC; 372c4762a1bSJed Brown options->useNearNullspace = PETSC_TRUE; 3739566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(options->dmType, DMPLEX, 256)); 374c4762a1bSJed Brown 375d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Linear Elasticity Problem Options", "DMPLEX"); 3769566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-deform_type", "Type of domain deformation", "ex17.c", deformTypes, NUM_DEFORM_TYPES, deformTypes[options->deform], &def, NULL)); 377482dd828SMatthew G. Knepley options->deform = (DeformType) def; 3789566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-sol_type", "Type of exact solution", "ex17.c", solutionTypes, NUM_SOLUTION_TYPES, solutionTypes[options->solType], &sol, NULL)); 379c4762a1bSJed Brown options->solType = (SolutionType) sol; 3809566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-near_nullspace", "Use the rigid body modes as an AMG near nullspace", "ex17.c", options->useNearNullspace, &options->useNearNullspace, NULL)); 3819566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-dm_type", "Convert DMPlex to another format", "ex17.c", DMList, options->dmType, options->dmType, 256, NULL)); 382d0609cedSBarry Smith PetscOptionsEnd(); 383c4762a1bSJed Brown PetscFunctionReturn(0); 384c4762a1bSJed Brown } 385c4762a1bSJed Brown 386d12cd81dSMatthew G. Knepley static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *ctx) 387d12cd81dSMatthew G. Knepley { 388d12cd81dSMatthew G. Knepley PetscBag bag; 389d12cd81dSMatthew G. Knepley Parameter *p; 390d12cd81dSMatthew G. Knepley 391d12cd81dSMatthew G. Knepley PetscFunctionBeginUser; 392d12cd81dSMatthew G. Knepley /* setup PETSc parameter bag */ 3939566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(ctx->bag,(void**)&p)); 3949566063dSJacob Faibussowitsch PetscCall(PetscBagSetName(ctx->bag,"par","Elastic Parameters")); 395d12cd81dSMatthew G. Knepley bag = ctx->bag; 3969566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterScalar(bag, &p->mu, 1.0, "mu", "Shear Modulus, Pa")); 3979566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterScalar(bag, &p->lambda, 1.0, "lambda", "Lame's first parameter, Pa")); 3989566063dSJacob Faibussowitsch PetscCall(PetscBagSetFromOptions(bag)); 399d12cd81dSMatthew G. Knepley { 400d12cd81dSMatthew G. Knepley PetscViewer viewer; 401d12cd81dSMatthew G. Knepley PetscViewerFormat format; 402d12cd81dSMatthew G. Knepley PetscBool flg; 403d12cd81dSMatthew G. Knepley 4049566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(comm, NULL, NULL, "-param_view", &viewer, &format, &flg)); 405d12cd81dSMatthew G. Knepley if (flg) { 4069566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, format)); 4079566063dSJacob Faibussowitsch PetscCall(PetscBagView(bag, viewer)); 4089566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 4099566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 4109566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 411d12cd81dSMatthew G. Knepley } 412d12cd81dSMatthew G. Knepley } 413d12cd81dSMatthew G. Knepley PetscFunctionReturn(0); 414d12cd81dSMatthew G. Knepley } 415d12cd81dSMatthew G. Knepley 416482dd828SMatthew G. Knepley static PetscErrorCode DMPlexDistortGeometry(DM dm) 417482dd828SMatthew G. Knepley { 418482dd828SMatthew G. Knepley DM cdm; 419482dd828SMatthew G. Knepley DMLabel label; 420482dd828SMatthew G. Knepley Vec coordinates; 421482dd828SMatthew G. Knepley PetscScalar *coords; 422482dd828SMatthew G. Knepley PetscReal mid = 0.5; 423482dd828SMatthew G. Knepley PetscInt cdim, d, vStart, vEnd, v; 424482dd828SMatthew G. Knepley 425482dd828SMatthew G. Knepley PetscFunctionBeginUser; 4269566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 4279566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim)); 4289566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 4299566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 4309566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 4319566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(coordinates, &coords)); 432482dd828SMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 433482dd828SMatthew G. Knepley PetscScalar *pcoords, shift; 434482dd828SMatthew G. Knepley PetscInt val; 435482dd828SMatthew G. Knepley 4369566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(label, v, &val)); 437482dd828SMatthew G. Knepley if (val >= 0) continue; 4389566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRef(cdm, v, coords, &pcoords)); 439482dd828SMatthew G. Knepley shift = 0.2 * PetscAbsScalar(pcoords[0] - mid); 440482dd828SMatthew G. Knepley shift = PetscRealPart(pcoords[0]) > mid ? shift : -shift; 441482dd828SMatthew G. Knepley for (d = 1; d < cdim; ++d) pcoords[d] += shift; 442482dd828SMatthew G. Knepley } 4439566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(coordinates, &coords)); 444482dd828SMatthew G. Knepley PetscFunctionReturn(0); 445482dd828SMatthew G. Knepley } 446482dd828SMatthew G. Knepley 447c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 448c4762a1bSJed Brown { 449c4762a1bSJed Brown PetscFunctionBeginUser; 4509566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 4519566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 4529566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 453482dd828SMatthew G. Knepley switch (user->deform) { 454482dd828SMatthew G. Knepley case DEFORM_NONE: break; 4559566063dSJacob Faibussowitsch case DEFORM_SHEAR: PetscCall(DMPlexShearGeometry(*dm, DM_X, NULL));break; 4569566063dSJacob Faibussowitsch case DEFORM_STEP: PetscCall(DMPlexDistortGeometry(*dm));break; 45763a3b9bcSJacob Faibussowitsch default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid deformation type: %s (%d)", deformTypes[PetscMin(user->deform, NUM_DEFORM_TYPES)], user->deform); 458482dd828SMatthew G. Knepley } 4599566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*dm, user)); 4609566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 461c4762a1bSJed Brown PetscFunctionReturn(0); 462c4762a1bSJed Brown } 463c4762a1bSJed Brown 464c4762a1bSJed Brown static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 465c4762a1bSJed Brown { 466c4762a1bSJed Brown PetscErrorCode (*exact)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 467d12cd81dSMatthew G. Knepley Parameter *param; 468482dd828SMatthew G. Knepley PetscDS ds; 46945480ffeSMatthew G. Knepley PetscWeakForm wf; 47045480ffeSMatthew G. Knepley DMLabel label; 47145480ffeSMatthew G. Knepley PetscInt id, bd; 472c4762a1bSJed Brown PetscInt dim; 473c4762a1bSJed Brown 474c4762a1bSJed Brown PetscFunctionBeginUser; 4759566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 4769566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf)); 4779566063dSJacob Faibussowitsch PetscCall(PetscDSGetSpatialDimension(ds, &dim)); 4789566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **) ¶m)); 479c4762a1bSJed Brown switch (user->solType) { 480482dd828SMatthew G. Knepley case SOL_MASS_QUADRATIC: 4819566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_mass_u, NULL)); 4829566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_mass_uu, NULL, NULL, NULL)); 4839566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, 0, 0, 1, f0_mass_quadratic_u, 0, NULL)); 484c4762a1bSJed Brown switch (dim) { 485c4762a1bSJed Brown case 2: exact = quadratic_2d_u;break; 486c4762a1bSJed Brown case 3: exact = quadratic_3d_u;break; 48763a3b9bcSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim); 488482dd828SMatthew G. Knepley } 489482dd828SMatthew G. Knepley break; 490482dd828SMatthew G. Knepley case SOL_VLAP_QUADRATIC: 4919566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_vlap_quadratic_u, f1_vlap_u)); 4929566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_vlap_uu)); 493482dd828SMatthew G. Knepley switch (dim) { 494482dd828SMatthew G. Knepley case 2: exact = quadratic_2d_u;break; 495482dd828SMatthew G. Knepley case 3: exact = quadratic_3d_u;break; 49663a3b9bcSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim); 497c4762a1bSJed Brown } 498c4762a1bSJed Brown break; 499c4762a1bSJed Brown case SOL_ELAS_QUADRATIC: 5009566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_elas_quadratic_u, f1_elas_u)); 5019566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu)); 502c4762a1bSJed Brown switch (dim) { 503c4762a1bSJed Brown case 2: exact = quadratic_2d_u;break; 504c4762a1bSJed Brown case 3: exact = quadratic_3d_u;break; 50563a3b9bcSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim); 506c4762a1bSJed Brown } 507c4762a1bSJed Brown break; 508c4762a1bSJed Brown case SOL_VLAP_TRIG: 5099566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_vlap_trig_u, f1_vlap_u)); 5109566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_vlap_uu)); 511c4762a1bSJed Brown switch (dim) { 512c4762a1bSJed Brown case 2: exact = trig_2d_u;break; 513c4762a1bSJed Brown case 3: exact = trig_3d_u;break; 51463a3b9bcSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim); 515c4762a1bSJed Brown } 516c4762a1bSJed Brown break; 517c4762a1bSJed Brown case SOL_ELAS_TRIG: 5189566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_elas_trig_u, f1_elas_u)); 5199566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu)); 520c4762a1bSJed Brown switch (dim) { 521c4762a1bSJed Brown case 2: exact = trig_2d_u;break; 522c4762a1bSJed Brown case 3: exact = trig_3d_u;break; 52363a3b9bcSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim); 524c4762a1bSJed Brown } 525c4762a1bSJed Brown break; 526c4762a1bSJed Brown case SOL_ELAS_AXIAL_DISP: 5279566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, NULL, f1_elas_u)); 5289566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu)); 52945480ffeSMatthew G. Knepley id = dim == 3 ? 5 : 2; 5309566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 5319566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "right", label, 1, &id, 0, 0, NULL, (void (*)(void)) NULL, NULL, user, &bd)); 5329566063dSJacob Faibussowitsch PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 5339566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_elas_axial_disp_bd_u, 0, NULL)); 534c4762a1bSJed Brown exact = axial_disp_u; 535c4762a1bSJed Brown break; 536c4762a1bSJed Brown case SOL_ELAS_UNIFORM_STRAIN: 5379566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, NULL, f1_elas_u)); 5389566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu)); 539c4762a1bSJed Brown exact = uniform_strain_u; 540c4762a1bSJed Brown break; 541d12cd81dSMatthew G. Knepley case SOL_ELAS_GE: 5429566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, NULL, f1_elas_u)); 5439566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu)); 544d12cd81dSMatthew G. Knepley exact = zero; /* No exact solution available */ 545d12cd81dSMatthew G. Knepley break; 54663a3b9bcSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Invalid solution type: %s (%d)", solutionTypes[PetscMin(user->solType, NUM_SOLUTION_TYPES)], user->solType); 547c4762a1bSJed Brown } 5489566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, exact, user)); 5499566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 550c4762a1bSJed Brown if (user->solType == SOL_ELAS_AXIAL_DISP) { 551c4762a1bSJed Brown PetscInt cmp; 552c4762a1bSJed Brown 553c4762a1bSJed Brown id = dim == 3 ? 6 : 4; 554c4762a1bSJed Brown cmp = 0; 5559566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "left", label, 1, &id, 0, 1, &cmp, (void (*)(void)) zero, NULL, user, NULL)); 556c4762a1bSJed Brown cmp = dim == 3 ? 2 : 1; 557c4762a1bSJed Brown id = dim == 3 ? 1 : 1; 5589566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "bottom", label, 1, &id, 0, 1, &cmp, (void (*)(void)) zero, NULL, user, NULL)); 559c4762a1bSJed Brown if (dim == 3) { 560c4762a1bSJed Brown cmp = 1; 561c4762a1bSJed Brown id = 3; 5629566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "front", label, 1, &id, 0, 1, &cmp, (void (*)(void)) zero, NULL, user, NULL)); 563c4762a1bSJed Brown } 564d12cd81dSMatthew G. Knepley } else if (user->solType == SOL_ELAS_GE) { 565d12cd81dSMatthew G. Knepley PetscInt cmp; 566d12cd81dSMatthew G. Knepley 567d12cd81dSMatthew G. Knepley id = dim == 3 ? 6 : 4; 5689566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "left", label, 1, &id, 0, 0, NULL, (void (*)(void)) zero, NULL, user, NULL)); 569d12cd81dSMatthew G. Knepley id = dim == 3 ? 5 : 2; 570d12cd81dSMatthew G. Knepley cmp = 0; 5719566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "right", label, 1, &id, 0, 1, &cmp, (void (*)(void)) ge_shift, NULL, user, NULL)); 572c4762a1bSJed Brown } else { 573c4762a1bSJed Brown id = 1; 5749566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) exact, NULL, user, NULL)); 575c4762a1bSJed Brown } 576d12cd81dSMatthew G. Knepley /* Setup constants */ 577d12cd81dSMatthew G. Knepley { 578d12cd81dSMatthew G. Knepley PetscScalar constants[2]; 579d12cd81dSMatthew G. Knepley 580d12cd81dSMatthew G. Knepley constants[0] = param->mu; /* shear modulus, Pa */ 581d12cd81dSMatthew G. Knepley constants[1] = param->lambda; /* Lame's first parameter, Pa */ 5829566063dSJacob Faibussowitsch PetscCall(PetscDSSetConstants(ds, 2, constants)); 583d12cd81dSMatthew G. Knepley } 584c4762a1bSJed Brown PetscFunctionReturn(0); 585c4762a1bSJed Brown } 586c4762a1bSJed Brown 5878cda7954SMatthew G. Knepley static PetscErrorCode CreateElasticityNullSpace(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullspace) 588c4762a1bSJed Brown { 589c4762a1bSJed Brown PetscFunctionBegin; 5909566063dSJacob Faibussowitsch PetscCall(DMPlexCreateRigidBody(dm, origField, nullspace)); 591c4762a1bSJed Brown PetscFunctionReturn(0); 592c4762a1bSJed Brown } 593c4762a1bSJed Brown 59430602db0SMatthew G. Knepley PetscErrorCode SetupFE(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), void *ctx) 595c4762a1bSJed Brown { 596c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 597c4762a1bSJed Brown DM cdm = dm; 598c4762a1bSJed Brown PetscFE fe; 599c4762a1bSJed Brown char prefix[PETSC_MAX_PATH_LEN]; 60030602db0SMatthew G. Knepley DMPolytopeType ct; 60130602db0SMatthew G. Knepley PetscBool simplex; 60230602db0SMatthew G. Knepley PetscInt dim, cStart; 603c4762a1bSJed Brown 604c4762a1bSJed Brown PetscFunctionBegin; 605c4762a1bSJed Brown /* Create finite element */ 6069566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 6079566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 6089566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 60930602db0SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 6109566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name)); 6119566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, dim, simplex, name ? prefix : NULL, -1, &fe)); 6129566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe, name)); 613c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 6149566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe)); 6159566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 6169566063dSJacob Faibussowitsch PetscCall((*setup)(dm, user)); 617c4762a1bSJed Brown while (cdm) { 6189566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 6199566063dSJacob Faibussowitsch if (user->useNearNullspace) PetscCall(DMSetNearNullSpaceConstructor(cdm, 0, CreateElasticityNullSpace)); 620c4762a1bSJed Brown /* TODO: Check whether the boundary of coarse meshes is marked */ 6219566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 622c4762a1bSJed Brown } 6239566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 624c4762a1bSJed Brown PetscFunctionReturn(0); 625c4762a1bSJed Brown } 626c4762a1bSJed Brown 627c4762a1bSJed Brown int main(int argc, char **argv) 628c4762a1bSJed Brown { 629c4762a1bSJed Brown DM dm; /* Problem specification */ 630c4762a1bSJed Brown SNES snes; /* Nonlinear solver */ 631c4762a1bSJed Brown Vec u; /* Solutions */ 632c4762a1bSJed Brown AppCtx user; /* User-defined work context */ 633c4762a1bSJed Brown 6349566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 6359566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 6369566063dSJacob Faibussowitsch PetscCall(PetscBagCreate(PETSC_COMM_SELF, sizeof(Parameter), &user.bag)); 6379566063dSJacob Faibussowitsch PetscCall(SetupParameters(PETSC_COMM_WORLD, &user)); 638c4762a1bSJed Brown /* Primal system */ 6399566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 6409566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 6419566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 6429566063dSJacob Faibussowitsch PetscCall(SetupFE(dm, "displacement", SetupPrimalProblem, &user)); 6439566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 6449566063dSJacob Faibussowitsch PetscCall(VecSet(u, 0.0)); 6459566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) u, "displacement")); 6469566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 6479566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 6489566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u)); 6499566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 6509566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &u)); 6519566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-displacement_view")); 652c4762a1bSJed Brown /* Cleanup */ 6539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 6549566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 6559566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 6569566063dSJacob Faibussowitsch PetscCall(PetscBagDestroy(&user.bag)); 6579566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 658b122ec5aSJacob Faibussowitsch return 0; 659c4762a1bSJed Brown } 660c4762a1bSJed Brown 661c4762a1bSJed Brown /*TEST 662c4762a1bSJed Brown 66330602db0SMatthew G. Knepley testset: 66430602db0SMatthew G. Knepley args: -dm_plex_box_faces 1,1,1 66530602db0SMatthew G. Knepley 666c4762a1bSJed Brown test: 667c4762a1bSJed Brown suffix: 2d_p1_quad_vlap 668c4762a1bSJed Brown requires: triangle 669c4762a1bSJed Brown args: -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate 670c4762a1bSJed Brown test: 671c4762a1bSJed Brown suffix: 2d_p2_quad_vlap 672c4762a1bSJed Brown requires: triangle 673c4762a1bSJed Brown args: -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001 674c4762a1bSJed Brown test: 675c4762a1bSJed Brown suffix: 2d_p3_quad_vlap 676c4762a1bSJed Brown requires: triangle 677c4762a1bSJed Brown args: -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001 678c4762a1bSJed Brown test: 679c4762a1bSJed Brown suffix: 2d_q1_quad_vlap 68030602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate 681c4762a1bSJed Brown test: 682c4762a1bSJed Brown suffix: 2d_q2_quad_vlap 68330602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001 684c4762a1bSJed Brown test: 685c4762a1bSJed Brown suffix: 2d_q3_quad_vlap 686c4762a1bSJed Brown requires: !single 68730602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001 688c4762a1bSJed Brown test: 689c4762a1bSJed Brown suffix: 2d_p1_quad_elas 690c4762a1bSJed Brown requires: triangle 691c4762a1bSJed Brown args: -sol_type elas_quad -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate 692c4762a1bSJed Brown test: 693c4762a1bSJed Brown suffix: 2d_p2_quad_elas 694c4762a1bSJed Brown requires: triangle 695c4762a1bSJed Brown args: -sol_type elas_quad -displacement_petscspace_degree 2 -dmsnes_check .0001 696c4762a1bSJed Brown test: 697c4762a1bSJed Brown suffix: 2d_p3_quad_elas 698c4762a1bSJed Brown requires: triangle 699c4762a1bSJed Brown args: -sol_type elas_quad -displacement_petscspace_degree 3 -dmsnes_check .0001 700c4762a1bSJed Brown test: 701c4762a1bSJed Brown suffix: 2d_q1_quad_elas 70230602db0SMatthew 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 703c4762a1bSJed Brown test: 704c4762a1bSJed Brown suffix: 2d_q1_quad_elas_shear 705482dd828SMatthew 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 706c4762a1bSJed Brown test: 707c4762a1bSJed Brown suffix: 2d_q2_quad_elas 70830602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dmsnes_check .0001 709c4762a1bSJed Brown test: 710c4762a1bSJed Brown suffix: 2d_q2_quad_elas_shear 711482dd828SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_simplex 0 -deform_type shear -displacement_petscspace_degree 2 -dmsnes_check 712c4762a1bSJed Brown test: 713c4762a1bSJed Brown suffix: 2d_q3_quad_elas 71430602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dmsnes_check .0001 715c4762a1bSJed Brown test: 716c4762a1bSJed Brown suffix: 2d_q3_quad_elas_shear 717c4762a1bSJed Brown requires: !single 718482dd828SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_simplex 0 -deform_type shear -displacement_petscspace_degree 3 -dmsnes_check 719c4762a1bSJed Brown 720c4762a1bSJed Brown test: 721c4762a1bSJed Brown suffix: 3d_p1_quad_vlap 722c4762a1bSJed Brown requires: ctetgen 72330602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_refine 1 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate 724c4762a1bSJed Brown test: 725c4762a1bSJed Brown suffix: 3d_p2_quad_vlap 726c4762a1bSJed Brown requires: ctetgen 72730602db0SMatthew G. Knepley args: -dm_plex_dim 3 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001 728c4762a1bSJed Brown test: 729c4762a1bSJed Brown suffix: 3d_p3_quad_vlap 730c4762a1bSJed Brown requires: ctetgen 73130602db0SMatthew G. Knepley args: -dm_plex_dim 3 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001 732c4762a1bSJed Brown test: 733c4762a1bSJed Brown suffix: 3d_q1_quad_vlap 73430602db0SMatthew 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 735c4762a1bSJed Brown test: 736c4762a1bSJed Brown suffix: 3d_q2_quad_vlap 73730602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001 738c4762a1bSJed Brown test: 739c4762a1bSJed Brown suffix: 3d_q3_quad_vlap 74030602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001 741c4762a1bSJed Brown test: 742c4762a1bSJed Brown suffix: 3d_p1_quad_elas 743c4762a1bSJed Brown requires: ctetgen 74430602db0SMatthew 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 745c4762a1bSJed Brown test: 746c4762a1bSJed Brown suffix: 3d_p2_quad_elas 747c4762a1bSJed Brown requires: ctetgen 74830602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_dim 3 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001 749c4762a1bSJed Brown test: 750c4762a1bSJed Brown suffix: 3d_p3_quad_elas 751c4762a1bSJed Brown requires: ctetgen 75230602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_dim 3 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001 753c4762a1bSJed Brown test: 754c4762a1bSJed Brown suffix: 3d_q1_quad_elas 75530602db0SMatthew 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 756c4762a1bSJed Brown test: 757c4762a1bSJed Brown suffix: 3d_q2_quad_elas 75830602db0SMatthew 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 759c4762a1bSJed Brown test: 760c4762a1bSJed Brown suffix: 3d_q3_quad_elas 761c4762a1bSJed Brown requires: !single 76230602db0SMatthew 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 763c4762a1bSJed Brown 764c4762a1bSJed Brown test: 765c4762a1bSJed Brown suffix: 2d_p1_trig_vlap 766c4762a1bSJed Brown requires: triangle 767c4762a1bSJed Brown args: -sol_type vlap_trig -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 768c4762a1bSJed Brown test: 769c4762a1bSJed Brown suffix: 2d_p2_trig_vlap 770c4762a1bSJed Brown requires: triangle 771c4762a1bSJed Brown args: -sol_type vlap_trig -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 772c4762a1bSJed Brown test: 773c4762a1bSJed Brown suffix: 2d_p3_trig_vlap 774c4762a1bSJed Brown requires: triangle 775c4762a1bSJed Brown args: -sol_type vlap_trig -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 776c4762a1bSJed Brown test: 777c4762a1bSJed Brown suffix: 2d_q1_trig_vlap 77830602db0SMatthew 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 779c4762a1bSJed Brown test: 780c4762a1bSJed Brown suffix: 2d_q2_trig_vlap 78130602db0SMatthew 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 782c4762a1bSJed Brown test: 783c4762a1bSJed Brown suffix: 2d_q3_trig_vlap 78430602db0SMatthew 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 785c4762a1bSJed Brown test: 786c4762a1bSJed Brown suffix: 2d_p1_trig_elas 787c4762a1bSJed Brown requires: triangle 788c4762a1bSJed Brown args: -sol_type elas_trig -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 789c4762a1bSJed Brown test: 790c4762a1bSJed Brown suffix: 2d_p2_trig_elas 791c4762a1bSJed Brown requires: triangle 792c4762a1bSJed Brown args: -sol_type elas_trig -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 793c4762a1bSJed Brown test: 794c4762a1bSJed Brown suffix: 2d_p3_trig_elas 795c4762a1bSJed Brown requires: triangle 796c4762a1bSJed Brown args: -sol_type elas_trig -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 797c4762a1bSJed Brown test: 798c4762a1bSJed Brown suffix: 2d_q1_trig_elas 79930602db0SMatthew 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 800c4762a1bSJed Brown test: 801c4762a1bSJed Brown suffix: 2d_q1_trig_elas_shear 802482dd828SMatthew 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 803c4762a1bSJed Brown test: 804c4762a1bSJed Brown suffix: 2d_q2_trig_elas 80530602db0SMatthew 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 806c4762a1bSJed Brown test: 807c4762a1bSJed Brown suffix: 2d_q2_trig_elas_shear 808482dd828SMatthew 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 809c4762a1bSJed Brown test: 810c4762a1bSJed Brown suffix: 2d_q3_trig_elas 81130602db0SMatthew 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 812c4762a1bSJed Brown test: 813c4762a1bSJed Brown suffix: 2d_q3_trig_elas_shear 814482dd828SMatthew 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 815c4762a1bSJed Brown 816c4762a1bSJed Brown test: 817c4762a1bSJed Brown suffix: 3d_p1_trig_vlap 818c4762a1bSJed Brown requires: ctetgen 81930602db0SMatthew 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 820c4762a1bSJed Brown test: 821c4762a1bSJed Brown suffix: 3d_p2_trig_vlap 822c4762a1bSJed Brown requires: ctetgen 82330602db0SMatthew 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 824c4762a1bSJed Brown test: 825c4762a1bSJed Brown suffix: 3d_p3_trig_vlap 826c4762a1bSJed Brown requires: ctetgen 82730602db0SMatthew 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 828c4762a1bSJed Brown test: 829c4762a1bSJed Brown suffix: 3d_q1_trig_vlap 83030602db0SMatthew 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 831c4762a1bSJed Brown test: 832c4762a1bSJed Brown suffix: 3d_q2_trig_vlap 83330602db0SMatthew 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 834c4762a1bSJed Brown test: 835c4762a1bSJed Brown suffix: 3d_q3_trig_vlap 836c4762a1bSJed Brown requires: !__float128 83730602db0SMatthew 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 838c4762a1bSJed Brown test: 839c4762a1bSJed Brown suffix: 3d_p1_trig_elas 840c4762a1bSJed Brown requires: ctetgen 84130602db0SMatthew 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 842c4762a1bSJed Brown test: 843c4762a1bSJed Brown suffix: 3d_p2_trig_elas 844c4762a1bSJed Brown requires: ctetgen 84530602db0SMatthew 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 846c4762a1bSJed Brown test: 847c4762a1bSJed Brown suffix: 3d_p3_trig_elas 848c4762a1bSJed Brown requires: ctetgen 84930602db0SMatthew 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 850c4762a1bSJed Brown test: 851c4762a1bSJed Brown suffix: 3d_q1_trig_elas 85230602db0SMatthew 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 853c4762a1bSJed Brown test: 854c4762a1bSJed Brown suffix: 3d_q2_trig_elas 85530602db0SMatthew 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 856c4762a1bSJed Brown test: 857c4762a1bSJed Brown suffix: 3d_q3_trig_elas 858c4762a1bSJed Brown requires: !__float128 85930602db0SMatthew 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 860c4762a1bSJed Brown 861c4762a1bSJed Brown test: 862c4762a1bSJed Brown suffix: 2d_p1_axial_elas 863c4762a1bSJed Brown requires: triangle 864c4762a1bSJed Brown args: -sol_type elas_axial_disp -displacement_petscspace_degree 1 -dm_plex_separate_marker -dm_refine 2 -dmsnes_check .0001 -pc_type lu 865c4762a1bSJed Brown test: 866c4762a1bSJed Brown suffix: 2d_p2_axial_elas 867c4762a1bSJed Brown requires: triangle 868c4762a1bSJed Brown args: -sol_type elas_axial_disp -displacement_petscspace_degree 2 -dm_plex_separate_marker -dmsnes_check .0001 -pc_type lu 869c4762a1bSJed Brown test: 870c4762a1bSJed Brown suffix: 2d_p3_axial_elas 871c4762a1bSJed Brown requires: triangle 872c4762a1bSJed Brown args: -sol_type elas_axial_disp -displacement_petscspace_degree 3 -dm_plex_separate_marker -dmsnes_check .0001 -pc_type lu 873c4762a1bSJed Brown test: 874c4762a1bSJed Brown suffix: 2d_q1_axial_elas 87530602db0SMatthew 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 876c4762a1bSJed Brown test: 877c4762a1bSJed Brown suffix: 2d_q2_axial_elas 87830602db0SMatthew 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 879c4762a1bSJed Brown test: 880c4762a1bSJed Brown suffix: 2d_q3_axial_elas 88130602db0SMatthew 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 882c4762a1bSJed Brown 883c4762a1bSJed Brown test: 884c4762a1bSJed Brown suffix: 2d_p1_uniform_elas 885c4762a1bSJed Brown requires: triangle 886c4762a1bSJed Brown args: -sol_type elas_uniform_strain -displacement_petscspace_degree 1 -dm_refine 2 -dmsnes_check .0001 -pc_type lu 887c4762a1bSJed Brown test: 888c4762a1bSJed Brown suffix: 2d_p2_uniform_elas 889c4762a1bSJed Brown requires: triangle 890c4762a1bSJed Brown args: -sol_type elas_uniform_strain -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001 -pc_type lu 891c4762a1bSJed Brown test: 892c4762a1bSJed Brown suffix: 2d_p3_uniform_elas 893c4762a1bSJed Brown requires: triangle 894c4762a1bSJed Brown args: -sol_type elas_uniform_strain -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001 -pc_type lu 895c4762a1bSJed Brown test: 896c4762a1bSJed Brown suffix: 2d_q1_uniform_elas 89730602db0SMatthew 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 898c4762a1bSJed Brown test: 899c4762a1bSJed Brown suffix: 2d_q2_uniform_elas 900c4762a1bSJed Brown requires: !single 90130602db0SMatthew 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 902c4762a1bSJed Brown test: 903c4762a1bSJed Brown suffix: 2d_q3_uniform_elas 904c4762a1bSJed Brown requires: !single 90530602db0SMatthew 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 906482dd828SMatthew G. Knepley test: 907482dd828SMatthew G. Knepley suffix: 2d_p1_uniform_elas_step 908482dd828SMatthew G. Knepley requires: triangle 909482dd828SMatthew G. Knepley args: -sol_type elas_uniform_strain -deform_type step -displacement_petscspace_degree 1 -dm_refine 2 -dmsnes_check .0001 -pc_type lu 910482dd828SMatthew G. Knepley 911482dd828SMatthew G. Knepley testset: 912482dd828SMatthew 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 913482dd828SMatthew G. Knepley 914482dd828SMatthew G. Knepley test: 915482dd828SMatthew G. Knepley suffix: 2d_q1_uniform_elas_step 916482dd828SMatthew G. Knepley args: -sol_type elas_uniform_strain -dm_refine 2 917482dd828SMatthew G. Knepley test: 918482dd828SMatthew G. Knepley suffix: 2d_q1_quad_vlap_step 919482dd828SMatthew G. Knepley args: 920482dd828SMatthew G. Knepley test: 921482dd828SMatthew G. Knepley suffix: 2d_q2_quad_vlap_step 922482dd828SMatthew G. Knepley args: -displacement_petscspace_degree 2 923482dd828SMatthew G. Knepley test: 924482dd828SMatthew G. Knepley suffix: 2d_q1_quad_mass_step 925482dd828SMatthew G. Knepley args: -sol_type mass_quad 926c4762a1bSJed Brown 927d12cd81dSMatthew G. Knepley testset: 928908b9b43SStefano Zampini filter: grep -v "variant HERMITIAN" 929d12cd81dSMatthew 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 \ 930d12cd81dSMatthew G. Knepley -dm_plex_box_faces 5,5,2 -dm_plex_separate_marker -dm_refine 0 -petscpartitioner_type simple \ 931908b9b43SStefano Zampini -sol_type elas_ge 932908b9b43SStefano Zampini 933908b9b43SStefano Zampini test: 934908b9b43SStefano Zampini suffix: ge_q1_0 935908b9b43SStefano Zampini args: -displacement_petscspace_degree 1 \ 936d12cd81dSMatthew G. Knepley -snes_max_it 2 -snes_rtol 1.e-10 \ 937d12cd81dSMatthew G. Knepley -ksp_type cg -ksp_rtol 1.e-10 -ksp_max_it 100 -ksp_norm_type unpreconditioned \ 938d12cd81dSMatthew G. Knepley -pc_type gamg -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 \ 939d12cd81dSMatthew G. Knepley -pc_gamg_coarse_eq_limit 10 -pc_gamg_reuse_interpolation true \ 940*bae903cbSmarkadams4 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 \ 941d12cd81dSMatthew 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 \ 942d12cd81dSMatthew G. Knepley -matptap_via scalable 943d12cd81dSMatthew G. Knepley test: 944908b9b43SStefano Zampini nsize: 5 945908b9b43SStefano Zampini suffix: ge_q1_gdsw 946908b9b43SStefano 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 947d12cd81dSMatthew G. Knepley 948c4762a1bSJed Brown TEST*/ 949