1 static char help[] = "Check that a DM can accurately represent and interpolate functions of a given polynomial order\n\n"; 2 3 #include <petscdmplex.h> 4 #include <petscdm.h> 5 #include <petscdmda.h> 6 #include <petscfe.h> 7 #include <petscds.h> 8 #include <petscksp.h> 9 #include <petscsnes.h> 10 11 typedef struct { 12 /* Domain and mesh definition */ 13 PetscBool useDA; /* Flag DMDA tensor product mesh */ 14 PetscBool shearCoords; /* Flag for shear transform */ 15 PetscBool nonaffineCoords; /* Flag for non-affine transform */ 16 /* Element definition */ 17 PetscInt qorder; /* Order of the quadrature */ 18 PetscInt numComponents; /* Number of field components */ 19 PetscFE fe; /* The finite element */ 20 /* Testing space */ 21 PetscInt porder; /* Order of polynomials to test */ 22 PetscBool convergence; /* Test for order of convergence */ 23 PetscBool convRefine; /* Test for convergence using refinement, otherwise use coarsening */ 24 PetscBool constraints; /* Test local constraints */ 25 PetscBool tree; /* Test tree routines */ 26 PetscBool testFEjacobian; /* Test finite element Jacobian assembly */ 27 PetscBool testFVgrad; /* Test finite difference gradient routine */ 28 PetscBool testInjector; /* Test finite element injection routines */ 29 PetscInt treeCell; /* Cell to refine in tree test */ 30 PetscReal constants[3]; /* Constant values for each dimension */ 31 } AppCtx; 32 33 /* u = 1 */ 34 PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx) 35 { 36 AppCtx *user = (AppCtx *) ctx; 37 PetscInt d; 38 for (d = 0; d < dim; ++d) u[d] = user->constants[d]; 39 return 0; 40 } 41 PetscErrorCode constantDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) 42 { 43 PetscInt d; 44 for (d = 0; d < dim; ++d) u[d] = 0.0; 45 return 0; 46 } 47 48 /* u = x */ 49 PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx) 50 { 51 PetscInt d; 52 for (d = 0; d < dim; ++d) u[d] = coords[d]; 53 return 0; 54 } 55 PetscErrorCode linearDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) 56 { 57 PetscInt d, e; 58 for (d = 0; d < dim; ++d) { 59 u[d] = 0.0; 60 for (e = 0; e < dim; ++e) u[d] += (d == e ? 1.0 : 0.0) * n[e]; 61 } 62 return 0; 63 } 64 65 /* u = x^2 or u = (x^2, xy) or u = (xy, yz, zx) */ 66 PetscErrorCode quadratic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx) 67 { 68 if (dim > 2) {u[0] = coords[0]*coords[1]; u[1] = coords[1]*coords[2]; u[2] = coords[2]*coords[0];} 69 else if (dim > 1) {u[0] = coords[0]*coords[0]; u[1] = coords[0]*coords[1];} 70 else if (dim > 0) {u[0] = coords[0]*coords[0];} 71 return 0; 72 } 73 PetscErrorCode quadraticDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) 74 { 75 if (dim > 2) {u[0] = coords[1]*n[0] + coords[0]*n[1]; u[1] = coords[2]*n[1] + coords[1]*n[2]; u[2] = coords[2]*n[0] + coords[0]*n[2];} 76 else if (dim > 1) {u[0] = 2.0*coords[0]*n[0]; u[1] = coords[1]*n[0] + coords[0]*n[1];} 77 else if (dim > 0) {u[0] = 2.0*coords[0]*n[0];} 78 return 0; 79 } 80 81 /* u = x^3 or u = (x^3, x^2y) or u = (x^2y, y^2z, z^2x) */ 82 PetscErrorCode cubic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx) 83 { 84 if (dim > 2) {u[0] = coords[0]*coords[0]*coords[1]; u[1] = coords[1]*coords[1]*coords[2]; u[2] = coords[2]*coords[2]*coords[0];} 85 else if (dim > 1) {u[0] = coords[0]*coords[0]*coords[0]; u[1] = coords[0]*coords[0]*coords[1];} 86 else if (dim > 0) {u[0] = coords[0]*coords[0]*coords[0];} 87 return 0; 88 } 89 PetscErrorCode cubicDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) 90 { 91 if (dim > 2) {u[0] = 2.0*coords[0]*coords[1]*n[0] + coords[0]*coords[0]*n[1]; u[1] = 2.0*coords[1]*coords[2]*n[1] + coords[1]*coords[1]*n[2]; u[2] = 2.0*coords[2]*coords[0]*n[2] + coords[2]*coords[2]*n[0];} 92 else if (dim > 1) {u[0] = 3.0*coords[0]*coords[0]*n[0]; u[1] = 2.0*coords[0]*coords[1]*n[0] + coords[0]*coords[0]*n[1];} 93 else if (dim > 0) {u[0] = 3.0*coords[0]*coords[0]*n[0];} 94 return 0; 95 } 96 97 /* u = tanh(x) */ 98 PetscErrorCode trig(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx) 99 { 100 PetscInt d; 101 for (d = 0; d < dim; ++d) u[d] = PetscTanhReal(coords[d] - 0.5); 102 return 0; 103 } 104 PetscErrorCode trigDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) 105 { 106 PetscInt d; 107 for (d = 0; d < dim; ++d) u[d] = 1.0/PetscSqr(PetscCoshReal(coords[d] - 0.5)) * n[d]; 108 return 0; 109 } 110 111 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 112 { 113 PetscInt n = 3; 114 PetscErrorCode ierr; 115 116 PetscFunctionBeginUser; 117 options->useDA = PETSC_FALSE; 118 options->shearCoords = PETSC_FALSE; 119 options->nonaffineCoords = PETSC_FALSE; 120 options->qorder = 0; 121 options->numComponents = PETSC_DEFAULT; 122 options->porder = 0; 123 options->convergence = PETSC_FALSE; 124 options->convRefine = PETSC_TRUE; 125 options->constraints = PETSC_FALSE; 126 options->tree = PETSC_FALSE; 127 options->treeCell = 0; 128 options->testFEjacobian = PETSC_FALSE; 129 options->testFVgrad = PETSC_FALSE; 130 options->testInjector = PETSC_FALSE; 131 options->constants[0] = 1.0; 132 options->constants[1] = 2.0; 133 options->constants[2] = 3.0; 134 135 ierr = PetscOptionsBegin(comm, "", "Projection Test Options", "DMPlex");CHKERRQ(ierr); 136 CHKERRQ(PetscOptionsBool("-use_da", "Flag for DMDA mesh", "ex3.c", options->useDA, &options->useDA, NULL)); 137 CHKERRQ(PetscOptionsBool("-shear_coords", "Transform coordinates with a shear", "ex3.c", options->shearCoords, &options->shearCoords, NULL)); 138 CHKERRQ(PetscOptionsBool("-non_affine_coords", "Transform coordinates with a non-affine transform", "ex3.c", options->nonaffineCoords, &options->nonaffineCoords, NULL)); 139 CHKERRQ(PetscOptionsBoundedInt("-qorder", "The quadrature order", "ex3.c", options->qorder, &options->qorder, NULL,0)); 140 CHKERRQ(PetscOptionsBoundedInt("-num_comp", "The number of field components", "ex3.c", options->numComponents, &options->numComponents, NULL,PETSC_DEFAULT)); 141 CHKERRQ(PetscOptionsBoundedInt("-porder", "The order of polynomials to test", "ex3.c", options->porder, &options->porder, NULL,0)); 142 CHKERRQ(PetscOptionsBool("-convergence", "Check the convergence rate", "ex3.c", options->convergence, &options->convergence, NULL)); 143 CHKERRQ(PetscOptionsBool("-conv_refine", "Use refinement for the convergence rate", "ex3.c", options->convRefine, &options->convRefine, NULL)); 144 CHKERRQ(PetscOptionsBool("-constraints", "Test local constraints (serial only)", "ex3.c", options->constraints, &options->constraints, NULL)); 145 CHKERRQ(PetscOptionsBool("-tree", "Test tree routines", "ex3.c", options->tree, &options->tree, NULL)); 146 CHKERRQ(PetscOptionsBoundedInt("-tree_cell", "cell to refine in tree test", "ex3.c", options->treeCell, &options->treeCell, NULL,0)); 147 CHKERRQ(PetscOptionsBool("-test_fe_jacobian", "Test finite element Jacobian assembly", "ex3.c", options->testFEjacobian, &options->testFEjacobian, NULL)); 148 CHKERRQ(PetscOptionsBool("-test_fv_grad", "Test finite volume gradient reconstruction", "ex3.c", options->testFVgrad, &options->testFVgrad, NULL)); 149 CHKERRQ(PetscOptionsBool("-test_injector","Test finite element injection", "ex3.c", options->testInjector, &options->testInjector,NULL)); 150 CHKERRQ(PetscOptionsRealArray("-constants","Set the constant values", "ex3.c", options->constants, &n,NULL)); 151 ierr = PetscOptionsEnd();CHKERRQ(ierr); 152 153 PetscFunctionReturn(0); 154 } 155 156 static PetscErrorCode TransformCoordinates(DM dm, AppCtx *user) 157 { 158 PetscSection coordSection; 159 Vec coordinates; 160 PetscScalar *coords; 161 PetscInt vStart, vEnd, v; 162 163 PetscFunctionBeginUser; 164 if (user->nonaffineCoords) { 165 /* x' = r^(1/p) (x/r), y' = r^(1/p) (y/r), z' = z */ 166 CHKERRQ(DMGetCoordinateSection(dm, &coordSection)); 167 CHKERRQ(DMGetCoordinatesLocal(dm, &coordinates)); 168 CHKERRQ(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 169 CHKERRQ(VecGetArray(coordinates, &coords)); 170 for (v = vStart; v < vEnd; ++v) { 171 PetscInt dof, off; 172 PetscReal p = 4.0, r; 173 174 CHKERRQ(PetscSectionGetDof(coordSection, v, &dof)); 175 CHKERRQ(PetscSectionGetOffset(coordSection, v, &off)); 176 switch (dof) { 177 case 2: 178 r = PetscSqr(PetscRealPart(coords[off+0])) + PetscSqr(PetscRealPart(coords[off+1])); 179 coords[off+0] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p)/(2*p))*coords[off+0]; 180 coords[off+1] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p)/(2*p))*coords[off+1]; 181 break; 182 case 3: 183 r = PetscSqr(PetscRealPart(coords[off+0])) + PetscSqr(PetscRealPart(coords[off+1])); 184 coords[off+0] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p)/(2*p))*coords[off+0]; 185 coords[off+1] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p)/(2*p))*coords[off+1]; 186 coords[off+2] = coords[off+2]; 187 break; 188 } 189 } 190 CHKERRQ(VecRestoreArray(coordinates, &coords)); 191 } 192 if (user->shearCoords) { 193 /* x' = x + m y + m z, y' = y + m z, z' = z */ 194 CHKERRQ(DMGetCoordinateSection(dm, &coordSection)); 195 CHKERRQ(DMGetCoordinatesLocal(dm, &coordinates)); 196 CHKERRQ(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 197 CHKERRQ(VecGetArray(coordinates, &coords)); 198 for (v = vStart; v < vEnd; ++v) { 199 PetscInt dof, off; 200 PetscReal m = 1.0; 201 202 CHKERRQ(PetscSectionGetDof(coordSection, v, &dof)); 203 CHKERRQ(PetscSectionGetOffset(coordSection, v, &off)); 204 switch (dof) { 205 case 2: 206 coords[off+0] = coords[off+0] + m*coords[off+1]; 207 coords[off+1] = coords[off+1]; 208 break; 209 case 3: 210 coords[off+0] = coords[off+0] + m*coords[off+1] + m*coords[off+2]; 211 coords[off+1] = coords[off+1] + m*coords[off+2]; 212 coords[off+2] = coords[off+2]; 213 break; 214 } 215 } 216 CHKERRQ(VecRestoreArray(coordinates, &coords)); 217 } 218 PetscFunctionReturn(0); 219 } 220 221 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 222 { 223 PetscInt dim = 2; 224 PetscBool simplex; 225 226 PetscFunctionBeginUser; 227 if (user->useDA) { 228 switch (dim) { 229 case 2: 230 CHKERRQ(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 2, 2, PETSC_DETERMINE, PETSC_DETERMINE, 1, 1, NULL, NULL, dm)); 231 CHKERRQ(DMSetFromOptions(*dm)); 232 CHKERRQ(DMSetUp(*dm)); 233 CHKERRQ(DMDASetVertexCoordinates(*dm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 234 break; 235 default: 236 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot create structured mesh of dimension %d", dim); 237 } 238 CHKERRQ(PetscObjectSetName((PetscObject) *dm, "Hexahedral Mesh")); 239 } else { 240 CHKERRQ(DMCreate(comm, dm)); 241 CHKERRQ(DMSetType(*dm, DMPLEX)); 242 CHKERRQ(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 243 CHKERRQ(DMSetFromOptions(*dm)); 244 245 CHKERRQ(DMGetDimension(*dm, &dim)); 246 CHKERRQ(DMPlexIsSimplex(*dm, &simplex)); 247 CHKERRMPI(MPI_Bcast(&simplex, 1, MPIU_BOOL, 0, comm)); 248 if (user->tree) { 249 DM refTree, ncdm = NULL; 250 251 CHKERRQ(DMPlexCreateDefaultReferenceTree(comm,dim,simplex,&refTree)); 252 CHKERRQ(DMViewFromOptions(refTree,NULL,"-reftree_dm_view")); 253 CHKERRQ(DMPlexSetReferenceTree(*dm,refTree)); 254 CHKERRQ(DMDestroy(&refTree)); 255 CHKERRQ(DMPlexTreeRefineCell(*dm,user->treeCell,&ncdm)); 256 if (ncdm) { 257 CHKERRQ(DMDestroy(dm)); 258 *dm = ncdm; 259 CHKERRQ(DMPlexSetRefinementUniform(*dm, PETSC_FALSE)); 260 } 261 CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) *dm, "tree_")); 262 CHKERRQ(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 263 CHKERRQ(DMSetFromOptions(*dm)); 264 CHKERRQ(DMViewFromOptions(*dm,NULL,"-dm_view")); 265 } else { 266 CHKERRQ(DMPlexSetRefinementUniform(*dm, PETSC_TRUE)); 267 } 268 CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) *dm, "dist_")); 269 CHKERRQ(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 270 CHKERRQ(DMSetFromOptions(*dm)); 271 CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) *dm, NULL)); 272 if (simplex) CHKERRQ(PetscObjectSetName((PetscObject) *dm, "Simplicial Mesh")); 273 else CHKERRQ(PetscObjectSetName((PetscObject) *dm, "Hexahedral Mesh")); 274 } 275 CHKERRQ(DMSetFromOptions(*dm)); 276 CHKERRQ(TransformCoordinates(*dm, user)); 277 CHKERRQ(DMViewFromOptions(*dm,NULL,"-dm_view")); 278 PetscFunctionReturn(0); 279 } 280 281 static void simple_mass(PetscInt dim, PetscInt Nf, PetscInt NfAux, 282 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 283 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 284 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 285 { 286 PetscInt d, e; 287 for (d = 0, e = 0; d < dim; d++, e+=dim+1) { 288 g0[e] = 1.; 289 } 290 } 291 292 /* < \nabla v, 1/2(\nabla u + {\nabla u}^T) > */ 293 static void symmetric_gradient_inner_product(PetscInt dim, PetscInt Nf, PetscInt NfAux, 294 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 295 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 296 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar C[]) 297 { 298 PetscInt compI, compJ, d, e; 299 300 for (compI = 0; compI < dim; ++compI) { 301 for (compJ = 0; compJ < dim; ++compJ) { 302 for (d = 0; d < dim; ++d) { 303 for (e = 0; e < dim; e++) { 304 if (d == e && d == compI && d == compJ) { 305 C[((compI*dim+compJ)*dim+d)*dim+e] = 1.0; 306 } else if ((d == compJ && e == compI) || (d == e && compI == compJ)) { 307 C[((compI*dim+compJ)*dim+d)*dim+e] = 0.5; 308 } else { 309 C[((compI*dim+compJ)*dim+d)*dim+e] = 0.0; 310 } 311 } 312 } 313 } 314 } 315 } 316 317 static PetscErrorCode SetupSection(DM dm, AppCtx *user) 318 { 319 PetscFunctionBeginUser; 320 if (user->constraints) { 321 /* test local constraints */ 322 DM coordDM; 323 PetscInt fStart, fEnd, f, vStart, vEnd, v; 324 PetscInt edgesx = 2, vertsx; 325 PetscInt edgesy = 2, vertsy; 326 PetscMPIInt size; 327 PetscInt numConst; 328 PetscSection aSec; 329 PetscInt *anchors; 330 PetscInt offset; 331 IS aIS; 332 MPI_Comm comm = PetscObjectComm((PetscObject)dm); 333 334 CHKERRMPI(MPI_Comm_size(comm,&size)); 335 PetscCheckFalse(size > 1,comm,PETSC_ERR_SUP,"Local constraint test can only be performed in serial"); 336 337 /* we are going to test constraints by using them to enforce periodicity 338 * in one direction, and comparing to the existing method of enforcing 339 * periodicity */ 340 341 /* first create the coordinate section so that it does not clone the 342 * constraints */ 343 CHKERRQ(DMGetCoordinateDM(dm,&coordDM)); 344 345 /* create the constrained-to-anchor section */ 346 CHKERRQ(DMPlexGetDepthStratum(dm,0,&vStart,&vEnd)); 347 CHKERRQ(DMPlexGetDepthStratum(dm,1,&fStart,&fEnd)); 348 CHKERRQ(PetscSectionCreate(PETSC_COMM_SELF,&aSec)); 349 CHKERRQ(PetscSectionSetChart(aSec,PetscMin(fStart,vStart),PetscMax(fEnd,vEnd))); 350 351 /* define the constraints */ 352 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-da_grid_x",&edgesx,NULL)); 353 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-da_grid_y",&edgesy,NULL)); 354 vertsx = edgesx + 1; 355 vertsy = edgesy + 1; 356 numConst = vertsy + edgesy; 357 CHKERRQ(PetscMalloc1(numConst,&anchors)); 358 offset = 0; 359 for (v = vStart + edgesx; v < vEnd; v+= vertsx) { 360 CHKERRQ(PetscSectionSetDof(aSec,v,1)); 361 anchors[offset++] = v - edgesx; 362 } 363 for (f = fStart + edgesx * vertsy + edgesx * edgesy; f < fEnd; f++) { 364 CHKERRQ(PetscSectionSetDof(aSec,f,1)); 365 anchors[offset++] = f - edgesx * edgesy; 366 } 367 CHKERRQ(PetscSectionSetUp(aSec)); 368 CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,numConst,anchors,PETSC_OWN_POINTER,&aIS)); 369 370 CHKERRQ(DMPlexSetAnchors(dm,aSec,aIS)); 371 CHKERRQ(PetscSectionDestroy(&aSec)); 372 CHKERRQ(ISDestroy(&aIS)); 373 } 374 CHKERRQ(DMSetNumFields(dm, 1)); 375 CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) user->fe)); 376 CHKERRQ(DMCreateDS(dm)); 377 if (user->constraints) { 378 /* test getting local constraint matrix that matches section */ 379 PetscSection aSec; 380 IS aIS; 381 382 CHKERRQ(DMPlexGetAnchors(dm,&aSec,&aIS)); 383 if (aSec) { 384 PetscDS ds; 385 PetscSection cSec, section; 386 PetscInt cStart, cEnd, c, numComp; 387 Mat cMat, mass; 388 Vec local; 389 const PetscInt *anchors; 390 391 CHKERRQ(DMGetLocalSection(dm,§ion)); 392 /* this creates the matrix and preallocates the matrix structure: we 393 * just have to fill in the values */ 394 CHKERRQ(DMGetDefaultConstraints(dm,&cSec,&cMat,NULL)); 395 CHKERRQ(PetscSectionGetChart(cSec,&cStart,&cEnd)); 396 CHKERRQ(ISGetIndices(aIS,&anchors)); 397 CHKERRQ(PetscFEGetNumComponents(user->fe, &numComp)); 398 for (c = cStart; c < cEnd; c++) { 399 PetscInt cDof; 400 401 /* is this point constrained? (does it have an anchor?) */ 402 CHKERRQ(PetscSectionGetDof(aSec,c,&cDof)); 403 if (cDof) { 404 PetscInt cOff, a, aDof, aOff, j; 405 PetscCheckFalse(cDof != 1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Found %d anchor points: should be just one",cDof); 406 407 /* find the anchor point */ 408 CHKERRQ(PetscSectionGetOffset(aSec,c,&cOff)); 409 a = anchors[cOff]; 410 411 /* find the constrained dofs (row in constraint matrix) */ 412 CHKERRQ(PetscSectionGetDof(cSec,c,&cDof)); 413 CHKERRQ(PetscSectionGetOffset(cSec,c,&cOff)); 414 415 /* find the anchor dofs (column in constraint matrix) */ 416 CHKERRQ(PetscSectionGetDof(section,a,&aDof)); 417 CHKERRQ(PetscSectionGetOffset(section,a,&aOff)); 418 419 PetscCheckFalse(cDof != aDof,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Point and anchor have different number of dofs: %d, %d",cDof,aDof); 420 PetscCheckFalse(cDof % numComp,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Point dofs not divisible by field components: %d, %d",cDof,numComp); 421 422 /* put in a simple equality constraint */ 423 for (j = 0; j < cDof; j++) { 424 CHKERRQ(MatSetValue(cMat,cOff+j,aOff+j,1.,INSERT_VALUES)); 425 } 426 } 427 } 428 CHKERRQ(MatAssemblyBegin(cMat,MAT_FINAL_ASSEMBLY)); 429 CHKERRQ(MatAssemblyEnd(cMat,MAT_FINAL_ASSEMBLY)); 430 CHKERRQ(ISRestoreIndices(aIS,&anchors)); 431 432 /* Now that we have constructed the constraint matrix, any FE matrix 433 * that we construct will apply the constraints during construction */ 434 435 CHKERRQ(DMCreateMatrix(dm,&mass)); 436 /* get a dummy local variable to serve as the solution */ 437 CHKERRQ(DMGetLocalVector(dm,&local)); 438 CHKERRQ(DMGetDS(dm,&ds)); 439 /* set the jacobian to be the mass matrix */ 440 CHKERRQ(PetscDSSetJacobian(ds, 0, 0, simple_mass, NULL, NULL, NULL)); 441 /* build the mass matrix */ 442 CHKERRQ(DMPlexSNESComputeJacobianFEM(dm,local,mass,mass,NULL)); 443 CHKERRQ(MatView(mass,PETSC_VIEWER_STDOUT_WORLD)); 444 CHKERRQ(MatDestroy(&mass)); 445 CHKERRQ(DMRestoreLocalVector(dm,&local)); 446 } 447 } 448 PetscFunctionReturn(0); 449 } 450 451 static PetscErrorCode TestFEJacobian(DM dm, AppCtx *user) 452 { 453 PetscFunctionBeginUser; 454 if (!user->useDA) { 455 Vec local; 456 const Vec *vecs; 457 Mat E; 458 MatNullSpace sp; 459 PetscBool isNullSpace, hasConst; 460 PetscInt dim, n, i; 461 Vec res = NULL, localX, localRes; 462 PetscDS ds; 463 464 CHKERRQ(DMGetDimension(dm, &dim)); 465 PetscCheckFalse(user->numComponents != dim,PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "The number of components %d must be equal to the dimension %d for this test", user->numComponents, dim); 466 CHKERRQ(DMGetDS(dm,&ds)); 467 CHKERRQ(PetscDSSetJacobian(ds,0,0,NULL,NULL,NULL,symmetric_gradient_inner_product)); 468 CHKERRQ(DMCreateMatrix(dm,&E)); 469 CHKERRQ(DMGetLocalVector(dm,&local)); 470 CHKERRQ(DMPlexSNESComputeJacobianFEM(dm,local,E,E,NULL)); 471 CHKERRQ(DMPlexCreateRigidBody(dm,0,&sp)); 472 CHKERRQ(MatNullSpaceGetVecs(sp,&hasConst,&n,&vecs)); 473 if (n) CHKERRQ(VecDuplicate(vecs[0],&res)); 474 CHKERRQ(DMCreateLocalVector(dm,&localX)); 475 CHKERRQ(DMCreateLocalVector(dm,&localRes)); 476 for (i = 0; i < n; i++) { /* also test via matrix-free Jacobian application */ 477 PetscReal resNorm; 478 479 CHKERRQ(VecSet(localRes,0.)); 480 CHKERRQ(VecSet(localX,0.)); 481 CHKERRQ(VecSet(local,0.)); 482 CHKERRQ(VecSet(res,0.)); 483 CHKERRQ(DMGlobalToLocalBegin(dm,vecs[i],INSERT_VALUES,localX)); 484 CHKERRQ(DMGlobalToLocalEnd(dm,vecs[i],INSERT_VALUES,localX)); 485 CHKERRQ(DMSNESComputeJacobianAction(dm,local,localX,localRes,NULL)); 486 CHKERRQ(DMLocalToGlobalBegin(dm,localRes,ADD_VALUES,res)); 487 CHKERRQ(DMLocalToGlobalEnd(dm,localRes,ADD_VALUES,res)); 488 CHKERRQ(VecNorm(res,NORM_2,&resNorm)); 489 if (resNorm > PETSC_SMALL) { 490 CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)dm),"Symmetric gradient action null space vector %D residual: %E\n",i,resNorm)); 491 } 492 } 493 CHKERRQ(VecDestroy(&localRes)); 494 CHKERRQ(VecDestroy(&localX)); 495 CHKERRQ(VecDestroy(&res)); 496 CHKERRQ(MatNullSpaceTest(sp,E,&isNullSpace)); 497 if (isNullSpace) { 498 CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)dm),"Symmetric gradient null space: PASS\n")); 499 } else { 500 CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)dm),"Symmetric gradient null space: FAIL\n")); 501 } 502 CHKERRQ(MatNullSpaceDestroy(&sp)); 503 CHKERRQ(MatDestroy(&E)); 504 CHKERRQ(DMRestoreLocalVector(dm,&local)); 505 } 506 PetscFunctionReturn(0); 507 } 508 509 static PetscErrorCode TestInjector(DM dm, AppCtx *user) 510 { 511 DM refTree; 512 PetscMPIInt rank; 513 514 PetscFunctionBegin; 515 CHKERRQ(DMPlexGetReferenceTree(dm,&refTree)); 516 if (refTree) { 517 Mat inj; 518 519 CHKERRQ(DMPlexComputeInjectorReferenceTree(refTree,&inj)); 520 CHKERRQ(PetscObjectSetName((PetscObject)inj,"Reference Tree Injector")); 521 CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 522 if (rank == 0) { 523 CHKERRQ(MatView(inj,PETSC_VIEWER_STDOUT_SELF)); 524 } 525 CHKERRQ(MatDestroy(&inj)); 526 } 527 PetscFunctionReturn(0); 528 } 529 530 static PetscErrorCode TestFVGrad(DM dm, AppCtx *user) 531 { 532 MPI_Comm comm; 533 DM dmRedist, dmfv, dmgrad, dmCell, refTree; 534 PetscFV fv; 535 PetscInt dim, nvecs, v, cStart, cEnd, cEndInterior; 536 PetscMPIInt size; 537 Vec cellgeom, grad, locGrad; 538 const PetscScalar *cgeom; 539 PetscReal allVecMaxDiff = 0., fvTol = 100. * PETSC_MACHINE_EPSILON; 540 541 PetscFunctionBeginUser; 542 comm = PetscObjectComm((PetscObject)dm); 543 CHKERRQ(DMGetDimension(dm, &dim)); 544 /* duplicate DM, give dup. a FV discretization */ 545 CHKERRQ(DMSetBasicAdjacency(dm,PETSC_TRUE,PETSC_FALSE)); 546 CHKERRMPI(MPI_Comm_size(comm,&size)); 547 dmRedist = NULL; 548 if (size > 1) { 549 CHKERRQ(DMPlexDistributeOverlap(dm,1,NULL,&dmRedist)); 550 } 551 if (!dmRedist) { 552 dmRedist = dm; 553 CHKERRQ(PetscObjectReference((PetscObject)dmRedist)); 554 } 555 CHKERRQ(PetscFVCreate(comm,&fv)); 556 CHKERRQ(PetscFVSetType(fv,PETSCFVLEASTSQUARES)); 557 CHKERRQ(PetscFVSetNumComponents(fv,user->numComponents)); 558 CHKERRQ(PetscFVSetSpatialDimension(fv,dim)); 559 CHKERRQ(PetscFVSetFromOptions(fv)); 560 CHKERRQ(PetscFVSetUp(fv)); 561 CHKERRQ(DMPlexConstructGhostCells(dmRedist,NULL,NULL,&dmfv)); 562 CHKERRQ(DMDestroy(&dmRedist)); 563 CHKERRQ(DMSetNumFields(dmfv,1)); 564 CHKERRQ(DMSetField(dmfv, 0, NULL, (PetscObject) fv)); 565 CHKERRQ(DMCreateDS(dmfv)); 566 CHKERRQ(DMPlexGetReferenceTree(dm,&refTree)); 567 if (refTree) CHKERRQ(DMCopyDisc(dmfv,refTree)); 568 CHKERRQ(DMPlexGetGradientDM(dmfv, fv, &dmgrad)); 569 CHKERRQ(DMPlexGetHeightStratum(dmfv,0,&cStart,&cEnd)); 570 nvecs = dim * (dim+1) / 2; 571 CHKERRQ(DMPlexGetGeometryFVM(dmfv,NULL,&cellgeom,NULL)); 572 CHKERRQ(VecGetDM(cellgeom,&dmCell)); 573 CHKERRQ(VecGetArrayRead(cellgeom,&cgeom)); 574 CHKERRQ(DMGetGlobalVector(dmgrad,&grad)); 575 CHKERRQ(DMGetLocalVector(dmgrad,&locGrad)); 576 CHKERRQ(DMPlexGetGhostCellStratum(dmgrad,&cEndInterior,NULL)); 577 cEndInterior = (cEndInterior < 0) ? cEnd: cEndInterior; 578 for (v = 0; v < nvecs; v++) { 579 Vec locX; 580 PetscInt c; 581 PetscScalar trueGrad[3][3] = {{0.}}; 582 const PetscScalar *gradArray; 583 PetscReal maxDiff, maxDiffGlob; 584 585 CHKERRQ(DMGetLocalVector(dmfv,&locX)); 586 /* get the local projection of the rigid body mode */ 587 for (c = cStart; c < cEnd; c++) { 588 PetscFVCellGeom *cg; 589 PetscScalar cx[3] = {0.,0.,0.}; 590 591 CHKERRQ(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); 592 if (v < dim) { 593 cx[v] = 1.; 594 } else { 595 PetscInt w = v - dim; 596 597 cx[(w + 1) % dim] = cg->centroid[(w + 2) % dim]; 598 cx[(w + 2) % dim] = -cg->centroid[(w + 1) % dim]; 599 } 600 CHKERRQ(DMPlexVecSetClosure(dmfv,NULL,locX,c,cx,INSERT_ALL_VALUES)); 601 } 602 /* TODO: this isn't in any header */ 603 CHKERRQ(DMPlexReconstructGradientsFVM(dmfv,locX,grad)); 604 CHKERRQ(DMGlobalToLocalBegin(dmgrad,grad,INSERT_VALUES,locGrad)); 605 CHKERRQ(DMGlobalToLocalEnd(dmgrad,grad,INSERT_VALUES,locGrad)); 606 CHKERRQ(VecGetArrayRead(locGrad,&gradArray)); 607 /* compare computed gradient to exact gradient */ 608 if (v >= dim) { 609 PetscInt w = v - dim; 610 611 trueGrad[(w + 1) % dim][(w + 2) % dim] = 1.; 612 trueGrad[(w + 2) % dim][(w + 1) % dim] = -1.; 613 } 614 maxDiff = 0.; 615 for (c = cStart; c < cEndInterior; c++) { 616 PetscScalar *compGrad; 617 PetscInt i, j, k; 618 PetscReal FrobDiff = 0.; 619 620 CHKERRQ(DMPlexPointLocalRead(dmgrad, c, gradArray, &compGrad)); 621 622 for (i = 0, k = 0; i < dim; i++) { 623 for (j = 0; j < dim; j++, k++) { 624 PetscScalar diff = compGrad[k] - trueGrad[i][j]; 625 FrobDiff += PetscRealPart(diff * PetscConj(diff)); 626 } 627 } 628 FrobDiff = PetscSqrtReal(FrobDiff); 629 maxDiff = PetscMax(maxDiff,FrobDiff); 630 } 631 CHKERRMPI(MPI_Allreduce(&maxDiff,&maxDiffGlob,1,MPIU_REAL,MPIU_MAX,comm)); 632 allVecMaxDiff = PetscMax(allVecMaxDiff,maxDiffGlob); 633 CHKERRQ(VecRestoreArrayRead(locGrad,&gradArray)); 634 CHKERRQ(DMRestoreLocalVector(dmfv,&locX)); 635 } 636 if (allVecMaxDiff < fvTol) { 637 CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)dm),"Finite volume gradient reconstruction: PASS\n")); 638 } else { 639 CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)dm),"Finite volume gradient reconstruction: FAIL at tolerance %g with max difference %g\n",fvTol,allVecMaxDiff)); 640 } 641 CHKERRQ(DMRestoreLocalVector(dmgrad,&locGrad)); 642 CHKERRQ(DMRestoreGlobalVector(dmgrad,&grad)); 643 CHKERRQ(VecRestoreArrayRead(cellgeom,&cgeom)); 644 CHKERRQ(DMDestroy(&dmfv)); 645 CHKERRQ(PetscFVDestroy(&fv)); 646 PetscFunctionReturn(0); 647 } 648 649 static PetscErrorCode ComputeError(DM dm, PetscErrorCode (**exactFuncs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), 650 PetscErrorCode (**exactFuncDers)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), 651 void **exactCtxs, PetscReal *error, PetscReal *errorDer, AppCtx *user) 652 { 653 Vec u; 654 PetscReal n[3] = {1.0, 1.0, 1.0}; 655 656 PetscFunctionBeginUser; 657 CHKERRQ(DMGetGlobalVector(dm, &u)); 658 /* Project function into FE function space */ 659 CHKERRQ(DMProjectFunction(dm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, u)); 660 CHKERRQ(VecViewFromOptions(u, NULL, "-projection_view")); 661 /* Compare approximation to exact in L_2 */ 662 CHKERRQ(DMComputeL2Diff(dm, 0.0, exactFuncs, exactCtxs, u, error)); 663 CHKERRQ(DMComputeL2GradientDiff(dm, 0.0, exactFuncDers, exactCtxs, u, n, errorDer)); 664 CHKERRQ(DMRestoreGlobalVector(dm, &u)); 665 PetscFunctionReturn(0); 666 } 667 668 static PetscErrorCode CheckFunctions(DM dm, PetscInt order, AppCtx *user) 669 { 670 PetscErrorCode (*exactFuncs[1]) (PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 671 PetscErrorCode (*exactFuncDers[1]) (PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx); 672 void *exactCtxs[3]; 673 MPI_Comm comm; 674 PetscReal error, errorDer, tol = PETSC_SMALL; 675 676 PetscFunctionBeginUser; 677 exactCtxs[0] = user; 678 exactCtxs[1] = user; 679 exactCtxs[2] = user; 680 CHKERRQ(PetscObjectGetComm((PetscObject)dm, &comm)); 681 /* Setup functions to approximate */ 682 switch (order) { 683 case 0: 684 exactFuncs[0] = constant; 685 exactFuncDers[0] = constantDer; 686 break; 687 case 1: 688 exactFuncs[0] = linear; 689 exactFuncDers[0] = linearDer; 690 break; 691 case 2: 692 exactFuncs[0] = quadratic; 693 exactFuncDers[0] = quadraticDer; 694 break; 695 case 3: 696 exactFuncs[0] = cubic; 697 exactFuncDers[0] = cubicDer; 698 break; 699 default: 700 SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Could not determine functions to test for order %d", order); 701 } 702 CHKERRQ(ComputeError(dm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user)); 703 /* Report result */ 704 if (error > tol) CHKERRQ(PetscPrintf(comm, "Function tests FAIL for order %D at tolerance %g error %g\n", order, (double)tol,(double) error)); 705 else CHKERRQ(PetscPrintf(comm, "Function tests pass for order %D at tolerance %g\n", order, (double)tol)); 706 if (errorDer > tol) CHKERRQ(PetscPrintf(comm, "Function tests FAIL for order %D derivatives at tolerance %g error %g\n", order, (double)tol, (double)errorDer)); 707 else CHKERRQ(PetscPrintf(comm, "Function tests pass for order %D derivatives at tolerance %g\n", order, (double)tol)); 708 PetscFunctionReturn(0); 709 } 710 711 static PetscErrorCode CheckInterpolation(DM dm, PetscBool checkRestrict, PetscInt order, AppCtx *user) 712 { 713 PetscErrorCode (*exactFuncs[1]) (PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx); 714 PetscErrorCode (*exactFuncDers[1]) (PetscInt, PetscReal, const PetscReal x[], const PetscReal n[], PetscInt, PetscScalar *u, void *ctx); 715 PetscReal n[3] = {1.0, 1.0, 1.0}; 716 void *exactCtxs[3]; 717 DM rdm, idm, fdm; 718 Mat Interp; 719 Vec iu, fu, scaling; 720 MPI_Comm comm; 721 PetscInt dim; 722 PetscReal error, errorDer, tol = PETSC_SMALL; 723 724 PetscFunctionBeginUser; 725 exactCtxs[0] = user; 726 exactCtxs[1] = user; 727 exactCtxs[2] = user; 728 CHKERRQ(PetscObjectGetComm((PetscObject)dm,&comm)); 729 CHKERRQ(DMGetDimension(dm, &dim)); 730 CHKERRQ(DMRefine(dm, comm, &rdm)); 731 CHKERRQ(DMSetCoarseDM(rdm, dm)); 732 CHKERRQ(DMPlexSetRegularRefinement(rdm, user->convRefine)); 733 if (user->tree) { 734 DM refTree; 735 CHKERRQ(DMPlexGetReferenceTree(dm,&refTree)); 736 CHKERRQ(DMPlexSetReferenceTree(rdm,refTree)); 737 } 738 if (user->useDA) CHKERRQ(DMDASetVertexCoordinates(rdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 739 CHKERRQ(SetupSection(rdm, user)); 740 /* Setup functions to approximate */ 741 switch (order) { 742 case 0: 743 exactFuncs[0] = constant; 744 exactFuncDers[0] = constantDer; 745 break; 746 case 1: 747 exactFuncs[0] = linear; 748 exactFuncDers[0] = linearDer; 749 break; 750 case 2: 751 exactFuncs[0] = quadratic; 752 exactFuncDers[0] = quadraticDer; 753 break; 754 case 3: 755 exactFuncs[0] = cubic; 756 exactFuncDers[0] = cubicDer; 757 break; 758 default: 759 SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Could not determine functions to test for dimension %D order %D", dim, order); 760 } 761 idm = checkRestrict ? rdm : dm; 762 fdm = checkRestrict ? dm : rdm; 763 CHKERRQ(DMGetGlobalVector(idm, &iu)); 764 CHKERRQ(DMGetGlobalVector(fdm, &fu)); 765 CHKERRQ(DMSetApplicationContext(dm, user)); 766 CHKERRQ(DMSetApplicationContext(rdm, user)); 767 CHKERRQ(DMCreateInterpolation(dm, rdm, &Interp, &scaling)); 768 /* Project function into initial FE function space */ 769 CHKERRQ(DMProjectFunction(idm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, iu)); 770 /* Interpolate function into final FE function space */ 771 if (checkRestrict) {CHKERRQ(MatRestrict(Interp, iu, fu));CHKERRQ(VecPointwiseMult(fu, scaling, fu));} 772 else CHKERRQ(MatInterpolate(Interp, iu, fu)); 773 /* Compare approximation to exact in L_2 */ 774 CHKERRQ(DMComputeL2Diff(fdm, 0.0, exactFuncs, exactCtxs, fu, &error)); 775 CHKERRQ(DMComputeL2GradientDiff(fdm, 0.0, exactFuncDers, exactCtxs, fu, n, &errorDer)); 776 /* Report result */ 777 if (error > tol) CHKERRQ(PetscPrintf(comm, "Interpolation tests FAIL for order %D at tolerance %g error %g\n", order, (double)tol, (double)error)); 778 else CHKERRQ(PetscPrintf(comm, "Interpolation tests pass for order %D at tolerance %g\n", order, (double)tol)); 779 if (errorDer > tol) CHKERRQ(PetscPrintf(comm, "Interpolation tests FAIL for order %D derivatives at tolerance %g error %g\n", order, (double)tol, (double)errorDer)); 780 else CHKERRQ(PetscPrintf(comm, "Interpolation tests pass for order %D derivatives at tolerance %g\n", order, (double)tol)); 781 CHKERRQ(DMRestoreGlobalVector(idm, &iu)); 782 CHKERRQ(DMRestoreGlobalVector(fdm, &fu)); 783 CHKERRQ(MatDestroy(&Interp)); 784 CHKERRQ(VecDestroy(&scaling)); 785 CHKERRQ(DMDestroy(&rdm)); 786 PetscFunctionReturn(0); 787 } 788 789 static PetscErrorCode CheckConvergence(DM dm, PetscInt Nr, AppCtx *user) 790 { 791 DM odm = dm, rdm = NULL, cdm = NULL; 792 PetscErrorCode (*exactFuncs[1]) (PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) = {trig}; 793 PetscErrorCode (*exactFuncDers[1]) (PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) = {trigDer}; 794 void *exactCtxs[3]; 795 PetscInt r, c, cStart, cEnd; 796 PetscReal errorOld, errorDerOld, error, errorDer, rel, len, lenOld; 797 double p; 798 799 PetscFunctionBeginUser; 800 if (!user->convergence) PetscFunctionReturn(0); 801 exactCtxs[0] = user; 802 exactCtxs[1] = user; 803 exactCtxs[2] = user; 804 CHKERRQ(PetscObjectReference((PetscObject) odm)); 805 if (!user->convRefine) { 806 for (r = 0; r < Nr; ++r) { 807 CHKERRQ(DMRefine(odm, PetscObjectComm((PetscObject) dm), &rdm)); 808 CHKERRQ(DMDestroy(&odm)); 809 odm = rdm; 810 } 811 CHKERRQ(SetupSection(odm, user)); 812 } 813 CHKERRQ(ComputeError(odm, exactFuncs, exactFuncDers, exactCtxs, &errorOld, &errorDerOld, user)); 814 if (user->convRefine) { 815 for (r = 0; r < Nr; ++r) { 816 CHKERRQ(DMRefine(odm, PetscObjectComm((PetscObject) dm), &rdm)); 817 if (user->useDA) CHKERRQ(DMDASetVertexCoordinates(rdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 818 CHKERRQ(SetupSection(rdm, user)); 819 CHKERRQ(ComputeError(rdm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user)); 820 p = PetscLog2Real(errorOld/error); 821 CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) dm), "Function convergence rate at refinement %D: %.2f\n", r, (double)p)); 822 p = PetscLog2Real(errorDerOld/errorDer); 823 CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) dm), "Derivative convergence rate at refinement %D: %.2f\n", r, (double)p)); 824 CHKERRQ(DMDestroy(&odm)); 825 odm = rdm; 826 errorOld = error; 827 errorDerOld = errorDer; 828 } 829 } else { 830 /* CHKERRQ(ComputeLongestEdge(dm, &lenOld)); */ 831 CHKERRQ(DMPlexGetHeightStratum(odm, 0, &cStart, &cEnd)); 832 lenOld = cEnd - cStart; 833 for (c = 0; c < Nr; ++c) { 834 CHKERRQ(DMCoarsen(odm, PetscObjectComm((PetscObject) dm), &cdm)); 835 if (user->useDA) CHKERRQ(DMDASetVertexCoordinates(cdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 836 CHKERRQ(SetupSection(cdm, user)); 837 CHKERRQ(ComputeError(cdm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user)); 838 /* CHKERRQ(ComputeLongestEdge(cdm, &len)); */ 839 CHKERRQ(DMPlexGetHeightStratum(cdm, 0, &cStart, &cEnd)); 840 len = cEnd - cStart; 841 rel = error/errorOld; 842 p = PetscLogReal(rel) / PetscLogReal(lenOld / len); 843 CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) dm), "Function convergence rate at coarsening %D: %.2f\n", c, (double)p)); 844 rel = errorDer/errorDerOld; 845 p = PetscLogReal(rel) / PetscLogReal(lenOld / len); 846 CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) dm), "Derivative convergence rate at coarsening %D: %.2f\n", c, (double)p)); 847 CHKERRQ(DMDestroy(&odm)); 848 odm = cdm; 849 errorOld = error; 850 errorDerOld = errorDer; 851 lenOld = len; 852 } 853 } 854 CHKERRQ(DMDestroy(&odm)); 855 PetscFunctionReturn(0); 856 } 857 858 int main(int argc, char **argv) 859 { 860 DM dm; 861 AppCtx user; /* user-defined work context */ 862 PetscInt dim = 2; 863 PetscBool simplex = PETSC_FALSE; 864 PetscErrorCode ierr; 865 866 ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 867 CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user)); 868 CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 869 if (!user.useDA) { 870 CHKERRQ(DMGetDimension(dm, &dim)); 871 CHKERRQ(DMPlexIsSimplex(dm, &simplex)); 872 } 873 CHKERRQ(DMPlexMetricSetFromOptions(dm)); 874 user.numComponents = user.numComponents < 0 ? dim : user.numComponents; 875 CHKERRQ(PetscFECreateDefault(PETSC_COMM_WORLD, dim, user.numComponents, simplex, NULL, user.qorder, &user.fe)); 876 CHKERRQ(SetupSection(dm, &user)); 877 if (user.testFEjacobian) CHKERRQ(TestFEJacobian(dm, &user)); 878 if (user.testFVgrad) CHKERRQ(TestFVGrad(dm, &user)); 879 if (user.testInjector) CHKERRQ(TestInjector(dm, &user)); 880 CHKERRQ(CheckFunctions(dm, user.porder, &user)); 881 { 882 PetscDualSpace dsp; 883 PetscInt k; 884 885 CHKERRQ(PetscFEGetDualSpace(user.fe, &dsp)); 886 CHKERRQ(PetscDualSpaceGetDeRahm(dsp, &k)); 887 if (dim == 2 && simplex == PETSC_TRUE && user.tree == PETSC_FALSE && k == 0) { 888 CHKERRQ(CheckInterpolation(dm, PETSC_FALSE, user.porder, &user)); 889 CHKERRQ(CheckInterpolation(dm, PETSC_TRUE, user.porder, &user)); 890 } 891 } 892 CHKERRQ(CheckConvergence(dm, 3, &user)); 893 CHKERRQ(PetscFEDestroy(&user.fe)); 894 CHKERRQ(DMDestroy(&dm)); 895 ierr = PetscFinalize(); 896 return ierr; 897 } 898 899 /*TEST 900 901 test: 902 suffix: 1 903 requires: triangle 904 905 # 2D P_1 on a triangle 906 test: 907 suffix: p1_2d_0 908 requires: triangle 909 args: -petscspace_degree 1 -qorder 1 -convergence 910 test: 911 suffix: p1_2d_1 912 requires: triangle 913 args: -petscspace_degree 1 -qorder 1 -porder 1 914 test: 915 suffix: p1_2d_2 916 requires: triangle 917 args: -petscspace_degree 1 -qorder 1 -porder 2 918 test: 919 suffix: p1_2d_3 920 requires: triangle mmg 921 args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0 922 test: 923 suffix: p1_2d_4 924 requires: triangle mmg 925 args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0 926 test: 927 suffix: p1_2d_5 928 requires: triangle mmg 929 args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0 930 931 # 3D P_1 on a tetrahedron 932 test: 933 suffix: p1_3d_0 934 requires: ctetgen 935 args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -convergence 936 test: 937 suffix: p1_3d_1 938 requires: ctetgen 939 args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -porder 1 940 test: 941 suffix: p1_3d_2 942 requires: ctetgen 943 args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -porder 2 944 test: 945 suffix: p1_3d_3 946 requires: ctetgen mmg 947 args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0 948 test: 949 suffix: p1_3d_4 950 requires: ctetgen mmg 951 args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0 952 test: 953 suffix: p1_3d_5 954 requires: ctetgen mmg 955 args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0 956 957 # 2D P_2 on a triangle 958 test: 959 suffix: p2_2d_0 960 requires: triangle 961 args: -petscspace_degree 2 -qorder 2 -convergence 962 test: 963 suffix: p2_2d_1 964 requires: triangle 965 args: -petscspace_degree 2 -qorder 2 -porder 1 966 test: 967 suffix: p2_2d_2 968 requires: triangle 969 args: -petscspace_degree 2 -qorder 2 -porder 2 970 test: 971 suffix: p2_2d_3 972 requires: triangle mmg 973 args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -convergence -conv_refine 0 974 test: 975 suffix: p2_2d_4 976 requires: triangle mmg 977 args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 1 -conv_refine 0 978 test: 979 suffix: p2_2d_5 980 requires: triangle mmg 981 args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 2 -conv_refine 0 982 983 # 3D P_2 on a tetrahedron 984 test: 985 suffix: p2_3d_0 986 requires: ctetgen 987 args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -convergence 988 test: 989 suffix: p2_3d_1 990 requires: ctetgen 991 args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -porder 1 992 test: 993 suffix: p2_3d_2 994 requires: ctetgen 995 args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -porder 2 996 test: 997 suffix: p2_3d_3 998 requires: ctetgen mmg 999 args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -convergence -conv_refine 0 1000 test: 1001 suffix: p2_3d_4 1002 requires: ctetgen mmg 1003 args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 1 -conv_refine 0 1004 test: 1005 suffix: p2_3d_5 1006 requires: ctetgen mmg 1007 args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 2 -conv_refine 0 1008 1009 # 2D Q_1 on a quadrilaterial DA 1010 test: 1011 suffix: q1_2d_da_0 1012 requires: broken 1013 args: -use_da 1 -petscspace_degree 1 -qorder 1 -convergence 1014 test: 1015 suffix: q1_2d_da_1 1016 requires: broken 1017 args: -use_da 1 -petscspace_degree 1 -qorder 1 -porder 1 1018 test: 1019 suffix: q1_2d_da_2 1020 requires: broken 1021 args: -use_da 1 -petscspace_degree 1 -qorder 1 -porder 2 1022 1023 # 2D Q_1 on a quadrilaterial Plex 1024 test: 1025 suffix: q1_2d_plex_0 1026 args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -convergence 1027 test: 1028 suffix: q1_2d_plex_1 1029 args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 1 1030 test: 1031 suffix: q1_2d_plex_2 1032 args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 2 1033 test: 1034 suffix: q1_2d_plex_3 1035 args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 1 -shear_coords 1036 test: 1037 suffix: q1_2d_plex_4 1038 args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 2 -shear_coords 1039 test: 1040 suffix: q1_2d_plex_5 1041 args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 0 -non_affine_coords -convergence 1042 test: 1043 suffix: q1_2d_plex_6 1044 args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 1 -non_affine_coords -convergence 1045 test: 1046 suffix: q1_2d_plex_7 1047 args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 2 -non_affine_coords -convergence 1048 1049 # 2D Q_2 on a quadrilaterial 1050 test: 1051 suffix: q2_2d_plex_0 1052 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -convergence 1053 test: 1054 suffix: q2_2d_plex_1 1055 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 1056 test: 1057 suffix: q2_2d_plex_2 1058 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 1059 test: 1060 suffix: q2_2d_plex_3 1061 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 -shear_coords 1062 test: 1063 suffix: q2_2d_plex_4 1064 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 -shear_coords 1065 test: 1066 suffix: q2_2d_plex_5 1067 args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 0 -non_affine_coords -convergence 1068 test: 1069 suffix: q2_2d_plex_6 1070 args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 1 -non_affine_coords -convergence 1071 test: 1072 suffix: q2_2d_plex_7 1073 args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 2 -non_affine_coords -convergence 1074 1075 # 2D P_3 on a triangle 1076 test: 1077 suffix: p3_2d_0 1078 requires: triangle !single 1079 args: -petscspace_degree 3 -qorder 3 -convergence 1080 test: 1081 suffix: p3_2d_1 1082 requires: triangle !single 1083 args: -petscspace_degree 3 -qorder 3 -porder 1 1084 test: 1085 suffix: p3_2d_2 1086 requires: triangle !single 1087 args: -petscspace_degree 3 -qorder 3 -porder 2 1088 test: 1089 suffix: p3_2d_3 1090 requires: triangle !single 1091 args: -petscspace_degree 3 -qorder 3 -porder 3 1092 test: 1093 suffix: p3_2d_4 1094 requires: triangle mmg 1095 args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -convergence -conv_refine 0 1096 test: 1097 suffix: p3_2d_5 1098 requires: triangle mmg 1099 args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -porder 1 -conv_refine 0 1100 test: 1101 suffix: p3_2d_6 1102 requires: triangle mmg 1103 args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -porder 3 -conv_refine 0 1104 1105 # 2D Q_3 on a quadrilaterial 1106 test: 1107 suffix: q3_2d_0 1108 requires: !single 1109 args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -convergence 1110 test: 1111 suffix: q3_2d_1 1112 requires: !single 1113 args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 1 1114 test: 1115 suffix: q3_2d_2 1116 requires: !single 1117 args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 2 1118 test: 1119 suffix: q3_2d_3 1120 requires: !single 1121 args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 3 1122 1123 # 2D P_1disc on a triangle/quadrilateral 1124 test: 1125 suffix: p1d_2d_0 1126 requires: triangle 1127 args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -convergence 1128 test: 1129 suffix: p1d_2d_1 1130 requires: triangle 1131 args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 1 1132 test: 1133 suffix: p1d_2d_2 1134 requires: triangle 1135 args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 2 1136 test: 1137 suffix: p1d_2d_3 1138 requires: triangle 1139 args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -convergence 1140 filter: sed -e "s/convergence rate at refinement 0: 2/convergence rate at refinement 0: 1.9/g" 1141 test: 1142 suffix: p1d_2d_4 1143 requires: triangle 1144 args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 1 1145 test: 1146 suffix: p1d_2d_5 1147 requires: triangle 1148 args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 2 1149 1150 # 2D BDM_1 on a triangle 1151 test: 1152 suffix: bdm1_2d_0 1153 requires: triangle 1154 args: -petscspace_degree 1 -petscdualspace_type bdm \ 1155 -num_comp 2 -qorder 1 -convergence 1156 test: 1157 suffix: bdm1_2d_1 1158 requires: triangle 1159 args: -petscspace_degree 1 -petscdualspace_type bdm \ 1160 -num_comp 2 -qorder 1 -porder 1 1161 test: 1162 suffix: bdm1_2d_2 1163 requires: triangle 1164 args: -petscspace_degree 1 -petscdualspace_type bdm \ 1165 -num_comp 2 -qorder 1 -porder 2 1166 1167 # 2D BDM_1 on a quadrilateral 1168 test: 1169 suffix: bdm1q_2d_0 1170 requires: triangle 1171 args: -petscspace_degree 1 -petscdualspace_type bdm \ 1172 -petscdualspace_lagrange_tensor 1 \ 1173 -dm_plex_simplex 0 -num_comp 2 -qorder 1 -convergence 1174 test: 1175 suffix: bdm1q_2d_1 1176 requires: triangle 1177 args: -petscspace_degree 1 -petscdualspace_type bdm \ 1178 -petscdualspace_lagrange_tensor 1 \ 1179 -dm_plex_simplex 0 -num_comp 2 -qorder 1 -porder 1 1180 test: 1181 suffix: bdm1q_2d_2 1182 requires: triangle 1183 args: -petscspace_degree 1 -petscdualspace_type bdm \ 1184 -petscdualspace_lagrange_tensor 1 \ 1185 -dm_plex_simplex 0 -num_comp 2 -qorder 1 -porder 2 1186 1187 # Test high order quadrature 1188 test: 1189 suffix: p1_quad_2 1190 requires: triangle 1191 args: -petscspace_degree 1 -qorder 2 -porder 1 1192 test: 1193 suffix: p1_quad_5 1194 requires: triangle 1195 args: -petscspace_degree 1 -qorder 5 -porder 1 1196 test: 1197 suffix: p2_quad_3 1198 requires: triangle 1199 args: -petscspace_degree 2 -qorder 3 -porder 2 1200 test: 1201 suffix: p2_quad_5 1202 requires: triangle 1203 args: -petscspace_degree 2 -qorder 5 -porder 2 1204 test: 1205 suffix: q1_quad_2 1206 args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 2 -porder 1 1207 test: 1208 suffix: q1_quad_5 1209 args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 5 -porder 1 1210 test: 1211 suffix: q2_quad_3 1212 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 3 -porder 1 1213 test: 1214 suffix: q2_quad_5 1215 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 5 -porder 1 1216 1217 # Nonconforming tests 1218 test: 1219 suffix: constraints 1220 args: -dm_coord_space 0 -dm_plex_simplex 0 -petscspace_type tensor -petscspace_degree 1 -qorder 0 -constraints 1221 test: 1222 suffix: nonconforming_tensor_2 1223 nsize: 4 1224 args: -dist_dm_distribute -test_fe_jacobian -test_injector -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_max_projection_height 1 -petscspace_type tensor -petscspace_degree 2 -qorder 2 -dm_view ascii::ASCII_INFO_DETAIL 1225 test: 1226 suffix: nonconforming_tensor_3 1227 nsize: 4 1228 args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_plex_max_projection_height 2 -petscspace_type tensor -petscspace_degree 1 -qorder 1 -dm_view ascii::ASCII_INFO_DETAIL 1229 test: 1230 suffix: nonconforming_tensor_2_fv 1231 nsize: 4 1232 args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -dm_plex_simplex 0 -num_comp 2 1233 test: 1234 suffix: nonconforming_tensor_3_fv 1235 nsize: 4 1236 args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -num_comp 3 1237 test: 1238 suffix: nonconforming_tensor_2_hi 1239 requires: !single 1240 nsize: 4 1241 args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_max_projection_height 1 -petscspace_type tensor -petscspace_degree 4 -qorder 4 1242 test: 1243 suffix: nonconforming_tensor_3_hi 1244 requires: !single skip 1245 nsize: 4 1246 args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_plex_max_projection_height 2 -petscspace_type tensor -petscspace_degree 4 -qorder 4 1247 test: 1248 suffix: nonconforming_simplex_2 1249 requires: triangle 1250 nsize: 4 1251 args: -dist_dm_distribute -test_fe_jacobian -test_injector -petscpartitioner_type simple -tree -dm_plex_max_projection_height 1 -petscspace_degree 2 -qorder 2 -dm_view ascii::ASCII_INFO_DETAIL 1252 test: 1253 suffix: nonconforming_simplex_2_hi 1254 requires: triangle !single 1255 nsize: 4 1256 args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_max_projection_height 1 -petscspace_degree 4 -qorder 4 1257 test: 1258 suffix: nonconforming_simplex_2_fv 1259 requires: triangle 1260 nsize: 4 1261 args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -num_comp 2 1262 test: 1263 suffix: nonconforming_simplex_3 1264 requires: ctetgen 1265 nsize: 4 1266 args: -dist_dm_distribute -test_fe_jacobian -test_injector -petscpartitioner_type simple -tree -dm_plex_dim 3 -dm_plex_max_projection_height 2 -petscspace_degree 2 -qorder 2 -dm_view ascii::ASCII_INFO_DETAIL 1267 test: 1268 suffix: nonconforming_simplex_3_hi 1269 requires: ctetgen skip 1270 nsize: 4 1271 args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_dim 3 -dm_plex_max_projection_height 2 -petscspace_degree 4 -qorder 4 1272 test: 1273 suffix: nonconforming_simplex_3_fv 1274 requires: ctetgen 1275 nsize: 4 1276 args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -dm_plex_dim 3 -num_comp 3 1277 1278 # 3D WXY on a triangular prism 1279 test: 1280 suffix: wxy_0 1281 args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -qorder 2 -porder 0 \ 1282 -petscspace_type sum \ 1283 -petscspace_variables 3 \ 1284 -petscspace_components 3 \ 1285 -petscspace_sum_spaces 2 \ 1286 -petscspace_sum_concatenate false \ 1287 -sumcomp_0_petscspace_variables 3 \ 1288 -sumcomp_0_petscspace_components 3 \ 1289 -sumcomp_0_petscspace_degree 1 \ 1290 -sumcomp_1_petscspace_variables 3 \ 1291 -sumcomp_1_petscspace_components 3 \ 1292 -sumcomp_1_petscspace_type wxy \ 1293 -petscdualspace_refcell triangular_prism \ 1294 -petscdualspace_form_degree 0 \ 1295 -petscdualspace_order 1 \ 1296 -petscdualspace_components 3 1297 1298 TEST*/ 1299 1300 /* 1301 # 2D Q_2 on a quadrilaterial Plex 1302 test: 1303 suffix: q2_2d_plex_0 1304 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -convergence 1305 test: 1306 suffix: q2_2d_plex_1 1307 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 1308 test: 1309 suffix: q2_2d_plex_2 1310 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 1311 test: 1312 suffix: q2_2d_plex_3 1313 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 -shear_coords 1314 test: 1315 suffix: q2_2d_plex_4 1316 args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 -shear_coords 1317 test: 1318 suffix: q2_2d_plex_5 1319 args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 0 -non_affine_coords 1320 test: 1321 suffix: q2_2d_plex_6 1322 args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 1 -non_affine_coords 1323 test: 1324 suffix: q2_2d_plex_7 1325 args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 2 -non_affine_coords 1326 1327 test: 1328 suffix: p1d_2d_6 1329 requires: mmg 1330 args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0 1331 test: 1332 suffix: p1d_2d_7 1333 requires: mmg 1334 args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0 1335 test: 1336 suffix: p1d_2d_8 1337 requires: mmg 1338 args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0 1339 */ 1340