1c4762a1bSJed Brown static char help[] = "Linear elasticity in 2d and 3d with finite elements.\n\ 2c4762a1bSJed Brown We solve the elasticity problem in a rectangular\n\ 3c4762a1bSJed Brown domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\ 4c4762a1bSJed Brown This example supports automatic convergence estimation\n\ 5c4762a1bSJed Brown and eventually adaptivity.\n\n\n"; 6c4762a1bSJed Brown 7c4762a1bSJed Brown /* 8c4762a1bSJed Brown https://en.wikipedia.org/wiki/Linear_elasticity 9c4762a1bSJed Brown */ 10c4762a1bSJed Brown 11c4762a1bSJed Brown #include <petscdmplex.h> 12c4762a1bSJed Brown #include <petscsnes.h> 13c4762a1bSJed Brown #include <petscds.h> 14c4762a1bSJed Brown #include <petscconvest.h> 15c4762a1bSJed Brown 16c4762a1bSJed Brown typedef enum {SOL_VLAP_QUADRATIC, SOL_ELAS_QUADRATIC, SOL_VLAP_TRIG, SOL_ELAS_TRIG, SOL_ELAS_AXIAL_DISP, SOL_ELAS_UNIFORM_STRAIN, NUM_SOLUTION_TYPES} SolutionType; 17c4762a1bSJed Brown const char *solutionTypes[NUM_SOLUTION_TYPES+1] = {"vlap_quad", "elas_quad", "vlap_trig", "elas_trig", "elas_axial_disp", "elas_uniform_strain", "unknown"}; 18c4762a1bSJed Brown 19c4762a1bSJed Brown typedef struct { 20c4762a1bSJed Brown /* Domain and mesh definition */ 21c4762a1bSJed Brown char dmType[256]; /* DM type for the solve */ 22c4762a1bSJed Brown PetscBool shear; /* Shear the domain */ 23c4762a1bSJed Brown /* Problem definition */ 24c4762a1bSJed Brown SolutionType solType; /* Type of exact solution */ 25c4762a1bSJed Brown /* Solver definition */ 26c4762a1bSJed Brown PetscBool useNearNullspace; /* Use the rigid body modes as a near nullspace for AMG */ 27c4762a1bSJed Brown } AppCtx; 28c4762a1bSJed Brown 29c4762a1bSJed Brown static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 30c4762a1bSJed Brown { 31c4762a1bSJed Brown PetscInt d; 32c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[d] = 0.0; 33c4762a1bSJed Brown return 0; 34c4762a1bSJed Brown } 35c4762a1bSJed Brown 36c4762a1bSJed Brown static PetscErrorCode quadratic_2d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 37c4762a1bSJed Brown { 38c4762a1bSJed Brown u[0] = x[0]*x[0]; 39c4762a1bSJed Brown u[1] = x[1]*x[1] - 2.0*x[0]*x[1]; 40c4762a1bSJed Brown return 0; 41c4762a1bSJed Brown } 42c4762a1bSJed Brown 43c4762a1bSJed Brown static PetscErrorCode quadratic_3d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 44c4762a1bSJed Brown { 45c4762a1bSJed Brown u[0] = x[0]*x[0]; 46c4762a1bSJed Brown u[1] = x[1]*x[1] - 2.0*x[0]*x[1]; 47c4762a1bSJed Brown u[2] = x[2]*x[2] - 2.0*x[1]*x[2]; 48c4762a1bSJed Brown return 0; 49c4762a1bSJed Brown } 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* 52c4762a1bSJed Brown u = x^2 53c4762a1bSJed Brown v = y^2 - 2xy 54c4762a1bSJed Brown Delta <u,v> - f = <2, 2> - <2, 2> 55c4762a1bSJed Brown 56c4762a1bSJed Brown u = x^2 57c4762a1bSJed Brown v = y^2 - 2xy 58c4762a1bSJed Brown w = z^2 - 2yz 59c4762a1bSJed Brown Delta <u,v,w> - f = <2, 2, 2> - <2, 2, 2> 60c4762a1bSJed Brown */ 61c4762a1bSJed Brown static void f0_vlap_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 62c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 63c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 64c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 65c4762a1bSJed Brown { 66c4762a1bSJed Brown PetscInt d; 67c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[d] += 2.0; 68c4762a1bSJed Brown } 69c4762a1bSJed Brown 70c4762a1bSJed Brown /* 71c4762a1bSJed Brown u = x^2 72c4762a1bSJed Brown v = y^2 - 2xy 73c4762a1bSJed Brown \varepsilon = / 2x -y \ 74c4762a1bSJed Brown \ -y 2y - 2x / 75c4762a1bSJed Brown Tr(\varepsilon) = div u = 2y 76c4762a1bSJed Brown div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij} 77c4762a1bSJed Brown = \lambda \partial_j (2y) + 2\mu < 2-1, 2 > 78c4762a1bSJed Brown = \lambda < 0, 2 > + \mu < 2, 4 > 79c4762a1bSJed Brown 80c4762a1bSJed Brown u = x^2 81c4762a1bSJed Brown v = y^2 - 2xy 82c4762a1bSJed Brown w = z^2 - 2yz 83c4762a1bSJed Brown \varepsilon = / 2x -y 0 \ 84c4762a1bSJed Brown | -y 2y - 2x -z | 85c4762a1bSJed Brown \ 0 -z 2z - 2y/ 86c4762a1bSJed Brown Tr(\varepsilon) = div u = 2z 87c4762a1bSJed Brown div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij} 88c4762a1bSJed Brown = \lambda \partial_j (2z) + 2\mu < 2-1, 2-1, 2 > 89c4762a1bSJed Brown = \lambda < 0, 0, 2 > + \mu < 2, 2, 4 > 90c4762a1bSJed Brown */ 91c4762a1bSJed Brown static void f0_elas_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 92c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 93c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 94c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 95c4762a1bSJed Brown { 96c4762a1bSJed Brown const PetscReal mu = 1.0; 97c4762a1bSJed Brown const PetscReal lambda = 1.0; 98c4762a1bSJed Brown PetscInt d; 99c4762a1bSJed Brown 100c4762a1bSJed Brown for (d = 0; d < dim-1; ++d) f0[d] += 2.0*mu; 101c4762a1bSJed Brown f0[dim-1] += 2.0*lambda + 4.0*mu; 102c4762a1bSJed Brown } 103c4762a1bSJed Brown 104c4762a1bSJed Brown static PetscErrorCode trig_2d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 105c4762a1bSJed Brown { 106c4762a1bSJed Brown u[0] = PetscSinReal(2.0*PETSC_PI*x[0]); 107c4762a1bSJed Brown u[1] = PetscSinReal(2.0*PETSC_PI*x[1]) - 2.0*x[0]*x[1]; 108c4762a1bSJed Brown return 0; 109c4762a1bSJed Brown } 110c4762a1bSJed Brown 111c4762a1bSJed Brown static PetscErrorCode trig_3d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 112c4762a1bSJed Brown { 113c4762a1bSJed Brown u[0] = PetscSinReal(2.0*PETSC_PI*x[0]); 114c4762a1bSJed Brown u[1] = PetscSinReal(2.0*PETSC_PI*x[1]) - 2.0*x[0]*x[1]; 115c4762a1bSJed Brown u[2] = PetscSinReal(2.0*PETSC_PI*x[2]) - 2.0*x[1]*x[2]; 116c4762a1bSJed Brown return 0; 117c4762a1bSJed Brown } 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* 120c4762a1bSJed Brown u = sin(2 pi x) 121c4762a1bSJed Brown v = sin(2 pi y) - 2xy 122c4762a1bSJed 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)> 123c4762a1bSJed Brown 124c4762a1bSJed Brown u = sin(2 pi x) 125c4762a1bSJed Brown v = sin(2 pi y) - 2xy 126c4762a1bSJed Brown w = sin(2 pi z) - 2yz 127c4762a1bSJed 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)> 128c4762a1bSJed Brown */ 129c4762a1bSJed Brown static void f0_vlap_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 130c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 131c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 132c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 133c4762a1bSJed Brown { 134c4762a1bSJed Brown PetscInt d; 135c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[d] += -4.0*PetscSqr(PETSC_PI)*PetscSinReal(2.0*PETSC_PI*x[d]); 136c4762a1bSJed Brown } 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* 139c4762a1bSJed Brown u = sin(2 pi x) 140c4762a1bSJed Brown v = sin(2 pi y) - 2xy 141c4762a1bSJed Brown \varepsilon = / 2 pi cos(2 pi x) -y \ 142c4762a1bSJed Brown \ -y 2 pi cos(2 pi y) - 2x / 143c4762a1bSJed Brown Tr(\varepsilon) = div u = 2 pi (cos(2 pi x) + cos(2 pi y)) - 2 x 144c4762a1bSJed Brown div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij} 145c4762a1bSJed 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) > 146c4762a1bSJed 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) > 147c4762a1bSJed Brown 148c4762a1bSJed Brown u = sin(2 pi x) 149c4762a1bSJed Brown v = sin(2 pi y) - 2xy 150c4762a1bSJed Brown w = sin(2 pi z) - 2yz 151c4762a1bSJed Brown \varepsilon = / 2 pi cos(2 pi x) -y 0 \ 152c4762a1bSJed Brown | -y 2 pi cos(2 pi y) - 2x -z | 153c4762a1bSJed Brown \ 0 -z 2 pi cos(2 pi z) - 2y / 154c4762a1bSJed Brown Tr(\varepsilon) = div u = 2 pi (cos(2 pi x) + cos(2 pi y) + cos(2 pi z)) - 2 x - 2 y 155c4762a1bSJed Brown div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij} 156c4762a1bSJed 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) > 157c4762a1bSJed 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) > 158c4762a1bSJed Brown */ 159c4762a1bSJed Brown static void f0_elas_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 160c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 161c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 162c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 163c4762a1bSJed Brown { 164c4762a1bSJed Brown const PetscReal mu = 1.0; 165c4762a1bSJed Brown const PetscReal lambda = 1.0; 166c4762a1bSJed Brown const PetscReal fact = 4.0*PetscSqr(PETSC_PI); 167c4762a1bSJed Brown PetscInt d; 168c4762a1bSJed Brown 169c4762a1bSJed 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); 170c4762a1bSJed Brown } 171c4762a1bSJed Brown 172c4762a1bSJed Brown static PetscErrorCode axial_disp_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 173c4762a1bSJed Brown { 174c4762a1bSJed Brown const PetscReal mu = 1.0; 175c4762a1bSJed Brown const PetscReal lambda = 1.0; 176c4762a1bSJed Brown const PetscReal N = 1.0; 177c4762a1bSJed Brown PetscInt d; 178c4762a1bSJed Brown 179c4762a1bSJed 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]; 180c4762a1bSJed Brown u[1] = -0.25*lambda/mu/(lambda+mu)*N*x[1]; 181c4762a1bSJed Brown for (d = 2; d < dim; ++d) u[d] = 0.0; 182c4762a1bSJed Brown return 0; 183c4762a1bSJed Brown } 184c4762a1bSJed Brown 185c4762a1bSJed Brown /* 186c4762a1bSJed Brown We will pull/push on the right side of a block of linearly elastic material. The uniform traction conditions on the 187c4762a1bSJed Brown right side of the box will result in a uniform strain along x and y. The Neumann BC is given by 188c4762a1bSJed Brown 189c4762a1bSJed Brown n_i \sigma_{ij} = t_i 190c4762a1bSJed Brown 191c4762a1bSJed Brown u = (1/(2\mu) - 1) x 192c4762a1bSJed Brown v = -y 193c4762a1bSJed Brown f = 0 194c4762a1bSJed Brown t = <4\mu/\lambda (\lambda + \mu), 0> 195c4762a1bSJed Brown \varepsilon = / 1/(2\mu) - 1 0 \ 196c4762a1bSJed Brown \ 0 -1 / 197c4762a1bSJed Brown Tr(\varepsilon) = div u = 1/(2\mu) - 2 198c4762a1bSJed Brown div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij} 199c4762a1bSJed Brown = \lambda \partial_j (1/(2\mu) - 2) + 2\mu < 0, 0 > 200c4762a1bSJed Brown = \lambda < 0, 0 > + \mu < 0, 0 > = 0 201c4762a1bSJed Brown NBC = <1,0> . <4\mu/\lambda (\lambda + \mu), 0> = 4\mu/\lambda (\lambda + \mu) 202c4762a1bSJed Brown 203c4762a1bSJed Brown u = x - 1/2 204c4762a1bSJed Brown v = 0 205c4762a1bSJed Brown w = 0 206c4762a1bSJed Brown \varepsilon = / x 0 0 \ 207c4762a1bSJed Brown | 0 0 0 | 208c4762a1bSJed Brown \ 0 0 0 / 209c4762a1bSJed Brown Tr(\varepsilon) = div u = x 210c4762a1bSJed Brown div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij} 211c4762a1bSJed Brown = \lambda \partial_j x + 2\mu < 1, 0, 0 > 212c4762a1bSJed Brown = \lambda < 1, 0, 0 > + \mu < 2, 0, 0 > 213c4762a1bSJed Brown */ 214c4762a1bSJed Brown static void f0_elas_axial_disp_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 215c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 216c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 217c4762a1bSJed Brown PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 218c4762a1bSJed Brown { 219c4762a1bSJed Brown const PetscReal N = -1.0; 220c4762a1bSJed Brown 221c4762a1bSJed Brown f0[0] = N; 222c4762a1bSJed Brown } 223c4762a1bSJed Brown 224c4762a1bSJed Brown static PetscErrorCode uniform_strain_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 225c4762a1bSJed Brown { 226c4762a1bSJed Brown const PetscReal eps_xx = 0.1; 227c4762a1bSJed Brown const PetscReal eps_xy = 0.3; 228c4762a1bSJed Brown const PetscReal eps_yy = 0.25; 229c4762a1bSJed Brown PetscInt d; 230c4762a1bSJed Brown 231c4762a1bSJed Brown u[0] = eps_xx*x[0] + eps_xy*x[1]; 232c4762a1bSJed Brown u[1] = eps_xy*x[0] + eps_yy*x[1]; 233c4762a1bSJed Brown for (d = 2; d < dim; ++d) u[d] = 0.0; 234c4762a1bSJed Brown return 0; 235c4762a1bSJed Brown } 236c4762a1bSJed Brown 237c4762a1bSJed Brown static void f1_vlap_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 238c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 239c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 240c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 241c4762a1bSJed Brown { 242c4762a1bSJed Brown const PetscInt Nc = dim; 243c4762a1bSJed Brown PetscInt c, d; 244c4762a1bSJed Brown 245c4762a1bSJed Brown for (c = 0; c < Nc; ++c) for (d = 0; d < dim; ++d) f1[c*dim+d] += u_x[c*dim+d]; 246c4762a1bSJed Brown } 247c4762a1bSJed Brown 248c4762a1bSJed Brown static void f1_elas_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, 249c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 250c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 251c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 252c4762a1bSJed Brown { 253c4762a1bSJed Brown const PetscInt Nc = dim; 254c4762a1bSJed Brown const PetscReal mu = 1.0; 255c4762a1bSJed Brown const PetscReal lambda = 1.0; 256c4762a1bSJed Brown PetscInt c, d; 257c4762a1bSJed Brown 258c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 259c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 260c4762a1bSJed Brown f1[c*dim+d] += mu*(u_x[c*dim+d] + u_x[d*dim+c]); 261c4762a1bSJed Brown f1[c*dim+c] += lambda*u_x[d*dim+d]; 262c4762a1bSJed Brown } 263c4762a1bSJed Brown } 264c4762a1bSJed Brown } 265c4762a1bSJed Brown 266c4762a1bSJed Brown static void g3_vlap_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 267c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 268c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 269c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 270c4762a1bSJed Brown { 271c4762a1bSJed Brown const PetscInt Nc = dim; 272c4762a1bSJed Brown PetscInt c, d; 273c4762a1bSJed Brown 274c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 275c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 276c4762a1bSJed Brown g3[((c*Nc + c)*dim + d)*dim + d] = 1.0; 277c4762a1bSJed Brown } 278c4762a1bSJed Brown } 279c4762a1bSJed Brown } 280c4762a1bSJed Brown 281c4762a1bSJed Brown /* 282c4762a1bSJed Brown \partial_df \phi_fc g_{fc,gc,df,dg} \partial_dg \phi_gc 283c4762a1bSJed Brown 284c4762a1bSJed Brown \partial_df \phi_fc \lambda \delta_{fc,df} \sum_gc \partial_dg \phi_gc \delta_{gc,dg} 285c4762a1bSJed Brown = \partial_fc \phi_fc \sum_gc \partial_gc \phi_gc 286c4762a1bSJed Brown */ 287c4762a1bSJed Brown static void g3_elas_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, 288c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 289c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 290c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 291c4762a1bSJed Brown { 292c4762a1bSJed Brown const PetscInt Nc = dim; 293c4762a1bSJed Brown const PetscReal mu = 1.0; 294c4762a1bSJed Brown const PetscReal lambda = 1.0; 295c4762a1bSJed Brown PetscInt c, d; 296c4762a1bSJed Brown 297c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 298c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 299c4762a1bSJed Brown g3[((c*Nc + c)*dim + d)*dim + d] += mu; 300c4762a1bSJed Brown g3[((c*Nc + d)*dim + d)*dim + c] += mu; 301c4762a1bSJed Brown g3[((c*Nc + d)*dim + c)*dim + d] += lambda; 302c4762a1bSJed Brown } 303c4762a1bSJed Brown } 304c4762a1bSJed Brown } 305c4762a1bSJed Brown 306c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 307c4762a1bSJed Brown { 308*30602db0SMatthew G. Knepley PetscInt sol; 309c4762a1bSJed Brown PetscErrorCode ierr; 310c4762a1bSJed Brown 311c4762a1bSJed Brown PetscFunctionBeginUser; 312c4762a1bSJed Brown options->shear = PETSC_FALSE; 313c4762a1bSJed Brown options->solType = SOL_VLAP_QUADRATIC; 314c4762a1bSJed Brown options->useNearNullspace = PETSC_TRUE; 315c4762a1bSJed Brown ierr = PetscStrncpy(options->dmType, DMPLEX, 256);CHKERRQ(ierr); 316c4762a1bSJed Brown 317c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Linear Elasticity Problem Options", "DMPLEX");CHKERRQ(ierr); 318c4762a1bSJed Brown ierr = PetscOptionsBool("-shear", "Shear the domain", "ex17.c", options->shear, &options->shear, NULL);CHKERRQ(ierr); 319c4762a1bSJed Brown sol = options->solType; 320c4762a1bSJed Brown ierr = PetscOptionsEList("-sol_type", "Type of exact solution", "ex17.c", solutionTypes, NUM_SOLUTION_TYPES, solutionTypes[options->solType], &sol, NULL);CHKERRQ(ierr); 321c4762a1bSJed Brown options->solType = (SolutionType) sol; 322c4762a1bSJed Brown ierr = PetscOptionsBool("-near_nullspace", "Use the rigid body modes as an AMG near nullspace", "ex17.c", options->useNearNullspace, &options->useNearNullspace, NULL);CHKERRQ(ierr); 323c4762a1bSJed Brown ierr = PetscOptionsFList("-dm_type", "Convert DMPlex to another format", "ex17.c", DMList, options->dmType, options->dmType, 256, NULL);CHKERRQ(ierr); 324c4762a1bSJed Brown ierr = PetscOptionsEnd(); 325c4762a1bSJed Brown PetscFunctionReturn(0); 326c4762a1bSJed Brown } 327c4762a1bSJed Brown 328c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 329c4762a1bSJed Brown { 330c4762a1bSJed Brown PetscErrorCode ierr; 331c4762a1bSJed Brown 332c4762a1bSJed Brown PetscFunctionBeginUser; 333*30602db0SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 334*30602db0SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 335c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 336*30602db0SMatthew G. Knepley if (user->shear) {ierr = DMPlexShearGeometry(*dm, DM_X, NULL);CHKERRQ(ierr);} 337*30602db0SMatthew G. Knepley ierr = DMSetApplicationContext(*dm, user);CHKERRQ(ierr); 338c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 339c4762a1bSJed Brown PetscFunctionReturn(0); 340c4762a1bSJed Brown } 341c4762a1bSJed Brown 342c4762a1bSJed Brown static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 343c4762a1bSJed Brown { 344c4762a1bSJed Brown PetscErrorCode (*exact)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 345c4762a1bSJed Brown PetscDS prob; 34645480ffeSMatthew G. Knepley PetscWeakForm wf; 34745480ffeSMatthew G. Knepley DMLabel label; 34845480ffeSMatthew G. Knepley PetscInt id, bd; 349c4762a1bSJed Brown PetscInt dim; 350c4762a1bSJed Brown PetscErrorCode ierr; 351c4762a1bSJed Brown 352c4762a1bSJed Brown PetscFunctionBeginUser; 353c4762a1bSJed Brown ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 354c4762a1bSJed Brown ierr = PetscDSGetSpatialDimension(prob, &dim);CHKERRQ(ierr); 355c4762a1bSJed Brown switch (user->solType) { 356c4762a1bSJed Brown case SOL_VLAP_QUADRATIC: 357c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 0, f0_vlap_quadratic_u, f1_vlap_u);CHKERRQ(ierr); 358c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_vlap_uu);CHKERRQ(ierr); 359c4762a1bSJed Brown switch (dim) { 360c4762a1bSJed Brown case 2: exact = quadratic_2d_u;break; 361c4762a1bSJed Brown case 3: exact = quadratic_3d_u;break; 362c4762a1bSJed Brown default: SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Invalid dimension: %D", dim); 363c4762a1bSJed Brown } 364c4762a1bSJed Brown break; 365c4762a1bSJed Brown case SOL_ELAS_QUADRATIC: 366c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 0, f0_elas_quadratic_u, f1_elas_u);CHKERRQ(ierr); 367c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_elas_uu);CHKERRQ(ierr); 368c4762a1bSJed Brown switch (dim) { 369c4762a1bSJed Brown case 2: exact = quadratic_2d_u;break; 370c4762a1bSJed Brown case 3: exact = quadratic_3d_u;break; 371c4762a1bSJed Brown default: SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Invalid dimension: %D", dim); 372c4762a1bSJed Brown } 373c4762a1bSJed Brown break; 374c4762a1bSJed Brown case SOL_VLAP_TRIG: 375c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 0, f0_vlap_trig_u, f1_vlap_u);CHKERRQ(ierr); 376c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_vlap_uu);CHKERRQ(ierr); 377c4762a1bSJed Brown switch (dim) { 378c4762a1bSJed Brown case 2: exact = trig_2d_u;break; 379c4762a1bSJed Brown case 3: exact = trig_3d_u;break; 380c4762a1bSJed Brown default: SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Invalid dimension: %D", dim); 381c4762a1bSJed Brown } 382c4762a1bSJed Brown break; 383c4762a1bSJed Brown case SOL_ELAS_TRIG: 384c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 0, f0_elas_trig_u, f1_elas_u);CHKERRQ(ierr); 385c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_elas_uu);CHKERRQ(ierr); 386c4762a1bSJed Brown switch (dim) { 387c4762a1bSJed Brown case 2: exact = trig_2d_u;break; 388c4762a1bSJed Brown case 3: exact = trig_3d_u;break; 389c4762a1bSJed Brown default: SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Invalid dimension: %D", dim); 390c4762a1bSJed Brown } 391c4762a1bSJed Brown break; 392c4762a1bSJed Brown case SOL_ELAS_AXIAL_DISP: 393c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 0, NULL, f1_elas_u);CHKERRQ(ierr); 394c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_elas_uu);CHKERRQ(ierr); 39545480ffeSMatthew G. Knepley id = dim == 3 ? 5 : 2; 39645480ffeSMatthew G. Knepley ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr); 39745480ffeSMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_NATURAL, "right", label, 1, &id, 0, 0, NULL, (void (*)(void)) NULL, NULL, user, &bd);CHKERRQ(ierr); 39845480ffeSMatthew G. Knepley ierr = PetscDSGetBoundary(prob, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 39945480ffeSMatthew G. Knepley ierr = PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, f0_elas_axial_disp_bd_u, 0, NULL);CHKERRQ(ierr); 400c4762a1bSJed Brown exact = axial_disp_u; 401c4762a1bSJed Brown break; 402c4762a1bSJed Brown case SOL_ELAS_UNIFORM_STRAIN: 403c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 0, NULL, f1_elas_u);CHKERRQ(ierr); 404c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_elas_uu);CHKERRQ(ierr); 405c4762a1bSJed Brown exact = uniform_strain_u; 406c4762a1bSJed Brown break; 407c4762a1bSJed Brown default: SETERRQ2(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Invalid solution type: %s (%D)", solutionTypes[PetscMin(user->solType, NUM_SOLUTION_TYPES)], user->solType); 408c4762a1bSJed Brown } 409c4762a1bSJed Brown ierr = PetscDSSetExactSolution(prob, 0, exact, user);CHKERRQ(ierr); 41045480ffeSMatthew G. Knepley ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr); 411c4762a1bSJed Brown if (user->solType == SOL_ELAS_AXIAL_DISP) { 412c4762a1bSJed Brown PetscInt cmp; 413c4762a1bSJed Brown 414c4762a1bSJed Brown id = dim == 3 ? 6 : 4; 415c4762a1bSJed Brown cmp = 0; 41645480ffeSMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "left", label, 1, &id, 0, 1, &cmp, (void (*)(void)) zero, NULL, user, NULL);CHKERRQ(ierr); 417c4762a1bSJed Brown cmp = dim == 3 ? 2 : 1; 418c4762a1bSJed Brown id = dim == 3 ? 1 : 1; 41945480ffeSMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "bottom", label, 1, &id, 0, 1, &cmp, (void (*)(void)) zero, NULL, user, NULL);CHKERRQ(ierr); 420c4762a1bSJed Brown if (dim == 3) { 421c4762a1bSJed Brown cmp = 1; 422c4762a1bSJed Brown id = 3; 42345480ffeSMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "front", label, 1, &id, 0, 1, &cmp, (void (*)(void)) zero, NULL, user, NULL);CHKERRQ(ierr); 424c4762a1bSJed Brown } 425c4762a1bSJed Brown } else { 426c4762a1bSJed Brown id = 1; 42745480ffeSMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) exact, NULL, user, NULL);CHKERRQ(ierr); 428c4762a1bSJed Brown } 429c4762a1bSJed Brown PetscFunctionReturn(0); 430c4762a1bSJed Brown } 431c4762a1bSJed Brown 4328cda7954SMatthew G. Knepley static PetscErrorCode CreateElasticityNullSpace(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullspace) 433c4762a1bSJed Brown { 434c4762a1bSJed Brown PetscErrorCode ierr; 435c4762a1bSJed Brown 436c4762a1bSJed Brown PetscFunctionBegin; 4378cda7954SMatthew G. Knepley ierr = DMPlexCreateRigidBody(dm, origField, nullspace);CHKERRQ(ierr); 438c4762a1bSJed Brown PetscFunctionReturn(0); 439c4762a1bSJed Brown } 440c4762a1bSJed Brown 441*30602db0SMatthew G. Knepley PetscErrorCode SetupFE(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), void *ctx) 442c4762a1bSJed Brown { 443c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 444c4762a1bSJed Brown DM cdm = dm; 445c4762a1bSJed Brown PetscFE fe; 446c4762a1bSJed Brown char prefix[PETSC_MAX_PATH_LEN]; 447*30602db0SMatthew G. Knepley DMPolytopeType ct; 448*30602db0SMatthew G. Knepley PetscBool simplex; 449*30602db0SMatthew G. Knepley PetscInt dim, cStart; 450c4762a1bSJed Brown PetscErrorCode ierr; 451c4762a1bSJed Brown 452c4762a1bSJed Brown PetscFunctionBegin; 453c4762a1bSJed Brown /* Create finite element */ 454c4762a1bSJed Brown ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 455*30602db0SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr); 456*30602db0SMatthew G. Knepley ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr); 457*30602db0SMatthew G. Knepley simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE; 458c4762a1bSJed Brown ierr = PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name);CHKERRQ(ierr); 459*30602db0SMatthew G. Knepley ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, dim, simplex, name ? prefix : NULL, -1, &fe);CHKERRQ(ierr); 460c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr); 461c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 462c4762a1bSJed Brown ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 463c4762a1bSJed Brown ierr = DMCreateDS(dm);CHKERRQ(ierr); 464c4762a1bSJed Brown ierr = (*setup)(dm, user);CHKERRQ(ierr); 465c4762a1bSJed Brown while (cdm) { 466c4762a1bSJed Brown ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr); 467c4762a1bSJed Brown if (user->useNearNullspace) {ierr = DMSetNearNullSpaceConstructor(cdm, 0, CreateElasticityNullSpace);CHKERRQ(ierr);} 468c4762a1bSJed Brown /* TODO: Check whether the boundary of coarse meshes is marked */ 469c4762a1bSJed Brown ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 470c4762a1bSJed Brown } 471c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 472c4762a1bSJed Brown PetscFunctionReturn(0); 473c4762a1bSJed Brown } 474c4762a1bSJed Brown 475c4762a1bSJed Brown int main(int argc, char **argv) 476c4762a1bSJed Brown { 477c4762a1bSJed Brown DM dm; /* Problem specification */ 478c4762a1bSJed Brown SNES snes; /* Nonlinear solver */ 479c4762a1bSJed Brown Vec u; /* Solutions */ 480c4762a1bSJed Brown AppCtx user; /* User-defined work context */ 481c4762a1bSJed Brown PetscErrorCode ierr; 482c4762a1bSJed Brown 483c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 484c4762a1bSJed Brown ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 485c4762a1bSJed Brown /* Primal system */ 486c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr); 487c4762a1bSJed Brown ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 488c4762a1bSJed Brown ierr = SNESSetDM(snes, dm);CHKERRQ(ierr); 489*30602db0SMatthew G. Knepley ierr = SetupFE(dm, "displacement", SetupPrimalProblem, &user);CHKERRQ(ierr); 490c4762a1bSJed Brown ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 491c4762a1bSJed Brown ierr = VecSet(u, 0.0);CHKERRQ(ierr); 492c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) u, "displacement");CHKERRQ(ierr); 493c4762a1bSJed Brown ierr = DMPlexSetSNESLocalFEM(dm, &user, &user, &user);CHKERRQ(ierr); 494c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 495348a1646SMatthew G. Knepley ierr = DMSNESCheckFromOptions(snes, u);CHKERRQ(ierr); 496c4762a1bSJed Brown ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr); 497c4762a1bSJed Brown ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr); 498c4762a1bSJed Brown ierr = VecViewFromOptions(u, NULL, "-displacement_view");CHKERRQ(ierr); 499c4762a1bSJed Brown /* Cleanup */ 500c4762a1bSJed Brown ierr = VecDestroy(&u);CHKERRQ(ierr); 501c4762a1bSJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 502c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 503c4762a1bSJed Brown ierr = PetscFinalize(); 504c4762a1bSJed Brown return ierr; 505c4762a1bSJed Brown } 506c4762a1bSJed Brown 507c4762a1bSJed Brown /*TEST 508c4762a1bSJed Brown 509*30602db0SMatthew G. Knepley testset: 510*30602db0SMatthew G. Knepley args: -dm_plex_box_faces 1,1,1 511*30602db0SMatthew G. Knepley 512c4762a1bSJed Brown test: 513c4762a1bSJed Brown suffix: 2d_p1_quad_vlap 514c4762a1bSJed Brown requires: triangle 515c4762a1bSJed Brown args: -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate 516c4762a1bSJed Brown test: 517c4762a1bSJed Brown suffix: 2d_p2_quad_vlap 518c4762a1bSJed Brown requires: triangle 519c4762a1bSJed Brown args: -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001 520c4762a1bSJed Brown test: 521c4762a1bSJed Brown suffix: 2d_p3_quad_vlap 522c4762a1bSJed Brown requires: triangle 523c4762a1bSJed Brown args: -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001 524c4762a1bSJed Brown test: 525c4762a1bSJed Brown suffix: 2d_q1_quad_vlap 526*30602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate 527c4762a1bSJed Brown test: 528c4762a1bSJed Brown suffix: 2d_q2_quad_vlap 529*30602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001 530c4762a1bSJed Brown test: 531c4762a1bSJed Brown suffix: 2d_q3_quad_vlap 532c4762a1bSJed Brown requires: !single 533*30602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001 534c4762a1bSJed Brown test: 535c4762a1bSJed Brown suffix: 2d_p1_quad_elas 536c4762a1bSJed Brown requires: triangle 537c4762a1bSJed Brown args: -sol_type elas_quad -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate 538c4762a1bSJed Brown test: 539c4762a1bSJed Brown suffix: 2d_p2_quad_elas 540c4762a1bSJed Brown requires: triangle 541c4762a1bSJed Brown args: -sol_type elas_quad -displacement_petscspace_degree 2 -dmsnes_check .0001 542c4762a1bSJed Brown test: 543c4762a1bSJed Brown suffix: 2d_p3_quad_elas 544c4762a1bSJed Brown requires: triangle 545c4762a1bSJed Brown args: -sol_type elas_quad -displacement_petscspace_degree 3 -dmsnes_check .0001 546c4762a1bSJed Brown test: 547c4762a1bSJed Brown suffix: 2d_q1_quad_elas 548*30602db0SMatthew 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 549c4762a1bSJed Brown test: 550c4762a1bSJed Brown suffix: 2d_q1_quad_elas_shear 551*30602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_simplex 0 -shear -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 552c4762a1bSJed Brown test: 553c4762a1bSJed Brown suffix: 2d_q2_quad_elas 554*30602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dmsnes_check .0001 555c4762a1bSJed Brown test: 556c4762a1bSJed Brown suffix: 2d_q2_quad_elas_shear 557*30602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_simplex 0 -shear -displacement_petscspace_degree 2 -dmsnes_check 558c4762a1bSJed Brown test: 559c4762a1bSJed Brown suffix: 2d_q3_quad_elas 560*30602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dmsnes_check .0001 561c4762a1bSJed Brown test: 562c4762a1bSJed Brown suffix: 2d_q3_quad_elas_shear 563c4762a1bSJed Brown requires: !single 564*30602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_simplex 0 -shear -displacement_petscspace_degree 3 -dmsnes_check 565c4762a1bSJed Brown 566c4762a1bSJed Brown test: 567c4762a1bSJed Brown suffix: 3d_p1_quad_vlap 568c4762a1bSJed Brown requires: ctetgen 569*30602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_refine 1 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate 570c4762a1bSJed Brown test: 571c4762a1bSJed Brown suffix: 3d_p2_quad_vlap 572c4762a1bSJed Brown requires: ctetgen 573*30602db0SMatthew G. Knepley args: -dm_plex_dim 3 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001 574c4762a1bSJed Brown test: 575c4762a1bSJed Brown suffix: 3d_p3_quad_vlap 576c4762a1bSJed Brown requires: ctetgen 577*30602db0SMatthew G. Knepley args: -dm_plex_dim 3 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001 578c4762a1bSJed Brown test: 579c4762a1bSJed Brown suffix: 3d_q1_quad_vlap 580*30602db0SMatthew 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 581c4762a1bSJed Brown test: 582c4762a1bSJed Brown suffix: 3d_q2_quad_vlap 583*30602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001 584c4762a1bSJed Brown test: 585c4762a1bSJed Brown suffix: 3d_q3_quad_vlap 586*30602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001 587c4762a1bSJed Brown test: 588c4762a1bSJed Brown suffix: 3d_p1_quad_elas 589c4762a1bSJed Brown requires: ctetgen 590*30602db0SMatthew 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 591c4762a1bSJed Brown test: 592c4762a1bSJed Brown suffix: 3d_p2_quad_elas 593c4762a1bSJed Brown requires: ctetgen 594*30602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_dim 3 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001 595c4762a1bSJed Brown test: 596c4762a1bSJed Brown suffix: 3d_p3_quad_elas 597c4762a1bSJed Brown requires: ctetgen 598*30602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_dim 3 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001 599c4762a1bSJed Brown test: 600c4762a1bSJed Brown suffix: 3d_q1_quad_elas 601*30602db0SMatthew 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 602c4762a1bSJed Brown test: 603c4762a1bSJed Brown suffix: 3d_q2_quad_elas 604*30602db0SMatthew 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 605c4762a1bSJed Brown test: 606c4762a1bSJed Brown suffix: 3d_q3_quad_elas 607c4762a1bSJed Brown requires: !single 608*30602db0SMatthew 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 609c4762a1bSJed Brown 610c4762a1bSJed Brown test: 611c4762a1bSJed Brown suffix: 2d_p1_trig_vlap 612c4762a1bSJed Brown requires: triangle 613c4762a1bSJed Brown args: -sol_type vlap_trig -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 614c4762a1bSJed Brown test: 615c4762a1bSJed Brown suffix: 2d_p2_trig_vlap 616c4762a1bSJed Brown requires: triangle 617c4762a1bSJed Brown args: -sol_type vlap_trig -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 618c4762a1bSJed Brown test: 619c4762a1bSJed Brown suffix: 2d_p3_trig_vlap 620c4762a1bSJed Brown requires: triangle 621c4762a1bSJed Brown args: -sol_type vlap_trig -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 622c4762a1bSJed Brown test: 623c4762a1bSJed Brown suffix: 2d_q1_trig_vlap 624*30602db0SMatthew 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 625c4762a1bSJed Brown test: 626c4762a1bSJed Brown suffix: 2d_q2_trig_vlap 627*30602db0SMatthew 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 628c4762a1bSJed Brown test: 629c4762a1bSJed Brown suffix: 2d_q3_trig_vlap 630*30602db0SMatthew 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 631c4762a1bSJed Brown test: 632c4762a1bSJed Brown suffix: 2d_p1_trig_elas 633c4762a1bSJed Brown requires: triangle 634c4762a1bSJed Brown args: -sol_type elas_trig -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 635c4762a1bSJed Brown test: 636c4762a1bSJed Brown suffix: 2d_p2_trig_elas 637c4762a1bSJed Brown requires: triangle 638c4762a1bSJed Brown args: -sol_type elas_trig -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 639c4762a1bSJed Brown test: 640c4762a1bSJed Brown suffix: 2d_p3_trig_elas 641c4762a1bSJed Brown requires: triangle 642c4762a1bSJed Brown args: -sol_type elas_trig -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 643c4762a1bSJed Brown test: 644c4762a1bSJed Brown suffix: 2d_q1_trig_elas 645*30602db0SMatthew 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 646c4762a1bSJed Brown test: 647c4762a1bSJed Brown suffix: 2d_q1_trig_elas_shear 648*30602db0SMatthew G. Knepley args: -sol_type elas_trig -dm_plex_simplex 0 -shear -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 649c4762a1bSJed Brown test: 650c4762a1bSJed Brown suffix: 2d_q2_trig_elas 651*30602db0SMatthew 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 652c4762a1bSJed Brown test: 653c4762a1bSJed Brown suffix: 2d_q2_trig_elas_shear 654*30602db0SMatthew G. Knepley args: -sol_type elas_trig -dm_plex_simplex 0 -shear -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 655c4762a1bSJed Brown test: 656c4762a1bSJed Brown suffix: 2d_q3_trig_elas 657*30602db0SMatthew 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 658c4762a1bSJed Brown test: 659c4762a1bSJed Brown suffix: 2d_q3_trig_elas_shear 660*30602db0SMatthew G. Knepley args: -sol_type elas_trig -dm_plex_simplex 0 -shear -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 661c4762a1bSJed Brown 662c4762a1bSJed Brown test: 663c4762a1bSJed Brown suffix: 3d_p1_trig_vlap 664c4762a1bSJed Brown requires: ctetgen 665*30602db0SMatthew 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 666c4762a1bSJed Brown test: 667c4762a1bSJed Brown suffix: 3d_p2_trig_vlap 668c4762a1bSJed Brown requires: ctetgen 669*30602db0SMatthew 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 670c4762a1bSJed Brown test: 671c4762a1bSJed Brown suffix: 3d_p3_trig_vlap 672c4762a1bSJed Brown requires: ctetgen 673*30602db0SMatthew 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 674c4762a1bSJed Brown test: 675c4762a1bSJed Brown suffix: 3d_q1_trig_vlap 676*30602db0SMatthew 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 677c4762a1bSJed Brown test: 678c4762a1bSJed Brown suffix: 3d_q2_trig_vlap 679*30602db0SMatthew 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 680c4762a1bSJed Brown test: 681c4762a1bSJed Brown suffix: 3d_q3_trig_vlap 682c4762a1bSJed Brown requires: !__float128 683*30602db0SMatthew 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 684c4762a1bSJed Brown test: 685c4762a1bSJed Brown suffix: 3d_p1_trig_elas 686c4762a1bSJed Brown requires: ctetgen 687*30602db0SMatthew 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 688c4762a1bSJed Brown test: 689c4762a1bSJed Brown suffix: 3d_p2_trig_elas 690c4762a1bSJed Brown requires: ctetgen 691*30602db0SMatthew 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 692c4762a1bSJed Brown test: 693c4762a1bSJed Brown suffix: 3d_p3_trig_elas 694c4762a1bSJed Brown requires: ctetgen 695*30602db0SMatthew 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 696c4762a1bSJed Brown test: 697c4762a1bSJed Brown suffix: 3d_q1_trig_elas 698*30602db0SMatthew 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 699c4762a1bSJed Brown test: 700c4762a1bSJed Brown suffix: 3d_q2_trig_elas 701*30602db0SMatthew 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 702c4762a1bSJed Brown test: 703c4762a1bSJed Brown suffix: 3d_q3_trig_elas 704c4762a1bSJed Brown requires: !__float128 705*30602db0SMatthew 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 706c4762a1bSJed Brown 707c4762a1bSJed Brown test: 708c4762a1bSJed Brown suffix: 2d_p1_axial_elas 709c4762a1bSJed Brown requires: triangle 710c4762a1bSJed Brown args: -sol_type elas_axial_disp -displacement_petscspace_degree 1 -dm_plex_separate_marker -dm_refine 2 -dmsnes_check .0001 -pc_type lu 711c4762a1bSJed Brown test: 712c4762a1bSJed Brown suffix: 2d_p2_axial_elas 713c4762a1bSJed Brown requires: triangle 714c4762a1bSJed Brown args: -sol_type elas_axial_disp -displacement_petscspace_degree 2 -dm_plex_separate_marker -dmsnes_check .0001 -pc_type lu 715c4762a1bSJed Brown test: 716c4762a1bSJed Brown suffix: 2d_p3_axial_elas 717c4762a1bSJed Brown requires: triangle 718c4762a1bSJed Brown args: -sol_type elas_axial_disp -displacement_petscspace_degree 3 -dm_plex_separate_marker -dmsnes_check .0001 -pc_type lu 719c4762a1bSJed Brown test: 720c4762a1bSJed Brown suffix: 2d_q1_axial_elas 721*30602db0SMatthew 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 722c4762a1bSJed Brown test: 723c4762a1bSJed Brown suffix: 2d_q2_axial_elas 724*30602db0SMatthew 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 725c4762a1bSJed Brown test: 726c4762a1bSJed Brown suffix: 2d_q3_axial_elas 727*30602db0SMatthew 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 728c4762a1bSJed Brown 729c4762a1bSJed Brown test: 730c4762a1bSJed Brown suffix: 2d_p1_uniform_elas 731c4762a1bSJed Brown requires: triangle 732c4762a1bSJed Brown args: -sol_type elas_uniform_strain -displacement_petscspace_degree 1 -dm_refine 2 -dmsnes_check .0001 -pc_type lu 733c4762a1bSJed Brown test: 734c4762a1bSJed Brown suffix: 2d_p2_uniform_elas 735c4762a1bSJed Brown requires: triangle 736c4762a1bSJed Brown args: -sol_type elas_uniform_strain -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001 -pc_type lu 737c4762a1bSJed Brown test: 738c4762a1bSJed Brown suffix: 2d_p3_uniform_elas 739c4762a1bSJed Brown requires: triangle 740c4762a1bSJed Brown args: -sol_type elas_uniform_strain -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001 -pc_type lu 741c4762a1bSJed Brown test: 742c4762a1bSJed Brown suffix: 2d_q1_uniform_elas 743*30602db0SMatthew 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 744c4762a1bSJed Brown test: 745c4762a1bSJed Brown suffix: 2d_q2_uniform_elas 746c4762a1bSJed Brown requires: !single 747*30602db0SMatthew 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 748c4762a1bSJed Brown test: 749c4762a1bSJed Brown suffix: 2d_q3_uniform_elas 750c4762a1bSJed Brown requires: !single 751*30602db0SMatthew 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 752c4762a1bSJed Brown 753c4762a1bSJed Brown TEST*/ 754