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