1c4762a1bSJed Brown static char help[] = "Nonlinear elasticity problem in 3d with simplicial finite elements.\n\ 2c4762a1bSJed Brown We solve a nonlinear elasticity problem, modelled as an incompressible Neo-Hookean solid, \n\ 3c4762a1bSJed Brown with pressure loading in a rectangular domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n"; 4c4762a1bSJed Brown 5c4762a1bSJed Brown /* 6c4762a1bSJed Brown Nonlinear elasticity problem, which we discretize using the finite 7c4762a1bSJed Brown element method on an unstructured mesh. This uses both Dirichlet boundary conditions (fixed faces) 8c4762a1bSJed Brown and nonlinear Neumann boundary conditions (pressure loading). 9c4762a1bSJed Brown The Lagrangian density (modulo boundary conditions) for this problem is given by 10c4762a1bSJed Brown \begin{equation} 11c4762a1bSJed Brown \frac{\mu}{2} (\mathrm{Tr}{C}-3) + J p + \frac{\kappa}{2} (J-1). 12c4762a1bSJed Brown \end{equation} 13c4762a1bSJed Brown 14c4762a1bSJed Brown Discretization: 15c4762a1bSJed Brown 16c4762a1bSJed Brown We use PetscFE to generate a tabulation of the finite element basis functions 17c4762a1bSJed Brown at quadrature points. We can currently generate an arbitrary order Lagrange 18c4762a1bSJed Brown element. 19c4762a1bSJed Brown 20c4762a1bSJed Brown Field Data: 21c4762a1bSJed Brown 22c4762a1bSJed Brown DMPLEX data is organized by point, and the closure operation just stacks up the 23c4762a1bSJed Brown data from each sieve point in the closure. Thus, for a P_2-P_1 Stokes element, we 24c4762a1bSJed Brown have 25c4762a1bSJed Brown 26c4762a1bSJed Brown cl{e} = {f e_0 e_1 e_2 v_0 v_1 v_2} 27c4762a1bSJed Brown x = [u_{e_0} v_{e_0} u_{e_1} v_{e_1} u_{e_2} v_{e_2} u_{v_0} v_{v_0} p_{v_0} u_{v_1} v_{v_1} p_{v_1} u_{v_2} v_{v_2} p_{v_2}] 28c4762a1bSJed Brown 29c4762a1bSJed Brown The problem here is that we would like to loop over each field separately for 30c4762a1bSJed Brown integration. Therefore, the closure visitor in DMPlexVecGetClosure() reorders 31c4762a1bSJed Brown the data so that each field is contiguous 32c4762a1bSJed Brown 33c4762a1bSJed Brown x' = [u_{e_0} v_{e_0} u_{e_1} v_{e_1} u_{e_2} v_{e_2} u_{v_0} v_{v_0} u_{v_1} v_{v_1} u_{v_2} v_{v_2} p_{v_0} p_{v_1} p_{v_2}] 34c4762a1bSJed Brown 35c4762a1bSJed Brown Likewise, DMPlexVecSetClosure() takes data partitioned by field, and correctly 36c4762a1bSJed Brown puts it into the Sieve ordering. 37c4762a1bSJed Brown 38c4762a1bSJed Brown */ 39c4762a1bSJed Brown 40c4762a1bSJed Brown #include <petscdmplex.h> 41c4762a1bSJed Brown #include <petscsnes.h> 42c4762a1bSJed Brown #include <petscds.h> 43c4762a1bSJed Brown 44c4762a1bSJed Brown typedef enum {RUN_FULL, RUN_TEST} RunType; 45c4762a1bSJed Brown 46c4762a1bSJed Brown typedef struct { 47c4762a1bSJed Brown RunType runType; /* Whether to run tests, or solve the full problem */ 48c4762a1bSJed Brown PetscReal mu; /* The shear modulus */ 49c4762a1bSJed Brown PetscReal p_wall; /* The wall pressure */ 50c4762a1bSJed Brown } AppCtx; 51c4762a1bSJed Brown 52c4762a1bSJed Brown #if 0 539fbee547SJacob Faibussowitsch static inline void Det2D(PetscReal *detJ, const PetscReal J[]) 54c4762a1bSJed Brown { 55c4762a1bSJed Brown *detJ = J[0]*J[3] - J[1]*J[2]; 56c4762a1bSJed Brown } 57c4762a1bSJed Brown #endif 58c4762a1bSJed Brown 599fbee547SJacob Faibussowitsch static inline void Det3D(PetscReal *detJ, const PetscScalar J[]) 60c4762a1bSJed Brown { 61c4762a1bSJed Brown *detJ = PetscRealPart(J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) + 62c4762a1bSJed Brown J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) + 63c4762a1bSJed Brown J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0])); 64c4762a1bSJed Brown } 65c4762a1bSJed Brown 66c4762a1bSJed Brown #if 0 679fbee547SJacob Faibussowitsch static inline void Cof2D(PetscReal C[], const PetscReal A[]) 68c4762a1bSJed Brown { 69c4762a1bSJed Brown C[0] = A[3]; 70c4762a1bSJed Brown C[1] = -A[2]; 71c4762a1bSJed Brown C[2] = -A[1]; 72c4762a1bSJed Brown C[3] = A[0]; 73c4762a1bSJed Brown } 74c4762a1bSJed Brown #endif 75c4762a1bSJed Brown 769fbee547SJacob Faibussowitsch static inline void Cof3D(PetscReal C[], const PetscScalar A[]) 77c4762a1bSJed Brown { 78c4762a1bSJed Brown C[0*3+0] = PetscRealPart(A[1*3+1]*A[2*3+2] - A[1*3+2]*A[2*3+1]); 79c4762a1bSJed Brown C[0*3+1] = PetscRealPart(A[1*3+2]*A[2*3+0] - A[1*3+0]*A[2*3+2]); 80c4762a1bSJed Brown C[0*3+2] = PetscRealPart(A[1*3+0]*A[2*3+1] - A[1*3+1]*A[2*3+0]); 81c4762a1bSJed Brown C[1*3+0] = PetscRealPart(A[0*3+2]*A[2*3+1] - A[0*3+1]*A[2*3+2]); 82c4762a1bSJed Brown C[1*3+1] = PetscRealPart(A[0*3+0]*A[2*3+2] - A[0*3+2]*A[2*3+0]); 83c4762a1bSJed Brown C[1*3+2] = PetscRealPart(A[0*3+1]*A[2*3+0] - A[0*3+0]*A[2*3+1]); 84c4762a1bSJed Brown C[2*3+0] = PetscRealPart(A[0*3+1]*A[1*3+2] - A[0*3+2]*A[1*3+1]); 85c4762a1bSJed Brown C[2*3+1] = PetscRealPart(A[0*3+2]*A[1*3+0] - A[0*3+0]*A[1*3+2]); 86c4762a1bSJed Brown C[2*3+2] = PetscRealPart(A[0*3+0]*A[1*3+1] - A[0*3+1]*A[1*3+0]); 87c4762a1bSJed Brown } 88c4762a1bSJed Brown 89c4762a1bSJed Brown PetscErrorCode zero_scalar(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 90c4762a1bSJed Brown { 91c4762a1bSJed Brown u[0] = 0.0; 92c4762a1bSJed Brown return 0; 93c4762a1bSJed Brown } 94c4762a1bSJed Brown 95c4762a1bSJed Brown PetscErrorCode zero_vector(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 96c4762a1bSJed Brown { 97c4762a1bSJed Brown const PetscInt Ncomp = dim; 98c4762a1bSJed Brown 99c4762a1bSJed Brown PetscInt comp; 100c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) u[comp] = 0.0; 101c4762a1bSJed Brown return 0; 102c4762a1bSJed Brown } 103c4762a1bSJed Brown 104c4762a1bSJed Brown PetscErrorCode coordinates(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 105c4762a1bSJed Brown { 106c4762a1bSJed Brown const PetscInt Ncomp = dim; 107c4762a1bSJed Brown 108c4762a1bSJed Brown PetscInt comp; 109c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) u[comp] = x[comp]; 110c4762a1bSJed Brown return 0; 111c4762a1bSJed Brown } 112c4762a1bSJed Brown 113c4762a1bSJed Brown PetscErrorCode elasticityMaterial(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 114c4762a1bSJed Brown { 115c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 116c4762a1bSJed Brown u[0] = user->mu; 117c4762a1bSJed Brown return 0; 118c4762a1bSJed Brown } 119c4762a1bSJed Brown 120c4762a1bSJed Brown PetscErrorCode wallPressure(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 121c4762a1bSJed Brown { 122c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 123c4762a1bSJed Brown u[0] = user->p_wall; 124c4762a1bSJed Brown return 0; 125c4762a1bSJed Brown } 126c4762a1bSJed Brown 127c4762a1bSJed Brown void f1_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 128c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 129c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 130c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 131c4762a1bSJed Brown { 132c4762a1bSJed Brown const PetscInt Ncomp = dim; 133c4762a1bSJed Brown const PetscReal mu = PetscRealPart(a[0]), kappa = 3.0; 134c4762a1bSJed Brown PetscReal cofu_x[9/*Ncomp*dim*/], detu_x, p = PetscRealPart(u[Ncomp]); 135c4762a1bSJed Brown PetscInt comp, d; 136c4762a1bSJed Brown 137c4762a1bSJed Brown Cof3D(cofu_x, u_x); 138c4762a1bSJed Brown Det3D(&detu_x, u_x); 139c4762a1bSJed Brown p += kappa * (detu_x - 1.0); 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* f1 is the first Piola-Kirchhoff tensor */ 142c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) { 143c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 144c4762a1bSJed Brown f1[comp*dim+d] = mu * u_x[comp*dim+d] + p * cofu_x[comp*dim+d]; 145c4762a1bSJed Brown } 146c4762a1bSJed Brown } 147c4762a1bSJed Brown } 148c4762a1bSJed Brown 149c4762a1bSJed Brown void g3_uu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 150c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 151c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 152c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 153c4762a1bSJed Brown { 154c4762a1bSJed Brown const PetscInt Ncomp = dim; 155c4762a1bSJed Brown const PetscReal mu = PetscRealPart(a[0]), kappa = 3.0; 156c4762a1bSJed Brown PetscReal cofu_x[9/*Ncomp*dim*/], detu_x, pp, pm, p = PetscRealPart(u[Ncomp]); 157c4762a1bSJed Brown PetscInt compI, compJ, d1, d2; 158c4762a1bSJed Brown 159c4762a1bSJed Brown Cof3D(cofu_x, u_x); 160c4762a1bSJed Brown Det3D(&detu_x, u_x); 161c4762a1bSJed Brown p += kappa * (detu_x - 1.0); 162c4762a1bSJed Brown pp = p/detu_x + kappa; 163c4762a1bSJed Brown pm = p/detu_x; 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* g3 is the first elasticity tensor, i.e. A_i^I_j^J = S^{IJ} g_{ij} + C^{KILJ} F^k_K F^l_L g_{kj} g_{li} */ 166c4762a1bSJed Brown for (compI = 0; compI < Ncomp; ++compI) { 167c4762a1bSJed Brown for (compJ = 0; compJ < Ncomp; ++compJ) { 168c4762a1bSJed Brown const PetscReal G = (compI == compJ) ? 1.0 : 0.0; 169c4762a1bSJed Brown 170c4762a1bSJed Brown for (d1 = 0; d1 < dim; ++d1) { 171c4762a1bSJed Brown for (d2 = 0; d2 < dim; ++d2) { 172c4762a1bSJed Brown const PetscReal g = (d1 == d2) ? 1.0 : 0.0; 173c4762a1bSJed Brown 174c4762a1bSJed Brown g3[((compI*Ncomp+compJ)*dim+d1)*dim+d2] = g * G * mu + pp * cofu_x[compI*dim+d1] * cofu_x[compJ*dim+d2] - pm * cofu_x[compI*dim+d2] * cofu_x[compJ*dim+d1]; 175c4762a1bSJed Brown } 176c4762a1bSJed Brown } 177c4762a1bSJed Brown } 178c4762a1bSJed Brown } 179c4762a1bSJed Brown } 180c4762a1bSJed Brown 181c4762a1bSJed Brown void f0_bd_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 182c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 183c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 184c4762a1bSJed Brown PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 185c4762a1bSJed Brown { 186c4762a1bSJed Brown const PetscInt Ncomp = dim; 187c4762a1bSJed Brown const PetscScalar p = a[aOff[1]]; 188c4762a1bSJed Brown PetscReal cofu_x[9/*Ncomp*dim*/]; 189c4762a1bSJed Brown PetscInt comp, d; 190c4762a1bSJed Brown 191c4762a1bSJed Brown Cof3D(cofu_x, u_x); 192c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) { 193c4762a1bSJed Brown for (d = 0, f0[comp] = 0.0; d < dim; ++d) f0[comp] += cofu_x[comp*dim+d] * n[d]; 194c4762a1bSJed Brown f0[comp] *= p; 195c4762a1bSJed Brown } 196c4762a1bSJed Brown } 197c4762a1bSJed Brown 198c4762a1bSJed Brown void g1_bd_uu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 199c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 200c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 201c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 202c4762a1bSJed Brown { 203c4762a1bSJed Brown const PetscInt Ncomp = dim; 204c4762a1bSJed Brown PetscScalar p = a[aOff[1]]; 205c4762a1bSJed Brown PetscReal cofu_x[9/*Ncomp*dim*/], m[3/*Ncomp*/], detu_x; 206c4762a1bSJed Brown PetscInt comp, compI, compJ, d; 207c4762a1bSJed Brown 208c4762a1bSJed Brown Cof3D(cofu_x, u_x); 209c4762a1bSJed Brown Det3D(&detu_x, u_x); 210c4762a1bSJed Brown p /= detu_x; 211c4762a1bSJed Brown 212c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) for (d = 0, m[comp] = 0.0; d < dim; ++d) m[comp] += cofu_x[comp*dim+d] * n[d]; 213c4762a1bSJed Brown for (compI = 0; compI < Ncomp; ++compI) { 214c4762a1bSJed Brown for (compJ = 0; compJ < Ncomp; ++compJ) { 215c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 216c4762a1bSJed Brown g1[(compI*Ncomp+compJ)*dim+d] = p * (m[compI] * cofu_x[compJ*dim+d] - cofu_x[compI*dim+d] * m[compJ]); 217c4762a1bSJed Brown } 218c4762a1bSJed Brown } 219c4762a1bSJed Brown } 220c4762a1bSJed Brown } 221c4762a1bSJed Brown 222c4762a1bSJed Brown void f0_p_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 223c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 224c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 225c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 226c4762a1bSJed Brown { 227c4762a1bSJed Brown PetscReal detu_x; 228c4762a1bSJed Brown Det3D(&detu_x, u_x); 229c4762a1bSJed Brown f0[0] = detu_x - 1.0; 230c4762a1bSJed Brown } 231c4762a1bSJed Brown 232c4762a1bSJed Brown void g1_pu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 233c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 234c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 235c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 236c4762a1bSJed Brown { 237c4762a1bSJed Brown PetscReal cofu_x[9/*Ncomp*dim*/]; 238c4762a1bSJed Brown PetscInt compI, d; 239c4762a1bSJed Brown 240c4762a1bSJed Brown Cof3D(cofu_x, u_x); 241c4762a1bSJed Brown for (compI = 0; compI < dim; ++compI) 242c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[compI*dim+d] = cofu_x[compI*dim+d]; 243c4762a1bSJed Brown } 244c4762a1bSJed Brown 245c4762a1bSJed Brown void g2_up_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 246c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 247c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 248c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 249c4762a1bSJed Brown { 250c4762a1bSJed Brown PetscReal cofu_x[9/*Ncomp*dim*/]; 251c4762a1bSJed Brown PetscInt compI, d; 252c4762a1bSJed Brown 253c4762a1bSJed Brown Cof3D(cofu_x, u_x); 254c4762a1bSJed Brown for (compI = 0; compI < dim; ++compI) 255c4762a1bSJed Brown for (d = 0; d < dim; ++d) g2[compI*dim+d] = cofu_x[compI*dim+d]; 256c4762a1bSJed Brown } 257c4762a1bSJed Brown 258c4762a1bSJed Brown PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 259c4762a1bSJed Brown { 260c4762a1bSJed Brown const char *runTypes[2] = {"full", "test"}; 261c4762a1bSJed Brown PetscInt run; 262c4762a1bSJed Brown 263c4762a1bSJed Brown PetscFunctionBeginUser; 264c4762a1bSJed Brown options->runType = RUN_FULL; 265c4762a1bSJed Brown options->mu = 1.0; 266c4762a1bSJed Brown options->p_wall = 0.4; 267d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Nonlinear elasticity problem options", "DMPLEX"); 268c4762a1bSJed Brown run = options->runType; 2699566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-run_type", "The run type", "ex77.c", runTypes, 2, runTypes[options->runType], &run, NULL)); 270c4762a1bSJed Brown options->runType = (RunType) run; 2719566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-shear_modulus", "The shear modulus", "ex77.c", options->mu, &options->mu, NULL)); 2729566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-wall_pressure", "The wall pressure", "ex77.c", options->p_wall, &options->p_wall, NULL)); 273d0609cedSBarry Smith PetscOptionsEnd(); 274c4762a1bSJed Brown PetscFunctionReturn(0); 275c4762a1bSJed Brown } 276c4762a1bSJed Brown 277c4762a1bSJed Brown PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 278c4762a1bSJed Brown { 279c4762a1bSJed Brown PetscFunctionBeginUser; 28030602db0SMatthew G. Knepley /* TODO The P1 coordinate space gives wrong results when compared to the affine version. Track this down */ 28130602db0SMatthew G. Knepley if (0) { 2829566063dSJacob Faibussowitsch PetscCall(DMPlexCreateBoxMesh(comm, 3, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm)); 28330602db0SMatthew G. Knepley } else { 2849566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 2859566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 28630602db0SMatthew G. Knepley } 2879566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 288c4762a1bSJed Brown /* Label the faces (bit of a hack here, until it is properly implemented for simplices) */ 2899566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view")); 290c4762a1bSJed Brown { 291c4762a1bSJed Brown DM cdm; 292c4762a1bSJed Brown DMLabel label; 293c4762a1bSJed Brown IS is; 29430602db0SMatthew G. Knepley PetscInt d, dim, b, f, Nf; 295c4762a1bSJed Brown const PetscInt *faces; 296c4762a1bSJed Brown PetscInt csize; 297c4762a1bSJed Brown PetscScalar *coords = NULL; 298c4762a1bSJed Brown PetscSection cs; 299c4762a1bSJed Brown Vec coordinates ; 300c4762a1bSJed Brown 3019566063dSJacob Faibussowitsch PetscCall(DMGetDimension(*dm, &dim)); 3029566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(*dm, "boundary")); 3039566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dm, "boundary", &label)); 3049566063dSJacob Faibussowitsch PetscCall(DMPlexMarkBoundaryFaces(*dm, 1, label)); 3059566063dSJacob Faibussowitsch PetscCall(DMGetStratumIS(*dm, "boundary", 1, &is)); 306c4762a1bSJed Brown if (is) { 307c4762a1bSJed Brown PetscReal faceCoord; 308c4762a1bSJed Brown PetscInt v; 309c4762a1bSJed Brown 3109566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &Nf)); 3119566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &faces)); 312c4762a1bSJed Brown 3139566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(*dm, &coordinates)); 3149566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(*dm, &cdm)); 3159566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(cdm, &cs)); 316c4762a1bSJed Brown 317c4762a1bSJed Brown /* Check for each boundary face if any component of its centroid is either 0.0 or 1.0 */ 318c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 3199566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(cdm, cs, coordinates, faces[f], &csize, &coords)); 320c4762a1bSJed Brown /* Calculate mean coordinate vector */ 321c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 322c4762a1bSJed Brown const PetscInt Nv = csize/dim; 323c4762a1bSJed Brown faceCoord = 0.0; 324c4762a1bSJed Brown for (v = 0; v < Nv; ++v) faceCoord += PetscRealPart(coords[v*dim+d]); 325c4762a1bSJed Brown faceCoord /= Nv; 326c4762a1bSJed Brown for (b = 0; b < 2; ++b) { 327c4762a1bSJed Brown if (PetscAbs(faceCoord - b*1.0) < PETSC_SMALL) { 3289566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, "Faces", faces[f], d*2+b+1)); 329c4762a1bSJed Brown } 330c4762a1bSJed Brown } 331c4762a1bSJed Brown } 3329566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(cdm, cs, coordinates, faces[f], &csize, &coords)); 333c4762a1bSJed Brown } 3349566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &faces)); 335c4762a1bSJed Brown } 3369566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 337c4762a1bSJed Brown } 3389566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 339c4762a1bSJed Brown PetscFunctionReturn(0); 340c4762a1bSJed Brown } 341c4762a1bSJed Brown 342c4762a1bSJed Brown PetscErrorCode SetupProblem(DM dm, PetscInt dim, AppCtx *user) 343c4762a1bSJed Brown { 34445480ffeSMatthew G. Knepley PetscDS ds; 34545480ffeSMatthew G. Knepley PetscWeakForm wf; 34645480ffeSMatthew G. Knepley DMLabel label; 34745480ffeSMatthew G. Knepley PetscInt bd; 348c4762a1bSJed Brown 349c4762a1bSJed Brown PetscFunctionBeginUser; 3509566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 3519566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, NULL, f1_u_3d)); 3529566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_p_3d, NULL)); 3539566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu_3d)); 3549566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up_3d, NULL)); 3559566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu_3d, NULL, NULL)); 356c4762a1bSJed Brown 3579566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Faces", &label)); 3589566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "pressure", label, 0, NULL, 0, 0, NULL, NULL, NULL, user, &bd)); 3599566063dSJacob Faibussowitsch PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 3609566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, 0, f0_bd_u_3d, 0, NULL)); 3619566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexBdJacobian(wf, label, 1, 0, 0, 0, 0, NULL, 0, g1_bd_uu_3d, 0, NULL, 0, NULL)); 362c4762a1bSJed Brown 3639566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "fixed", label, 0, NULL, 0, 0, NULL, (void (*)(void)) coordinates, NULL, user, NULL)); 364c4762a1bSJed Brown PetscFunctionReturn(0); 365c4762a1bSJed Brown } 366c4762a1bSJed Brown 367c4762a1bSJed Brown PetscErrorCode SetupMaterial(DM dm, DM dmAux, AppCtx *user) 368c4762a1bSJed Brown { 369c4762a1bSJed Brown PetscErrorCode (*matFuncs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx) = {elasticityMaterial, wallPressure}; 370c4762a1bSJed Brown Vec A; 371c4762a1bSJed Brown void *ctxs[2]; 372c4762a1bSJed Brown 373c4762a1bSJed Brown PetscFunctionBegin; 374c4762a1bSJed Brown ctxs[0] = user; ctxs[1] = user; 3759566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dmAux, &A)); 3769566063dSJacob Faibussowitsch PetscCall(DMProjectFunctionLocal(dmAux, 0.0, matFuncs, ctxs, INSERT_ALL_VALUES, A)); 3779566063dSJacob Faibussowitsch PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, A)); 3789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&A)); 379c4762a1bSJed Brown PetscFunctionReturn(0); 380c4762a1bSJed Brown } 381c4762a1bSJed Brown 382c4762a1bSJed Brown PetscErrorCode SetupNearNullSpace(DM dm, AppCtx *user) 383c4762a1bSJed Brown { 384c4762a1bSJed Brown /* Set up the near null space (a.k.a. rigid body modes) that will be used by the multigrid preconditioner */ 385c4762a1bSJed Brown DM subdm; 386c4762a1bSJed Brown MatNullSpace nearNullSpace; 387c4762a1bSJed Brown PetscInt fields = 0; 388c4762a1bSJed Brown PetscObject deformation; 389c4762a1bSJed Brown 390c4762a1bSJed Brown PetscFunctionBeginUser; 3919566063dSJacob Faibussowitsch PetscCall(DMCreateSubDM(dm, 1, &fields, NULL, &subdm)); 3929566063dSJacob Faibussowitsch PetscCall(DMPlexCreateRigidBody(subdm, 0, &nearNullSpace)); 3939566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, 0, NULL, &deformation)); 3949566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose(deformation, "nearnullspace", (PetscObject) nearNullSpace)); 3959566063dSJacob Faibussowitsch PetscCall(DMDestroy(&subdm)); 3969566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nearNullSpace)); 397c4762a1bSJed Brown PetscFunctionReturn(0); 398c4762a1bSJed Brown } 399c4762a1bSJed Brown 400c4762a1bSJed Brown static PetscErrorCode SetupAuxDM(DM dm, PetscInt NfAux, PetscFE feAux[], AppCtx *user) 401c4762a1bSJed Brown { 402c4762a1bSJed Brown DM dmAux, coordDM; 403c4762a1bSJed Brown PetscInt f; 404c4762a1bSJed Brown 405c4762a1bSJed Brown PetscFunctionBegin; 406c4762a1bSJed Brown /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */ 4079566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &coordDM)); 408c4762a1bSJed Brown if (!feAux) PetscFunctionReturn(0); 4099566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmAux)); 4109566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateDM(dmAux, coordDM)); 4119566063dSJacob Faibussowitsch for (f = 0; f < NfAux; ++f) PetscCall(DMSetField(dmAux, f, NULL, (PetscObject) feAux[f])); 4129566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dmAux)); 4139566063dSJacob Faibussowitsch PetscCall(SetupMaterial(dm, dmAux, user)); 4149566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmAux)); 415c4762a1bSJed Brown PetscFunctionReturn(0); 416c4762a1bSJed Brown } 417c4762a1bSJed Brown 418c4762a1bSJed Brown PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 419c4762a1bSJed Brown { 420c4762a1bSJed Brown DM cdm = dm; 421c4762a1bSJed Brown PetscFE fe[2], feAux[2]; 42230602db0SMatthew G. Knepley PetscBool simplex; 42330602db0SMatthew G. Knepley PetscInt dim; 424c4762a1bSJed Brown MPI_Comm comm; 425c4762a1bSJed Brown 426c4762a1bSJed Brown PetscFunctionBeginUser; 4279566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) dm, &comm)); 4289566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 4299566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 430c4762a1bSJed Brown /* Create finite element */ 4319566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "def_", PETSC_DEFAULT, &fe[0])); 4329566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe[0], "deformation")); 4339566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1])); 4349566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 435c4762a1bSJed Brown 4369566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) fe[1], "pressure")); 437c4762a1bSJed Brown 4389566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "elastMat_", PETSC_DEFAULT, &feAux[0])); 4399566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) feAux[0], "elasticityMaterial")); 4409566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], feAux[0])); 441c4762a1bSJed Brown /* It is not yet possible to define a field on a submesh (e.g. a boundary), so we will use a normal finite element */ 4429566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "wall_pres_", PETSC_DEFAULT, &feAux[1])); 4439566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) feAux[1], "wall_pressure")); 4449566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], feAux[1])); 445c4762a1bSJed Brown 446c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 4479566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe[0])); 4489566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 1, NULL, (PetscObject) fe[1])); 4499566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 4509566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, dim, user)); 451c4762a1bSJed Brown while (cdm) { 4529566063dSJacob Faibussowitsch PetscCall(SetupAuxDM(cdm, 2, feAux, user)); 4539566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 4549566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 455c4762a1bSJed Brown } 4569566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[0])); 4579566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[1])); 4589566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&feAux[0])); 4599566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&feAux[1])); 460c4762a1bSJed Brown PetscFunctionReturn(0); 461c4762a1bSJed Brown } 462c4762a1bSJed Brown 463c4762a1bSJed Brown int main(int argc, char **argv) 464c4762a1bSJed Brown { 465c4762a1bSJed Brown SNES snes; /* nonlinear solver */ 466c4762a1bSJed Brown DM dm; /* problem definition */ 467c4762a1bSJed Brown Vec u,r; /* solution, residual vectors */ 468c4762a1bSJed Brown Mat A,J; /* Jacobian matrix */ 469c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 470c4762a1bSJed Brown PetscInt its; /* iterations for convergence */ 471c4762a1bSJed Brown 472*327415f7SBarry Smith PetscFunctionBeginUser; 4739566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 4749566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 4759566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 4769566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 4779566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 4789566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dm, &user)); 479c4762a1bSJed Brown 4809566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, &user)); 4819566063dSJacob Faibussowitsch PetscCall(DMPlexCreateClosureIndex(dm, NULL)); 4829566063dSJacob Faibussowitsch PetscCall(SetupNearNullSpace(dm, &user)); 483c4762a1bSJed Brown 4849566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 4859566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r)); 486c4762a1bSJed Brown 4879566063dSJacob Faibussowitsch PetscCall(DMSetMatType(dm,MATAIJ)); 4889566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &J)); 489c4762a1bSJed Brown A = J; 490c4762a1bSJed Brown 4919566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm,&user,&user,&user)); 4929566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, A, J, NULL, NULL)); 493c4762a1bSJed Brown 4949566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 495c4762a1bSJed Brown 496c4762a1bSJed Brown { 497c4762a1bSJed Brown PetscErrorCode (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void* ctx); 498c4762a1bSJed Brown initialGuess[0] = coordinates; initialGuess[1] = zero_scalar; 4999566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u)); 500c4762a1bSJed Brown } 501c4762a1bSJed Brown 502c4762a1bSJed Brown if (user.runType == RUN_FULL) { 5039566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 5049566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 50563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its)); 506c4762a1bSJed Brown } else { 507c4762a1bSJed Brown PetscReal res = 0.0; 508c4762a1bSJed Brown 509c4762a1bSJed Brown /* Check initial guess */ 5109566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n")); 5119566063dSJacob Faibussowitsch PetscCall(VecView(u, PETSC_VIEWER_STDOUT_WORLD)); 512c4762a1bSJed Brown /* Check residual */ 5139566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, u, r)); 5149566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n")); 5159566063dSJacob Faibussowitsch PetscCall(VecChop(r, 1.0e-10)); 5169566063dSJacob Faibussowitsch PetscCall(VecView(r, PETSC_VIEWER_STDOUT_WORLD)); 5179566063dSJacob Faibussowitsch PetscCall(VecNorm(r, NORM_2, &res)); 5189566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", (double)res)); 519c4762a1bSJed Brown /* Check Jacobian */ 520c4762a1bSJed Brown { 521c4762a1bSJed Brown Vec b; 522c4762a1bSJed Brown 5239566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, u, A, A)); 5249566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &b)); 5259566063dSJacob Faibussowitsch PetscCall(VecSet(r, 0.0)); 5269566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, r, b)); 5279566063dSJacob Faibussowitsch PetscCall(MatMult(A, u, r)); 5289566063dSJacob Faibussowitsch PetscCall(VecAXPY(r, 1.0, b)); 5299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 5309566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Au - b = Au + F(0)\n")); 5319566063dSJacob Faibussowitsch PetscCall(VecChop(r, 1.0e-10)); 5329566063dSJacob Faibussowitsch PetscCall(VecView(r, PETSC_VIEWER_STDOUT_WORLD)); 5339566063dSJacob Faibussowitsch PetscCall(VecNorm(r, NORM_2, &res)); 5349566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", (double)res)); 535c4762a1bSJed Brown } 536c4762a1bSJed Brown } 5379566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 538c4762a1bSJed Brown 5399566063dSJacob Faibussowitsch if (A != J) PetscCall(MatDestroy(&A)); 5409566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 5419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 5429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 5439566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 5449566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 5459566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 546b122ec5aSJacob Faibussowitsch return 0; 547c4762a1bSJed Brown } 548c4762a1bSJed Brown 549c4762a1bSJed Brown /*TEST 550c4762a1bSJed Brown 551c4762a1bSJed Brown build: 552c4762a1bSJed Brown requires: !complex 553c4762a1bSJed Brown 55430602db0SMatthew G. Knepley testset: 55530602db0SMatthew G. Knepley requires: ctetgen 55630602db0SMatthew G. Knepley args: -run_type full -dm_plex_dim 3 \ 55730602db0SMatthew G. Knepley -def_petscspace_degree 2 -pres_petscspace_degree 1 -elastMat_petscspace_degree 0 -wall_pres_petscspace_degree 0 \ 55830602db0SMatthew G. Knepley -snes_rtol 1e-05 -snes_monitor_short -snes_converged_reason \ 55930602db0SMatthew G. Knepley -ksp_type fgmres -ksp_rtol 1e-10 -ksp_monitor_short -ksp_converged_reason \ 56030602db0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper \ 56130602db0SMatthew G. Knepley -fieldsplit_deformation_ksp_type preonly -fieldsplit_deformation_pc_type lu \ 56230602db0SMatthew G. Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 56330602db0SMatthew G. Knepley 564c4762a1bSJed Brown test: 565c4762a1bSJed Brown suffix: 0 56630602db0SMatthew G. Knepley requires: !single 56730602db0SMatthew G. Knepley args: -dm_refine 2 \ 56830602db0SMatthew G. Knepley -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4 569c4762a1bSJed Brown 570c4762a1bSJed Brown test: 571c4762a1bSJed Brown suffix: 1 57230602db0SMatthew G. Knepley requires: superlu_dist 573c4762a1bSJed Brown nsize: 2 574e600fa54SMatthew G. Knepley args: -dm_refine 0 -petscpartitioner_type simple \ 57530602db0SMatthew G. Knepley -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4 576c4762a1bSJed Brown timeoutfactor: 2 577c4762a1bSJed Brown 578c4762a1bSJed Brown test: 579c4762a1bSJed Brown suffix: 4 58030602db0SMatthew G. Knepley requires: superlu_dist 581c4762a1bSJed Brown nsize: 2 582e600fa54SMatthew G. Knepley args: -dm_refine 0 -petscpartitioner_type simple \ 58330602db0SMatthew G. Knepley -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4 584c4762a1bSJed Brown output_file: output/ex77_1.out 585c4762a1bSJed Brown 586c4762a1bSJed Brown test: 58730602db0SMatthew G. Knepley suffix: 2 58830602db0SMatthew G. Knepley requires: !single 58930602db0SMatthew G. Knepley args: -dm_refine 2 \ 59030602db0SMatthew G. Knepley -bc_fixed 3,4,5,6 -bc_pressure 2 -wall_pressure 1.0 591c4762a1bSJed Brown 592c4762a1bSJed Brown test: 593c4762a1bSJed Brown suffix: 2_par 59430602db0SMatthew G. Knepley requires: superlu_dist !single 595c4762a1bSJed Brown nsize: 4 596e600fa54SMatthew G. Knepley args: -dm_refine 2 -petscpartitioner_type simple \ 59730602db0SMatthew G. Knepley -bc_fixed 3,4,5,6 -bc_pressure 2 -wall_pressure 1.0 598c4762a1bSJed Brown output_file: output/ex77_2.out 599c4762a1bSJed Brown 600c4762a1bSJed Brown TEST*/ 601