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