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