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