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 449371c9d4SSatish Balay typedef enum { 459371c9d4SSatish Balay RUN_FULL, 469371c9d4SSatish Balay RUN_TEST 479371c9d4SSatish Balay } RunType; 48c4762a1bSJed Brown 49c4762a1bSJed Brown typedef struct { 50c4762a1bSJed Brown RunType runType; /* Whether to run tests, or solve the full problem */ 51c4762a1bSJed Brown PetscReal mu; /* The shear modulus */ 52c4762a1bSJed Brown PetscReal p_wall; /* The wall pressure */ 53c4762a1bSJed Brown } AppCtx; 54c4762a1bSJed Brown 55c4762a1bSJed Brown #if 0 569fbee547SJacob Faibussowitsch static inline void Det2D(PetscReal *detJ, const PetscReal J[]) 57c4762a1bSJed Brown { 58c4762a1bSJed Brown *detJ = J[0]*J[3] - J[1]*J[2]; 59c4762a1bSJed Brown } 60c4762a1bSJed Brown #endif 61c4762a1bSJed Brown 629371c9d4SSatish Balay static inline void Det3D(PetscReal *detJ, const PetscScalar J[]) { 639371c9d4SSatish Balay *detJ = PetscRealPart(J[0 * 3 + 0] * (J[1 * 3 + 1] * J[2 * 3 + 2] - J[1 * 3 + 2] * J[2 * 3 + 1]) + J[0 * 3 + 1] * (J[1 * 3 + 2] * J[2 * 3 + 0] - J[1 * 3 + 0] * J[2 * 3 + 2]) + 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 769371c9d4SSatish Balay static inline void Cof3D(PetscReal C[], const PetscScalar A[]) { 77c4762a1bSJed Brown C[0 * 3 + 0] = PetscRealPart(A[1 * 3 + 1] * A[2 * 3 + 2] - A[1 * 3 + 2] * A[2 * 3 + 1]); 78c4762a1bSJed Brown C[0 * 3 + 1] = PetscRealPart(A[1 * 3 + 2] * A[2 * 3 + 0] - A[1 * 3 + 0] * A[2 * 3 + 2]); 79c4762a1bSJed Brown C[0 * 3 + 2] = PetscRealPart(A[1 * 3 + 0] * A[2 * 3 + 1] - A[1 * 3 + 1] * A[2 * 3 + 0]); 80c4762a1bSJed Brown C[1 * 3 + 0] = PetscRealPart(A[0 * 3 + 2] * A[2 * 3 + 1] - A[0 * 3 + 1] * A[2 * 3 + 2]); 81c4762a1bSJed Brown C[1 * 3 + 1] = PetscRealPart(A[0 * 3 + 0] * A[2 * 3 + 2] - A[0 * 3 + 2] * A[2 * 3 + 0]); 82c4762a1bSJed Brown C[1 * 3 + 2] = PetscRealPart(A[0 * 3 + 1] * A[2 * 3 + 0] - A[0 * 3 + 0] * A[2 * 3 + 1]); 83c4762a1bSJed Brown C[2 * 3 + 0] = PetscRealPart(A[0 * 3 + 1] * A[1 * 3 + 2] - A[0 * 3 + 2] * A[1 * 3 + 1]); 84c4762a1bSJed Brown C[2 * 3 + 1] = PetscRealPart(A[0 * 3 + 2] * A[1 * 3 + 0] - A[0 * 3 + 0] * A[1 * 3 + 2]); 85c4762a1bSJed Brown C[2 * 3 + 2] = PetscRealPart(A[0 * 3 + 0] * A[1 * 3 + 1] - A[0 * 3 + 1] * A[1 * 3 + 0]); 86c4762a1bSJed Brown } 87c4762a1bSJed Brown 889371c9d4SSatish Balay PetscErrorCode zero_scalar(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 89c4762a1bSJed Brown u[0] = 0.0; 90c4762a1bSJed Brown return 0; 91c4762a1bSJed Brown } 92c4762a1bSJed Brown 939371c9d4SSatish Balay PetscErrorCode zero_vector(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 94c4762a1bSJed Brown const PetscInt Ncomp = dim; 95c4762a1bSJed Brown 96c4762a1bSJed Brown PetscInt comp; 97c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) u[comp] = 0.0; 98c4762a1bSJed Brown return 0; 99c4762a1bSJed Brown } 100c4762a1bSJed Brown 1019371c9d4SSatish Balay PetscErrorCode coordinates(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 102c4762a1bSJed Brown const PetscInt Ncomp = dim; 103c4762a1bSJed Brown 104c4762a1bSJed Brown PetscInt comp; 105c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) u[comp] = x[comp]; 106c4762a1bSJed Brown return 0; 107c4762a1bSJed Brown } 108c4762a1bSJed Brown 1099371c9d4SSatish Balay PetscErrorCode elasticityMaterial(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 110c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 111c4762a1bSJed Brown u[0] = user->mu; 112c4762a1bSJed Brown return 0; 113c4762a1bSJed Brown } 114c4762a1bSJed Brown 1159371c9d4SSatish Balay PetscErrorCode wallPressure(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 116c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 117c4762a1bSJed Brown u[0] = user->p_wall; 118c4762a1bSJed Brown return 0; 119c4762a1bSJed Brown } 120c4762a1bSJed Brown 1219371c9d4SSatish Balay void f1_u_3d(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[]) { 122c4762a1bSJed Brown const PetscInt Ncomp = dim; 123c4762a1bSJed Brown const PetscReal mu = PetscRealPart(a[0]), kappa = 3.0; 124c4762a1bSJed Brown PetscReal cofu_x[9 /*Ncomp*dim*/], detu_x, p = PetscRealPart(u[Ncomp]); 125c4762a1bSJed Brown PetscInt comp, d; 126c4762a1bSJed Brown 127c4762a1bSJed Brown Cof3D(cofu_x, u_x); 128c4762a1bSJed Brown Det3D(&detu_x, u_x); 129c4762a1bSJed Brown p += kappa * (detu_x - 1.0); 130c4762a1bSJed Brown 131c4762a1bSJed Brown /* f1 is the first Piola-Kirchhoff tensor */ 132c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) { 1339371c9d4SSatish Balay for (d = 0; d < dim; ++d) { f1[comp * dim + d] = mu * u_x[comp * dim + d] + p * cofu_x[comp * dim + d]; } 134c4762a1bSJed Brown } 135c4762a1bSJed Brown } 136c4762a1bSJed Brown 1379371c9d4SSatish Balay void g3_uu_3d(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[]) { 138c4762a1bSJed Brown const PetscInt Ncomp = dim; 139c4762a1bSJed Brown const PetscReal mu = PetscRealPart(a[0]), kappa = 3.0; 140c4762a1bSJed Brown PetscReal cofu_x[9 /*Ncomp*dim*/], detu_x, pp, pm, p = PetscRealPart(u[Ncomp]); 141c4762a1bSJed Brown PetscInt compI, compJ, d1, d2; 142c4762a1bSJed Brown 143c4762a1bSJed Brown Cof3D(cofu_x, u_x); 144c4762a1bSJed Brown Det3D(&detu_x, u_x); 145c4762a1bSJed Brown p += kappa * (detu_x - 1.0); 146c4762a1bSJed Brown pp = p / detu_x + kappa; 147c4762a1bSJed Brown pm = p / detu_x; 148c4762a1bSJed Brown 149c4762a1bSJed 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} */ 150c4762a1bSJed Brown for (compI = 0; compI < Ncomp; ++compI) { 151c4762a1bSJed Brown for (compJ = 0; compJ < Ncomp; ++compJ) { 152c4762a1bSJed Brown const PetscReal G = (compI == compJ) ? 1.0 : 0.0; 153c4762a1bSJed Brown 154c4762a1bSJed Brown for (d1 = 0; d1 < dim; ++d1) { 155c4762a1bSJed Brown for (d2 = 0; d2 < dim; ++d2) { 156c4762a1bSJed Brown const PetscReal g = (d1 == d2) ? 1.0 : 0.0; 157c4762a1bSJed Brown 158c4762a1bSJed 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]; 159c4762a1bSJed Brown } 160c4762a1bSJed Brown } 161c4762a1bSJed Brown } 162c4762a1bSJed Brown } 163c4762a1bSJed Brown } 164c4762a1bSJed Brown 1659371c9d4SSatish Balay void f0_bd_u_3d(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[]) { 166c4762a1bSJed Brown const PetscInt Ncomp = dim; 167c4762a1bSJed Brown const PetscScalar p = a[aOff[1]]; 168c4762a1bSJed Brown PetscReal cofu_x[9 /*Ncomp*dim*/]; 169c4762a1bSJed Brown PetscInt comp, d; 170c4762a1bSJed Brown 171c4762a1bSJed Brown Cof3D(cofu_x, u_x); 172c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) { 173c4762a1bSJed Brown for (d = 0, f0[comp] = 0.0; d < dim; ++d) f0[comp] += cofu_x[comp * dim + d] * n[d]; 174c4762a1bSJed Brown f0[comp] *= p; 175c4762a1bSJed Brown } 176c4762a1bSJed Brown } 177c4762a1bSJed Brown 1789371c9d4SSatish Balay void g1_bd_uu_3d(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[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) { 179c4762a1bSJed Brown const PetscInt Ncomp = dim; 180c4762a1bSJed Brown PetscScalar p = a[aOff[1]]; 181c4762a1bSJed Brown PetscReal cofu_x[9 /*Ncomp*dim*/], m[3 /*Ncomp*/], detu_x; 182c4762a1bSJed Brown PetscInt comp, compI, compJ, d; 183c4762a1bSJed Brown 184c4762a1bSJed Brown Cof3D(cofu_x, u_x); 185c4762a1bSJed Brown Det3D(&detu_x, u_x); 186c4762a1bSJed Brown p /= detu_x; 187c4762a1bSJed Brown 1889371c9d4SSatish Balay for (comp = 0; comp < Ncomp; ++comp) 1899371c9d4SSatish Balay for (d = 0, m[comp] = 0.0; d < dim; ++d) m[comp] += cofu_x[comp * dim + d] * n[d]; 190c4762a1bSJed Brown for (compI = 0; compI < Ncomp; ++compI) { 191c4762a1bSJed Brown for (compJ = 0; compJ < Ncomp; ++compJ) { 1929371c9d4SSatish Balay for (d = 0; d < dim; ++d) { g1[(compI * Ncomp + compJ) * dim + d] = p * (m[compI] * cofu_x[compJ * dim + d] - cofu_x[compI * dim + d] * m[compJ]); } 193c4762a1bSJed Brown } 194c4762a1bSJed Brown } 195c4762a1bSJed Brown } 196c4762a1bSJed Brown 1979371c9d4SSatish Balay void f0_p_3d(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[]) { 198c4762a1bSJed Brown PetscReal detu_x; 199c4762a1bSJed Brown Det3D(&detu_x, u_x); 200c4762a1bSJed Brown f0[0] = detu_x - 1.0; 201c4762a1bSJed Brown } 202c4762a1bSJed Brown 2039371c9d4SSatish Balay void g1_pu_3d(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 g1[]) { 204c4762a1bSJed Brown PetscReal cofu_x[9 /*Ncomp*dim*/]; 205c4762a1bSJed Brown PetscInt compI, d; 206c4762a1bSJed Brown 207c4762a1bSJed Brown Cof3D(cofu_x, u_x); 208c4762a1bSJed Brown for (compI = 0; compI < dim; ++compI) 209c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[compI * dim + d] = cofu_x[compI * dim + d]; 210c4762a1bSJed Brown } 211c4762a1bSJed Brown 2129371c9d4SSatish Balay void g2_up_3d(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 g2[]) { 213c4762a1bSJed Brown PetscReal cofu_x[9 /*Ncomp*dim*/]; 214c4762a1bSJed Brown PetscInt compI, d; 215c4762a1bSJed Brown 216c4762a1bSJed Brown Cof3D(cofu_x, u_x); 217c4762a1bSJed Brown for (compI = 0; compI < dim; ++compI) 218c4762a1bSJed Brown for (d = 0; d < dim; ++d) g2[compI * dim + d] = cofu_x[compI * dim + d]; 219c4762a1bSJed Brown } 220c4762a1bSJed Brown 2219371c9d4SSatish Balay PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) { 222c4762a1bSJed Brown const char *runTypes[2] = {"full", "test"}; 223c4762a1bSJed Brown PetscInt run; 224c4762a1bSJed Brown 225c4762a1bSJed Brown PetscFunctionBeginUser; 226c4762a1bSJed Brown options->runType = RUN_FULL; 227c4762a1bSJed Brown options->mu = 1.0; 228c4762a1bSJed Brown options->p_wall = 0.4; 229d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Nonlinear elasticity problem options", "DMPLEX"); 230c4762a1bSJed Brown run = options->runType; 2319566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-run_type", "The run type", "ex77.c", runTypes, 2, runTypes[options->runType], &run, NULL)); 232c4762a1bSJed Brown options->runType = (RunType)run; 2339566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-shear_modulus", "The shear modulus", "ex77.c", options->mu, &options->mu, NULL)); 2349566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-wall_pressure", "The wall pressure", "ex77.c", options->p_wall, &options->p_wall, NULL)); 235d0609cedSBarry Smith PetscOptionsEnd(); 236c4762a1bSJed Brown PetscFunctionReturn(0); 237c4762a1bSJed Brown } 238c4762a1bSJed Brown 2399371c9d4SSatish Balay PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) { 240c4762a1bSJed Brown PetscFunctionBeginUser; 24130602db0SMatthew G. Knepley /* TODO The P1 coordinate space gives wrong results when compared to the affine version. Track this down */ 24230602db0SMatthew G. Knepley if (0) { 2439566063dSJacob Faibussowitsch PetscCall(DMPlexCreateBoxMesh(comm, 3, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm)); 24430602db0SMatthew G. Knepley } else { 2459566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 2469566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 24730602db0SMatthew G. Knepley } 2489566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 249c4762a1bSJed Brown /* Label the faces (bit of a hack here, until it is properly implemented for simplices) */ 2509566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view")); 251c4762a1bSJed Brown { 252c4762a1bSJed Brown DM cdm; 253c4762a1bSJed Brown DMLabel label; 254c4762a1bSJed Brown IS is; 25530602db0SMatthew G. Knepley PetscInt d, dim, b, f, Nf; 256c4762a1bSJed Brown const PetscInt *faces; 257c4762a1bSJed Brown PetscInt csize; 258c4762a1bSJed Brown PetscScalar *coords = NULL; 259c4762a1bSJed Brown PetscSection cs; 260c4762a1bSJed Brown Vec coordinates; 261c4762a1bSJed Brown 2629566063dSJacob Faibussowitsch PetscCall(DMGetDimension(*dm, &dim)); 2639566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(*dm, "boundary")); 2649566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dm, "boundary", &label)); 2659566063dSJacob Faibussowitsch PetscCall(DMPlexMarkBoundaryFaces(*dm, 1, label)); 2669566063dSJacob Faibussowitsch PetscCall(DMGetStratumIS(*dm, "boundary", 1, &is)); 267c4762a1bSJed Brown if (is) { 268c4762a1bSJed Brown PetscReal faceCoord; 269c4762a1bSJed Brown PetscInt v; 270c4762a1bSJed Brown 2719566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &Nf)); 2729566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &faces)); 273c4762a1bSJed Brown 2749566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(*dm, &coordinates)); 2759566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(*dm, &cdm)); 2769566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(cdm, &cs)); 277c4762a1bSJed Brown 278c4762a1bSJed Brown /* Check for each boundary face if any component of its centroid is either 0.0 or 1.0 */ 279c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 2809566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(cdm, cs, coordinates, faces[f], &csize, &coords)); 281c4762a1bSJed Brown /* Calculate mean coordinate vector */ 282c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 283c4762a1bSJed Brown const PetscInt Nv = csize / dim; 284c4762a1bSJed Brown faceCoord = 0.0; 285c4762a1bSJed Brown for (v = 0; v < Nv; ++v) faceCoord += PetscRealPart(coords[v * dim + d]); 286c4762a1bSJed Brown faceCoord /= Nv; 287c4762a1bSJed Brown for (b = 0; b < 2; ++b) { 288*48a46eb9SPierre Jolivet if (PetscAbs(faceCoord - b * 1.0) < PETSC_SMALL) PetscCall(DMSetLabelValue(*dm, "Faces", faces[f], d * 2 + b + 1)); 289c4762a1bSJed Brown } 290c4762a1bSJed Brown } 2919566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(cdm, cs, coordinates, faces[f], &csize, &coords)); 292c4762a1bSJed Brown } 2939566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &faces)); 294c4762a1bSJed Brown } 2959566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 296c4762a1bSJed Brown } 2979566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 298c4762a1bSJed Brown PetscFunctionReturn(0); 299c4762a1bSJed Brown } 300c4762a1bSJed Brown 3019371c9d4SSatish Balay PetscErrorCode SetupProblem(DM dm, PetscInt dim, AppCtx *user) { 30245480ffeSMatthew G. Knepley PetscDS ds; 30345480ffeSMatthew G. Knepley PetscWeakForm wf; 30445480ffeSMatthew G. Knepley DMLabel label; 30545480ffeSMatthew G. Knepley PetscInt bd; 306c4762a1bSJed Brown 307c4762a1bSJed Brown PetscFunctionBeginUser; 3089566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 3099566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, NULL, f1_u_3d)); 3109566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_p_3d, NULL)); 3119566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu_3d)); 3129566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up_3d, NULL)); 3139566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu_3d, NULL, NULL)); 314c4762a1bSJed Brown 3159566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Faces", &label)); 3169566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "pressure", label, 0, NULL, 0, 0, NULL, NULL, NULL, user, &bd)); 3179566063dSJacob Faibussowitsch PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 3189566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, 0, f0_bd_u_3d, 0, NULL)); 3199566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexBdJacobian(wf, label, 1, 0, 0, 0, 0, NULL, 0, g1_bd_uu_3d, 0, NULL, 0, NULL)); 320c4762a1bSJed Brown 3219566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "fixed", label, 0, NULL, 0, 0, NULL, (void (*)(void))coordinates, NULL, user, NULL)); 322c4762a1bSJed Brown PetscFunctionReturn(0); 323c4762a1bSJed Brown } 324c4762a1bSJed Brown 3259371c9d4SSatish Balay PetscErrorCode SetupMaterial(DM dm, DM dmAux, AppCtx *user) { 326c4762a1bSJed Brown PetscErrorCode (*matFuncs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx) = {elasticityMaterial, wallPressure}; 327c4762a1bSJed Brown Vec A; 328c4762a1bSJed Brown void *ctxs[2]; 329c4762a1bSJed Brown 330c4762a1bSJed Brown PetscFunctionBegin; 3319371c9d4SSatish Balay ctxs[0] = user; 3329371c9d4SSatish Balay ctxs[1] = user; 3339566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dmAux, &A)); 3349566063dSJacob Faibussowitsch PetscCall(DMProjectFunctionLocal(dmAux, 0.0, matFuncs, ctxs, INSERT_ALL_VALUES, A)); 3359566063dSJacob Faibussowitsch PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, A)); 3369566063dSJacob Faibussowitsch PetscCall(VecDestroy(&A)); 337c4762a1bSJed Brown PetscFunctionReturn(0); 338c4762a1bSJed Brown } 339c4762a1bSJed Brown 3409371c9d4SSatish Balay PetscErrorCode SetupNearNullSpace(DM dm, AppCtx *user) { 341c4762a1bSJed Brown /* Set up the near null space (a.k.a. rigid body modes) that will be used by the multigrid preconditioner */ 342c4762a1bSJed Brown DM subdm; 343c4762a1bSJed Brown MatNullSpace nearNullSpace; 344c4762a1bSJed Brown PetscInt fields = 0; 345c4762a1bSJed Brown PetscObject deformation; 346c4762a1bSJed Brown 347c4762a1bSJed Brown PetscFunctionBeginUser; 3489566063dSJacob Faibussowitsch PetscCall(DMCreateSubDM(dm, 1, &fields, NULL, &subdm)); 3499566063dSJacob Faibussowitsch PetscCall(DMPlexCreateRigidBody(subdm, 0, &nearNullSpace)); 3509566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, 0, NULL, &deformation)); 3519566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose(deformation, "nearnullspace", (PetscObject)nearNullSpace)); 3529566063dSJacob Faibussowitsch PetscCall(DMDestroy(&subdm)); 3539566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nearNullSpace)); 354c4762a1bSJed Brown PetscFunctionReturn(0); 355c4762a1bSJed Brown } 356c4762a1bSJed Brown 3579371c9d4SSatish Balay static PetscErrorCode SetupAuxDM(DM dm, PetscInt NfAux, PetscFE feAux[], AppCtx *user) { 358c4762a1bSJed Brown DM dmAux, coordDM; 359c4762a1bSJed Brown PetscInt f; 360c4762a1bSJed Brown 361c4762a1bSJed Brown PetscFunctionBegin; 362c4762a1bSJed Brown /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */ 3639566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &coordDM)); 364c4762a1bSJed Brown if (!feAux) PetscFunctionReturn(0); 3659566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmAux)); 3669566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateDM(dmAux, coordDM)); 3679566063dSJacob Faibussowitsch for (f = 0; f < NfAux; ++f) PetscCall(DMSetField(dmAux, f, NULL, (PetscObject)feAux[f])); 3689566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dmAux)); 3699566063dSJacob Faibussowitsch PetscCall(SetupMaterial(dm, dmAux, user)); 3709566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmAux)); 371c4762a1bSJed Brown PetscFunctionReturn(0); 372c4762a1bSJed Brown } 373c4762a1bSJed Brown 3749371c9d4SSatish Balay PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) { 375c4762a1bSJed Brown DM cdm = dm; 376c4762a1bSJed Brown PetscFE fe[2], feAux[2]; 37730602db0SMatthew G. Knepley PetscBool simplex; 37830602db0SMatthew G. Knepley PetscInt dim; 379c4762a1bSJed Brown MPI_Comm comm; 380c4762a1bSJed Brown 381c4762a1bSJed Brown PetscFunctionBeginUser; 3829566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 3839566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 3849566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 385c4762a1bSJed Brown /* Create finite element */ 3869566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "def_", PETSC_DEFAULT, &fe[0])); 3879566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[0], "deformation")); 3889566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1])); 3899566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 390c4762a1bSJed Brown 3919566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[1], "pressure")); 392c4762a1bSJed Brown 3939566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "elastMat_", PETSC_DEFAULT, &feAux[0])); 3949566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)feAux[0], "elasticityMaterial")); 3959566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], feAux[0])); 396c4762a1bSJed 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 */ 3979566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "wall_pres_", PETSC_DEFAULT, &feAux[1])); 3989566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)feAux[1], "wall_pressure")); 3999566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], feAux[1])); 400c4762a1bSJed Brown 401c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 4029566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe[0])); 4039566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe[1])); 4049566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 4059566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, dim, user)); 406c4762a1bSJed Brown while (cdm) { 4079566063dSJacob Faibussowitsch PetscCall(SetupAuxDM(cdm, 2, feAux, user)); 4089566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 4099566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 410c4762a1bSJed Brown } 4119566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[0])); 4129566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[1])); 4139566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&feAux[0])); 4149566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&feAux[1])); 415c4762a1bSJed Brown PetscFunctionReturn(0); 416c4762a1bSJed Brown } 417c4762a1bSJed Brown 4189371c9d4SSatish Balay int main(int argc, char **argv) { 419c4762a1bSJed Brown SNES snes; /* nonlinear solver */ 420c4762a1bSJed Brown DM dm; /* problem definition */ 421c4762a1bSJed Brown Vec u, r; /* solution, residual vectors */ 422c4762a1bSJed Brown Mat A, J; /* Jacobian matrix */ 423c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 424c4762a1bSJed Brown PetscInt its; /* iterations for convergence */ 425c4762a1bSJed Brown 426327415f7SBarry Smith PetscFunctionBeginUser; 4279566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 4289566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 4299566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 4309566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 4319566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 4329566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dm, &user)); 433c4762a1bSJed Brown 4349566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, &user)); 4359566063dSJacob Faibussowitsch PetscCall(DMPlexCreateClosureIndex(dm, NULL)); 4369566063dSJacob Faibussowitsch PetscCall(SetupNearNullSpace(dm, &user)); 437c4762a1bSJed Brown 4389566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u)); 4399566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r)); 440c4762a1bSJed Brown 4419566063dSJacob Faibussowitsch PetscCall(DMSetMatType(dm, MATAIJ)); 4429566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &J)); 443c4762a1bSJed Brown A = J; 444c4762a1bSJed Brown 4459566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 4469566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, A, J, NULL, NULL)); 447c4762a1bSJed Brown 4489566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 449c4762a1bSJed Brown 450c4762a1bSJed Brown { 451c4762a1bSJed Brown PetscErrorCode (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 4529371c9d4SSatish Balay initialGuess[0] = coordinates; 4539371c9d4SSatish Balay initialGuess[1] = zero_scalar; 4549566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u)); 455c4762a1bSJed Brown } 456c4762a1bSJed Brown 457c4762a1bSJed Brown if (user.runType == RUN_FULL) { 4589566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 4599566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 46063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its)); 461c4762a1bSJed Brown } else { 462c4762a1bSJed Brown PetscReal res = 0.0; 463c4762a1bSJed Brown 464c4762a1bSJed Brown /* Check initial guess */ 4659566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n")); 4669566063dSJacob Faibussowitsch PetscCall(VecView(u, PETSC_VIEWER_STDOUT_WORLD)); 467c4762a1bSJed Brown /* Check residual */ 4689566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, u, r)); 4699566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n")); 4709566063dSJacob Faibussowitsch PetscCall(VecChop(r, 1.0e-10)); 4719566063dSJacob Faibussowitsch PetscCall(VecView(r, PETSC_VIEWER_STDOUT_WORLD)); 4729566063dSJacob Faibussowitsch PetscCall(VecNorm(r, NORM_2, &res)); 4739566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", (double)res)); 474c4762a1bSJed Brown /* Check Jacobian */ 475c4762a1bSJed Brown { 476c4762a1bSJed Brown Vec b; 477c4762a1bSJed Brown 4789566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, u, A, A)); 4799566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &b)); 4809566063dSJacob Faibussowitsch PetscCall(VecSet(r, 0.0)); 4819566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, r, b)); 4829566063dSJacob Faibussowitsch PetscCall(MatMult(A, u, r)); 4839566063dSJacob Faibussowitsch PetscCall(VecAXPY(r, 1.0, b)); 4849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 4859566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Au - b = Au + F(0)\n")); 4869566063dSJacob Faibussowitsch PetscCall(VecChop(r, 1.0e-10)); 4879566063dSJacob Faibussowitsch PetscCall(VecView(r, PETSC_VIEWER_STDOUT_WORLD)); 4889566063dSJacob Faibussowitsch PetscCall(VecNorm(r, NORM_2, &res)); 4899566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", (double)res)); 490c4762a1bSJed Brown } 491c4762a1bSJed Brown } 4929566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view")); 493c4762a1bSJed Brown 4949566063dSJacob Faibussowitsch if (A != J) PetscCall(MatDestroy(&A)); 4959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 4969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 4979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 4989566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 4999566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 5009566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 501b122ec5aSJacob Faibussowitsch return 0; 502c4762a1bSJed Brown } 503c4762a1bSJed Brown 504c4762a1bSJed Brown /*TEST 505c4762a1bSJed Brown 506c4762a1bSJed Brown build: 507c4762a1bSJed Brown requires: !complex 508c4762a1bSJed Brown 50930602db0SMatthew G. Knepley testset: 51030602db0SMatthew G. Knepley requires: ctetgen 51130602db0SMatthew G. Knepley args: -run_type full -dm_plex_dim 3 \ 51230602db0SMatthew G. Knepley -def_petscspace_degree 2 -pres_petscspace_degree 1 -elastMat_petscspace_degree 0 -wall_pres_petscspace_degree 0 \ 51330602db0SMatthew G. Knepley -snes_rtol 1e-05 -snes_monitor_short -snes_converged_reason \ 51430602db0SMatthew G. Knepley -ksp_type fgmres -ksp_rtol 1e-10 -ksp_monitor_short -ksp_converged_reason \ 51530602db0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper \ 51630602db0SMatthew G. Knepley -fieldsplit_deformation_ksp_type preonly -fieldsplit_deformation_pc_type lu \ 51730602db0SMatthew G. Knepley -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 51830602db0SMatthew G. Knepley 519c4762a1bSJed Brown test: 520c4762a1bSJed Brown suffix: 0 52130602db0SMatthew G. Knepley requires: !single 52230602db0SMatthew G. Knepley args: -dm_refine 2 \ 52330602db0SMatthew G. Knepley -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4 524c4762a1bSJed Brown 525c4762a1bSJed Brown test: 526c4762a1bSJed Brown suffix: 1 52730602db0SMatthew G. Knepley requires: superlu_dist 528c4762a1bSJed Brown nsize: 2 529e600fa54SMatthew G. Knepley args: -dm_refine 0 -petscpartitioner_type simple \ 53030602db0SMatthew G. Knepley -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4 531c4762a1bSJed Brown timeoutfactor: 2 532c4762a1bSJed Brown 533c4762a1bSJed Brown test: 534c4762a1bSJed Brown suffix: 4 53530602db0SMatthew G. Knepley requires: superlu_dist 536c4762a1bSJed Brown nsize: 2 537e600fa54SMatthew G. Knepley args: -dm_refine 0 -petscpartitioner_type simple \ 53830602db0SMatthew G. Knepley -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4 539c4762a1bSJed Brown output_file: output/ex77_1.out 540c4762a1bSJed Brown 541c4762a1bSJed Brown test: 54230602db0SMatthew G. Knepley suffix: 2 54330602db0SMatthew G. Knepley requires: !single 54430602db0SMatthew G. Knepley args: -dm_refine 2 \ 54530602db0SMatthew G. Knepley -bc_fixed 3,4,5,6 -bc_pressure 2 -wall_pressure 1.0 546c4762a1bSJed Brown 547c4762a1bSJed Brown test: 548c4762a1bSJed Brown suffix: 2_par 54930602db0SMatthew G. Knepley requires: superlu_dist !single 550c4762a1bSJed Brown nsize: 4 551e600fa54SMatthew G. Knepley args: -dm_refine 2 -petscpartitioner_type simple \ 55230602db0SMatthew G. Knepley -bc_fixed 3,4,5,6 -bc_pressure 2 -wall_pressure 1.0 553c4762a1bSJed Brown output_file: output/ex77_2.out 554c4762a1bSJed Brown 555c4762a1bSJed Brown TEST*/ 556