1c4762a1bSJed Brown static char help[] = "Linear elasticity in 2d and 3d with finite elements.\n\ 2c4762a1bSJed Brown We solve the elasticity problem in a rectangular\n\ 3c4762a1bSJed Brown domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\ 4c4762a1bSJed Brown This example supports automatic convergence estimation\n\ 5c4762a1bSJed Brown and eventually adaptivity.\n\n\n"; 6c4762a1bSJed Brown 7c4762a1bSJed Brown /* 8c4762a1bSJed Brown https://en.wikipedia.org/wiki/Linear_elasticity 9d12cd81dSMatthew G. Knepley 10d12cd81dSMatthew G. Knepley Converting elastic constants: 11d12cd81dSMatthew G. Knepley lambda = E nu / ((1 + nu) (1 - 2 nu)) 12d12cd81dSMatthew G. Knepley mu = E / (2 (1 + nu)) 13c4762a1bSJed Brown */ 14c4762a1bSJed Brown 15c4762a1bSJed Brown #include <petscdmplex.h> 16c4762a1bSJed Brown #include <petscsnes.h> 17c4762a1bSJed Brown #include <petscds.h> 18d12cd81dSMatthew G. Knepley #include <petscbag.h> 19c4762a1bSJed Brown #include <petscconvest.h> 20c4762a1bSJed Brown 219371c9d4SSatish Balay typedef enum { 229371c9d4SSatish Balay SOL_VLAP_QUADRATIC, 239371c9d4SSatish Balay SOL_ELAS_QUADRATIC, 249371c9d4SSatish Balay SOL_VLAP_TRIG, 259371c9d4SSatish Balay SOL_ELAS_TRIG, 269371c9d4SSatish Balay SOL_ELAS_AXIAL_DISP, 279371c9d4SSatish Balay SOL_ELAS_UNIFORM_STRAIN, 289371c9d4SSatish Balay SOL_ELAS_GE, 299371c9d4SSatish Balay SOL_MASS_QUADRATIC, 309371c9d4SSatish Balay NUM_SOLUTION_TYPES 319371c9d4SSatish Balay } SolutionType; 32d12cd81dSMatthew G. Knepley const char *solutionTypes[NUM_SOLUTION_TYPES + 1] = {"vlap_quad", "elas_quad", "vlap_trig", "elas_trig", "elas_axial_disp", "elas_uniform_strain", "elas_ge", "mass_quad", "unknown"}; 33482dd828SMatthew G. Knepley 349371c9d4SSatish Balay typedef enum { 359371c9d4SSatish Balay DEFORM_NONE, 369371c9d4SSatish Balay DEFORM_SHEAR, 379371c9d4SSatish Balay DEFORM_STEP, 389371c9d4SSatish Balay NUM_DEFORM_TYPES 399371c9d4SSatish Balay } DeformType; 40482dd828SMatthew G. Knepley const char *deformTypes[NUM_DEFORM_TYPES + 1] = {"none", "shear", "step", "unknown"}; 41c4762a1bSJed Brown 42c4762a1bSJed Brown typedef struct { 43d12cd81dSMatthew G. Knepley PetscScalar mu; /* shear modulus */ 44d12cd81dSMatthew G. Knepley PetscScalar lambda; /* Lame's first parameter */ 45f9ea58a4SMatthew G. Knepley PetscScalar N; /* Tension force on right wall */ 46d12cd81dSMatthew G. Knepley } Parameter; 47d12cd81dSMatthew G. Knepley 48d12cd81dSMatthew G. Knepley typedef struct { 49c4762a1bSJed Brown /* Domain and mesh definition */ 50c4762a1bSJed Brown char dmType[256]; /* DM type for the solve */ 51482dd828SMatthew G. Knepley DeformType deform; /* Domain deformation type */ 52c4762a1bSJed Brown /* Problem definition */ 53c4762a1bSJed Brown SolutionType solType; /* Type of exact solution */ 54d12cd81dSMatthew G. Knepley PetscBag bag; /* Problem parameters */ 55c4762a1bSJed Brown /* Solver definition */ 56c4762a1bSJed Brown PetscBool useNearNullspace; /* Use the rigid body modes as a near nullspace for AMG */ 57c4762a1bSJed Brown } AppCtx; 58c4762a1bSJed Brown 59d71ae5a4SJacob Faibussowitsch static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 60d71ae5a4SJacob Faibussowitsch { 61c4762a1bSJed Brown PetscInt d; 62c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[d] = 0.0; 633ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 64c4762a1bSJed Brown } 65c4762a1bSJed Brown 66d71ae5a4SJacob Faibussowitsch static PetscErrorCode ge_shift(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 67d71ae5a4SJacob Faibussowitsch { 68d12cd81dSMatthew G. Knepley PetscInt d; 69d12cd81dSMatthew G. Knepley u[0] = 0.1; 70d12cd81dSMatthew G. Knepley for (d = 1; d < dim; ++d) u[d] = 0.0; 713ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 72d12cd81dSMatthew G. Knepley } 73d12cd81dSMatthew G. Knepley 74d71ae5a4SJacob Faibussowitsch static PetscErrorCode quadratic_2d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 75d71ae5a4SJacob Faibussowitsch { 76c4762a1bSJed Brown u[0] = x[0] * x[0]; 77c4762a1bSJed Brown u[1] = x[1] * x[1] - 2.0 * x[0] * x[1]; 783ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 79c4762a1bSJed Brown } 80c4762a1bSJed Brown 81d71ae5a4SJacob Faibussowitsch static PetscErrorCode quadratic_3d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 82d71ae5a4SJacob Faibussowitsch { 83c4762a1bSJed Brown u[0] = x[0] * x[0]; 84c4762a1bSJed Brown u[1] = x[1] * x[1] - 2.0 * x[0] * x[1]; 85c4762a1bSJed Brown u[2] = x[2] * x[2] - 2.0 * x[1] * x[2]; 863ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 87c4762a1bSJed Brown } 88c4762a1bSJed Brown 89c4762a1bSJed Brown /* 90c4762a1bSJed Brown u = x^2 91c4762a1bSJed Brown v = y^2 - 2xy 92c4762a1bSJed Brown Delta <u,v> - f = <2, 2> - <2, 2> 93c4762a1bSJed Brown 94c4762a1bSJed Brown u = x^2 95c4762a1bSJed Brown v = y^2 - 2xy 96c4762a1bSJed Brown w = z^2 - 2yz 97c4762a1bSJed Brown Delta <u,v,w> - f = <2, 2, 2> - <2, 2, 2> 98c4762a1bSJed Brown */ 99d71ae5a4SJacob Faibussowitsch static void f0_vlap_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 100d71ae5a4SJacob Faibussowitsch { 101c4762a1bSJed Brown PetscInt d; 102c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[d] += 2.0; 103c4762a1bSJed Brown } 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* 106c4762a1bSJed Brown u = x^2 107c4762a1bSJed Brown v = y^2 - 2xy 108c4762a1bSJed Brown \varepsilon = / 2x -y \ 109c4762a1bSJed Brown \ -y 2y - 2x / 110c4762a1bSJed Brown Tr(\varepsilon) = div u = 2y 111c4762a1bSJed Brown div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij} 112c4762a1bSJed Brown = \lambda \partial_j (2y) + 2\mu < 2-1, 2 > 113c4762a1bSJed Brown = \lambda < 0, 2 > + \mu < 2, 4 > 114c4762a1bSJed Brown 115c4762a1bSJed Brown u = x^2 116c4762a1bSJed Brown v = y^2 - 2xy 117c4762a1bSJed Brown w = z^2 - 2yz 118c4762a1bSJed Brown \varepsilon = / 2x -y 0 \ 119c4762a1bSJed Brown | -y 2y - 2x -z | 120c4762a1bSJed Brown \ 0 -z 2z - 2y/ 121c4762a1bSJed Brown Tr(\varepsilon) = div u = 2z 122c4762a1bSJed Brown div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij} 123c4762a1bSJed Brown = \lambda \partial_j (2z) + 2\mu < 2-1, 2-1, 2 > 124c4762a1bSJed Brown = \lambda < 0, 0, 2 > + \mu < 2, 2, 4 > 125c4762a1bSJed Brown */ 126d71ae5a4SJacob Faibussowitsch static void f0_elas_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 127d71ae5a4SJacob Faibussowitsch { 128f9ea58a4SMatthew G. Knepley const PetscReal mu = PetscRealPart(constants[0]); 129f9ea58a4SMatthew G. Knepley const PetscReal lambda = PetscRealPart(constants[1]); 130c4762a1bSJed Brown 131f9ea58a4SMatthew G. Knepley for (PetscInt d = 0; d < dim - 1; ++d) f0[d] += 2.0 * mu; 132c4762a1bSJed Brown f0[dim - 1] += 2.0 * lambda + 4.0 * mu; 133c4762a1bSJed Brown } 134c4762a1bSJed Brown 135d71ae5a4SJacob Faibussowitsch static void f0_mass_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 136d71ae5a4SJacob Faibussowitsch { 137482dd828SMatthew G. Knepley if (dim == 2) { 138482dd828SMatthew G. Knepley f0[0] -= x[0] * x[0]; 139482dd828SMatthew G. Knepley f0[1] -= x[1] * x[1] - 2.0 * x[0] * x[1]; 140482dd828SMatthew G. Knepley } else { 141482dd828SMatthew G. Knepley f0[0] -= x[0] * x[0]; 142482dd828SMatthew G. Knepley f0[1] -= x[1] * x[1] - 2.0 * x[0] * x[1]; 143482dd828SMatthew G. Knepley f0[2] -= x[2] * x[2] - 2.0 * x[1] * x[2]; 144482dd828SMatthew G. Knepley } 145482dd828SMatthew G. Knepley } 146482dd828SMatthew G. Knepley 147d71ae5a4SJacob Faibussowitsch static PetscErrorCode trig_2d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 148d71ae5a4SJacob Faibussowitsch { 149c4762a1bSJed Brown u[0] = PetscSinReal(2.0 * PETSC_PI * x[0]); 150c4762a1bSJed Brown u[1] = PetscSinReal(2.0 * PETSC_PI * x[1]) - 2.0 * x[0] * x[1]; 1513ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 152c4762a1bSJed Brown } 153c4762a1bSJed Brown 154d71ae5a4SJacob Faibussowitsch static PetscErrorCode trig_3d_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 155d71ae5a4SJacob Faibussowitsch { 156c4762a1bSJed Brown u[0] = PetscSinReal(2.0 * PETSC_PI * x[0]); 157c4762a1bSJed Brown u[1] = PetscSinReal(2.0 * PETSC_PI * x[1]) - 2.0 * x[0] * x[1]; 158c4762a1bSJed Brown u[2] = PetscSinReal(2.0 * PETSC_PI * x[2]) - 2.0 * x[1] * x[2]; 1593ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 160c4762a1bSJed Brown } 161c4762a1bSJed Brown 162c4762a1bSJed Brown /* 163c4762a1bSJed Brown u = sin(2 pi x) 164c4762a1bSJed Brown v = sin(2 pi y) - 2xy 165c4762a1bSJed 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)> 166c4762a1bSJed Brown 167c4762a1bSJed Brown u = sin(2 pi x) 168c4762a1bSJed Brown v = sin(2 pi y) - 2xy 169c4762a1bSJed Brown w = sin(2 pi z) - 2yz 170c4762a1bSJed 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)> 171c4762a1bSJed Brown */ 172d71ae5a4SJacob Faibussowitsch static void f0_vlap_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 173d71ae5a4SJacob Faibussowitsch { 174c4762a1bSJed Brown PetscInt d; 175c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[d] += -4.0 * PetscSqr(PETSC_PI) * PetscSinReal(2.0 * PETSC_PI * x[d]); 176c4762a1bSJed Brown } 177c4762a1bSJed Brown 178c4762a1bSJed Brown /* 179c4762a1bSJed Brown u = sin(2 pi x) 180c4762a1bSJed Brown v = sin(2 pi y) - 2xy 181c4762a1bSJed Brown \varepsilon = / 2 pi cos(2 pi x) -y \ 182c4762a1bSJed Brown \ -y 2 pi cos(2 pi y) - 2x / 183c4762a1bSJed Brown Tr(\varepsilon) = div u = 2 pi (cos(2 pi x) + cos(2 pi y)) - 2 x 184c4762a1bSJed Brown div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij} 185c4762a1bSJed 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) > 186c4762a1bSJed 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) > 187c4762a1bSJed Brown 188c4762a1bSJed Brown u = sin(2 pi x) 189c4762a1bSJed Brown v = sin(2 pi y) - 2xy 190c4762a1bSJed Brown w = sin(2 pi z) - 2yz 191c4762a1bSJed Brown \varepsilon = / 2 pi cos(2 pi x) -y 0 \ 192c4762a1bSJed Brown | -y 2 pi cos(2 pi y) - 2x -z | 193c4762a1bSJed Brown \ 0 -z 2 pi cos(2 pi z) - 2y / 194c4762a1bSJed Brown Tr(\varepsilon) = div u = 2 pi (cos(2 pi x) + cos(2 pi y) + cos(2 pi z)) - 2 x - 2 y 195c4762a1bSJed Brown div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij} 196c4762a1bSJed 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) > 197c4762a1bSJed 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) > 198c4762a1bSJed Brown */ 199d71ae5a4SJacob Faibussowitsch static void f0_elas_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 200d71ae5a4SJacob Faibussowitsch { 201f9ea58a4SMatthew G. Knepley const PetscReal mu = PetscRealPart(constants[0]); 202f9ea58a4SMatthew G. Knepley const PetscReal lambda = PetscRealPart(constants[1]); 203c4762a1bSJed Brown const PetscReal fact = 4.0 * PetscSqr(PETSC_PI); 204c4762a1bSJed Brown 205f9ea58a4SMatthew G. Knepley for (PetscInt 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); 206c4762a1bSJed Brown } 207c4762a1bSJed Brown 208d71ae5a4SJacob Faibussowitsch static PetscErrorCode axial_disp_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 209d71ae5a4SJacob Faibussowitsch { 210f9ea58a4SMatthew G. Knepley AppCtx *user = (AppCtx *)ctx; 211f9ea58a4SMatthew G. Knepley Parameter *param; 212f9ea58a4SMatthew G. Knepley 213f9ea58a4SMatthew G. Knepley PetscCall(PetscBagGetData(user->bag, (void **)¶m)); 214f9ea58a4SMatthew G. Knepley { 215f9ea58a4SMatthew G. Knepley const PetscReal mu = PetscRealPart(param->mu); 216f9ea58a4SMatthew G. Knepley const PetscReal lambda = PetscRealPart(param->lambda); 217f9ea58a4SMatthew G. Knepley const PetscReal N = PetscRealPart(param->N); 218c4762a1bSJed Brown 219c4762a1bSJed 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]; 220f9ea58a4SMatthew G. Knepley u[1] = 0.25 * lambda / mu / (lambda + mu) * N * x[1]; 221f9ea58a4SMatthew G. Knepley for (PetscInt d = 2; d < dim; ++d) u[d] = 0.0; 222f9ea58a4SMatthew G. Knepley } 2233ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 224c4762a1bSJed Brown } 225c4762a1bSJed Brown 226c4762a1bSJed Brown /* 227c4762a1bSJed Brown We will pull/push on the right side of a block of linearly elastic material. The uniform traction conditions on the 228c4762a1bSJed Brown right side of the box will result in a uniform strain along x and y. The Neumann BC is given by 229c4762a1bSJed Brown 230c4762a1bSJed Brown n_i \sigma_{ij} = t_i 231c4762a1bSJed Brown 232c4762a1bSJed Brown u = (1/(2\mu) - 1) x 233c4762a1bSJed Brown v = -y 234c4762a1bSJed Brown f = 0 235c4762a1bSJed Brown t = <4\mu/\lambda (\lambda + \mu), 0> 236c4762a1bSJed Brown \varepsilon = / 1/(2\mu) - 1 0 \ 237c4762a1bSJed Brown \ 0 -1 / 238c4762a1bSJed Brown Tr(\varepsilon) = div u = 1/(2\mu) - 2 239c4762a1bSJed Brown div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij} 240c4762a1bSJed Brown = \lambda \partial_j (1/(2\mu) - 2) + 2\mu < 0, 0 > 241c4762a1bSJed Brown = \lambda < 0, 0 > + \mu < 0, 0 > = 0 242c4762a1bSJed Brown NBC = <1,0> . <4\mu/\lambda (\lambda + \mu), 0> = 4\mu/\lambda (\lambda + \mu) 243c4762a1bSJed Brown 244c4762a1bSJed Brown u = x - 1/2 245c4762a1bSJed Brown v = 0 246c4762a1bSJed Brown w = 0 247c4762a1bSJed Brown \varepsilon = / x 0 0 \ 248c4762a1bSJed Brown | 0 0 0 | 249c4762a1bSJed Brown \ 0 0 0 / 250c4762a1bSJed Brown Tr(\varepsilon) = div u = x 251c4762a1bSJed Brown div \sigma = \partial_i \lambda \delta_{ij} \varepsilon_{kk} + \partial_i 2\mu\varepsilon_{ij} 252c4762a1bSJed Brown = \lambda \partial_j x + 2\mu < 1, 0, 0 > 253c4762a1bSJed Brown = \lambda < 1, 0, 0 > + \mu < 2, 0, 0 > 254c4762a1bSJed Brown */ 255d71ae5a4SJacob Faibussowitsch static void f0_elas_axial_disp_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 256d71ae5a4SJacob Faibussowitsch { 257f9ea58a4SMatthew G. Knepley const PetscReal N = PetscRealPart(constants[2]); 258c4762a1bSJed Brown 259c4762a1bSJed Brown f0[0] = N; 260c4762a1bSJed Brown } 261c4762a1bSJed Brown 262d71ae5a4SJacob Faibussowitsch static PetscErrorCode uniform_strain_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 263d71ae5a4SJacob Faibussowitsch { 264c4762a1bSJed Brown const PetscReal eps_xx = 0.1; 265c4762a1bSJed Brown const PetscReal eps_xy = 0.3; 266c4762a1bSJed Brown const PetscReal eps_yy = 0.25; 267c4762a1bSJed Brown PetscInt d; 268c4762a1bSJed Brown 269c4762a1bSJed Brown u[0] = eps_xx * x[0] + eps_xy * x[1]; 270c4762a1bSJed Brown u[1] = eps_xy * x[0] + eps_yy * x[1]; 271c4762a1bSJed Brown for (d = 2; d < dim; ++d) u[d] = 0.0; 2723ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 273c4762a1bSJed Brown } 274c4762a1bSJed Brown 275d71ae5a4SJacob Faibussowitsch static void f0_mass_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 276d71ae5a4SJacob Faibussowitsch { 277482dd828SMatthew G. Knepley const PetscInt Nc = dim; 278482dd828SMatthew G. Knepley PetscInt c; 279482dd828SMatthew G. Knepley 280482dd828SMatthew G. Knepley for (c = 0; c < Nc; ++c) f0[c] = u[c]; 281482dd828SMatthew G. Knepley } 282482dd828SMatthew G. Knepley 283d71ae5a4SJacob Faibussowitsch static void f1_vlap_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 284d71ae5a4SJacob Faibussowitsch { 285c4762a1bSJed Brown const PetscInt Nc = dim; 286c4762a1bSJed Brown PetscInt c, d; 287c4762a1bSJed Brown 2889371c9d4SSatish Balay for (c = 0; c < Nc; ++c) 2899371c9d4SSatish Balay for (d = 0; d < dim; ++d) f1[c * dim + d] += u_x[c * dim + d]; 290c4762a1bSJed Brown } 291c4762a1bSJed Brown 292d71ae5a4SJacob Faibussowitsch static void f1_elas_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 293d71ae5a4SJacob Faibussowitsch { 294f9ea58a4SMatthew G. Knepley const PetscReal mu = PetscRealPart(constants[0]); 295f9ea58a4SMatthew G. Knepley const PetscReal lambda = PetscRealPart(constants[1]); 296c4762a1bSJed Brown const PetscInt Nc = dim; 297c4762a1bSJed Brown 298f9ea58a4SMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) { 299f9ea58a4SMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) { 300c4762a1bSJed Brown f1[c * dim + d] += mu * (u_x[c * dim + d] + u_x[d * dim + c]); 301c4762a1bSJed Brown f1[c * dim + c] += lambda * u_x[d * dim + d]; 302c4762a1bSJed Brown } 303c4762a1bSJed Brown } 304c4762a1bSJed Brown } 305c4762a1bSJed Brown 306d71ae5a4SJacob Faibussowitsch static void g0_mass_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 307d71ae5a4SJacob Faibussowitsch { 308482dd828SMatthew G. Knepley const PetscInt Nc = dim; 309482dd828SMatthew G. Knepley PetscInt c; 310482dd828SMatthew G. Knepley 311482dd828SMatthew G. Knepley for (c = 0; c < Nc; ++c) g0[c * Nc + c] = 1.0; 312482dd828SMatthew G. Knepley } 313482dd828SMatthew G. Knepley 314d71ae5a4SJacob Faibussowitsch static void g3_vlap_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 315d71ae5a4SJacob Faibussowitsch { 316c4762a1bSJed Brown const PetscInt Nc = dim; 317c4762a1bSJed Brown PetscInt c, d; 318c4762a1bSJed Brown 319c4762a1bSJed Brown for (c = 0; c < Nc; ++c) { 320ad540459SPierre Jolivet for (d = 0; d < dim; ++d) g3[((c * Nc + c) * dim + d) * dim + d] = 1.0; 321c4762a1bSJed Brown } 322c4762a1bSJed Brown } 323c4762a1bSJed Brown 324c4762a1bSJed Brown /* 325c4762a1bSJed Brown \partial_df \phi_fc g_{fc,gc,df,dg} \partial_dg \phi_gc 326c4762a1bSJed Brown 327c4762a1bSJed Brown \partial_df \phi_fc \lambda \delta_{fc,df} \sum_gc \partial_dg \phi_gc \delta_{gc,dg} 328c4762a1bSJed Brown = \partial_fc \phi_fc \sum_gc \partial_gc \phi_gc 329c4762a1bSJed Brown */ 330d71ae5a4SJacob Faibussowitsch static void g3_elas_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 331d71ae5a4SJacob Faibussowitsch { 332f9ea58a4SMatthew G. Knepley const PetscReal mu = PetscRealPart(constants[0]); 333f9ea58a4SMatthew G. Knepley const PetscReal lambda = PetscRealPart(constants[1]); 334c4762a1bSJed Brown const PetscInt Nc = dim; 335c4762a1bSJed Brown 336f9ea58a4SMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) { 337f9ea58a4SMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) { 338c4762a1bSJed Brown g3[((c * Nc + c) * dim + d) * dim + d] += mu; 339c4762a1bSJed Brown g3[((c * Nc + d) * dim + d) * dim + c] += mu; 340c4762a1bSJed Brown g3[((c * Nc + d) * dim + c) * dim + d] += lambda; 341c4762a1bSJed Brown } 342c4762a1bSJed Brown } 343c4762a1bSJed Brown } 344c4762a1bSJed Brown 345d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 346d71ae5a4SJacob Faibussowitsch { 347482dd828SMatthew G. Knepley PetscInt sol = 0, def = 0; 348c4762a1bSJed Brown 349c4762a1bSJed Brown PetscFunctionBeginUser; 350482dd828SMatthew G. Knepley options->deform = DEFORM_NONE; 351c4762a1bSJed Brown options->solType = SOL_VLAP_QUADRATIC; 352c4762a1bSJed Brown options->useNearNullspace = PETSC_TRUE; 3539566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(options->dmType, DMPLEX, 256)); 354c4762a1bSJed Brown 355d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Linear Elasticity Problem Options", "DMPLEX"); 3569566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-deform_type", "Type of domain deformation", "ex17.c", deformTypes, NUM_DEFORM_TYPES, deformTypes[options->deform], &def, NULL)); 357482dd828SMatthew G. Knepley options->deform = (DeformType)def; 3589566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-sol_type", "Type of exact solution", "ex17.c", solutionTypes, NUM_SOLUTION_TYPES, solutionTypes[options->solType], &sol, NULL)); 359c4762a1bSJed Brown options->solType = (SolutionType)sol; 3609566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-near_nullspace", "Use the rigid body modes as an AMG near nullspace", "ex17.c", options->useNearNullspace, &options->useNearNullspace, NULL)); 3619566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-dm_type", "Convert DMPlex to another format", "ex17.c", DMList, options->dmType, options->dmType, 256, NULL)); 362d0609cedSBarry Smith PetscOptionsEnd(); 3633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 364c4762a1bSJed Brown } 365c4762a1bSJed Brown 366d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *ctx) 367d71ae5a4SJacob Faibussowitsch { 368d12cd81dSMatthew G. Knepley PetscBag bag; 369d12cd81dSMatthew G. Knepley Parameter *p; 370d12cd81dSMatthew G. Knepley 371d12cd81dSMatthew G. Knepley PetscFunctionBeginUser; 372d12cd81dSMatthew G. Knepley /* setup PETSc parameter bag */ 3739566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(ctx->bag, (void **)&p)); 3749566063dSJacob Faibussowitsch PetscCall(PetscBagSetName(ctx->bag, "par", "Elastic Parameters")); 375d12cd81dSMatthew G. Knepley bag = ctx->bag; 3769566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterScalar(bag, &p->mu, 1.0, "mu", "Shear Modulus, Pa")); 3779566063dSJacob Faibussowitsch PetscCall(PetscBagRegisterScalar(bag, &p->lambda, 1.0, "lambda", "Lame's first parameter, Pa")); 378f9ea58a4SMatthew G. Knepley PetscCall(PetscBagRegisterScalar(bag, &p->N, -1.0, "N", "Tension on right wall, Pa")); 3799566063dSJacob Faibussowitsch PetscCall(PetscBagSetFromOptions(bag)); 380d12cd81dSMatthew G. Knepley { 381d12cd81dSMatthew G. Knepley PetscViewer viewer; 382d12cd81dSMatthew G. Knepley PetscViewerFormat format; 383d12cd81dSMatthew G. Knepley PetscBool flg; 384d12cd81dSMatthew G. Knepley 385648c30bcSBarry Smith PetscCall(PetscOptionsCreateViewer(comm, NULL, NULL, "-param_view", &viewer, &format, &flg)); 386d12cd81dSMatthew G. Knepley if (flg) { 3879566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, format)); 3889566063dSJacob Faibussowitsch PetscCall(PetscBagView(bag, viewer)); 3899566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 3909566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 391648c30bcSBarry Smith PetscCall(PetscViewerDestroy(&viewer)); 392d12cd81dSMatthew G. Knepley } 393d12cd81dSMatthew G. Knepley } 3943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 395d12cd81dSMatthew G. Knepley } 396d12cd81dSMatthew G. Knepley 397d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexDistortGeometry(DM dm) 398d71ae5a4SJacob Faibussowitsch { 399482dd828SMatthew G. Knepley DM cdm; 400482dd828SMatthew G. Knepley DMLabel label; 401482dd828SMatthew G. Knepley Vec coordinates; 402482dd828SMatthew G. Knepley PetscScalar *coords; 403482dd828SMatthew G. Knepley PetscReal mid = 0.5; 404482dd828SMatthew G. Knepley PetscInt cdim, d, vStart, vEnd, v; 405482dd828SMatthew G. Knepley 406482dd828SMatthew G. Knepley PetscFunctionBeginUser; 4079566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 4089566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim)); 4099566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 4109566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 4119566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 4129566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(coordinates, &coords)); 413482dd828SMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 414482dd828SMatthew G. Knepley PetscScalar *pcoords, shift; 415482dd828SMatthew G. Knepley PetscInt val; 416482dd828SMatthew G. Knepley 4179566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(label, v, &val)); 418482dd828SMatthew G. Knepley if (val >= 0) continue; 4199566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRef(cdm, v, coords, &pcoords)); 420482dd828SMatthew G. Knepley shift = 0.2 * PetscAbsScalar(pcoords[0] - mid); 421482dd828SMatthew G. Knepley shift = PetscRealPart(pcoords[0]) > mid ? shift : -shift; 422482dd828SMatthew G. Knepley for (d = 1; d < cdim; ++d) pcoords[d] += shift; 423482dd828SMatthew G. Knepley } 4249566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(coordinates, &coords)); 4253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 426482dd828SMatthew G. Knepley } 427482dd828SMatthew G. Knepley 428d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 429d71ae5a4SJacob Faibussowitsch { 430c4762a1bSJed Brown PetscFunctionBeginUser; 4319566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 4329566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 4339566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 434482dd828SMatthew G. Knepley switch (user->deform) { 435d71ae5a4SJacob Faibussowitsch case DEFORM_NONE: 436d71ae5a4SJacob Faibussowitsch break; 437d71ae5a4SJacob Faibussowitsch case DEFORM_SHEAR: 438d71ae5a4SJacob Faibussowitsch PetscCall(DMPlexShearGeometry(*dm, DM_X, NULL)); 439d71ae5a4SJacob Faibussowitsch break; 440d71ae5a4SJacob Faibussowitsch case DEFORM_STEP: 441d71ae5a4SJacob Faibussowitsch PetscCall(DMPlexDistortGeometry(*dm)); 442d71ae5a4SJacob Faibussowitsch break; 443d71ae5a4SJacob Faibussowitsch default: 444d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid deformation type: %s (%d)", deformTypes[PetscMin(user->deform, NUM_DEFORM_TYPES)], user->deform); 445482dd828SMatthew G. Knepley } 4469566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*dm, user)); 4479566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 4483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 449c4762a1bSJed Brown } 450c4762a1bSJed Brown 451d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 452d71ae5a4SJacob Faibussowitsch { 453c4762a1bSJed Brown PetscErrorCode (*exact)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *); 454d12cd81dSMatthew G. Knepley Parameter *param; 455482dd828SMatthew G. Knepley PetscDS ds; 45645480ffeSMatthew G. Knepley PetscWeakForm wf; 45745480ffeSMatthew G. Knepley DMLabel label; 45845480ffeSMatthew G. Knepley PetscInt id, bd; 459c4762a1bSJed Brown PetscInt dim; 460c4762a1bSJed Brown 461c4762a1bSJed Brown PetscFunctionBeginUser; 4629566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 4639566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf)); 4649566063dSJacob Faibussowitsch PetscCall(PetscDSGetSpatialDimension(ds, &dim)); 4659566063dSJacob Faibussowitsch PetscCall(PetscBagGetData(user->bag, (void **)¶m)); 466c4762a1bSJed Brown switch (user->solType) { 467482dd828SMatthew G. Knepley case SOL_MASS_QUADRATIC: 4689566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_mass_u, NULL)); 4699566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_mass_uu, NULL, NULL, NULL)); 4709566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, 0, 0, 1, f0_mass_quadratic_u, 0, NULL)); 471c4762a1bSJed Brown switch (dim) { 472d71ae5a4SJacob Faibussowitsch case 2: 473d71ae5a4SJacob Faibussowitsch exact = quadratic_2d_u; 474d71ae5a4SJacob Faibussowitsch break; 475d71ae5a4SJacob Faibussowitsch case 3: 476d71ae5a4SJacob Faibussowitsch exact = quadratic_3d_u; 477d71ae5a4SJacob Faibussowitsch break; 478d71ae5a4SJacob Faibussowitsch default: 479d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim); 480482dd828SMatthew G. Knepley } 481482dd828SMatthew G. Knepley break; 482482dd828SMatthew G. Knepley case SOL_VLAP_QUADRATIC: 4839566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_vlap_quadratic_u, f1_vlap_u)); 4849566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_vlap_uu)); 485482dd828SMatthew G. Knepley switch (dim) { 486d71ae5a4SJacob Faibussowitsch case 2: 487d71ae5a4SJacob Faibussowitsch exact = quadratic_2d_u; 488d71ae5a4SJacob Faibussowitsch break; 489d71ae5a4SJacob Faibussowitsch case 3: 490d71ae5a4SJacob Faibussowitsch exact = quadratic_3d_u; 491d71ae5a4SJacob Faibussowitsch break; 492d71ae5a4SJacob Faibussowitsch default: 493d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim); 494c4762a1bSJed Brown } 495c4762a1bSJed Brown break; 496c4762a1bSJed Brown case SOL_ELAS_QUADRATIC: 4979566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_elas_quadratic_u, f1_elas_u)); 4989566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu)); 499c4762a1bSJed Brown switch (dim) { 500d71ae5a4SJacob Faibussowitsch case 2: 501d71ae5a4SJacob Faibussowitsch exact = quadratic_2d_u; 502d71ae5a4SJacob Faibussowitsch break; 503d71ae5a4SJacob Faibussowitsch case 3: 504d71ae5a4SJacob Faibussowitsch exact = quadratic_3d_u; 505d71ae5a4SJacob Faibussowitsch break; 506d71ae5a4SJacob Faibussowitsch default: 507d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim); 508c4762a1bSJed Brown } 509c4762a1bSJed Brown break; 510c4762a1bSJed Brown case SOL_VLAP_TRIG: 5119566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_vlap_trig_u, f1_vlap_u)); 5129566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_vlap_uu)); 513c4762a1bSJed Brown switch (dim) { 514d71ae5a4SJacob Faibussowitsch case 2: 515d71ae5a4SJacob Faibussowitsch exact = trig_2d_u; 516d71ae5a4SJacob Faibussowitsch break; 517d71ae5a4SJacob Faibussowitsch case 3: 518d71ae5a4SJacob Faibussowitsch exact = trig_3d_u; 519d71ae5a4SJacob Faibussowitsch break; 520d71ae5a4SJacob Faibussowitsch default: 521d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim); 522c4762a1bSJed Brown } 523c4762a1bSJed Brown break; 524c4762a1bSJed Brown case SOL_ELAS_TRIG: 5259566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_elas_trig_u, f1_elas_u)); 5269566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu)); 527c4762a1bSJed Brown switch (dim) { 528d71ae5a4SJacob Faibussowitsch case 2: 529d71ae5a4SJacob Faibussowitsch exact = trig_2d_u; 530d71ae5a4SJacob Faibussowitsch break; 531d71ae5a4SJacob Faibussowitsch case 3: 532d71ae5a4SJacob Faibussowitsch exact = trig_3d_u; 533d71ae5a4SJacob Faibussowitsch break; 534d71ae5a4SJacob Faibussowitsch default: 535d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid dimension: %" PetscInt_FMT, dim); 536c4762a1bSJed Brown } 537c4762a1bSJed Brown break; 538c4762a1bSJed Brown case SOL_ELAS_AXIAL_DISP: 5399566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, NULL, f1_elas_u)); 5409566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu)); 54145480ffeSMatthew G. Knepley id = dim == 3 ? 5 : 2; 5429566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 5439566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "right", label, 1, &id, 0, 0, NULL, (void (*)(void))NULL, NULL, user, &bd)); 5449566063dSJacob Faibussowitsch PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 5459566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_elas_axial_disp_bd_u, 0, NULL)); 546c4762a1bSJed Brown exact = axial_disp_u; 547c4762a1bSJed Brown break; 548c4762a1bSJed Brown case SOL_ELAS_UNIFORM_STRAIN: 5499566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, NULL, f1_elas_u)); 5509566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu)); 551c4762a1bSJed Brown exact = uniform_strain_u; 552c4762a1bSJed Brown break; 553d12cd81dSMatthew G. Knepley case SOL_ELAS_GE: 5549566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, NULL, f1_elas_u)); 5559566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_elas_uu)); 556d12cd81dSMatthew G. Knepley exact = zero; /* No exact solution available */ 557d12cd81dSMatthew G. Knepley break; 558d71ae5a4SJacob Faibussowitsch default: 559d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Invalid solution type: %s (%d)", solutionTypes[PetscMin(user->solType, NUM_SOLUTION_TYPES)], user->solType); 560c4762a1bSJed Brown } 5619566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, exact, user)); 5629566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 563c4762a1bSJed Brown if (user->solType == SOL_ELAS_AXIAL_DISP) { 564c4762a1bSJed Brown PetscInt cmp; 565c4762a1bSJed Brown 566c4762a1bSJed Brown id = dim == 3 ? 6 : 4; 567c4762a1bSJed Brown cmp = 0; 5689566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "left", label, 1, &id, 0, 1, &cmp, (void (*)(void))zero, NULL, user, NULL)); 569c4762a1bSJed Brown cmp = dim == 3 ? 2 : 1; 570c4762a1bSJed Brown id = dim == 3 ? 1 : 1; 5719566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "bottom", label, 1, &id, 0, 1, &cmp, (void (*)(void))zero, NULL, user, NULL)); 572c4762a1bSJed Brown if (dim == 3) { 573c4762a1bSJed Brown cmp = 1; 574c4762a1bSJed Brown id = 3; 5759566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "front", label, 1, &id, 0, 1, &cmp, (void (*)(void))zero, NULL, user, NULL)); 576c4762a1bSJed Brown } 577d12cd81dSMatthew G. Knepley } else if (user->solType == SOL_ELAS_GE) { 578d12cd81dSMatthew G. Knepley PetscInt cmp; 579d12cd81dSMatthew G. Knepley 580d12cd81dSMatthew G. Knepley id = dim == 3 ? 6 : 4; 5819566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "left", label, 1, &id, 0, 0, NULL, (void (*)(void))zero, NULL, user, NULL)); 582d12cd81dSMatthew G. Knepley id = dim == 3 ? 5 : 2; 583d12cd81dSMatthew G. Knepley cmp = 0; 5849566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "right", label, 1, &id, 0, 1, &cmp, (void (*)(void))ge_shift, NULL, user, NULL)); 585c4762a1bSJed Brown } else { 586c4762a1bSJed Brown id = 1; 5879566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))exact, NULL, user, NULL)); 588c4762a1bSJed Brown } 589d12cd81dSMatthew G. Knepley /* Setup constants */ 590d12cd81dSMatthew G. Knepley { 591f9ea58a4SMatthew G. Knepley PetscScalar constants[3]; 592d12cd81dSMatthew G. Knepley 593d12cd81dSMatthew G. Knepley constants[0] = param->mu; /* shear modulus, Pa */ 594d12cd81dSMatthew G. Knepley constants[1] = param->lambda; /* Lame's first parameter, Pa */ 595f9ea58a4SMatthew G. Knepley constants[2] = param->N; /* Tension on right wall, Pa */ 596f9ea58a4SMatthew G. Knepley PetscCall(PetscDSSetConstants(ds, 3, constants)); 597d12cd81dSMatthew G. Knepley } 5983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 599c4762a1bSJed Brown } 600c4762a1bSJed Brown 601d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateElasticityNullSpace(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullspace) 602d71ae5a4SJacob Faibussowitsch { 603c4762a1bSJed Brown PetscFunctionBegin; 6049566063dSJacob Faibussowitsch PetscCall(DMPlexCreateRigidBody(dm, origField, nullspace)); 6053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 606c4762a1bSJed Brown } 607c4762a1bSJed Brown 608d71ae5a4SJacob Faibussowitsch PetscErrorCode SetupFE(DM dm, const char name[], PetscErrorCode (*setup)(DM, AppCtx *), void *ctx) 609d71ae5a4SJacob Faibussowitsch { 610c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 611c4762a1bSJed Brown DM cdm = dm; 612c4762a1bSJed Brown PetscFE fe; 613c4762a1bSJed Brown char prefix[PETSC_MAX_PATH_LEN]; 61430602db0SMatthew G. Knepley DMPolytopeType ct; 61530602db0SMatthew G. Knepley PetscInt dim, cStart; 616c4762a1bSJed Brown 617c4762a1bSJed Brown PetscFunctionBegin; 618c4762a1bSJed Brown /* Create finite element */ 6199566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 6209566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 6219566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 6229566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name)); 623f9ea58a4SMatthew G. Knepley PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, dim, ct, name ? prefix : NULL, -1, &fe)); 6249566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, name)); 625c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 6269566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 6279566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 6289566063dSJacob Faibussowitsch PetscCall((*setup)(dm, user)); 629c4762a1bSJed Brown while (cdm) { 6309566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 6319566063dSJacob Faibussowitsch if (user->useNearNullspace) PetscCall(DMSetNearNullSpaceConstructor(cdm, 0, CreateElasticityNullSpace)); 6329566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 633c4762a1bSJed Brown } 6349566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 6353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 636c4762a1bSJed Brown } 637c4762a1bSJed Brown 638d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 639d71ae5a4SJacob Faibussowitsch { 640c4762a1bSJed Brown DM dm; /* Problem specification */ 641c4762a1bSJed Brown SNES snes; /* Nonlinear solver */ 642c4762a1bSJed Brown Vec u; /* Solutions */ 643c4762a1bSJed Brown AppCtx user; /* User-defined work context */ 644c4762a1bSJed Brown 645327415f7SBarry Smith PetscFunctionBeginUser; 6469566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 6479566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 6489566063dSJacob Faibussowitsch PetscCall(PetscBagCreate(PETSC_COMM_SELF, sizeof(Parameter), &user.bag)); 6499566063dSJacob Faibussowitsch PetscCall(SetupParameters(PETSC_COMM_WORLD, &user)); 650c4762a1bSJed Brown /* Primal system */ 6519566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 6529566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 6539566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 6549566063dSJacob Faibussowitsch PetscCall(SetupFE(dm, "displacement", SetupPrimalProblem, &user)); 6559566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 6569566063dSJacob Faibussowitsch PetscCall(VecSet(u, 0.0)); 6579566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "displacement")); 6586493148fSStefano Zampini PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user)); 6599566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 6609566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u)); 6619566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 6629566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes, &u)); 6639566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-displacement_view")); 664c4762a1bSJed Brown /* Cleanup */ 6659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 6669566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 6679566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 6689566063dSJacob Faibussowitsch PetscCall(PetscBagDestroy(&user.bag)); 6699566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 670b122ec5aSJacob Faibussowitsch return 0; 671c4762a1bSJed Brown } 672c4762a1bSJed Brown 673c4762a1bSJed Brown /*TEST 674c4762a1bSJed Brown 67530602db0SMatthew G. Knepley testset: 67630602db0SMatthew G. Knepley args: -dm_plex_box_faces 1,1,1 67730602db0SMatthew G. Knepley 678c4762a1bSJed Brown test: 679c4762a1bSJed Brown suffix: 2d_p1_quad_vlap 680c4762a1bSJed Brown requires: triangle 681c4762a1bSJed Brown args: -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate 682c4762a1bSJed Brown test: 683c4762a1bSJed Brown suffix: 2d_p2_quad_vlap 684c4762a1bSJed Brown requires: triangle 685c4762a1bSJed Brown args: -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001 686c4762a1bSJed Brown test: 687c4762a1bSJed Brown suffix: 2d_p3_quad_vlap 688c4762a1bSJed Brown requires: triangle 689c4762a1bSJed Brown args: -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001 690c4762a1bSJed Brown test: 691c4762a1bSJed Brown suffix: 2d_q1_quad_vlap 69230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate 693c4762a1bSJed Brown test: 694c4762a1bSJed Brown suffix: 2d_q2_quad_vlap 69530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001 696c4762a1bSJed Brown test: 697c4762a1bSJed Brown suffix: 2d_q3_quad_vlap 698c4762a1bSJed Brown requires: !single 69930602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001 700c4762a1bSJed Brown test: 701c4762a1bSJed Brown suffix: 2d_p1_quad_elas 702c4762a1bSJed Brown requires: triangle 703c4762a1bSJed Brown args: -sol_type elas_quad -displacement_petscspace_degree 1 -dm_refine 2 -convest_num_refine 3 -snes_convergence_estimate 704c4762a1bSJed Brown test: 705c4762a1bSJed Brown suffix: 2d_p2_quad_elas 706c4762a1bSJed Brown requires: triangle 707c4762a1bSJed Brown args: -sol_type elas_quad -displacement_petscspace_degree 2 -dmsnes_check .0001 708c4762a1bSJed Brown test: 709c4762a1bSJed Brown suffix: 2d_p3_quad_elas 710c4762a1bSJed Brown requires: triangle 711c4762a1bSJed Brown args: -sol_type elas_quad -displacement_petscspace_degree 3 -dmsnes_check .0001 712c4762a1bSJed Brown test: 713c4762a1bSJed Brown suffix: 2d_q1_quad_elas 71430602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_simplex 0 -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 715c4762a1bSJed Brown test: 716c4762a1bSJed Brown suffix: 2d_q1_quad_elas_shear 717482dd828SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_simplex 0 -deform_type shear -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 718c4762a1bSJed Brown test: 719c4762a1bSJed Brown suffix: 2d_q2_quad_elas 72030602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dmsnes_check .0001 721c4762a1bSJed Brown test: 722c4762a1bSJed Brown suffix: 2d_q2_quad_elas_shear 723482dd828SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_simplex 0 -deform_type shear -displacement_petscspace_degree 2 -dmsnes_check 724c4762a1bSJed Brown test: 725c4762a1bSJed Brown suffix: 2d_q3_quad_elas 72630602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dmsnes_check .0001 727c4762a1bSJed Brown test: 728c4762a1bSJed Brown suffix: 2d_q3_quad_elas_shear 729c4762a1bSJed Brown requires: !single 730482dd828SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_simplex 0 -deform_type shear -displacement_petscspace_degree 3 -dmsnes_check 731c4762a1bSJed Brown 732c4762a1bSJed Brown test: 733c4762a1bSJed Brown suffix: 3d_p1_quad_vlap 734c4762a1bSJed Brown requires: ctetgen 73530602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_refine 1 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate 736c4762a1bSJed Brown test: 737c4762a1bSJed Brown suffix: 3d_p2_quad_vlap 738c4762a1bSJed Brown requires: ctetgen 73930602db0SMatthew G. Knepley args: -dm_plex_dim 3 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001 740c4762a1bSJed Brown test: 741c4762a1bSJed Brown suffix: 3d_p3_quad_vlap 742c4762a1bSJed Brown requires: ctetgen 74330602db0SMatthew G. Knepley args: -dm_plex_dim 3 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001 744c4762a1bSJed Brown test: 745c4762a1bSJed Brown suffix: 3d_q1_quad_vlap 74630602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_plex_simplex 0 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate 747c4762a1bSJed Brown test: 748c4762a1bSJed Brown suffix: 3d_q2_quad_vlap 74930602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001 750c4762a1bSJed Brown test: 751c4762a1bSJed Brown suffix: 3d_q3_quad_vlap 75230602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001 753c4762a1bSJed Brown test: 754c4762a1bSJed Brown suffix: 3d_p1_quad_elas 755c4762a1bSJed Brown requires: ctetgen 75630602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_dim 3 -dm_refine 1 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate 757c4762a1bSJed Brown test: 758c4762a1bSJed Brown suffix: 3d_p2_quad_elas 759c4762a1bSJed Brown requires: ctetgen 76030602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_dim 3 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001 761c4762a1bSJed Brown test: 762c4762a1bSJed Brown suffix: 3d_p3_quad_elas 763c4762a1bSJed Brown requires: ctetgen 76430602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_dim 3 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001 765c4762a1bSJed Brown test: 766c4762a1bSJed Brown suffix: 3d_q1_quad_elas 76730602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_plex_simplex 0 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate 768c4762a1bSJed Brown test: 769c4762a1bSJed Brown suffix: 3d_q2_quad_elas 77030602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 1 -dmsnes_check .0001 771c4762a1bSJed Brown test: 772c4762a1bSJed Brown suffix: 3d_q3_quad_elas 773c4762a1bSJed Brown requires: !single 77430602db0SMatthew G. Knepley args: -sol_type elas_quad -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 0 -dmsnes_check .0001 775c4762a1bSJed Brown 776c4762a1bSJed Brown test: 777c4762a1bSJed Brown suffix: 2d_p1_trig_vlap 778c4762a1bSJed Brown requires: triangle 779c4762a1bSJed Brown args: -sol_type vlap_trig -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 780c4762a1bSJed Brown test: 781c4762a1bSJed Brown suffix: 2d_p2_trig_vlap 782c4762a1bSJed Brown requires: triangle 783c4762a1bSJed Brown args: -sol_type vlap_trig -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 784c4762a1bSJed Brown test: 785c4762a1bSJed Brown suffix: 2d_p3_trig_vlap 786c4762a1bSJed Brown requires: triangle 787c4762a1bSJed Brown args: -sol_type vlap_trig -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 788c4762a1bSJed Brown test: 789c4762a1bSJed Brown suffix: 2d_q1_trig_vlap 79030602db0SMatthew G. Knepley args: -sol_type vlap_trig -dm_plex_simplex 0 -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 791c4762a1bSJed Brown test: 792c4762a1bSJed Brown suffix: 2d_q2_trig_vlap 79330602db0SMatthew G. Knepley args: -sol_type vlap_trig -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 794c4762a1bSJed Brown test: 795c4762a1bSJed Brown suffix: 2d_q3_trig_vlap 79630602db0SMatthew G. Knepley args: -sol_type vlap_trig -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 797c4762a1bSJed Brown test: 798c4762a1bSJed Brown suffix: 2d_p1_trig_elas 799c4762a1bSJed Brown requires: triangle 800c4762a1bSJed Brown args: -sol_type elas_trig -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 801c4762a1bSJed Brown test: 802c4762a1bSJed Brown suffix: 2d_p2_trig_elas 803c4762a1bSJed Brown requires: triangle 804c4762a1bSJed Brown args: -sol_type elas_trig -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 805c4762a1bSJed Brown test: 806c4762a1bSJed Brown suffix: 2d_p3_trig_elas 807c4762a1bSJed Brown requires: triangle 808c4762a1bSJed Brown args: -sol_type elas_trig -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 809c4762a1bSJed Brown test: 810c4762a1bSJed Brown suffix: 2d_q1_trig_elas 81130602db0SMatthew G. Knepley args: -sol_type elas_trig -dm_plex_simplex 0 -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 812c4762a1bSJed Brown test: 813c4762a1bSJed Brown suffix: 2d_q1_trig_elas_shear 814482dd828SMatthew G. Knepley args: -sol_type elas_trig -dm_plex_simplex 0 -deform_type shear -displacement_petscspace_degree 1 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 815c4762a1bSJed Brown test: 816c4762a1bSJed Brown suffix: 2d_q2_trig_elas 81730602db0SMatthew G. Knepley args: -sol_type elas_trig -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 818c4762a1bSJed Brown test: 819c4762a1bSJed Brown suffix: 2d_q2_trig_elas_shear 820482dd828SMatthew G. Knepley args: -sol_type elas_trig -dm_plex_simplex 0 -deform_type shear -displacement_petscspace_degree 2 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 821c4762a1bSJed Brown test: 822c4762a1bSJed Brown suffix: 2d_q3_trig_elas 82330602db0SMatthew G. Knepley args: -sol_type elas_trig -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 824c4762a1bSJed Brown test: 825c4762a1bSJed Brown suffix: 2d_q3_trig_elas_shear 826482dd828SMatthew G. Knepley args: -sol_type elas_trig -dm_plex_simplex 0 -deform_type shear -displacement_petscspace_degree 3 -dm_refine 1 -convest_num_refine 3 -snes_convergence_estimate 827c4762a1bSJed Brown 828c4762a1bSJed Brown test: 829c4762a1bSJed Brown suffix: 3d_p1_trig_vlap 830c4762a1bSJed Brown requires: ctetgen 83130602db0SMatthew G. Knepley args: -sol_type vlap_trig -dm_plex_dim 3 -dm_refine 1 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate 832c4762a1bSJed Brown test: 833c4762a1bSJed Brown suffix: 3d_p2_trig_vlap 834c4762a1bSJed Brown requires: ctetgen 83530602db0SMatthew G. Knepley args: -sol_type vlap_trig -dm_plex_dim 3 -displacement_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate 836c4762a1bSJed Brown test: 837c4762a1bSJed Brown suffix: 3d_p3_trig_vlap 838c4762a1bSJed Brown requires: ctetgen 83930602db0SMatthew G. Knepley args: -sol_type vlap_trig -dm_plex_dim 3 -displacement_petscspace_degree 3 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate 840c4762a1bSJed Brown test: 841c4762a1bSJed Brown suffix: 3d_q1_trig_vlap 84230602db0SMatthew G. Knepley args: -sol_type vlap_trig -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_plex_simplex 0 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate 843c4762a1bSJed Brown test: 844c4762a1bSJed Brown suffix: 3d_q2_trig_vlap 84530602db0SMatthew G. Knepley args: -sol_type vlap_trig -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate 846c4762a1bSJed Brown test: 847c4762a1bSJed Brown suffix: 3d_q3_trig_vlap 848c4762a1bSJed Brown requires: !__float128 84930602db0SMatthew G. Knepley args: -sol_type vlap_trig -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate 850c4762a1bSJed Brown test: 851c4762a1bSJed Brown suffix: 3d_p1_trig_elas 852c4762a1bSJed Brown requires: ctetgen 85330602db0SMatthew G. Knepley args: -sol_type elas_trig -dm_plex_dim 3 -dm_refine 1 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate 854c4762a1bSJed Brown test: 855c4762a1bSJed Brown suffix: 3d_p2_trig_elas 856c4762a1bSJed Brown requires: ctetgen 85730602db0SMatthew G. Knepley args: -sol_type elas_trig -dm_plex_dim 3 -displacement_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate 858c4762a1bSJed Brown test: 859c4762a1bSJed Brown suffix: 3d_p3_trig_elas 860c4762a1bSJed Brown requires: ctetgen 86130602db0SMatthew G. Knepley args: -sol_type elas_trig -dm_plex_dim 3 -displacement_petscspace_degree 3 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate 862c4762a1bSJed Brown test: 863c4762a1bSJed Brown suffix: 3d_q1_trig_elas 86430602db0SMatthew G. Knepley args: -sol_type elas_trig -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_plex_simplex 0 -displacement_petscspace_degree 1 -convest_num_refine 2 -snes_convergence_estimate 865c4762a1bSJed Brown test: 866c4762a1bSJed Brown suffix: 3d_q2_trig_elas 86730602db0SMatthew G. Knepley args: -sol_type elas_trig -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate 868c4762a1bSJed Brown test: 869c4762a1bSJed Brown suffix: 3d_q3_trig_elas 870c4762a1bSJed Brown requires: !__float128 87130602db0SMatthew G. Knepley args: -sol_type elas_trig -dm_plex_dim 3 -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate 872c4762a1bSJed Brown 873c4762a1bSJed Brown test: 874c4762a1bSJed Brown suffix: 2d_p1_axial_elas 875c4762a1bSJed Brown requires: triangle 876c4762a1bSJed Brown args: -sol_type elas_axial_disp -displacement_petscspace_degree 1 -dm_plex_separate_marker -dm_refine 2 -dmsnes_check .0001 -pc_type lu 877c4762a1bSJed Brown test: 878c4762a1bSJed Brown suffix: 2d_p2_axial_elas 879c4762a1bSJed Brown requires: triangle 880c4762a1bSJed Brown args: -sol_type elas_axial_disp -displacement_petscspace_degree 2 -dm_plex_separate_marker -dmsnes_check .0001 -pc_type lu 881c4762a1bSJed Brown test: 882c4762a1bSJed Brown suffix: 2d_p3_axial_elas 883c4762a1bSJed Brown requires: triangle 884c4762a1bSJed Brown args: -sol_type elas_axial_disp -displacement_petscspace_degree 3 -dm_plex_separate_marker -dmsnes_check .0001 -pc_type lu 885c4762a1bSJed Brown test: 886c4762a1bSJed Brown suffix: 2d_q1_axial_elas 88730602db0SMatthew G. Knepley args: -sol_type elas_axial_disp -dm_plex_simplex 0 -displacement_petscspace_degree 1 -dm_plex_separate_marker -dm_refine 1 -dmsnes_check .0001 -pc_type lu 888c4762a1bSJed Brown test: 889c4762a1bSJed Brown suffix: 2d_q2_axial_elas 89030602db0SMatthew G. Knepley args: -sol_type elas_axial_disp -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_plex_separate_marker -dmsnes_check .0001 -pc_type lu 891c4762a1bSJed Brown test: 892c4762a1bSJed Brown suffix: 2d_q3_axial_elas 89330602db0SMatthew G. Knepley args: -sol_type elas_axial_disp -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_plex_separate_marker -dmsnes_check .0001 -pc_type lu 894c4762a1bSJed Brown 895c4762a1bSJed Brown test: 896c4762a1bSJed Brown suffix: 2d_p1_uniform_elas 897c4762a1bSJed Brown requires: triangle 898c4762a1bSJed Brown args: -sol_type elas_uniform_strain -displacement_petscspace_degree 1 -dm_refine 2 -dmsnes_check .0001 -pc_type lu 899c4762a1bSJed Brown test: 900c4762a1bSJed Brown suffix: 2d_p2_uniform_elas 901c4762a1bSJed Brown requires: triangle 902c4762a1bSJed Brown args: -sol_type elas_uniform_strain -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001 -pc_type lu 903c4762a1bSJed Brown test: 904c4762a1bSJed Brown suffix: 2d_p3_uniform_elas 905c4762a1bSJed Brown requires: triangle 906c4762a1bSJed Brown args: -sol_type elas_uniform_strain -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001 -pc_type lu 907c4762a1bSJed Brown test: 908c4762a1bSJed Brown suffix: 2d_q1_uniform_elas 90930602db0SMatthew G. Knepley args: -sol_type elas_uniform_strain -dm_plex_simplex 0 -displacement_petscspace_degree 1 -dm_refine 2 -dmsnes_check .0001 -pc_type lu 910c4762a1bSJed Brown test: 911c4762a1bSJed Brown suffix: 2d_q2_uniform_elas 912c4762a1bSJed Brown requires: !single 91330602db0SMatthew G. Knepley args: -sol_type elas_uniform_strain -dm_plex_simplex 0 -displacement_petscspace_degree 2 -dm_refine 2 -dmsnes_check .0001 -pc_type lu 914c4762a1bSJed Brown test: 915c4762a1bSJed Brown suffix: 2d_q3_uniform_elas 916c4762a1bSJed Brown requires: !single 91730602db0SMatthew G. Knepley args: -sol_type elas_uniform_strain -dm_plex_simplex 0 -displacement_petscspace_degree 3 -dm_refine 2 -dmsnes_check .0001 -pc_type lu 918482dd828SMatthew G. Knepley test: 919482dd828SMatthew G. Knepley suffix: 2d_p1_uniform_elas_step 920482dd828SMatthew G. Knepley requires: triangle 921482dd828SMatthew G. Knepley args: -sol_type elas_uniform_strain -deform_type step -displacement_petscspace_degree 1 -dm_refine 2 -dmsnes_check .0001 -pc_type lu 922482dd828SMatthew G. Knepley 923482dd828SMatthew G. Knepley testset: 924482dd828SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 3,3 -deform_type step -displacement_petscspace_degree 1 -dmsnes_check .0001 -pc_type lu 925482dd828SMatthew G. Knepley 926482dd828SMatthew G. Knepley test: 927482dd828SMatthew G. Knepley suffix: 2d_q1_uniform_elas_step 928482dd828SMatthew G. Knepley args: -sol_type elas_uniform_strain -dm_refine 2 929482dd828SMatthew G. Knepley test: 930482dd828SMatthew G. Knepley suffix: 2d_q1_quad_vlap_step 931482dd828SMatthew G. Knepley args: 932482dd828SMatthew G. Knepley test: 933482dd828SMatthew G. Knepley suffix: 2d_q2_quad_vlap_step 934482dd828SMatthew G. Knepley args: -displacement_petscspace_degree 2 935482dd828SMatthew G. Knepley test: 936482dd828SMatthew G. Knepley suffix: 2d_q1_quad_mass_step 937482dd828SMatthew G. Knepley args: -sol_type mass_quad 938c4762a1bSJed Brown 939d12cd81dSMatthew G. Knepley testset: 940908b9b43SStefano Zampini filter: grep -v "variant HERMITIAN" 941d12cd81dSMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower -5,-5,-0.25 -dm_plex_box_upper 5,5,0.25 \ 942d12cd81dSMatthew G. Knepley -dm_plex_box_faces 5,5,2 -dm_plex_separate_marker -dm_refine 0 -petscpartitioner_type simple \ 943908b9b43SStefano Zampini -sol_type elas_ge 944908b9b43SStefano Zampini 945908b9b43SStefano Zampini test: 946908b9b43SStefano Zampini suffix: ge_q1_0 947908b9b43SStefano Zampini args: -displacement_petscspace_degree 1 \ 948d12cd81dSMatthew G. Knepley -snes_max_it 2 -snes_rtol 1.e-10 \ 949d12cd81dSMatthew G. Knepley -ksp_type cg -ksp_rtol 1.e-10 -ksp_max_it 100 -ksp_norm_type unpreconditioned \ 950d12cd81dSMatthew G. Knepley -pc_type gamg -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 \ 951d12cd81dSMatthew G. Knepley -pc_gamg_coarse_eq_limit 10 -pc_gamg_reuse_interpolation true \ 952bae903cbSmarkadams4 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 \ 953d12cd81dSMatthew G. Knepley -mg_levels_ksp_max_it 2 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.1 -mg_levels_pc_type jacobi \ 954d12cd81dSMatthew G. Knepley -matptap_via scalable 955*e0008caeSPierre Jolivet output_file: output/empty.out 956d12cd81dSMatthew G. Knepley test: 957f9ea58a4SMatthew G. Knepley suffix: ge_q1_gmg 958f9ea58a4SMatthew G. Knepley args: -displacement_petscspace_degree 1 \ 959f9ea58a4SMatthew G. Knepley -dm_plex_box_faces 2,2 -dm_refine_hierarchy 3 \ 960f9ea58a4SMatthew G. Knepley -snes_max_it 2 -snes_rtol 1.e-10 \ 961f9ea58a4SMatthew G. Knepley -ksp_type cg -ksp_rtol 1.e-10 -ksp_max_it 100 -ksp_norm_type unpreconditioned \ 962f9ea58a4SMatthew G. Knepley -pc_type mg -pc_mg_type full \ 963f9ea58a4SMatthew G. Knepley -mg_levels_ksp_max_it 4 -mg_levels_esteig_ksp_type cg \ 964f9ea58a4SMatthew G. Knepley -mg_levels_esteig_ksp_max_it 10 -mg_levels_ksp_chebyshev_esteig 0,0.1,0,1.1 \ 965f9ea58a4SMatthew G. Knepley -mg_levels_pc_type jacobi 966*e0008caeSPierre Jolivet output_file: output/empty.out 967f9ea58a4SMatthew G. Knepley test: 968908b9b43SStefano Zampini nsize: 5 969908b9b43SStefano Zampini suffix: ge_q1_gdsw 970908b9b43SStefano Zampini args: -snes_max_it 1 -ksp_type cg -ksp_norm_type natural -displacement_petscspace_degree 1 -snes_monitor_short -ksp_monitor_short -pc_type mg -pc_mg_adapt_interp_coarse_space gdsw -pc_mg_levels 2 -pc_mg_galerkin -mg_levels_pc_type bjacobi -mg_levels_esteig_ksp_type cg -mg_levels_sub_pc_type icc -mg_coarse_redundant_pc_type cholesky -ksp_view 971d12cd81dSMatthew G. Knepley 972c4762a1bSJed Brown TEST*/ 973