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