1 static char help[] = "Nonlinear elasticity problem in 3d with simplicial finite elements.\n\ 2 We solve a nonlinear elasticity problem, modelled as an incompressible Neo-Hookean solid, \n\ 3 with pressure loading in a rectangular domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n"; 4 5 /* 6 Nonlinear elasticity problem, which we discretize using the finite 7 element method on an unstructured mesh. This uses both Dirichlet boundary conditions (fixed faces) 8 and nonlinear Neumann boundary conditions (pressure loading). 9 The Lagrangian density (modulo boundary conditions) for this problem is given by 10 \begin{equation} 11 \frac{\mu}{2} (\mathrm{Tr}{C}-3) + J p + \frac{\kappa}{2} (J-1). 12 \end{equation} 13 14 Discretization: 15 16 We use PetscFE to generate a tabulation of the finite element basis functions 17 at quadrature points. We can currently generate an arbitrary order Lagrange 18 element. 19 20 Field Data: 21 22 DMPLEX data is organized by point, and the closure operation just stacks up the 23 data from each sieve point in the closure. Thus, for a P_2-P_1 Stokes element, we 24 have 25 26 cl{e} = {f e_0 e_1 e_2 v_0 v_1 v_2} 27 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 29 The problem here is that we would like to loop over each field separately for 30 integration. Therefore, the closure visitor in DMPlexVecGetClosure() reorders 31 the data so that each field is contiguous 32 33 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 35 Likewise, DMPlexVecSetClosure() takes data partitioned by field, and correctly 36 puts it into the Sieve ordering. 37 38 */ 39 40 #include <petscdmplex.h> 41 #include <petscsnes.h> 42 #include <petscds.h> 43 44 typedef enum {RUN_FULL, RUN_TEST} RunType; 45 46 typedef struct { 47 RunType runType; /* Whether to run tests, or solve the full problem */ 48 PetscReal mu; /* The shear modulus */ 49 PetscReal p_wall; /* The wall pressure */ 50 } AppCtx; 51 52 #if 0 53 PETSC_STATIC_INLINE void Det2D(PetscReal *detJ, const PetscReal J[]) 54 { 55 *detJ = J[0]*J[3] - J[1]*J[2]; 56 } 57 #endif 58 59 PETSC_STATIC_INLINE void Det3D(PetscReal *detJ, const PetscScalar J[]) 60 { 61 *detJ = PetscRealPart(J[0*3+0]*(J[1*3+1]*J[2*3+2] - J[1*3+2]*J[2*3+1]) + 62 J[0*3+1]*(J[1*3+2]*J[2*3+0] - J[1*3+0]*J[2*3+2]) + 63 J[0*3+2]*(J[1*3+0]*J[2*3+1] - J[1*3+1]*J[2*3+0])); 64 } 65 66 #if 0 67 PETSC_STATIC_INLINE void Cof2D(PetscReal C[], const PetscReal A[]) 68 { 69 C[0] = A[3]; 70 C[1] = -A[2]; 71 C[2] = -A[1]; 72 C[3] = A[0]; 73 } 74 #endif 75 76 PETSC_STATIC_INLINE void Cof3D(PetscReal C[], const PetscScalar A[]) 77 { 78 C[0*3+0] = PetscRealPart(A[1*3+1]*A[2*3+2] - A[1*3+2]*A[2*3+1]); 79 C[0*3+1] = PetscRealPart(A[1*3+2]*A[2*3+0] - A[1*3+0]*A[2*3+2]); 80 C[0*3+2] = PetscRealPart(A[1*3+0]*A[2*3+1] - A[1*3+1]*A[2*3+0]); 81 C[1*3+0] = PetscRealPart(A[0*3+2]*A[2*3+1] - A[0*3+1]*A[2*3+2]); 82 C[1*3+1] = PetscRealPart(A[0*3+0]*A[2*3+2] - A[0*3+2]*A[2*3+0]); 83 C[1*3+2] = PetscRealPart(A[0*3+1]*A[2*3+0] - A[0*3+0]*A[2*3+1]); 84 C[2*3+0] = PetscRealPart(A[0*3+1]*A[1*3+2] - A[0*3+2]*A[1*3+1]); 85 C[2*3+1] = PetscRealPart(A[0*3+2]*A[1*3+0] - A[0*3+0]*A[1*3+2]); 86 C[2*3+2] = PetscRealPart(A[0*3+0]*A[1*3+1] - A[0*3+1]*A[1*3+0]); 87 } 88 89 PetscErrorCode zero_scalar(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 90 { 91 u[0] = 0.0; 92 return 0; 93 } 94 95 PetscErrorCode zero_vector(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 96 { 97 const PetscInt Ncomp = dim; 98 99 PetscInt comp; 100 for (comp = 0; comp < Ncomp; ++comp) u[comp] = 0.0; 101 return 0; 102 } 103 104 PetscErrorCode coordinates(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 105 { 106 const PetscInt Ncomp = dim; 107 108 PetscInt comp; 109 for (comp = 0; comp < Ncomp; ++comp) u[comp] = x[comp]; 110 return 0; 111 } 112 113 PetscErrorCode elasticityMaterial(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 114 { 115 AppCtx *user = (AppCtx *) ctx; 116 u[0] = user->mu; 117 return 0; 118 } 119 120 PetscErrorCode wallPressure(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 121 { 122 AppCtx *user = (AppCtx *) ctx; 123 u[0] = user->p_wall; 124 return 0; 125 } 126 127 void f1_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 128 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 129 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 130 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 131 { 132 const PetscInt Ncomp = dim; 133 const PetscReal mu = PetscRealPart(a[0]), kappa = 3.0; 134 PetscReal cofu_x[9/*Ncomp*dim*/], detu_x, p = PetscRealPart(u[Ncomp]); 135 PetscInt comp, d; 136 137 Cof3D(cofu_x, u_x); 138 Det3D(&detu_x, u_x); 139 p += kappa * (detu_x - 1.0); 140 141 /* f1 is the first Piola-Kirchhoff tensor */ 142 for (comp = 0; comp < Ncomp; ++comp) { 143 for (d = 0; d < dim; ++d) { 144 f1[comp*dim+d] = mu * u_x[comp*dim+d] + p * cofu_x[comp*dim+d]; 145 } 146 } 147 } 148 149 void g3_uu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 150 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 151 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 152 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 153 { 154 const PetscInt Ncomp = dim; 155 const PetscReal mu = PetscRealPart(a[0]), kappa = 3.0; 156 PetscReal cofu_x[9/*Ncomp*dim*/], detu_x, pp, pm, p = PetscRealPart(u[Ncomp]); 157 PetscInt compI, compJ, d1, d2; 158 159 Cof3D(cofu_x, u_x); 160 Det3D(&detu_x, u_x); 161 p += kappa * (detu_x - 1.0); 162 pp = p/detu_x + kappa; 163 pm = p/detu_x; 164 165 /* 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} */ 166 for (compI = 0; compI < Ncomp; ++compI) { 167 for (compJ = 0; compJ < Ncomp; ++compJ) { 168 const PetscReal G = (compI == compJ) ? 1.0 : 0.0; 169 170 for (d1 = 0; d1 < dim; ++d1) { 171 for (d2 = 0; d2 < dim; ++d2) { 172 const PetscReal g = (d1 == d2) ? 1.0 : 0.0; 173 174 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]; 175 } 176 } 177 } 178 } 179 } 180 181 void f0_bd_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 182 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 183 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 184 PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 185 { 186 const PetscInt Ncomp = dim; 187 const PetscScalar p = a[aOff[1]]; 188 PetscReal cofu_x[9/*Ncomp*dim*/]; 189 PetscInt comp, d; 190 191 Cof3D(cofu_x, u_x); 192 for (comp = 0; comp < Ncomp; ++comp) { 193 for (d = 0, f0[comp] = 0.0; d < dim; ++d) f0[comp] += cofu_x[comp*dim+d] * n[d]; 194 f0[comp] *= p; 195 } 196 } 197 198 void g1_bd_uu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 199 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 200 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 201 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 202 { 203 const PetscInt Ncomp = dim; 204 PetscScalar p = a[aOff[1]]; 205 PetscReal cofu_x[9/*Ncomp*dim*/], m[3/*Ncomp*/], detu_x; 206 PetscInt comp, compI, compJ, d; 207 208 Cof3D(cofu_x, u_x); 209 Det3D(&detu_x, u_x); 210 p /= detu_x; 211 212 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]; 213 for (compI = 0; compI < Ncomp; ++compI) { 214 for (compJ = 0; compJ < Ncomp; ++compJ) { 215 for (d = 0; d < dim; ++d) { 216 g1[(compI*Ncomp+compJ)*dim+d] = p * (m[compI] * cofu_x[compJ*dim+d] - cofu_x[compI*dim+d] * m[compJ]); 217 } 218 } 219 } 220 } 221 222 void f0_p_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 223 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 224 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 225 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 226 { 227 PetscReal detu_x; 228 Det3D(&detu_x, u_x); 229 f0[0] = detu_x - 1.0; 230 } 231 232 void g1_pu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 233 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 234 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 235 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) 236 { 237 PetscReal cofu_x[9/*Ncomp*dim*/]; 238 PetscInt compI, d; 239 240 Cof3D(cofu_x, u_x); 241 for (compI = 0; compI < dim; ++compI) 242 for (d = 0; d < dim; ++d) g1[compI*dim+d] = cofu_x[compI*dim+d]; 243 } 244 245 void g2_up_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, 246 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 247 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 248 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) 249 { 250 PetscReal cofu_x[9/*Ncomp*dim*/]; 251 PetscInt compI, d; 252 253 Cof3D(cofu_x, u_x); 254 for (compI = 0; compI < dim; ++compI) 255 for (d = 0; d < dim; ++d) g2[compI*dim+d] = cofu_x[compI*dim+d]; 256 } 257 258 PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 259 { 260 const char *runTypes[2] = {"full", "test"}; 261 PetscInt run; 262 PetscErrorCode ierr; 263 264 PetscFunctionBeginUser; 265 options->runType = RUN_FULL; 266 options->mu = 1.0; 267 options->p_wall = 0.4; 268 269 ierr = PetscOptionsBegin(comm, "", "Nonlinear elasticity problem options", "DMPLEX");CHKERRQ(ierr); 270 run = options->runType; 271 ierr = PetscOptionsEList("-run_type", "The run type", "ex77.c", runTypes, 2, runTypes[options->runType], &run, NULL);CHKERRQ(ierr); 272 options->runType = (RunType) run; 273 ierr = PetscOptionsReal("-shear_modulus", "The shear modulus", "ex77.c", options->mu, &options->mu, NULL);CHKERRQ(ierr); 274 ierr = PetscOptionsReal("-wall_pressure", "The wall pressure", "ex77.c", options->p_wall, &options->p_wall, NULL);CHKERRQ(ierr); 275 ierr = PetscOptionsEnd(); 276 PetscFunctionReturn(0); 277 } 278 279 PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 280 { 281 PetscErrorCode ierr; 282 283 PetscFunctionBeginUser; 284 /* TODO The P1 coordinate space gives wrong results when compared to the affine version. Track this down */ 285 if (0) { 286 ierr = DMPlexCreateBoxMesh(comm, 3, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr); 287 } else { 288 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 289 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 290 } 291 ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 292 /* Label the faces (bit of a hack here, until it is properly implemented for simplices) */ 293 ierr = DMViewFromOptions(*dm, NULL, "-orig_dm_view");CHKERRQ(ierr); 294 { 295 DM cdm; 296 DMLabel label; 297 IS is; 298 PetscInt d, dim, b, f, Nf; 299 const PetscInt *faces; 300 PetscInt csize; 301 PetscScalar *coords = NULL; 302 PetscSection cs; 303 Vec coordinates ; 304 305 ierr = DMGetDimension(*dm, &dim);CHKERRQ(ierr); 306 ierr = DMCreateLabel(*dm, "boundary");CHKERRQ(ierr); 307 ierr = DMGetLabel(*dm, "boundary", &label);CHKERRQ(ierr); 308 ierr = DMPlexMarkBoundaryFaces(*dm, 1, label);CHKERRQ(ierr); 309 ierr = DMGetStratumIS(*dm, "boundary", 1, &is);CHKERRQ(ierr); 310 if (is) { 311 PetscReal faceCoord; 312 PetscInt v; 313 314 ierr = ISGetLocalSize(is, &Nf);CHKERRQ(ierr); 315 ierr = ISGetIndices(is, &faces);CHKERRQ(ierr); 316 317 ierr = DMGetCoordinatesLocal(*dm, &coordinates);CHKERRQ(ierr); 318 ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr); 319 ierr = DMGetLocalSection(cdm, &cs);CHKERRQ(ierr); 320 321 /* Check for each boundary face if any component of its centroid is either 0.0 or 1.0 */ 322 for (f = 0; f < Nf; ++f) { 323 ierr = DMPlexVecGetClosure(cdm, cs, coordinates, faces[f], &csize, &coords);CHKERRQ(ierr); 324 /* Calculate mean coordinate vector */ 325 for (d = 0; d < dim; ++d) { 326 const PetscInt Nv = csize/dim; 327 faceCoord = 0.0; 328 for (v = 0; v < Nv; ++v) faceCoord += PetscRealPart(coords[v*dim+d]); 329 faceCoord /= Nv; 330 for (b = 0; b < 2; ++b) { 331 if (PetscAbs(faceCoord - b*1.0) < PETSC_SMALL) { 332 ierr = DMSetLabelValue(*dm, "Faces", faces[f], d*2+b+1);CHKERRQ(ierr); 333 } 334 } 335 } 336 ierr = DMPlexVecRestoreClosure(cdm, cs, coordinates, faces[f], &csize, &coords);CHKERRQ(ierr); 337 } 338 ierr = ISRestoreIndices(is, &faces);CHKERRQ(ierr); 339 } 340 ierr = ISDestroy(&is);CHKERRQ(ierr); 341 } 342 ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 343 PetscFunctionReturn(0); 344 } 345 346 PetscErrorCode SetupProblem(DM dm, PetscInt dim, AppCtx *user) 347 { 348 PetscDS ds; 349 PetscWeakForm wf; 350 DMLabel label; 351 PetscInt bd; 352 PetscErrorCode ierr; 353 354 PetscFunctionBeginUser; 355 ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 356 ierr = PetscDSSetResidual(ds, 0, NULL, f1_u_3d);CHKERRQ(ierr); 357 ierr = PetscDSSetResidual(ds, 1, f0_p_3d, NULL);CHKERRQ(ierr); 358 ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_uu_3d);CHKERRQ(ierr); 359 ierr = PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up_3d, NULL);CHKERRQ(ierr); 360 ierr = PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu_3d, NULL, NULL);CHKERRQ(ierr); 361 362 ierr = DMGetLabel(dm, "Faces", &label);CHKERRQ(ierr); 363 ierr = DMAddBoundary(dm, DM_BC_NATURAL, "pressure", label, 0, NULL, 0, 0, NULL, NULL, NULL, user, &bd);CHKERRQ(ierr); 364 ierr = PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 365 ierr = PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, 0, f0_bd_u_3d, 0, NULL);CHKERRQ(ierr); 366 ierr = PetscWeakFormSetIndexBdJacobian(wf, label, 1, 0, 0, 0, 0, NULL, 0, g1_bd_uu_3d, 0, NULL, 0, NULL);CHKERRQ(ierr); 367 368 ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "fixed", label, 0, NULL, 0, 0, NULL, (void (*)(void)) coordinates, NULL, user, NULL);CHKERRQ(ierr); 369 PetscFunctionReturn(0); 370 } 371 372 PetscErrorCode SetupMaterial(DM dm, DM dmAux, AppCtx *user) 373 { 374 PetscErrorCode (*matFuncs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx) = {elasticityMaterial, wallPressure}; 375 Vec A; 376 void *ctxs[2]; 377 PetscErrorCode ierr; 378 379 PetscFunctionBegin; 380 ctxs[0] = user; ctxs[1] = user; 381 ierr = DMCreateLocalVector(dmAux, &A);CHKERRQ(ierr); 382 ierr = DMProjectFunctionLocal(dmAux, 0.0, matFuncs, ctxs, INSERT_ALL_VALUES, A);CHKERRQ(ierr); 383 ierr = DMSetAuxiliaryVec(dm, NULL, 0, A);CHKERRQ(ierr); 384 ierr = VecDestroy(&A);CHKERRQ(ierr); 385 PetscFunctionReturn(0); 386 } 387 388 PetscErrorCode SetupNearNullSpace(DM dm, AppCtx *user) 389 { 390 /* Set up the near null space (a.k.a. rigid body modes) that will be used by the multigrid preconditioner */ 391 DM subdm; 392 MatNullSpace nearNullSpace; 393 PetscInt fields = 0; 394 PetscObject deformation; 395 PetscErrorCode ierr; 396 397 PetscFunctionBeginUser; 398 ierr = DMCreateSubDM(dm, 1, &fields, NULL, &subdm);CHKERRQ(ierr); 399 ierr = DMPlexCreateRigidBody(subdm, 0, &nearNullSpace);CHKERRQ(ierr); 400 ierr = DMGetField(dm, 0, NULL, &deformation);CHKERRQ(ierr); 401 ierr = PetscObjectCompose(deformation, "nearnullspace", (PetscObject) nearNullSpace);CHKERRQ(ierr); 402 ierr = DMDestroy(&subdm);CHKERRQ(ierr); 403 ierr = MatNullSpaceDestroy(&nearNullSpace);CHKERRQ(ierr); 404 PetscFunctionReturn(0); 405 } 406 407 static PetscErrorCode SetupAuxDM(DM dm, PetscInt NfAux, PetscFE feAux[], AppCtx *user) 408 { 409 DM dmAux, coordDM; 410 PetscInt f; 411 PetscErrorCode ierr; 412 413 PetscFunctionBegin; 414 /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */ 415 ierr = DMGetCoordinateDM(dm, &coordDM);CHKERRQ(ierr); 416 if (!feAux) PetscFunctionReturn(0); 417 ierr = DMClone(dm, &dmAux);CHKERRQ(ierr); 418 ierr = DMSetCoordinateDM(dmAux, coordDM);CHKERRQ(ierr); 419 for (f = 0; f < NfAux; ++f) {ierr = DMSetField(dmAux, f, NULL, (PetscObject) feAux[f]);CHKERRQ(ierr);} 420 ierr = DMCreateDS(dmAux);CHKERRQ(ierr); 421 ierr = SetupMaterial(dm, dmAux, user);CHKERRQ(ierr); 422 ierr = DMDestroy(&dmAux);CHKERRQ(ierr); 423 PetscFunctionReturn(0); 424 } 425 426 PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 427 { 428 DM cdm = dm; 429 PetscFE fe[2], feAux[2]; 430 PetscBool simplex; 431 PetscInt dim; 432 MPI_Comm comm; 433 PetscErrorCode ierr; 434 435 PetscFunctionBeginUser; 436 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 437 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 438 ierr = DMPlexIsSimplex(dm, &simplex);CHKERRQ(ierr); 439 /* Create finite element */ 440 ierr = PetscFECreateDefault(comm, dim, dim, simplex, "def_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr); 441 ierr = PetscObjectSetName((PetscObject) fe[0], "deformation");CHKERRQ(ierr); 442 ierr = PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr); 443 ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr); 444 445 ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr); 446 447 ierr = PetscFECreateDefault(comm, dim, 1, simplex, "elastMat_", PETSC_DEFAULT, &feAux[0]);CHKERRQ(ierr); 448 ierr = PetscObjectSetName((PetscObject) feAux[0], "elasticityMaterial");CHKERRQ(ierr); 449 ierr = PetscFECopyQuadrature(fe[0], feAux[0]);CHKERRQ(ierr); 450 /* It is not yet possible to define a field on a submesh (e.g. a boundary), so we will use a normal finite element */ 451 ierr = PetscFECreateDefault(comm, dim, 1, simplex, "wall_pres_", PETSC_DEFAULT, &feAux[1]);CHKERRQ(ierr); 452 ierr = PetscObjectSetName((PetscObject) feAux[1], "wall_pressure");CHKERRQ(ierr); 453 ierr = PetscFECopyQuadrature(fe[0], feAux[1]);CHKERRQ(ierr); 454 455 /* Set discretization and boundary conditions for each mesh */ 456 ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr); 457 ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr); 458 ierr = DMCreateDS(dm);CHKERRQ(ierr); 459 ierr = SetupProblem(dm, dim, user);CHKERRQ(ierr); 460 while (cdm) { 461 ierr = SetupAuxDM(cdm, 2, feAux, user);CHKERRQ(ierr); 462 ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr); 463 ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr); 464 } 465 ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr); 466 ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr); 467 ierr = PetscFEDestroy(&feAux[0]);CHKERRQ(ierr); 468 ierr = PetscFEDestroy(&feAux[1]);CHKERRQ(ierr); 469 PetscFunctionReturn(0); 470 } 471 472 int main(int argc, char **argv) 473 { 474 SNES snes; /* nonlinear solver */ 475 DM dm; /* problem definition */ 476 Vec u,r; /* solution, residual vectors */ 477 Mat A,J; /* Jacobian matrix */ 478 AppCtx user; /* user-defined work context */ 479 PetscInt its; /* iterations for convergence */ 480 PetscErrorCode ierr; 481 482 ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 483 ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 484 ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr); 485 ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 486 ierr = SNESSetDM(snes, dm);CHKERRQ(ierr); 487 ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr); 488 489 ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr); 490 ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr); 491 ierr = SetupNearNullSpace(dm, &user);CHKERRQ(ierr); 492 493 ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr); 494 ierr = VecDuplicate(u, &r);CHKERRQ(ierr); 495 496 ierr = DMSetMatType(dm,MATAIJ);CHKERRQ(ierr); 497 ierr = DMCreateMatrix(dm, &J);CHKERRQ(ierr); 498 A = J; 499 500 ierr = DMPlexSetSNESLocalFEM(dm,&user,&user,&user);CHKERRQ(ierr); 501 ierr = SNESSetJacobian(snes, A, J, NULL, NULL);CHKERRQ(ierr); 502 503 ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 504 505 { 506 PetscErrorCode (*initialGuess[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void* ctx); 507 initialGuess[0] = coordinates; initialGuess[1] = zero_scalar; 508 ierr = DMProjectFunction(dm, 0.0, initialGuess, NULL, INSERT_VALUES, u);CHKERRQ(ierr); 509 } 510 511 if (user.runType == RUN_FULL) { 512 ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr); 513 ierr = SNESGetIterationNumber(snes, &its);CHKERRQ(ierr); 514 ierr = PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %D\n", its);CHKERRQ(ierr); 515 } else { 516 PetscReal res = 0.0; 517 518 /* Check initial guess */ 519 ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial guess\n");CHKERRQ(ierr); 520 ierr = VecView(u, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 521 /* Check residual */ 522 ierr = SNESComputeFunction(snes, u, r);CHKERRQ(ierr); 523 ierr = PetscPrintf(PETSC_COMM_WORLD, "Initial Residual\n");CHKERRQ(ierr); 524 ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr); 525 ierr = VecView(r, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 526 ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr); 527 ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Residual: %g\n", (double)res);CHKERRQ(ierr); 528 /* Check Jacobian */ 529 { 530 Vec b; 531 532 ierr = SNESComputeJacobian(snes, u, A, A);CHKERRQ(ierr); 533 ierr = VecDuplicate(u, &b);CHKERRQ(ierr); 534 ierr = VecSet(r, 0.0);CHKERRQ(ierr); 535 ierr = SNESComputeFunction(snes, r, b);CHKERRQ(ierr); 536 ierr = MatMult(A, u, r);CHKERRQ(ierr); 537 ierr = VecAXPY(r, 1.0, b);CHKERRQ(ierr); 538 ierr = VecDestroy(&b);CHKERRQ(ierr); 539 ierr = PetscPrintf(PETSC_COMM_WORLD, "Au - b = Au + F(0)\n");CHKERRQ(ierr); 540 ierr = VecChop(r, 1.0e-10);CHKERRQ(ierr); 541 ierr = VecView(r, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 542 ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr); 543 ierr = PetscPrintf(PETSC_COMM_WORLD, "Linear L_2 Residual: %g\n", (double)res);CHKERRQ(ierr); 544 } 545 } 546 ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr); 547 548 if (A != J) {ierr = MatDestroy(&A);CHKERRQ(ierr);} 549 ierr = MatDestroy(&J);CHKERRQ(ierr); 550 ierr = VecDestroy(&u);CHKERRQ(ierr); 551 ierr = VecDestroy(&r);CHKERRQ(ierr); 552 ierr = SNESDestroy(&snes);CHKERRQ(ierr); 553 ierr = DMDestroy(&dm);CHKERRQ(ierr); 554 ierr = PetscFinalize(); 555 return ierr; 556 } 557 558 /*TEST 559 560 build: 561 requires: !complex 562 563 testset: 564 requires: ctetgen 565 args: -run_type full -dm_plex_dim 3 \ 566 -def_petscspace_degree 2 -pres_petscspace_degree 1 -elastMat_petscspace_degree 0 -wall_pres_petscspace_degree 0 \ 567 -snes_rtol 1e-05 -snes_monitor_short -snes_converged_reason \ 568 -ksp_type fgmres -ksp_rtol 1e-10 -ksp_monitor_short -ksp_converged_reason \ 569 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper \ 570 -fieldsplit_deformation_ksp_type preonly -fieldsplit_deformation_pc_type lu \ 571 -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi 572 573 test: 574 suffix: 0 575 requires: !single 576 args: -dm_refine 2 \ 577 -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4 578 579 test: 580 suffix: 1 581 requires: superlu_dist 582 nsize: 2 583 args: -dm_distribute -dm_refine 0 -petscpartitioner_type simple \ 584 -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4 585 timeoutfactor: 2 586 587 test: 588 suffix: 4 589 requires: superlu_dist 590 nsize: 2 591 args: -dm_distribute -dm_refine 0 -petscpartitioner_type simple \ 592 -bc_fixed 1 -bc_pressure 2 -wall_pressure 0.4 593 output_file: output/ex77_1.out 594 595 test: 596 suffix: 2 597 requires: !single 598 args: -dm_refine 2 \ 599 -bc_fixed 3,4,5,6 -bc_pressure 2 -wall_pressure 1.0 600 601 test: 602 suffix: 2_par 603 requires: superlu_dist !single 604 nsize: 4 605 args: -dm_distribute -dm_refine 2 -petscpartitioner_type simple \ 606 -bc_fixed 3,4,5,6 -bc_pressure 2 -wall_pressure 1.0 607 output_file: output/ex77_2.out 608 609 TEST*/ 610