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 PetscInt debug; /* The debugging level */ 48c4762a1bSJed Brown RunType runType; /* Whether to run tests, or solve the full problem */ 49c4762a1bSJed Brown PetscLogEvent createMeshEvent; 50c4762a1bSJed Brown PetscBool showInitial, showSolution; 51c4762a1bSJed Brown /* Domain and mesh definition */ 52c4762a1bSJed Brown PetscInt dim; /* The topological mesh dimension */ 53c4762a1bSJed Brown PetscBool interpolate; /* Generate intermediate mesh elements */ 54c4762a1bSJed Brown PetscBool simplex; /* Use simplices or tensor product cells */ 55c4762a1bSJed Brown PetscReal refinementLimit; /* The largest allowable cell volume */ 56c4762a1bSJed Brown PetscBool testPartition; /* Use a fixed partitioning for testing */ 57c4762a1bSJed Brown PetscReal mu; /* The shear modulus */ 58c4762a1bSJed Brown PetscReal p_wall; /* The wall pressure */ 59c4762a1bSJed Brown } AppCtx; 60c4762a1bSJed Brown 61c4762a1bSJed Brown #if 0 62c4762a1bSJed Brown PETSC_STATIC_INLINE void Det2D(PetscReal *detJ, const PetscReal J[]) 63c4762a1bSJed Brown { 64c4762a1bSJed Brown *detJ = J[0]*J[3] - J[1]*J[2]; 65c4762a1bSJed Brown } 66c4762a1bSJed Brown #endif 67c4762a1bSJed Brown 68c4762a1bSJed Brown PETSC_STATIC_INLINE void Det3D(PetscReal *detJ, const PetscScalar J[]) 69c4762a1bSJed Brown { 70c4762a1bSJed Brown *detJ = PetscRealPart(J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) + 71c4762a1bSJed Brown J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) + 72c4762a1bSJed Brown J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0])); 73c4762a1bSJed Brown } 74c4762a1bSJed Brown 75c4762a1bSJed Brown #if 0 76c4762a1bSJed Brown PETSC_STATIC_INLINE void Cof2D(PetscReal C[], const PetscReal A[]) 77c4762a1bSJed Brown { 78c4762a1bSJed Brown C[0] = A[3]; 79c4762a1bSJed Brown C[1] = -A[2]; 80c4762a1bSJed Brown C[2] = -A[1]; 81c4762a1bSJed Brown C[3] = A[0]; 82c4762a1bSJed Brown } 83c4762a1bSJed Brown #endif 84c4762a1bSJed Brown 85c4762a1bSJed Brown PETSC_STATIC_INLINE void Cof3D(PetscReal C[], const PetscScalar A[]) 86c4762a1bSJed Brown { 87c4762a1bSJed Brown C[0*3+0] = PetscRealPart(A[1*3+1]*A[2*3+2] - A[1*3+2]*A[2*3+1]); 88c4762a1bSJed Brown C[0*3+1] = PetscRealPart(A[1*3+2]*A[2*3+0] - A[1*3+0]*A[2*3+2]); 89c4762a1bSJed Brown C[0*3+2] = PetscRealPart(A[1*3+0]*A[2*3+1] - A[1*3+1]*A[2*3+0]); 90c4762a1bSJed Brown C[1*3+0] = PetscRealPart(A[0*3+2]*A[2*3+1] - A[0*3+1]*A[2*3+2]); 91c4762a1bSJed Brown C[1*3+1] = PetscRealPart(A[0*3+0]*A[2*3+2] - A[0*3+2]*A[2*3+0]); 92c4762a1bSJed Brown C[1*3+2] = PetscRealPart(A[0*3+1]*A[2*3+0] - A[0*3+0]*A[2*3+1]); 93c4762a1bSJed Brown C[2*3+0] = PetscRealPart(A[0*3+1]*A[1*3+2] - A[0*3+2]*A[1*3+1]); 94c4762a1bSJed Brown C[2*3+1] = PetscRealPart(A[0*3+2]*A[1*3+0] - A[0*3+0]*A[1*3+2]); 95c4762a1bSJed Brown C[2*3+2] = PetscRealPart(A[0*3+0]*A[1*3+1] - A[0*3+1]*A[1*3+0]); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown 98c4762a1bSJed Brown PetscErrorCode zero_scalar(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 99c4762a1bSJed Brown { 100c4762a1bSJed Brown u[0] = 0.0; 101c4762a1bSJed Brown return 0; 102c4762a1bSJed Brown } 103c4762a1bSJed Brown 104c4762a1bSJed Brown PetscErrorCode zero_vector(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] = 0.0; 110c4762a1bSJed Brown return 0; 111c4762a1bSJed Brown } 112c4762a1bSJed Brown 113c4762a1bSJed Brown PetscErrorCode coordinates(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 114c4762a1bSJed Brown { 115c4762a1bSJed Brown const PetscInt Ncomp = dim; 116c4762a1bSJed Brown 117c4762a1bSJed Brown PetscInt comp; 118c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) u[comp] = x[comp]; 119c4762a1bSJed Brown return 0; 120c4762a1bSJed Brown } 121c4762a1bSJed Brown 122c4762a1bSJed Brown PetscErrorCode elasticityMaterial(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 123c4762a1bSJed Brown { 124c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 125c4762a1bSJed Brown u[0] = user->mu; 126c4762a1bSJed Brown return 0; 127c4762a1bSJed Brown } 128c4762a1bSJed Brown 129c4762a1bSJed Brown PetscErrorCode wallPressure(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 130c4762a1bSJed Brown { 131c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 132c4762a1bSJed Brown u[0] = user->p_wall; 133c4762a1bSJed Brown return 0; 134c4762a1bSJed Brown } 135c4762a1bSJed Brown 136c4762a1bSJed Brown void f1_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 137c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 138c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 139c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 140c4762a1bSJed Brown { 141c4762a1bSJed Brown const PetscInt Ncomp = dim; 142c4762a1bSJed Brown const PetscReal mu = PetscRealPart(a[0]), kappa = 3.0; 143c4762a1bSJed Brown PetscReal cofu_x[9/*Ncomp*dim*/], detu_x, p = PetscRealPart(u[Ncomp]); 144c4762a1bSJed Brown PetscInt comp, d; 145c4762a1bSJed Brown 146c4762a1bSJed Brown Cof3D(cofu_x, u_x); 147c4762a1bSJed Brown Det3D(&detu_x, u_x); 148c4762a1bSJed Brown p += kappa * (detu_x - 1.0); 149c4762a1bSJed Brown 150c4762a1bSJed Brown /* f1 is the first Piola-Kirchhoff tensor */ 151c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) { 152c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 153c4762a1bSJed Brown f1[comp*dim+d] = mu * u_x[comp*dim+d] + p * cofu_x[comp*dim+d]; 154c4762a1bSJed Brown } 155c4762a1bSJed Brown } 156c4762a1bSJed Brown } 157c4762a1bSJed Brown 158c4762a1bSJed Brown void g3_uu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 159c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 160c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 161c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 162c4762a1bSJed Brown { 163c4762a1bSJed Brown const PetscInt Ncomp = dim; 164c4762a1bSJed Brown const PetscReal mu = PetscRealPart(a[0]), kappa = 3.0; 165c4762a1bSJed Brown PetscReal cofu_x[9/*Ncomp*dim*/], detu_x, pp, pm, p = PetscRealPart(u[Ncomp]); 166c4762a1bSJed Brown PetscInt compI, compJ, d1, d2; 167c4762a1bSJed Brown 168c4762a1bSJed Brown Cof3D(cofu_x, u_x); 169c4762a1bSJed Brown Det3D(&detu_x, u_x); 170c4762a1bSJed Brown p += kappa * (detu_x - 1.0); 171c4762a1bSJed Brown pp = p/detu_x + kappa; 172c4762a1bSJed Brown pm = p/detu_x; 173c4762a1bSJed Brown 174c4762a1bSJed 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} */ 175c4762a1bSJed Brown for (compI = 0; compI < Ncomp; ++compI) { 176c4762a1bSJed Brown for (compJ = 0; compJ < Ncomp; ++compJ) { 177c4762a1bSJed Brown const PetscReal G = (compI == compJ) ? 1.0 : 0.0; 178c4762a1bSJed Brown 179c4762a1bSJed Brown for (d1 = 0; d1 < dim; ++d1) { 180c4762a1bSJed Brown for (d2 = 0; d2 < dim; ++d2) { 181c4762a1bSJed Brown const PetscReal g = (d1 == d2) ? 1.0 : 0.0; 182c4762a1bSJed Brown 183c4762a1bSJed 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]; 184c4762a1bSJed Brown } 185c4762a1bSJed Brown } 186c4762a1bSJed Brown } 187c4762a1bSJed Brown } 188c4762a1bSJed Brown } 189c4762a1bSJed Brown 190c4762a1bSJed Brown void f0_bd_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 191c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 192c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 193c4762a1bSJed Brown PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 194c4762a1bSJed Brown { 195c4762a1bSJed Brown const PetscInt Ncomp = dim; 196c4762a1bSJed Brown const PetscScalar p = a[aOff[1]]; 197c4762a1bSJed Brown PetscReal cofu_x[9/*Ncomp*dim*/]; 198c4762a1bSJed Brown PetscInt comp, d; 199c4762a1bSJed Brown 200c4762a1bSJed Brown Cof3D(cofu_x, u_x); 201c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) { 202c4762a1bSJed Brown for (d = 0, f0[comp] = 0.0; d < dim; ++d) f0[comp] += cofu_x[comp*dim+d] * n[d]; 203c4762a1bSJed Brown f0[comp] *= p; 204c4762a1bSJed Brown } 205c4762a1bSJed Brown } 206c4762a1bSJed Brown 207c4762a1bSJed Brown void g1_bd_uu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 208c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 209c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 210c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 211c4762a1bSJed Brown { 212c4762a1bSJed Brown const PetscInt Ncomp = dim; 213c4762a1bSJed Brown PetscScalar p = a[aOff[1]]; 214c4762a1bSJed Brown PetscReal cofu_x[9/*Ncomp*dim*/], m[3/*Ncomp*/], detu_x; 215c4762a1bSJed Brown PetscInt comp, compI, compJ, d; 216c4762a1bSJed Brown 217c4762a1bSJed Brown Cof3D(cofu_x, u_x); 218c4762a1bSJed Brown Det3D(&detu_x, u_x); 219c4762a1bSJed Brown p /= detu_x; 220c4762a1bSJed Brown 221c4762a1bSJed 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]; 222c4762a1bSJed Brown for (compI = 0; compI < Ncomp; ++compI) { 223c4762a1bSJed Brown for (compJ = 0; compJ < Ncomp; ++compJ) { 224c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 225c4762a1bSJed Brown g1[(compI*Ncomp+compJ)*dim+d] = p * (m[compI] * cofu_x[compJ*dim+d] - cofu_x[compI*dim+d] * m[compJ]); 226c4762a1bSJed Brown } 227c4762a1bSJed Brown } 228c4762a1bSJed Brown } 229c4762a1bSJed Brown } 230c4762a1bSJed Brown 231c4762a1bSJed Brown void f0_p_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 232c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 233c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 234c4762a1bSJed Brown PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 235c4762a1bSJed Brown { 236c4762a1bSJed Brown PetscReal detu_x; 237c4762a1bSJed Brown Det3D(&detu_x, u_x); 238c4762a1bSJed Brown f0[0] = detu_x - 1.0; 239c4762a1bSJed Brown } 240c4762a1bSJed Brown 241c4762a1bSJed Brown void g1_pu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 242c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 243c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 244c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 245c4762a1bSJed Brown { 246c4762a1bSJed Brown PetscReal cofu_x[9/*Ncomp*dim*/]; 247c4762a1bSJed Brown PetscInt compI, d; 248c4762a1bSJed Brown 249c4762a1bSJed Brown Cof3D(cofu_x, u_x); 250c4762a1bSJed Brown for (compI = 0; compI < dim; ++compI) 251c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[compI*dim+d] = cofu_x[compI*dim+d]; 252c4762a1bSJed Brown } 253c4762a1bSJed Brown 254c4762a1bSJed Brown void g2_up_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 255c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 256c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 257c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 258c4762a1bSJed Brown { 259c4762a1bSJed Brown PetscReal cofu_x[9/*Ncomp*dim*/]; 260c4762a1bSJed Brown PetscInt compI, d; 261c4762a1bSJed Brown 262c4762a1bSJed Brown Cof3D(cofu_x, u_x); 263c4762a1bSJed Brown for (compI = 0; compI < dim; ++compI) 264c4762a1bSJed Brown for (d = 0; d < dim; ++d) g2[compI*dim+d] = cofu_x[compI*dim+d]; 265c4762a1bSJed Brown } 266c4762a1bSJed Brown 267c4762a1bSJed Brown PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 268c4762a1bSJed Brown { 269c4762a1bSJed Brown const char *runTypes[2] = {"full", "test"}; 270c4762a1bSJed Brown PetscInt run; 271c4762a1bSJed Brown PetscErrorCode ierr; 272c4762a1bSJed Brown 273c4762a1bSJed Brown PetscFunctionBeginUser; 274c4762a1bSJed Brown options->debug = 0; 275c4762a1bSJed Brown options->runType = RUN_FULL; 276c4762a1bSJed Brown options->dim = 3; 277c4762a1bSJed Brown options->interpolate = PETSC_FALSE; 278c4762a1bSJed Brown options->simplex = PETSC_TRUE; 279c4762a1bSJed Brown options->refinementLimit = 0.0; 280c4762a1bSJed Brown options->mu = 1.0; 281c4762a1bSJed Brown options->p_wall = 0.4; 282c4762a1bSJed Brown options->testPartition = PETSC_FALSE; 283c4762a1bSJed Brown options->showInitial = PETSC_FALSE; 284c4762a1bSJed Brown options->showSolution = PETSC_TRUE; 285c4762a1bSJed Brown 286c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Nonlinear elasticity problem options", "DMPLEX");CHKERRQ(ierr); 287c4762a1bSJed Brown ierr = PetscOptionsInt("-debug", "The debugging level", "ex77.c", options->debug, &options->debug, NULL);CHKERRQ(ierr); 288c4762a1bSJed Brown run = options->runType; 289c4762a1bSJed Brown ierr = PetscOptionsEList("-run_type", "The run type", "ex77.c", runTypes, 2, runTypes[options->runType], &run, NULL);CHKERRQ(ierr); 290c4762a1bSJed Brown 291c4762a1bSJed Brown options->runType = (RunType) run; 292c4762a1bSJed Brown 293c4762a1bSJed Brown ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex77.c", options->dim, &options->dim, NULL);CHKERRQ(ierr); 294c4762a1bSJed Brown ierr = PetscOptionsBool("-interpolate", "Generate intermediate mesh elements", "ex77.c", options->interpolate, &options->interpolate, NULL);CHKERRQ(ierr); 295c4762a1bSJed Brown ierr = PetscOptionsBool("-simplex", "Use simplices or tensor product cells", "ex77.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr); 296c4762a1bSJed Brown ierr = PetscOptionsReal("-refinement_limit", "The largest allowable cell volume", "ex77.c", options->refinementLimit, &options->refinementLimit, NULL);CHKERRQ(ierr); 297c4762a1bSJed Brown ierr = PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex77.c", options->testPartition, &options->testPartition, NULL);CHKERRQ(ierr); 298c4762a1bSJed Brown ierr = PetscOptionsReal("-shear_modulus", "The shear modulus", "ex77.c", options->mu, &options->mu, NULL);CHKERRQ(ierr); 299c4762a1bSJed Brown ierr = PetscOptionsReal("-wall_pressure", "The wall pressure", "ex77.c", options->p_wall, &options->p_wall, NULL);CHKERRQ(ierr); 300c4762a1bSJed Brown 301c4762a1bSJed Brown ierr = PetscOptionsBool("-show_initial", "Output the initial guess for verification", "ex77.c", options->showInitial, &options->showInitial, NULL);CHKERRQ(ierr); 302c4762a1bSJed Brown ierr = PetscOptionsBool("-show_solution", "Output the solution for verification", "ex77.c", options->showSolution, &options->showSolution, NULL);CHKERRQ(ierr); 303c4762a1bSJed Brown ierr = PetscOptionsEnd(); 304c4762a1bSJed Brown 305c4762a1bSJed Brown ierr = PetscLogEventRegister("CreateMesh", DM_CLASSID, &options->createMeshEvent);CHKERRQ(ierr); 306c4762a1bSJed Brown PetscFunctionReturn(0); 307c4762a1bSJed Brown } 308c4762a1bSJed Brown 309c4762a1bSJed Brown PetscErrorCode DMVecViewLocal(DM dm, Vec v, PetscViewer viewer) 310c4762a1bSJed Brown { 311c4762a1bSJed Brown Vec lv; 312c4762a1bSJed Brown PetscInt p; 313c4762a1bSJed Brown PetscMPIInt rank, size; 314c4762a1bSJed Brown PetscErrorCode ierr; 315c4762a1bSJed Brown 316c4762a1bSJed Brown PetscFunctionBeginUser; 317c4762a1bSJed Brown ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 318c4762a1bSJed Brown ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size);CHKERRQ(ierr); 319c4762a1bSJed Brown ierr = DMGetLocalVector(dm, &lv);CHKERRQ(ierr); 320c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(dm, v, INSERT_VALUES, lv);CHKERRQ(ierr); 321c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(dm, v, INSERT_VALUES, lv);CHKERRQ(ierr); 322c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Local function\n");CHKERRQ(ierr); 323c4762a1bSJed Brown for (p = 0; p < size; ++p) { 324c4762a1bSJed Brown if (p == rank) {ierr = VecView(lv, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);} 325c4762a1bSJed Brown ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr); 326c4762a1bSJed Brown } 327c4762a1bSJed Brown ierr = DMRestoreLocalVector(dm, &lv);CHKERRQ(ierr); 328c4762a1bSJed Brown PetscFunctionReturn(0); 329c4762a1bSJed Brown } 330c4762a1bSJed Brown 331c4762a1bSJed Brown PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 332c4762a1bSJed Brown { 333c4762a1bSJed Brown PetscInt dim = user->dim; 334c4762a1bSJed Brown PetscBool interpolate = user->interpolate; 335c4762a1bSJed Brown PetscReal refinementLimit = user->refinementLimit; 336c4762a1bSJed Brown const PetscInt cells[3] = {3, 3, 3}; 337c4762a1bSJed Brown PetscErrorCode ierr; 338c4762a1bSJed Brown 339c4762a1bSJed Brown PetscFunctionBeginUser; 340c4762a1bSJed Brown ierr = PetscLogEventBegin(user->createMeshEvent,0,0,0,0);CHKERRQ(ierr); 341c4762a1bSJed Brown ierr = DMPlexCreateBoxMesh(comm, dim, user->simplex, user->simplex ? NULL : cells, NULL, NULL, NULL, interpolate, dm);CHKERRQ(ierr); 342c4762a1bSJed Brown /* Label the faces (bit of a hack here, until it is properly implemented for simplices) */ 343c4762a1bSJed Brown { 344c4762a1bSJed Brown DM cdm; 345c4762a1bSJed Brown DMLabel label; 346c4762a1bSJed Brown IS is; 347c4762a1bSJed Brown PetscInt d, dim = user->dim, b, f, Nf; 348c4762a1bSJed Brown const PetscInt *faces; 349c4762a1bSJed Brown PetscInt csize; 350c4762a1bSJed Brown PetscScalar *coords = NULL; 351c4762a1bSJed Brown PetscSection cs; 352c4762a1bSJed Brown Vec coordinates ; 353c4762a1bSJed Brown 354c4762a1bSJed Brown ierr = DMCreateLabel(*dm, "boundary");CHKERRQ(ierr); 355c4762a1bSJed Brown ierr = DMGetLabel(*dm, "boundary", &label);CHKERRQ(ierr); 356c4762a1bSJed Brown ierr = DMPlexMarkBoundaryFaces(*dm, 1, label);CHKERRQ(ierr); 357c4762a1bSJed Brown ierr = DMGetStratumIS(*dm, "boundary", 1, &is);CHKERRQ(ierr); 358c4762a1bSJed Brown if (is) { 359c4762a1bSJed Brown PetscReal faceCoord; 360c4762a1bSJed Brown PetscInt v; 361c4762a1bSJed Brown 362c4762a1bSJed Brown ierr = ISGetLocalSize(is, &Nf);CHKERRQ(ierr); 363c4762a1bSJed Brown ierr = ISGetIndices(is, &faces);CHKERRQ(ierr); 364c4762a1bSJed Brown 365c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(*dm, &coordinates);CHKERRQ(ierr); 366c4762a1bSJed Brown ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr); 367c4762a1bSJed Brown ierr = DMGetLocalSection(cdm, &cs);CHKERRQ(ierr); 368c4762a1bSJed Brown 369c4762a1bSJed Brown /* Check for each boundary face if any component of its centroid is either 0.0 or 1.0 */ 370c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 371c4762a1bSJed Brown ierr = DMPlexVecGetClosure(cdm, cs, coordinates, faces[f], &csize, &coords);CHKERRQ(ierr); 372c4762a1bSJed Brown /* Calculate mean coordinate vector */ 373c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 374c4762a1bSJed Brown const PetscInt Nv = csize/dim; 375c4762a1bSJed Brown faceCoord = 0.0; 376c4762a1bSJed Brown for (v = 0; v < Nv; ++v) faceCoord += PetscRealPart(coords[v*dim+d]); 377c4762a1bSJed Brown faceCoord /= Nv; 378c4762a1bSJed Brown for (b = 0; b < 2; ++b) { 379c4762a1bSJed Brown if (PetscAbs(faceCoord - b*1.0) < PETSC_SMALL) { 380c4762a1bSJed Brown ierr = DMSetLabelValue(*dm, "Faces", faces[f], d*2+b+1);CHKERRQ(ierr); 381c4762a1bSJed Brown } 382c4762a1bSJed Brown } 383c4762a1bSJed Brown } 384c4762a1bSJed Brown ierr = DMPlexVecRestoreClosure(cdm, cs, coordinates, faces[f], &csize, &coords);CHKERRQ(ierr); 385c4762a1bSJed Brown } 386c4762a1bSJed Brown ierr = ISRestoreIndices(is, &faces);CHKERRQ(ierr); 387c4762a1bSJed Brown } 388c4762a1bSJed Brown ierr = ISDestroy(&is);CHKERRQ(ierr); 389c4762a1bSJed Brown } 390c4762a1bSJed Brown { 391c4762a1bSJed Brown DM refinedMesh = NULL; 392c4762a1bSJed Brown DM distributedMesh = NULL; 393c4762a1bSJed Brown PetscPartitioner part; 394c4762a1bSJed Brown 395c4762a1bSJed Brown ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr); 396c4762a1bSJed Brown 397c4762a1bSJed Brown /* Refine mesh using a volume constraint */ 398c4762a1bSJed Brown ierr = DMPlexSetRefinementLimit(*dm, refinementLimit);CHKERRQ(ierr); 399c4762a1bSJed Brown if (user->simplex) {ierr = DMRefine(*dm, comm, &refinedMesh);CHKERRQ(ierr);} 400c4762a1bSJed Brown if (refinedMesh) { 401c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 402c4762a1bSJed Brown *dm = refinedMesh; 403c4762a1bSJed Brown } 404c4762a1bSJed Brown /* Setup test partitioning */ 405c4762a1bSJed Brown if (user->testPartition) { 406c4762a1bSJed Brown PetscInt triSizes_n2[2] = {4, 4}; 407c4762a1bSJed Brown PetscInt triPoints_n2[8] = {3, 5, 6, 7, 0, 1, 2, 4}; 408c4762a1bSJed Brown PetscInt triSizes_n3[3] = {2, 3, 3}; 409c4762a1bSJed Brown PetscInt triPoints_n3[8] = {3, 5, 1, 6, 7, 0, 2, 4}; 410c4762a1bSJed Brown PetscInt triSizes_n5[5] = {1, 2, 2, 1, 2}; 411c4762a1bSJed Brown PetscInt triPoints_n5[8] = {3, 5, 6, 4, 7, 0, 1, 2}; 412c4762a1bSJed Brown PetscInt triSizes_ref_n2[2] = {8, 8}; 413c4762a1bSJed Brown PetscInt triPoints_ref_n2[16] = {1, 5, 6, 7, 10, 11, 14, 15, 0, 2, 3, 4, 8, 9, 12, 13}; 414c4762a1bSJed Brown PetscInt triSizes_ref_n3[3] = {5, 6, 5}; 415c4762a1bSJed Brown PetscInt triPoints_ref_n3[16] = {1, 7, 10, 14, 15, 2, 6, 8, 11, 12, 13, 0, 3, 4, 5, 9}; 416c4762a1bSJed Brown PetscInt triSizes_ref_n5[5] = {3, 4, 3, 3, 3}; 417c4762a1bSJed Brown PetscInt triPoints_ref_n5[16] = {1, 7, 10, 2, 11, 13, 14, 5, 6, 15, 0, 8, 9, 3, 4, 12}; 418c4762a1bSJed Brown PetscInt tetSizes_n2[2] = {3, 3}; 419c4762a1bSJed Brown PetscInt tetPoints_n2[6] = {1, 2, 3, 0, 4, 5}; 420c4762a1bSJed Brown const PetscInt *sizes = NULL; 421c4762a1bSJed Brown const PetscInt *points = NULL; 422c4762a1bSJed Brown PetscInt cEnd; 423c4762a1bSJed Brown PetscMPIInt rank, size; 424c4762a1bSJed Brown 425c4762a1bSJed Brown ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 426c4762a1bSJed Brown ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 427c4762a1bSJed Brown ierr = DMPlexGetHeightStratum(*dm, 0, NULL, &cEnd);CHKERRQ(ierr); 428c4762a1bSJed Brown if (!rank) { 429c4762a1bSJed Brown if (dim == 2 && user->simplex && size == 2 && cEnd == 8) { 430c4762a1bSJed Brown sizes = triSizes_n2; points = triPoints_n2; 431c4762a1bSJed Brown } else if (dim == 2 && user->simplex && size == 3 && cEnd == 8) { 432c4762a1bSJed Brown sizes = triSizes_n3; points = triPoints_n3; 433c4762a1bSJed Brown } else if (dim == 2 && user->simplex && size == 5 && cEnd == 8) { 434c4762a1bSJed Brown sizes = triSizes_n5; points = triPoints_n5; 435c4762a1bSJed Brown } else if (dim == 2 && user->simplex && size == 2 && cEnd == 16) { 436c4762a1bSJed Brown sizes = triSizes_ref_n2; points = triPoints_ref_n2; 437c4762a1bSJed Brown } else if (dim == 2 && user->simplex && size == 3 && cEnd == 16) { 438c4762a1bSJed Brown sizes = triSizes_ref_n3; points = triPoints_ref_n3; 439c4762a1bSJed Brown } else if (dim == 2 && user->simplex && size == 5 && cEnd == 16) { 440c4762a1bSJed Brown sizes = triSizes_ref_n5; points = triPoints_ref_n5; 441c4762a1bSJed Brown } else if (dim == 3 && user->simplex && size == 2 && cEnd == 6) { 442c4762a1bSJed Brown sizes = tetSizes_n2; points = tetPoints_n2; 443c4762a1bSJed Brown } else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "No stored partition matching run parameters"); 444c4762a1bSJed Brown } 445c4762a1bSJed Brown ierr = PetscPartitionerSetType(part, PETSCPARTITIONERSHELL);CHKERRQ(ierr); 446c4762a1bSJed Brown ierr = PetscPartitionerShellSetPartition(part, size, sizes, points);CHKERRQ(ierr); 447c4762a1bSJed Brown } else { 448c4762a1bSJed Brown ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 449c4762a1bSJed Brown } 450c4762a1bSJed Brown /* Distribute mesh over processes */ 451c4762a1bSJed Brown ierr = DMPlexDistribute(*dm, 0, NULL, &distributedMesh);CHKERRQ(ierr); 452c4762a1bSJed Brown if (distributedMesh) { 453c4762a1bSJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 454c4762a1bSJed Brown *dm = distributedMesh; 455c4762a1bSJed Brown } 456c4762a1bSJed Brown } 457c4762a1bSJed Brown ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 458c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 459c4762a1bSJed Brown 460c4762a1bSJed Brown ierr = PetscLogEventEnd(user->createMeshEvent,0,0,0,0);CHKERRQ(ierr); 461c4762a1bSJed Brown PetscFunctionReturn(0); 462c4762a1bSJed Brown } 463c4762a1bSJed Brown 464c4762a1bSJed Brown PetscErrorCode SetupProblem(DM dm, PetscInt dim, AppCtx *user) 465c4762a1bSJed Brown { 466c4762a1bSJed Brown PetscDS prob; 467c4762a1bSJed Brown PetscErrorCode ierr; 468c4762a1bSJed Brown 469c4762a1bSJed Brown PetscFunctionBeginUser; 470c4762a1bSJed Brown ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 471c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 0, NULL, f1_u_3d);CHKERRQ(ierr); 472c4762a1bSJed Brown ierr = PetscDSSetResidual(prob, 1, f0_p_3d, NULL);CHKERRQ(ierr); 473c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu_3d);CHKERRQ(ierr); 474c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_up_3d, NULL);CHKERRQ(ierr); 475c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 1, 0, NULL, g1_pu_3d, NULL, NULL);CHKERRQ(ierr); 476c4762a1bSJed Brown 477c4762a1bSJed Brown ierr = PetscDSSetBdResidual(prob, 0, f0_bd_u_3d, NULL);CHKERRQ(ierr); 478c4762a1bSJed Brown ierr = PetscDSSetBdJacobian(prob, 0, 0, NULL, g1_bd_uu_3d, NULL, NULL);CHKERRQ(ierr); 479c4762a1bSJed Brown 480*408cafa0SMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "fixed", "Faces", 0, 0, NULL, (void (*)(void)) coordinates, 0, NULL, user);CHKERRQ(ierr); 481*408cafa0SMatthew G. Knepley ierr = DMAddBoundary(dm, DM_BC_NATURAL, "pressure", "Faces", 0, 0, NULL, NULL, 0, NULL, user);CHKERRQ(ierr); 482c4762a1bSJed Brown PetscFunctionReturn(0); 483c4762a1bSJed Brown } 484c4762a1bSJed Brown 485c4762a1bSJed Brown PetscErrorCode SetupMaterial(DM dm, DM dmAux, AppCtx *user) 486c4762a1bSJed Brown { 487c4762a1bSJed Brown PetscErrorCode (*matFuncs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx) = {elasticityMaterial, wallPressure}; 488c4762a1bSJed Brown Vec A; 489c4762a1bSJed Brown void *ctxs[2]; 490c4762a1bSJed Brown PetscErrorCode ierr; 491c4762a1bSJed Brown 492c4762a1bSJed Brown PetscFunctionBegin; 493c4762a1bSJed Brown ctxs[0] = user; ctxs[1] = user; 494c4762a1bSJed Brown ierr = DMCreateLocalVector(dmAux, &A);CHKERRQ(ierr); 495c4762a1bSJed Brown ierr = DMProjectFunctionLocal(dmAux, 0.0, matFuncs, ctxs, INSERT_ALL_VALUES, A);CHKERRQ(ierr); 496c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dm, "A", (PetscObject) A);CHKERRQ(ierr); 497c4762a1bSJed Brown ierr = VecDestroy(&A);CHKERRQ(ierr); 498c4762a1bSJed Brown PetscFunctionReturn(0); 499c4762a1bSJed Brown } 500c4762a1bSJed Brown 501c4762a1bSJed Brown PetscErrorCode SetupNearNullSpace(DM dm, AppCtx *user) 502c4762a1bSJed Brown { 503c4762a1bSJed Brown /* Set up the near null space (a.k.a. rigid body modes) that will be used by the multigrid preconditioner */ 504c4762a1bSJed Brown DM subdm; 505c4762a1bSJed Brown MatNullSpace nearNullSpace; 506c4762a1bSJed Brown PetscInt fields = 0; 507c4762a1bSJed Brown PetscObject deformation; 508c4762a1bSJed Brown PetscErrorCode ierr; 509c4762a1bSJed Brown 510c4762a1bSJed Brown PetscFunctionBeginUser; 511c4762a1bSJed Brown ierr = DMCreateSubDM(dm, 1, &fields, NULL, &subdm);CHKERRQ(ierr); 512c4762a1bSJed Brown ierr = DMPlexCreateRigidBody(subdm, &nearNullSpace);CHKERRQ(ierr); 513c4762a1bSJed Brown ierr = DMGetField(dm, 0, NULL, &deformation);CHKERRQ(ierr); 514c4762a1bSJed Brown ierr = PetscObjectCompose(deformation, "nearnullspace", (PetscObject) nearNullSpace);CHKERRQ(ierr); 515c4762a1bSJed Brown ierr = DMDestroy(&subdm);CHKERRQ(ierr); 516c4762a1bSJed Brown ierr = MatNullSpaceDestroy(&nearNullSpace);CHKERRQ(ierr); 517c4762a1bSJed Brown PetscFunctionReturn(0); 518c4762a1bSJed Brown } 519c4762a1bSJed Brown 520c4762a1bSJed Brown static PetscErrorCode SetupAuxDM(DM dm, PetscInt NfAux, PetscFE feAux[], AppCtx *user) 521c4762a1bSJed Brown { 522c4762a1bSJed Brown DM dmAux, coordDM; 523c4762a1bSJed Brown PetscInt f; 524c4762a1bSJed Brown PetscErrorCode ierr; 525c4762a1bSJed Brown 526c4762a1bSJed Brown PetscFunctionBegin; 527c4762a1bSJed Brown /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */ 528c4762a1bSJed Brown ierr = DMGetCoordinateDM(dm, &coordDM);CHKERRQ(ierr); 529c4762a1bSJed Brown if (!feAux) PetscFunctionReturn(0); 530c4762a1bSJed Brown ierr = DMClone(dm, &dmAux);CHKERRQ(ierr); 531c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dm, "dmAux", (PetscObject) dmAux);CHKERRQ(ierr); 532c4762a1bSJed Brown ierr = DMSetCoordinateDM(dmAux, coordDM);CHKERRQ(ierr); 533c4762a1bSJed Brown for (f = 0; f < NfAux; ++f) {ierr = DMSetField(dmAux, f, NULL, (PetscObject) feAux[f]);CHKERRQ(ierr);} 534c4762a1bSJed Brown ierr = DMCreateDS(dmAux);CHKERRQ(ierr); 535c4762a1bSJed Brown ierr = SetupMaterial(dm, dmAux, user);CHKERRQ(ierr); 536c4762a1bSJed Brown ierr = DMDestroy(&dmAux);CHKERRQ(ierr); 537c4762a1bSJed Brown PetscFunctionReturn(0); 538c4762a1bSJed Brown } 539c4762a1bSJed Brown 540c4762a1bSJed Brown PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 541c4762a1bSJed Brown { 542c4762a1bSJed Brown DM cdm = dm; 543c4762a1bSJed Brown const PetscInt dim = user->dim; 544c4762a1bSJed Brown const PetscBool simplex = user->simplex; 545c4762a1bSJed Brown PetscFE fe[2], feAux[2]; 546c4762a1bSJed Brown MPI_Comm comm; 547c4762a1bSJed Brown PetscErrorCode ierr; 548c4762a1bSJed Brown 549c4762a1bSJed Brown PetscFunctionBeginUser; 550c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 551c4762a1bSJed Brown /* Create finite element */ 552c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, dim, simplex, "def_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr); 553c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fe[0], "deformation");CHKERRQ(ierr); 554c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr); 555c4762a1bSJed Brown ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr); 556c4762a1bSJed Brown 557c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr); 558c4762a1bSJed Brown 559c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, 1, simplex, "elastMat_", PETSC_DEFAULT, &feAux[0]);CHKERRQ(ierr); 560c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) feAux[0], "elasticityMaterial");CHKERRQ(ierr); 561c4762a1bSJed Brown ierr = PetscFECopyQuadrature(fe[0], feAux[0]);CHKERRQ(ierr); 562c4762a1bSJed 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 */ 563c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, 1, simplex, "wall_pres_", PETSC_DEFAULT, &feAux[1]);CHKERRQ(ierr); 564c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) feAux[1], "wall_pressure");CHKERRQ(ierr); 565c4762a1bSJed Brown ierr = PetscFECopyQuadrature(fe[0], feAux[1]);CHKERRQ(ierr); 566c4762a1bSJed Brown 567c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */ 568c4762a1bSJed Brown ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr); 569c4762a1bSJed Brown ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr); 570c4762a1bSJed Brown ierr = DMCreateDS(dm);CHKERRQ(ierr); 571c4762a1bSJed Brown ierr = SetupProblem(dm, dim, user);CHKERRQ(ierr); 572c4762a1bSJed Brown while (cdm) { 573c4762a1bSJed Brown ierr = SetupAuxDM(cdm, 2, feAux, user);CHKERRQ(ierr); 574*408cafa0SMatthew G. Knepley ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr); 575c4762a1bSJed Brown ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 576c4762a1bSJed Brown } 577c4762a1bSJed Brown ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr); 578c4762a1bSJed Brown ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr); 579c4762a1bSJed Brown ierr = PetscFEDestroy(&feAux[0]);CHKERRQ(ierr); 580c4762a1bSJed Brown ierr = PetscFEDestroy(&feAux[1]);CHKERRQ(ierr); 581c4762a1bSJed Brown PetscFunctionReturn(0); 582c4762a1bSJed Brown } 583c4762a1bSJed Brown 584c4762a1bSJed Brown 585c4762a1bSJed Brown int main(int argc, char **argv) 586c4762a1bSJed Brown { 587c4762a1bSJed Brown SNES snes; /* nonlinear solver */ 588c4762a1bSJed Brown DM dm; /* problem definition */ 589c4762a1bSJed Brown Vec u,r; /* solution, residual vectors */ 590c4762a1bSJed Brown Mat A,J; /* Jacobian matrix */ 591c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 592c4762a1bSJed Brown PetscInt its; /* iterations for convergence */ 593c4762a1bSJed Brown PetscErrorCode ierr; 594c4762a1bSJed Brown 595c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 596c4762a1bSJed Brown ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 597c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr); 598c4762a1bSJed Brown ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 599c4762a1bSJed Brown ierr = SNESSetDM(snes, dm);CHKERRQ(ierr); 600c4762a1bSJed Brown ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr); 601c4762a1bSJed Brown 602c4762a1bSJed Brown ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr); 603c4762a1bSJed Brown ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr); 604c4762a1bSJed Brown ierr = SetupNearNullSpace(dm, &user);CHKERRQ(ierr); 605c4762a1bSJed Brown 606c4762a1bSJed Brown ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 607c4762a1bSJed Brown ierr = VecDuplicate(u, &r);CHKERRQ(ierr); 608c4762a1bSJed Brown 609c4762a1bSJed Brown ierr = DMSetMatType(dm,MATAIJ);CHKERRQ(ierr); 610c4762a1bSJed Brown ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr); 611c4762a1bSJed Brown A = J; 612c4762a1bSJed Brown 613c4762a1bSJed Brown ierr = DMPlexSetSNESLocalFEM(dm,&user,&user,&user);CHKERRQ(ierr); 614c4762a1bSJed Brown ierr = SNESSetJacobian(snes, A, J, NULL, NULL);CHKERRQ(ierr); 615c4762a1bSJed Brown 616c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 617c4762a1bSJed Brown 618c4762a1bSJed Brown { 619c4762a1bSJed Brown PetscErrorCode (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void* ctx); 620c4762a1bSJed Brown initialGuess[0] = coordinates; initialGuess[1] = zero_scalar; 621c4762a1bSJed Brown ierr = DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u);CHKERRQ(ierr); 622c4762a1bSJed Brown } 623c4762a1bSJed Brown if (user.showInitial) {ierr = DMVecViewLocal(dm, u, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);} 624c4762a1bSJed Brown 625c4762a1bSJed Brown if (user.runType == RUN_FULL) { 626c4762a1bSJed Brown if (user.debug) { 627c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n");CHKERRQ(ierr); 628c4762a1bSJed Brown ierr = VecView(u, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 629c4762a1bSJed Brown } 630c4762a1bSJed Brown ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr); 631c4762a1bSJed Brown ierr = SNESGetIterationNumber(snes, &its);CHKERRQ(ierr); 632c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %D\n", its);CHKERRQ(ierr); 633c4762a1bSJed Brown if (user.showSolution) { 634c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Solution\n");CHKERRQ(ierr); 635c4762a1bSJed Brown ierr = VecChop(u, 3.0e-9);CHKERRQ(ierr); 636c4762a1bSJed Brown ierr = VecView(u, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 637c4762a1bSJed Brown } 638c4762a1bSJed Brown } else { 639c4762a1bSJed Brown PetscReal res = 0.0; 640c4762a1bSJed Brown 641c4762a1bSJed Brown /* Check initial guess */ 642c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n");CHKERRQ(ierr); 643c4762a1bSJed Brown ierr = VecView(u, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 644c4762a1bSJed Brown /* Check residual */ 645c4762a1bSJed Brown ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr); 646c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n");CHKERRQ(ierr); 647c4762a1bSJed Brown ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr); 648c4762a1bSJed Brown ierr = VecView(r, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 649c4762a1bSJed Brown ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr); 650c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", (double)res);CHKERRQ(ierr); 651c4762a1bSJed Brown /* Check Jacobian */ 652c4762a1bSJed Brown { 653c4762a1bSJed Brown Vec b; 654c4762a1bSJed Brown 655c4762a1bSJed Brown ierr = SNESComputeJacobian(snes, u, A, A);CHKERRQ(ierr); 656c4762a1bSJed Brown ierr = VecDuplicate(u, &b);CHKERRQ(ierr); 657c4762a1bSJed Brown ierr = VecSet(r, 0.0);CHKERRQ(ierr); 658c4762a1bSJed Brown ierr = SNESComputeFunction(snes, r, b);CHKERRQ(ierr); 659c4762a1bSJed Brown ierr = MatMult(A, u, r);CHKERRQ(ierr); 660c4762a1bSJed Brown ierr = VecAXPY(r, 1.0, b);CHKERRQ(ierr); 661c4762a1bSJed Brown ierr = VecDestroy(&b);CHKERRQ(ierr); 662c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Au - b = Au + F(0)\n");CHKERRQ(ierr); 663c4762a1bSJed Brown ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr); 664c4762a1bSJed Brown ierr = VecView(r, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 665c4762a1bSJed Brown ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr); 666c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", (double)res);CHKERRQ(ierr); 667c4762a1bSJed Brown } 668c4762a1bSJed Brown } 669c4762a1bSJed Brown ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr); 670c4762a1bSJed Brown 671c4762a1bSJed Brown if (A != J) {ierr = MatDestroy(&A);CHKERRQ(ierr);} 672c4762a1bSJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); 673c4762a1bSJed Brown ierr = VecDestroy(&u);CHKERRQ(ierr); 674c4762a1bSJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 675c4762a1bSJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 676c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 677c4762a1bSJed Brown ierr = PetscFinalize(); 678c4762a1bSJed Brown return ierr; 679c4762a1bSJed Brown } 680c4762a1bSJed Brown 681c4762a1bSJed Brown /*TEST 682c4762a1bSJed Brown 683c4762a1bSJed Brown build: 684c4762a1bSJed Brown requires: !complex 685c4762a1bSJed Brown 686c4762a1bSJed Brown test: 687c4762a1bSJed Brown suffix: 0 688c4762a1bSJed Brown requires: ctetgen !single 689c4762a1bSJed Brown args: -run_type full -dim 3 -dm_refine 2 -interpolate 1 -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4 -def_petscspace_degree 2 -pres_petscspace_degree 1 -elastMat_petscspace_degree 0 -wall_pres_petscspace_degree 0 -snes_rtol 1e-05 -ksp_type fgmres -ksp_rtol 1e-10 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper -fieldsplit_deformation_ksp_type preonly -fieldsplit_deformation_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi -snes_monitor_short -ksp_monitor_short -snes_converged_reason -ksp_converged_reason -show_solution 0 690c4762a1bSJed Brown 691c4762a1bSJed Brown test: 692c4762a1bSJed Brown suffix: 1 693c4762a1bSJed Brown requires: ctetgen superlu_dist 694c4762a1bSJed Brown nsize: 2 695c4762a1bSJed Brown args: -run_type full -dim 3 -dm_refine 0 -interpolate 1 -test_partition -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4 -def_petscspace_degree 2 -pres_petscspace_degree 1 -elastMat_petscspace_degree 0 -wall_pres_petscspace_degree 0 -snes_rtol 1e-05 -ksp_type fgmres -ksp_rtol 1e-10 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper -fieldsplit_deformation_ksp_type preonly -fieldsplit_deformation_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi -snes_monitor_short -ksp_monitor_short -snes_converged_reason -ksp_converged_reason -show_solution 0 696c4762a1bSJed Brown timeoutfactor: 2 697c4762a1bSJed Brown 698c4762a1bSJed Brown test: 699c4762a1bSJed Brown suffix: 2 700c4762a1bSJed Brown requires: ctetgen !single 701c4762a1bSJed Brown args: -run_type full -dim 3 -dm_refine 2 -interpolate 1 -bc_fixed 3,4,5,6 -bc_pressure 2 -wall_pressure 1.0 -def_petscspace_degree 2 -pres_petscspace_degree 1 -elastMat_petscspace_degree 0 -wall_pres_petscspace_degree 0 -snes_rtol 1e-05 -ksp_type fgmres -ksp_rtol 1e-10 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper -fieldsplit_deformation_ksp_type preonly -fieldsplit_deformation_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi -snes_monitor_short -ksp_monitor_short -snes_converged_reason -ksp_converged_reason -show_solution 0 702c4762a1bSJed Brown 703c4762a1bSJed Brown test: 704c4762a1bSJed Brown requires: ctetgen superlu_dist 705c4762a1bSJed Brown suffix: 4 706c4762a1bSJed Brown nsize: 2 707c4762a1bSJed Brown args: -run_type full -dim 3 -dm_refine 0 -interpolate 1 -test_partition -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4 -def_petscspace_degree 2 -pres_petscspace_degree 1 -elastMat_petscspace_degree 0 -wall_pres_petscspace_degree 0 -snes_rtol 1e-05 -ksp_type fgmres -ksp_rtol 1e-10 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper -fieldsplit_deformation_ksp_type preonly -fieldsplit_deformation_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi -snes_monitor_short -ksp_monitor_short -snes_converged_reason -ksp_converged_reason -show_solution 0 708c4762a1bSJed Brown output_file: output/ex77_1.out 709c4762a1bSJed Brown 710c4762a1bSJed Brown test: 711c4762a1bSJed Brown requires: ctetgen !single 712c4762a1bSJed Brown suffix: 3 713c4762a1bSJed Brown args: -run_type full -dim 3 -dm_refine 2 -interpolate 1 -bc_fixed 3,4,5,6 -bc_pressure 2 -wall_pressure 1.0 -def_petscspace_degree 2 -pres_petscspace_degree 1 -elastMat_petscspace_degree 0 -wall_pres_petscspace_degree 0 -snes_rtol 1e-05 -ksp_type fgmres -ksp_rtol 1e-10 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper -fieldsplit_deformation_ksp_type preonly -fieldsplit_deformation_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi -snes_monitor_short -ksp_monitor_short -snes_converged_reason -ksp_converged_reason -show_solution 0 714c4762a1bSJed Brown output_file: output/ex77_2.out 715c4762a1bSJed Brown 716c4762a1bSJed Brown #TODO: this example deadlocks for me when using ParMETIS 717c4762a1bSJed Brown test: 718c4762a1bSJed Brown requires: ctetgen superlu_dist !single 719c4762a1bSJed Brown suffix: 2_par 720c4762a1bSJed Brown nsize: 4 721c4762a1bSJed Brown args: -run_type full -dim 3 -dm_refine 2 -interpolate 1 -bc_fixed 3,4,5,6 -bc_pressure 2 -wall_pressure 1.0 -def_petscspace_degree 2 -pres_petscspace_degree 1 -elastMat_petscspace_degree 0 -wall_pres_petscspace_degree 0 -snes_rtol 1e-05 -ksp_type fgmres -ksp_rtol 1e-10 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper -fieldsplit_deformation_ksp_type preonly -fieldsplit_deformation_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi -snes_monitor_short -ksp_monitor_short -snes_converged_reason -ksp_converged_reason -show_solution 0 -petscpartitioner_type simple 722c4762a1bSJed Brown output_file: output/ex77_2.out 723c4762a1bSJed Brown 724c4762a1bSJed Brown TEST*/ 725