1c4762a1bSJed Brown static char help[] = "Check that a DM can accurately represent and interpolate functions of a given polynomial order\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown #include <petscdm.h> 5c4762a1bSJed Brown #include <petscdmda.h> 6c4762a1bSJed Brown #include <petscfe.h> 7c4762a1bSJed Brown #include <petscds.h> 8c4762a1bSJed Brown #include <petscksp.h> 9c4762a1bSJed Brown #include <petscsnes.h> 10c4762a1bSJed Brown 11c4762a1bSJed Brown typedef struct { 12c4762a1bSJed Brown /* Domain and mesh definition */ 13c4762a1bSJed Brown PetscBool useDA; /* Flag DMDA tensor product mesh */ 14c4762a1bSJed Brown PetscBool shearCoords; /* Flag for shear transform */ 15c4762a1bSJed Brown PetscBool nonaffineCoords; /* Flag for non-affine transform */ 16c4762a1bSJed Brown /* Element definition */ 17c4762a1bSJed Brown PetscInt qorder; /* Order of the quadrature */ 18c4762a1bSJed Brown PetscInt numComponents; /* Number of field components */ 19c4762a1bSJed Brown PetscFE fe; /* The finite element */ 20c4762a1bSJed Brown /* Testing space */ 21c4762a1bSJed Brown PetscInt porder; /* Order of polynomials to test */ 22c4762a1bSJed Brown PetscBool convergence; /* Test for order of convergence */ 23c4762a1bSJed Brown PetscBool convRefine; /* Test for convergence using refinement, otherwise use coarsening */ 24c4762a1bSJed Brown PetscBool constraints; /* Test local constraints */ 25c4762a1bSJed Brown PetscBool tree; /* Test tree routines */ 26c4762a1bSJed Brown PetscBool testFEjacobian; /* Test finite element Jacobian assembly */ 27c4762a1bSJed Brown PetscBool testFVgrad; /* Test finite difference gradient routine */ 28c4762a1bSJed Brown PetscBool testInjector; /* Test finite element injection routines */ 29c4762a1bSJed Brown PetscInt treeCell; /* Cell to refine in tree test */ 30c4762a1bSJed Brown PetscReal constants[3]; /* Constant values for each dimension */ 31c4762a1bSJed Brown } AppCtx; 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* u = 1 */ 34c4762a1bSJed Brown PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx) 35c4762a1bSJed Brown { 36c4762a1bSJed Brown AppCtx *user = (AppCtx *) ctx; 37c4762a1bSJed Brown PetscInt d; 3830602db0SMatthew G. Knepley for (d = 0; d < dim; ++d) u[d] = user->constants[d]; 39c4762a1bSJed Brown return 0; 40c4762a1bSJed Brown } 41c4762a1bSJed Brown PetscErrorCode constantDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) 42c4762a1bSJed Brown { 43c4762a1bSJed Brown PetscInt d; 4430602db0SMatthew G. Knepley for (d = 0; d < dim; ++d) u[d] = 0.0; 45c4762a1bSJed Brown return 0; 46c4762a1bSJed Brown } 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* u = x */ 49c4762a1bSJed Brown PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx) 50c4762a1bSJed Brown { 51c4762a1bSJed Brown PetscInt d; 52c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[d] = coords[d]; 53c4762a1bSJed Brown return 0; 54c4762a1bSJed Brown } 55c4762a1bSJed Brown PetscErrorCode linearDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) 56c4762a1bSJed Brown { 57c4762a1bSJed Brown PetscInt d, e; 58c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 59c4762a1bSJed Brown u[d] = 0.0; 60c4762a1bSJed Brown for (e = 0; e < dim; ++e) u[d] += (d == e ? 1.0 : 0.0) * n[e]; 61c4762a1bSJed Brown } 62c4762a1bSJed Brown return 0; 63c4762a1bSJed Brown } 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* u = x^2 or u = (x^2, xy) or u = (xy, yz, zx) */ 66c4762a1bSJed Brown PetscErrorCode quadratic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx) 67c4762a1bSJed Brown { 6830602db0SMatthew G. Knepley if (dim > 2) {u[0] = coords[0]*coords[1]; u[1] = coords[1]*coords[2]; u[2] = coords[2]*coords[0];} 6930602db0SMatthew G. Knepley else if (dim > 1) {u[0] = coords[0]*coords[0]; u[1] = coords[0]*coords[1];} 7030602db0SMatthew G. Knepley else if (dim > 0) {u[0] = coords[0]*coords[0];} 71c4762a1bSJed Brown return 0; 72c4762a1bSJed Brown } 73c4762a1bSJed Brown PetscErrorCode quadraticDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) 74c4762a1bSJed Brown { 7530602db0SMatthew G. Knepley 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];} 7630602db0SMatthew G. Knepley else if (dim > 1) {u[0] = 2.0*coords[0]*n[0]; u[1] = coords[1]*n[0] + coords[0]*n[1];} 7730602db0SMatthew G. Knepley else if (dim > 0) {u[0] = 2.0*coords[0]*n[0];} 78c4762a1bSJed Brown return 0; 79c4762a1bSJed Brown } 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* u = x^3 or u = (x^3, x^2y) or u = (x^2y, y^2z, z^2x) */ 82c4762a1bSJed Brown PetscErrorCode cubic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx) 83c4762a1bSJed Brown { 8430602db0SMatthew G. Knepley 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];} 8530602db0SMatthew G. Knepley else if (dim > 1) {u[0] = coords[0]*coords[0]*coords[0]; u[1] = coords[0]*coords[0]*coords[1];} 8630602db0SMatthew G. Knepley else if (dim > 0) {u[0] = coords[0]*coords[0]*coords[0];} 87c4762a1bSJed Brown return 0; 88c4762a1bSJed Brown } 89c4762a1bSJed Brown PetscErrorCode cubicDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) 90c4762a1bSJed Brown { 9130602db0SMatthew G. Knepley 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];} 9230602db0SMatthew G. Knepley 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];} 9330602db0SMatthew G. Knepley else if (dim > 0) {u[0] = 3.0*coords[0]*coords[0]*n[0];} 94c4762a1bSJed Brown return 0; 95c4762a1bSJed Brown } 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* u = tanh(x) */ 98c4762a1bSJed Brown PetscErrorCode trig(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nf, PetscScalar *u, void *ctx) 99c4762a1bSJed Brown { 100c4762a1bSJed Brown PetscInt d; 10130602db0SMatthew G. Knepley for (d = 0; d < dim; ++d) u[d] = PetscTanhReal(coords[d] - 0.5); 102c4762a1bSJed Brown return 0; 103c4762a1bSJed Brown } 104c4762a1bSJed Brown PetscErrorCode trigDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) 105c4762a1bSJed Brown { 106c4762a1bSJed Brown PetscInt d; 10730602db0SMatthew G. Knepley for (d = 0; d < dim; ++d) u[d] = 1.0/PetscSqr(PetscCoshReal(coords[d] - 0.5)) * n[d]; 108c4762a1bSJed Brown return 0; 109c4762a1bSJed Brown } 110c4762a1bSJed Brown 111c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 112c4762a1bSJed Brown { 113c4762a1bSJed Brown PetscInt n = 3; 114c4762a1bSJed Brown 115c4762a1bSJed Brown PetscFunctionBeginUser; 11630602db0SMatthew G. Knepley options->useDA = PETSC_FALSE; 117c4762a1bSJed Brown options->shearCoords = PETSC_FALSE; 118c4762a1bSJed Brown options->nonaffineCoords = PETSC_FALSE; 119c4762a1bSJed Brown options->qorder = 0; 120c4762a1bSJed Brown options->numComponents = PETSC_DEFAULT; 121c4762a1bSJed Brown options->porder = 0; 122c4762a1bSJed Brown options->convergence = PETSC_FALSE; 123c4762a1bSJed Brown options->convRefine = PETSC_TRUE; 124c4762a1bSJed Brown options->constraints = PETSC_FALSE; 125c4762a1bSJed Brown options->tree = PETSC_FALSE; 126c4762a1bSJed Brown options->treeCell = 0; 127c4762a1bSJed Brown options->testFEjacobian = PETSC_FALSE; 128c4762a1bSJed Brown options->testFVgrad = PETSC_FALSE; 129c4762a1bSJed Brown options->testInjector = PETSC_FALSE; 130c4762a1bSJed Brown options->constants[0] = 1.0; 131c4762a1bSJed Brown options->constants[1] = 2.0; 132c4762a1bSJed Brown options->constants[2] = 3.0; 133c4762a1bSJed Brown 134d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Projection Test Options", "DMPlex"); 1359566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_da", "Flag for DMDA mesh", "ex3.c", options->useDA, &options->useDA, NULL)); 1369566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-shear_coords", "Transform coordinates with a shear", "ex3.c", options->shearCoords, &options->shearCoords, NULL)); 1379566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-non_affine_coords", "Transform coordinates with a non-affine transform", "ex3.c", options->nonaffineCoords, &options->nonaffineCoords, NULL)); 1389566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-qorder", "The quadrature order", "ex3.c", options->qorder, &options->qorder, NULL,0)); 1399566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-num_comp", "The number of field components", "ex3.c", options->numComponents, &options->numComponents, NULL,PETSC_DEFAULT)); 1409566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-porder", "The order of polynomials to test", "ex3.c", options->porder, &options->porder, NULL,0)); 1419566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-convergence", "Check the convergence rate", "ex3.c", options->convergence, &options->convergence, NULL)); 1429566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-conv_refine", "Use refinement for the convergence rate", "ex3.c", options->convRefine, &options->convRefine, NULL)); 1439566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-constraints", "Test local constraints (serial only)", "ex3.c", options->constraints, &options->constraints, NULL)); 1449566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tree", "Test tree routines", "ex3.c", options->tree, &options->tree, NULL)); 1459566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-tree_cell", "cell to refine in tree test", "ex3.c", options->treeCell, &options->treeCell, NULL,0)); 1469566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_fe_jacobian", "Test finite element Jacobian assembly", "ex3.c", options->testFEjacobian, &options->testFEjacobian, NULL)); 1479566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_fv_grad", "Test finite volume gradient reconstruction", "ex3.c", options->testFVgrad, &options->testFVgrad, NULL)); 1489566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_injector","Test finite element injection", "ex3.c", options->testInjector, &options->testInjector,NULL)); 1499566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-constants","Set the constant values", "ex3.c", options->constants, &n,NULL)); 150d0609cedSBarry Smith PetscOptionsEnd(); 151c4762a1bSJed Brown 152c4762a1bSJed Brown PetscFunctionReturn(0); 153c4762a1bSJed Brown } 154c4762a1bSJed Brown 155c4762a1bSJed Brown static PetscErrorCode TransformCoordinates(DM dm, AppCtx *user) 156c4762a1bSJed Brown { 157c4762a1bSJed Brown PetscSection coordSection; 158c4762a1bSJed Brown Vec coordinates; 159c4762a1bSJed Brown PetscScalar *coords; 160c4762a1bSJed Brown PetscInt vStart, vEnd, v; 161c4762a1bSJed Brown 162c4762a1bSJed Brown PetscFunctionBeginUser; 163c4762a1bSJed Brown if (user->nonaffineCoords) { 164c4762a1bSJed Brown /* x' = r^(1/p) (x/r), y' = r^(1/p) (y/r), z' = z */ 1659566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 1669566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 1679566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 1689566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 169c4762a1bSJed Brown for (v = vStart; v < vEnd; ++v) { 170c4762a1bSJed Brown PetscInt dof, off; 171c4762a1bSJed Brown PetscReal p = 4.0, r; 172c4762a1bSJed Brown 1739566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(coordSection, v, &dof)); 1749566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, v, &off)); 175c4762a1bSJed Brown switch (dof) { 176c4762a1bSJed Brown case 2: 177c4762a1bSJed Brown r = PetscSqr(PetscRealPart(coords[off+0])) + PetscSqr(PetscRealPart(coords[off+1])); 178c4762a1bSJed Brown coords[off+0] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p)/(2*p))*coords[off+0]; 179c4762a1bSJed Brown coords[off+1] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p)/(2*p))*coords[off+1]; 180c4762a1bSJed Brown break; 181c4762a1bSJed Brown case 3: 182c4762a1bSJed Brown r = PetscSqr(PetscRealPart(coords[off+0])) + PetscSqr(PetscRealPart(coords[off+1])); 183c4762a1bSJed Brown coords[off+0] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p)/(2*p))*coords[off+0]; 184c4762a1bSJed Brown coords[off+1] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p)/(2*p))*coords[off+1]; 185c4762a1bSJed Brown coords[off+2] = coords[off+2]; 186c4762a1bSJed Brown break; 187c4762a1bSJed Brown } 188c4762a1bSJed Brown } 1899566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 190c4762a1bSJed Brown } 191c4762a1bSJed Brown if (user->shearCoords) { 192c4762a1bSJed Brown /* x' = x + m y + m z, y' = y + m z, z' = z */ 1939566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dm, &coordSection)); 1949566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 1959566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 1969566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 197c4762a1bSJed Brown for (v = vStart; v < vEnd; ++v) { 198c4762a1bSJed Brown PetscInt dof, off; 199c4762a1bSJed Brown PetscReal m = 1.0; 200c4762a1bSJed Brown 2019566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(coordSection, v, &dof)); 2029566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(coordSection, v, &off)); 203c4762a1bSJed Brown switch (dof) { 204c4762a1bSJed Brown case 2: 205c4762a1bSJed Brown coords[off+0] = coords[off+0] + m*coords[off+1]; 206c4762a1bSJed Brown coords[off+1] = coords[off+1]; 207c4762a1bSJed Brown break; 208c4762a1bSJed Brown case 3: 209c4762a1bSJed Brown coords[off+0] = coords[off+0] + m*coords[off+1] + m*coords[off+2]; 210c4762a1bSJed Brown coords[off+1] = coords[off+1] + m*coords[off+2]; 211c4762a1bSJed Brown coords[off+2] = coords[off+2]; 212c4762a1bSJed Brown break; 213c4762a1bSJed Brown } 214c4762a1bSJed Brown } 2159566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 216c4762a1bSJed Brown } 217c4762a1bSJed Brown PetscFunctionReturn(0); 218c4762a1bSJed Brown } 219c4762a1bSJed Brown 220c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 221c4762a1bSJed Brown { 22230602db0SMatthew G. Knepley PetscInt dim = 2; 22330602db0SMatthew G. Knepley PetscBool simplex; 224c4762a1bSJed Brown 225c4762a1bSJed Brown PetscFunctionBeginUser; 22630602db0SMatthew G. Knepley if (user->useDA) { 22730602db0SMatthew G. Knepley switch (dim) { 228c4762a1bSJed Brown case 2: 2299566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 2, 2, PETSC_DETERMINE, PETSC_DETERMINE, 1, 1, NULL, NULL, dm)); 2309566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 2319566063dSJacob Faibussowitsch PetscCall(DMSetUp(*dm)); 2329566063dSJacob Faibussowitsch PetscCall(DMDASetVertexCoordinates(*dm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 233c4762a1bSJed Brown break; 234c4762a1bSJed Brown default: 235*63a3b9bcSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot create structured mesh of dimension %" PetscInt_FMT, dim); 236c4762a1bSJed Brown } 2379566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) *dm, "Hexahedral Mesh")); 23830602db0SMatthew G. Knepley } else { 2399566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 2409566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 2419566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 2429566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 243c4762a1bSJed Brown 2449566063dSJacob Faibussowitsch PetscCall(DMGetDimension(*dm, &dim)); 2459566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(*dm, &simplex)); 2469566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&simplex, 1, MPIU_BOOL, 0, comm)); 247c4762a1bSJed Brown if (user->tree) { 24830602db0SMatthew G. Knepley DM refTree, ncdm = NULL; 249c4762a1bSJed Brown 2509566063dSJacob Faibussowitsch PetscCall(DMPlexCreateDefaultReferenceTree(comm,dim,simplex,&refTree)); 2519566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(refTree,NULL,"-reftree_dm_view")); 2529566063dSJacob Faibussowitsch PetscCall(DMPlexSetReferenceTree(*dm,refTree)); 2539566063dSJacob Faibussowitsch PetscCall(DMDestroy(&refTree)); 2549566063dSJacob Faibussowitsch PetscCall(DMPlexTreeRefineCell(*dm,user->treeCell,&ncdm)); 255c4762a1bSJed Brown if (ncdm) { 2569566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 257c4762a1bSJed Brown *dm = ncdm; 2589566063dSJacob Faibussowitsch PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_FALSE)); 259c4762a1bSJed Brown } 2609566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *dm, "tree_")); 2619566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 2629566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 2639566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm,NULL,"-dm_view")); 264c4762a1bSJed Brown } else { 2659566063dSJacob Faibussowitsch PetscCall(DMPlexSetRefinementUniform(*dm, PETSC_TRUE)); 266c4762a1bSJed Brown } 2679566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *dm, "dist_")); 2689566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 2699566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 2709566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *dm, NULL)); 2719566063dSJacob Faibussowitsch if (simplex) PetscCall(PetscObjectSetName((PetscObject) *dm, "Simplicial Mesh")); 2729566063dSJacob Faibussowitsch else PetscCall(PetscObjectSetName((PetscObject) *dm, "Hexahedral Mesh")); 273c4762a1bSJed Brown } 2749566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 2759566063dSJacob Faibussowitsch PetscCall(TransformCoordinates(*dm, user)); 2769566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm,NULL,"-dm_view")); 277c4762a1bSJed Brown PetscFunctionReturn(0); 278c4762a1bSJed Brown } 279c4762a1bSJed Brown 280c4762a1bSJed Brown static void simple_mass(PetscInt dim, PetscInt Nf, PetscInt NfAux, 281c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 282c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 283c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 284c4762a1bSJed Brown { 285c4762a1bSJed Brown PetscInt d, e; 286c4762a1bSJed Brown for (d = 0, e = 0; d < dim; d++, e+=dim+1) { 287c4762a1bSJed Brown g0[e] = 1.; 288c4762a1bSJed Brown } 289c4762a1bSJed Brown } 290c4762a1bSJed Brown 291c4762a1bSJed Brown /* < \nabla v, 1/2(\nabla u + {\nabla u}^T) > */ 292c4762a1bSJed Brown static void symmetric_gradient_inner_product(PetscInt dim, PetscInt Nf, PetscInt NfAux, 293c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 294c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 295c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar C[]) 296c4762a1bSJed Brown { 297c4762a1bSJed Brown PetscInt compI, compJ, d, e; 298c4762a1bSJed Brown 299c4762a1bSJed Brown for (compI = 0; compI < dim; ++compI) { 300c4762a1bSJed Brown for (compJ = 0; compJ < dim; ++compJ) { 301c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 302c4762a1bSJed Brown for (e = 0; e < dim; e++) { 303c4762a1bSJed Brown if (d == e && d == compI && d == compJ) { 304c4762a1bSJed Brown C[((compI*dim+compJ)*dim+d)*dim+e] = 1.0; 305c4762a1bSJed Brown } else if ((d == compJ && e == compI) || (d == e && compI == compJ)) { 306c4762a1bSJed Brown C[((compI*dim+compJ)*dim+d)*dim+e] = 0.5; 307c4762a1bSJed Brown } else { 308c4762a1bSJed Brown C[((compI*dim+compJ)*dim+d)*dim+e] = 0.0; 309c4762a1bSJed Brown } 310c4762a1bSJed Brown } 311c4762a1bSJed Brown } 312c4762a1bSJed Brown } 313c4762a1bSJed Brown } 314c4762a1bSJed Brown } 315c4762a1bSJed Brown 316c4762a1bSJed Brown static PetscErrorCode SetupSection(DM dm, AppCtx *user) 317c4762a1bSJed Brown { 318c4762a1bSJed Brown PetscFunctionBeginUser; 31930602db0SMatthew G. Knepley if (user->constraints) { 320c4762a1bSJed Brown /* test local constraints */ 321c4762a1bSJed Brown DM coordDM; 322c4762a1bSJed Brown PetscInt fStart, fEnd, f, vStart, vEnd, v; 323c4762a1bSJed Brown PetscInt edgesx = 2, vertsx; 324c4762a1bSJed Brown PetscInt edgesy = 2, vertsy; 325c4762a1bSJed Brown PetscMPIInt size; 326c4762a1bSJed Brown PetscInt numConst; 327c4762a1bSJed Brown PetscSection aSec; 328c4762a1bSJed Brown PetscInt *anchors; 329c4762a1bSJed Brown PetscInt offset; 330c4762a1bSJed Brown IS aIS; 331c4762a1bSJed Brown MPI_Comm comm = PetscObjectComm((PetscObject)dm); 332c4762a1bSJed Brown 3339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 33408401ef6SPierre Jolivet PetscCheck(size <= 1,comm,PETSC_ERR_SUP,"Local constraint test can only be performed in serial"); 335c4762a1bSJed Brown 336c4762a1bSJed Brown /* we are going to test constraints by using them to enforce periodicity 337c4762a1bSJed Brown * in one direction, and comparing to the existing method of enforcing 338c4762a1bSJed Brown * periodicity */ 339c4762a1bSJed Brown 340c4762a1bSJed Brown /* first create the coordinate section so that it does not clone the 341c4762a1bSJed Brown * constraints */ 3429566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm,&coordDM)); 343c4762a1bSJed Brown 344c4762a1bSJed Brown /* create the constrained-to-anchor section */ 3459566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm,0,&vStart,&vEnd)); 3469566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm,1,&fStart,&fEnd)); 3479566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PETSC_COMM_SELF,&aSec)); 3489566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(aSec,PetscMin(fStart,vStart),PetscMax(fEnd,vEnd))); 349c4762a1bSJed Brown 350c4762a1bSJed Brown /* define the constraints */ 3519566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-da_grid_x",&edgesx,NULL)); 3529566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-da_grid_y",&edgesy,NULL)); 353c4762a1bSJed Brown vertsx = edgesx + 1; 354c4762a1bSJed Brown vertsy = edgesy + 1; 355c4762a1bSJed Brown numConst = vertsy + edgesy; 3569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numConst,&anchors)); 357c4762a1bSJed Brown offset = 0; 358c4762a1bSJed Brown for (v = vStart + edgesx; v < vEnd; v+= vertsx) { 3599566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(aSec,v,1)); 360c4762a1bSJed Brown anchors[offset++] = v - edgesx; 361c4762a1bSJed Brown } 362c4762a1bSJed Brown for (f = fStart + edgesx * vertsy + edgesx * edgesy; f < fEnd; f++) { 3639566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(aSec,f,1)); 364c4762a1bSJed Brown anchors[offset++] = f - edgesx * edgesy; 365c4762a1bSJed Brown } 3669566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(aSec)); 3679566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,numConst,anchors,PETSC_OWN_POINTER,&aIS)); 368c4762a1bSJed Brown 3699566063dSJacob Faibussowitsch PetscCall(DMPlexSetAnchors(dm,aSec,aIS)); 3709566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&aSec)); 3719566063dSJacob Faibussowitsch PetscCall(ISDestroy(&aIS)); 372c4762a1bSJed Brown } 3739566063dSJacob Faibussowitsch PetscCall(DMSetNumFields(dm, 1)); 3749566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject) user->fe)); 3759566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 37630602db0SMatthew G. Knepley if (user->constraints) { 377c4762a1bSJed Brown /* test getting local constraint matrix that matches section */ 378c4762a1bSJed Brown PetscSection aSec; 379c4762a1bSJed Brown IS aIS; 380c4762a1bSJed Brown 3819566063dSJacob Faibussowitsch PetscCall(DMPlexGetAnchors(dm,&aSec,&aIS)); 382c4762a1bSJed Brown if (aSec) { 383c4762a1bSJed Brown PetscDS ds; 384c4762a1bSJed Brown PetscSection cSec, section; 385c4762a1bSJed Brown PetscInt cStart, cEnd, c, numComp; 386c4762a1bSJed Brown Mat cMat, mass; 387c4762a1bSJed Brown Vec local; 388c4762a1bSJed Brown const PetscInt *anchors; 389c4762a1bSJed Brown 3909566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm,§ion)); 391c4762a1bSJed Brown /* this creates the matrix and preallocates the matrix structure: we 392c4762a1bSJed Brown * just have to fill in the values */ 3939566063dSJacob Faibussowitsch PetscCall(DMGetDefaultConstraints(dm,&cSec,&cMat,NULL)); 3949566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(cSec,&cStart,&cEnd)); 3959566063dSJacob Faibussowitsch PetscCall(ISGetIndices(aIS,&anchors)); 3969566063dSJacob Faibussowitsch PetscCall(PetscFEGetNumComponents(user->fe, &numComp)); 397c4762a1bSJed Brown for (c = cStart; c < cEnd; c++) { 398c4762a1bSJed Brown PetscInt cDof; 399c4762a1bSJed Brown 400c4762a1bSJed Brown /* is this point constrained? (does it have an anchor?) */ 4019566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(aSec,c,&cDof)); 402c4762a1bSJed Brown if (cDof) { 403c4762a1bSJed Brown PetscInt cOff, a, aDof, aOff, j; 404*63a3b9bcSJacob Faibussowitsch PetscCheck(cDof == 1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Found %" PetscInt_FMT " anchor points: should be just one",cDof); 405c4762a1bSJed Brown 406c4762a1bSJed Brown /* find the anchor point */ 4079566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(aSec,c,&cOff)); 408c4762a1bSJed Brown a = anchors[cOff]; 409c4762a1bSJed Brown 410c4762a1bSJed Brown /* find the constrained dofs (row in constraint matrix) */ 4119566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(cSec,c,&cDof)); 4129566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(cSec,c,&cOff)); 413c4762a1bSJed Brown 414c4762a1bSJed Brown /* find the anchor dofs (column in constraint matrix) */ 4159566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section,a,&aDof)); 4169566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section,a,&aOff)); 417c4762a1bSJed Brown 418*63a3b9bcSJacob Faibussowitsch PetscCheck(cDof == aDof,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Point and anchor have different number of dofs: %" PetscInt_FMT ", %" PetscInt_FMT,cDof,aDof); 419*63a3b9bcSJacob Faibussowitsch PetscCheck(cDof % numComp == 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Point dofs not divisible by field components: %" PetscInt_FMT ", %" PetscInt_FMT,cDof,numComp); 420c4762a1bSJed Brown 421c4762a1bSJed Brown /* put in a simple equality constraint */ 422c4762a1bSJed Brown for (j = 0; j < cDof; j++) { 4239566063dSJacob Faibussowitsch PetscCall(MatSetValue(cMat,cOff+j,aOff+j,1.,INSERT_VALUES)); 424c4762a1bSJed Brown } 425c4762a1bSJed Brown } 426c4762a1bSJed Brown } 4279566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(cMat,MAT_FINAL_ASSEMBLY)); 4289566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(cMat,MAT_FINAL_ASSEMBLY)); 4299566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(aIS,&anchors)); 430c4762a1bSJed Brown 431c4762a1bSJed Brown /* Now that we have constructed the constraint matrix, any FE matrix 432c4762a1bSJed Brown * that we construct will apply the constraints during construction */ 433c4762a1bSJed Brown 4349566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm,&mass)); 435c4762a1bSJed Brown /* get a dummy local variable to serve as the solution */ 4369566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm,&local)); 4379566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm,&ds)); 438c4762a1bSJed Brown /* set the jacobian to be the mass matrix */ 4399566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, simple_mass, NULL, NULL, NULL)); 440c4762a1bSJed Brown /* build the mass matrix */ 4419566063dSJacob Faibussowitsch PetscCall(DMPlexSNESComputeJacobianFEM(dm,local,mass,mass,NULL)); 4429566063dSJacob Faibussowitsch PetscCall(MatView(mass,PETSC_VIEWER_STDOUT_WORLD)); 4439566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mass)); 4449566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm,&local)); 445c4762a1bSJed Brown } 446c4762a1bSJed Brown } 447c4762a1bSJed Brown PetscFunctionReturn(0); 448c4762a1bSJed Brown } 449c4762a1bSJed Brown 450c4762a1bSJed Brown static PetscErrorCode TestFEJacobian(DM dm, AppCtx *user) 451c4762a1bSJed Brown { 452c4762a1bSJed Brown PetscFunctionBeginUser; 45330602db0SMatthew G. Knepley if (!user->useDA) { 454c4762a1bSJed Brown Vec local; 455c4762a1bSJed Brown const Vec *vecs; 456c4762a1bSJed Brown Mat E; 457c4762a1bSJed Brown MatNullSpace sp; 458c4762a1bSJed Brown PetscBool isNullSpace, hasConst; 45930602db0SMatthew G. Knepley PetscInt dim, n, i; 460c4762a1bSJed Brown Vec res = NULL, localX, localRes; 461c4762a1bSJed Brown PetscDS ds; 462c4762a1bSJed Brown 4639566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 464*63a3b9bcSJacob Faibussowitsch 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); 4659566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm,&ds)); 4669566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds,0,0,NULL,NULL,NULL,symmetric_gradient_inner_product)); 4679566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm,&E)); 4689566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm,&local)); 4699566063dSJacob Faibussowitsch PetscCall(DMPlexSNESComputeJacobianFEM(dm,local,E,E,NULL)); 4709566063dSJacob Faibussowitsch PetscCall(DMPlexCreateRigidBody(dm,0,&sp)); 4719566063dSJacob Faibussowitsch PetscCall(MatNullSpaceGetVecs(sp,&hasConst,&n,&vecs)); 4729566063dSJacob Faibussowitsch if (n) PetscCall(VecDuplicate(vecs[0],&res)); 4739566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dm,&localX)); 4749566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dm,&localRes)); 475c4762a1bSJed Brown for (i = 0; i < n; i++) { /* also test via matrix-free Jacobian application */ 476c4762a1bSJed Brown PetscReal resNorm; 477c4762a1bSJed Brown 4789566063dSJacob Faibussowitsch PetscCall(VecSet(localRes,0.)); 4799566063dSJacob Faibussowitsch PetscCall(VecSet(localX,0.)); 4809566063dSJacob Faibussowitsch PetscCall(VecSet(local,0.)); 4819566063dSJacob Faibussowitsch PetscCall(VecSet(res,0.)); 4829566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm,vecs[i],INSERT_VALUES,localX)); 4839566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm,vecs[i],INSERT_VALUES,localX)); 4849566063dSJacob Faibussowitsch PetscCall(DMSNESComputeJacobianAction(dm,local,localX,localRes,NULL)); 4859566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dm,localRes,ADD_VALUES,res)); 4869566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dm,localRes,ADD_VALUES,res)); 4879566063dSJacob Faibussowitsch PetscCall(VecNorm(res,NORM_2,&resNorm)); 488c4762a1bSJed Brown if (resNorm > PETSC_SMALL) { 489*63a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm),"Symmetric gradient action null space vector %" PetscInt_FMT " residual: %E\n",i,(double)resNorm)); 490c4762a1bSJed Brown } 491c4762a1bSJed Brown } 4929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&localRes)); 4939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&localX)); 4949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&res)); 4959566063dSJacob Faibussowitsch PetscCall(MatNullSpaceTest(sp,E,&isNullSpace)); 496c4762a1bSJed Brown if (isNullSpace) { 4979566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm),"Symmetric gradient null space: PASS\n")); 498c4762a1bSJed Brown } else { 4999566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm),"Symmetric gradient null space: FAIL\n")); 500c4762a1bSJed Brown } 5019566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&sp)); 5029566063dSJacob Faibussowitsch PetscCall(MatDestroy(&E)); 5039566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm,&local)); 504c4762a1bSJed Brown } 505c4762a1bSJed Brown PetscFunctionReturn(0); 506c4762a1bSJed Brown } 507c4762a1bSJed Brown 508c4762a1bSJed Brown static PetscErrorCode TestInjector(DM dm, AppCtx *user) 509c4762a1bSJed Brown { 510c4762a1bSJed Brown DM refTree; 511c4762a1bSJed Brown PetscMPIInt rank; 512c4762a1bSJed Brown 513c4762a1bSJed Brown PetscFunctionBegin; 5149566063dSJacob Faibussowitsch PetscCall(DMPlexGetReferenceTree(dm,&refTree)); 515c4762a1bSJed Brown if (refTree) { 516c4762a1bSJed Brown Mat inj; 517c4762a1bSJed Brown 5189566063dSJacob Faibussowitsch PetscCall(DMPlexComputeInjectorReferenceTree(refTree,&inj)); 5199566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)inj,"Reference Tree Injector")); 5209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 521dd400576SPatrick Sanan if (rank == 0) { 5229566063dSJacob Faibussowitsch PetscCall(MatView(inj,PETSC_VIEWER_STDOUT_SELF)); 523c4762a1bSJed Brown } 5249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&inj)); 525c4762a1bSJed Brown } 526c4762a1bSJed Brown PetscFunctionReturn(0); 527c4762a1bSJed Brown } 528c4762a1bSJed Brown 529c4762a1bSJed Brown static PetscErrorCode TestFVGrad(DM dm, AppCtx *user) 530c4762a1bSJed Brown { 531c4762a1bSJed Brown MPI_Comm comm; 532c4762a1bSJed Brown DM dmRedist, dmfv, dmgrad, dmCell, refTree; 533c4762a1bSJed Brown PetscFV fv; 53430602db0SMatthew G. Knepley PetscInt dim, nvecs, v, cStart, cEnd, cEndInterior; 535c4762a1bSJed Brown PetscMPIInt size; 536c4762a1bSJed Brown Vec cellgeom, grad, locGrad; 537c4762a1bSJed Brown const PetscScalar *cgeom; 538c4762a1bSJed Brown PetscReal allVecMaxDiff = 0., fvTol = 100. * PETSC_MACHINE_EPSILON; 539c4762a1bSJed Brown 540c4762a1bSJed Brown PetscFunctionBeginUser; 541c4762a1bSJed Brown comm = PetscObjectComm((PetscObject)dm); 5429566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 543c4762a1bSJed Brown /* duplicate DM, give dup. a FV discretization */ 5449566063dSJacob Faibussowitsch PetscCall(DMSetBasicAdjacency(dm,PETSC_TRUE,PETSC_FALSE)); 5459566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 546c4762a1bSJed Brown dmRedist = NULL; 547c4762a1bSJed Brown if (size > 1) { 5489566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeOverlap(dm,1,NULL,&dmRedist)); 549c4762a1bSJed Brown } 550c4762a1bSJed Brown if (!dmRedist) { 551c4762a1bSJed Brown dmRedist = dm; 5529566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dmRedist)); 553c4762a1bSJed Brown } 5549566063dSJacob Faibussowitsch PetscCall(PetscFVCreate(comm,&fv)); 5559566063dSJacob Faibussowitsch PetscCall(PetscFVSetType(fv,PETSCFVLEASTSQUARES)); 5569566063dSJacob Faibussowitsch PetscCall(PetscFVSetNumComponents(fv,user->numComponents)); 5579566063dSJacob Faibussowitsch PetscCall(PetscFVSetSpatialDimension(fv,dim)); 5589566063dSJacob Faibussowitsch PetscCall(PetscFVSetFromOptions(fv)); 5599566063dSJacob Faibussowitsch PetscCall(PetscFVSetUp(fv)); 5609566063dSJacob Faibussowitsch PetscCall(DMPlexConstructGhostCells(dmRedist,NULL,NULL,&dmfv)); 5619566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmRedist)); 5629566063dSJacob Faibussowitsch PetscCall(DMSetNumFields(dmfv,1)); 5639566063dSJacob Faibussowitsch PetscCall(DMSetField(dmfv, 0, NULL, (PetscObject) fv)); 5649566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dmfv)); 5659566063dSJacob Faibussowitsch PetscCall(DMPlexGetReferenceTree(dm,&refTree)); 5669566063dSJacob Faibussowitsch if (refTree) PetscCall(DMCopyDisc(dmfv,refTree)); 5679566063dSJacob Faibussowitsch PetscCall(DMPlexGetGradientDM(dmfv, fv, &dmgrad)); 5689566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmfv,0,&cStart,&cEnd)); 56930602db0SMatthew G. Knepley nvecs = dim * (dim+1) / 2; 5709566063dSJacob Faibussowitsch PetscCall(DMPlexGetGeometryFVM(dmfv,NULL,&cellgeom,NULL)); 5719566063dSJacob Faibussowitsch PetscCall(VecGetDM(cellgeom,&dmCell)); 5729566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(cellgeom,&cgeom)); 5739566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmgrad,&grad)); 5749566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dmgrad,&locGrad)); 5759566063dSJacob Faibussowitsch PetscCall(DMPlexGetGhostCellStratum(dmgrad,&cEndInterior,NULL)); 576c4762a1bSJed Brown cEndInterior = (cEndInterior < 0) ? cEnd: cEndInterior; 577c4762a1bSJed Brown for (v = 0; v < nvecs; v++) { 578c4762a1bSJed Brown Vec locX; 579c4762a1bSJed Brown PetscInt c; 580c4762a1bSJed Brown PetscScalar trueGrad[3][3] = {{0.}}; 581c4762a1bSJed Brown const PetscScalar *gradArray; 582c4762a1bSJed Brown PetscReal maxDiff, maxDiffGlob; 583c4762a1bSJed Brown 5849566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dmfv,&locX)); 585c4762a1bSJed Brown /* get the local projection of the rigid body mode */ 586c4762a1bSJed Brown for (c = cStart; c < cEnd; c++) { 587c4762a1bSJed Brown PetscFVCellGeom *cg; 588c4762a1bSJed Brown PetscScalar cx[3] = {0.,0.,0.}; 589c4762a1bSJed Brown 5909566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); 59130602db0SMatthew G. Knepley if (v < dim) { 592c4762a1bSJed Brown cx[v] = 1.; 593c4762a1bSJed Brown } else { 59430602db0SMatthew G. Knepley PetscInt w = v - dim; 595c4762a1bSJed Brown 59630602db0SMatthew G. Knepley cx[(w + 1) % dim] = cg->centroid[(w + 2) % dim]; 59730602db0SMatthew G. Knepley cx[(w + 2) % dim] = -cg->centroid[(w + 1) % dim]; 598c4762a1bSJed Brown } 5999566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetClosure(dmfv,NULL,locX,c,cx,INSERT_ALL_VALUES)); 600c4762a1bSJed Brown } 601c4762a1bSJed Brown /* TODO: this isn't in any header */ 6029566063dSJacob Faibussowitsch PetscCall(DMPlexReconstructGradientsFVM(dmfv,locX,grad)); 6039566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dmgrad,grad,INSERT_VALUES,locGrad)); 6049566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dmgrad,grad,INSERT_VALUES,locGrad)); 6059566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(locGrad,&gradArray)); 606c4762a1bSJed Brown /* compare computed gradient to exact gradient */ 60730602db0SMatthew G. Knepley if (v >= dim) { 60830602db0SMatthew G. Knepley PetscInt w = v - dim; 609c4762a1bSJed Brown 61030602db0SMatthew G. Knepley trueGrad[(w + 1) % dim][(w + 2) % dim] = 1.; 61130602db0SMatthew G. Knepley trueGrad[(w + 2) % dim][(w + 1) % dim] = -1.; 612c4762a1bSJed Brown } 613c4762a1bSJed Brown maxDiff = 0.; 614c4762a1bSJed Brown for (c = cStart; c < cEndInterior; c++) { 615c4762a1bSJed Brown PetscScalar *compGrad; 616c4762a1bSJed Brown PetscInt i, j, k; 617c4762a1bSJed Brown PetscReal FrobDiff = 0.; 618c4762a1bSJed Brown 6199566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(dmgrad, c, gradArray, &compGrad)); 620c4762a1bSJed Brown 62130602db0SMatthew G. Knepley for (i = 0, k = 0; i < dim; i++) { 62230602db0SMatthew G. Knepley for (j = 0; j < dim; j++, k++) { 623c4762a1bSJed Brown PetscScalar diff = compGrad[k] - trueGrad[i][j]; 624c4762a1bSJed Brown FrobDiff += PetscRealPart(diff * PetscConj(diff)); 625c4762a1bSJed Brown } 626c4762a1bSJed Brown } 627c4762a1bSJed Brown FrobDiff = PetscSqrtReal(FrobDiff); 628c4762a1bSJed Brown maxDiff = PetscMax(maxDiff,FrobDiff); 629c4762a1bSJed Brown } 6309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&maxDiff,&maxDiffGlob,1,MPIU_REAL,MPIU_MAX,comm)); 631c4762a1bSJed Brown allVecMaxDiff = PetscMax(allVecMaxDiff,maxDiffGlob); 6329566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(locGrad,&gradArray)); 6339566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dmfv,&locX)); 634c4762a1bSJed Brown } 635c4762a1bSJed Brown if (allVecMaxDiff < fvTol) { 6369566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm),"Finite volume gradient reconstruction: PASS\n")); 637c4762a1bSJed Brown } else { 638*63a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm),"Finite volume gradient reconstruction: FAIL at tolerance %g with max difference %g\n",(double)fvTol,(double)allVecMaxDiff)); 639c4762a1bSJed Brown } 6409566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dmgrad,&locGrad)); 6419566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmgrad,&grad)); 6429566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(cellgeom,&cgeom)); 6439566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmfv)); 6449566063dSJacob Faibussowitsch PetscCall(PetscFVDestroy(&fv)); 645c4762a1bSJed Brown PetscFunctionReturn(0); 646c4762a1bSJed Brown } 647c4762a1bSJed Brown 648c4762a1bSJed Brown static PetscErrorCode ComputeError(DM dm, PetscErrorCode (**exactFuncs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), 649c4762a1bSJed Brown PetscErrorCode (**exactFuncDers)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), 650c4762a1bSJed Brown void **exactCtxs, PetscReal *error, PetscReal *errorDer, AppCtx *user) 651c4762a1bSJed Brown { 652c4762a1bSJed Brown Vec u; 653c4762a1bSJed Brown PetscReal n[3] = {1.0, 1.0, 1.0}; 654c4762a1bSJed Brown 655c4762a1bSJed Brown PetscFunctionBeginUser; 6569566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &u)); 657c4762a1bSJed Brown /* Project function into FE function space */ 6589566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, u)); 6599566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-projection_view")); 660c4762a1bSJed Brown /* Compare approximation to exact in L_2 */ 6619566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, exactCtxs, u, error)); 6629566063dSJacob Faibussowitsch PetscCall(DMComputeL2GradientDiff(dm, 0.0, exactFuncDers, exactCtxs, u, n, errorDer)); 6639566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &u)); 664c4762a1bSJed Brown PetscFunctionReturn(0); 665c4762a1bSJed Brown } 666c4762a1bSJed Brown 667c4762a1bSJed Brown static PetscErrorCode CheckFunctions(DM dm, PetscInt order, AppCtx *user) 668c4762a1bSJed Brown { 669c4762a1bSJed Brown PetscErrorCode (*exactFuncs[1]) (PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 670c4762a1bSJed Brown PetscErrorCode (*exactFuncDers[1]) (PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx); 671c4762a1bSJed Brown void *exactCtxs[3]; 672c4762a1bSJed Brown MPI_Comm comm; 673c4762a1bSJed Brown PetscReal error, errorDer, tol = PETSC_SMALL; 674c4762a1bSJed Brown 675c4762a1bSJed Brown PetscFunctionBeginUser; 676c4762a1bSJed Brown exactCtxs[0] = user; 677c4762a1bSJed Brown exactCtxs[1] = user; 678c4762a1bSJed Brown exactCtxs[2] = user; 6799566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 680c4762a1bSJed Brown /* Setup functions to approximate */ 681c4762a1bSJed Brown switch (order) { 682c4762a1bSJed Brown case 0: 683c4762a1bSJed Brown exactFuncs[0] = constant; 684c4762a1bSJed Brown exactFuncDers[0] = constantDer; 685c4762a1bSJed Brown break; 686c4762a1bSJed Brown case 1: 687c4762a1bSJed Brown exactFuncs[0] = linear; 688c4762a1bSJed Brown exactFuncDers[0] = linearDer; 689c4762a1bSJed Brown break; 690c4762a1bSJed Brown case 2: 691c4762a1bSJed Brown exactFuncs[0] = quadratic; 692c4762a1bSJed Brown exactFuncDers[0] = quadraticDer; 693c4762a1bSJed Brown break; 694c4762a1bSJed Brown case 3: 695c4762a1bSJed Brown exactFuncs[0] = cubic; 696c4762a1bSJed Brown exactFuncDers[0] = cubicDer; 697c4762a1bSJed Brown break; 698c4762a1bSJed Brown default: 699*63a3b9bcSJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Could not determine functions to test for order %" PetscInt_FMT, order); 700c4762a1bSJed Brown } 7019566063dSJacob Faibussowitsch PetscCall(ComputeError(dm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user)); 702c4762a1bSJed Brown /* Report result */ 703*63a3b9bcSJacob Faibussowitsch 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*63a3b9bcSJacob Faibussowitsch else PetscCall(PetscPrintf(comm, "Function tests pass for order %" PetscInt_FMT " at tolerance %g\n", order, (double)tol)); 705*63a3b9bcSJacob Faibussowitsch 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*63a3b9bcSJacob Faibussowitsch else PetscCall(PetscPrintf(comm, "Function tests pass for order %" PetscInt_FMT " derivatives at tolerance %g\n", order, (double)tol)); 707c4762a1bSJed Brown PetscFunctionReturn(0); 708c4762a1bSJed Brown } 709c4762a1bSJed Brown 710c4762a1bSJed Brown static PetscErrorCode CheckInterpolation(DM dm, PetscBool checkRestrict, PetscInt order, AppCtx *user) 711c4762a1bSJed Brown { 712c4762a1bSJed Brown PetscErrorCode (*exactFuncs[1]) (PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx); 713c4762a1bSJed Brown PetscErrorCode (*exactFuncDers[1]) (PetscInt, PetscReal, const PetscReal x[], const PetscReal n[], PetscInt, PetscScalar *u, void *ctx); 714c4762a1bSJed Brown PetscReal n[3] = {1.0, 1.0, 1.0}; 715c4762a1bSJed Brown void *exactCtxs[3]; 716c4762a1bSJed Brown DM rdm, idm, fdm; 717c4762a1bSJed Brown Mat Interp; 718c4762a1bSJed Brown Vec iu, fu, scaling; 719c4762a1bSJed Brown MPI_Comm comm; 72030602db0SMatthew G. Knepley PetscInt dim; 721c4762a1bSJed Brown PetscReal error, errorDer, tol = PETSC_SMALL; 722c4762a1bSJed Brown 723c4762a1bSJed Brown PetscFunctionBeginUser; 724c4762a1bSJed Brown exactCtxs[0] = user; 725c4762a1bSJed Brown exactCtxs[1] = user; 726c4762a1bSJed Brown exactCtxs[2] = user; 7279566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm,&comm)); 7289566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 7299566063dSJacob Faibussowitsch PetscCall(DMRefine(dm, comm, &rdm)); 7309566063dSJacob Faibussowitsch PetscCall(DMSetCoarseDM(rdm, dm)); 7319566063dSJacob Faibussowitsch PetscCall(DMPlexSetRegularRefinement(rdm, user->convRefine)); 73230602db0SMatthew G. Knepley if (user->tree) { 733c4762a1bSJed Brown DM refTree; 7349566063dSJacob Faibussowitsch PetscCall(DMPlexGetReferenceTree(dm,&refTree)); 7359566063dSJacob Faibussowitsch PetscCall(DMPlexSetReferenceTree(rdm,refTree)); 736c4762a1bSJed Brown } 7379566063dSJacob Faibussowitsch if (user->useDA) PetscCall(DMDASetVertexCoordinates(rdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 7389566063dSJacob Faibussowitsch PetscCall(SetupSection(rdm, user)); 739c4762a1bSJed Brown /* Setup functions to approximate */ 740c4762a1bSJed Brown switch (order) { 741c4762a1bSJed Brown case 0: 742c4762a1bSJed Brown exactFuncs[0] = constant; 743c4762a1bSJed Brown exactFuncDers[0] = constantDer; 744c4762a1bSJed Brown break; 745c4762a1bSJed Brown case 1: 746c4762a1bSJed Brown exactFuncs[0] = linear; 747c4762a1bSJed Brown exactFuncDers[0] = linearDer; 748c4762a1bSJed Brown break; 749c4762a1bSJed Brown case 2: 750c4762a1bSJed Brown exactFuncs[0] = quadratic; 751c4762a1bSJed Brown exactFuncDers[0] = quadraticDer; 752c4762a1bSJed Brown break; 753c4762a1bSJed Brown case 3: 754c4762a1bSJed Brown exactFuncs[0] = cubic; 755c4762a1bSJed Brown exactFuncDers[0] = cubicDer; 756c4762a1bSJed Brown break; 757c4762a1bSJed Brown default: 758*63a3b9bcSJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Could not determine functions to test for dimension %" PetscInt_FMT " order %" PetscInt_FMT, dim, order); 759c4762a1bSJed Brown } 760c4762a1bSJed Brown idm = checkRestrict ? rdm : dm; 761c4762a1bSJed Brown fdm = checkRestrict ? dm : rdm; 7629566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(idm, &iu)); 7639566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(fdm, &fu)); 7649566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dm, user)); 7659566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(rdm, user)); 7669566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dm, rdm, &Interp, &scaling)); 767c4762a1bSJed Brown /* Project function into initial FE function space */ 7689566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(idm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, iu)); 769c4762a1bSJed Brown /* Interpolate function into final FE function space */ 7709566063dSJacob Faibussowitsch if (checkRestrict) {PetscCall(MatRestrict(Interp, iu, fu));PetscCall(VecPointwiseMult(fu, scaling, fu));} 7719566063dSJacob Faibussowitsch else PetscCall(MatInterpolate(Interp, iu, fu)); 772c4762a1bSJed Brown /* Compare approximation to exact in L_2 */ 7739566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(fdm, 0.0, exactFuncs, exactCtxs, fu, &error)); 7749566063dSJacob Faibussowitsch PetscCall(DMComputeL2GradientDiff(fdm, 0.0, exactFuncDers, exactCtxs, fu, n, &errorDer)); 775c4762a1bSJed Brown /* Report result */ 776*63a3b9bcSJacob Faibussowitsch 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*63a3b9bcSJacob Faibussowitsch else PetscCall(PetscPrintf(comm, "Interpolation tests pass for order %" PetscInt_FMT " at tolerance %g\n", order, (double)tol)); 778*63a3b9bcSJacob Faibussowitsch 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*63a3b9bcSJacob Faibussowitsch else PetscCall(PetscPrintf(comm, "Interpolation tests pass for order %" PetscInt_FMT " derivatives at tolerance %g\n", order, (double)tol)); 7809566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(idm, &iu)); 7819566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(fdm, &fu)); 7829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Interp)); 7839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&scaling)); 7849566063dSJacob Faibussowitsch PetscCall(DMDestroy(&rdm)); 785c4762a1bSJed Brown PetscFunctionReturn(0); 786c4762a1bSJed Brown } 787c4762a1bSJed Brown 788c4762a1bSJed Brown static PetscErrorCode CheckConvergence(DM dm, PetscInt Nr, AppCtx *user) 789c4762a1bSJed Brown { 790c4762a1bSJed Brown DM odm = dm, rdm = NULL, cdm = NULL; 791c4762a1bSJed Brown PetscErrorCode (*exactFuncs[1]) (PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) = {trig}; 792c4762a1bSJed Brown PetscErrorCode (*exactFuncDers[1]) (PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) = {trigDer}; 793c4762a1bSJed Brown void *exactCtxs[3]; 794c4762a1bSJed Brown PetscInt r, c, cStart, cEnd; 795c4762a1bSJed Brown PetscReal errorOld, errorDerOld, error, errorDer, rel, len, lenOld; 796c4762a1bSJed Brown double p; 797c4762a1bSJed Brown 798c4762a1bSJed Brown PetscFunctionBeginUser; 799c4762a1bSJed Brown if (!user->convergence) PetscFunctionReturn(0); 800c4762a1bSJed Brown exactCtxs[0] = user; 801c4762a1bSJed Brown exactCtxs[1] = user; 802c4762a1bSJed Brown exactCtxs[2] = user; 8039566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject) odm)); 804c4762a1bSJed Brown if (!user->convRefine) { 805c4762a1bSJed Brown for (r = 0; r < Nr; ++r) { 8069566063dSJacob Faibussowitsch PetscCall(DMRefine(odm, PetscObjectComm((PetscObject) dm), &rdm)); 8079566063dSJacob Faibussowitsch PetscCall(DMDestroy(&odm)); 808c4762a1bSJed Brown odm = rdm; 809c4762a1bSJed Brown } 8109566063dSJacob Faibussowitsch PetscCall(SetupSection(odm, user)); 811c4762a1bSJed Brown } 8129566063dSJacob Faibussowitsch PetscCall(ComputeError(odm, exactFuncs, exactFuncDers, exactCtxs, &errorOld, &errorDerOld, user)); 813c4762a1bSJed Brown if (user->convRefine) { 814c4762a1bSJed Brown for (r = 0; r < Nr; ++r) { 8159566063dSJacob Faibussowitsch PetscCall(DMRefine(odm, PetscObjectComm((PetscObject) dm), &rdm)); 8169566063dSJacob Faibussowitsch if (user->useDA) PetscCall(DMDASetVertexCoordinates(rdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 8179566063dSJacob Faibussowitsch PetscCall(SetupSection(rdm, user)); 8189566063dSJacob Faibussowitsch PetscCall(ComputeError(rdm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user)); 819c4762a1bSJed Brown p = PetscLog2Real(errorOld/error); 820*63a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject) dm), "Function convergence rate at refinement %" PetscInt_FMT ": %.2f\n", r, (double)p)); 821c4762a1bSJed Brown p = PetscLog2Real(errorDerOld/errorDer); 822*63a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject) dm), "Derivative convergence rate at refinement %" PetscInt_FMT ": %.2f\n", r, (double)p)); 8239566063dSJacob Faibussowitsch PetscCall(DMDestroy(&odm)); 824c4762a1bSJed Brown odm = rdm; 825c4762a1bSJed Brown errorOld = error; 826c4762a1bSJed Brown errorDerOld = errorDer; 827c4762a1bSJed Brown } 828c4762a1bSJed Brown } else { 8299566063dSJacob Faibussowitsch /* PetscCall(ComputeLongestEdge(dm, &lenOld)); */ 8309566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(odm, 0, &cStart, &cEnd)); 831c4762a1bSJed Brown lenOld = cEnd - cStart; 832c4762a1bSJed Brown for (c = 0; c < Nr; ++c) { 8339566063dSJacob Faibussowitsch PetscCall(DMCoarsen(odm, PetscObjectComm((PetscObject) dm), &cdm)); 8349566063dSJacob Faibussowitsch if (user->useDA) PetscCall(DMDASetVertexCoordinates(cdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 8359566063dSJacob Faibussowitsch PetscCall(SetupSection(cdm, user)); 8369566063dSJacob Faibussowitsch PetscCall(ComputeError(cdm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user)); 8379566063dSJacob Faibussowitsch /* PetscCall(ComputeLongestEdge(cdm, &len)); */ 8389566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(cdm, 0, &cStart, &cEnd)); 839c4762a1bSJed Brown len = cEnd - cStart; 840c4762a1bSJed Brown rel = error/errorOld; 841c4762a1bSJed Brown p = PetscLogReal(rel) / PetscLogReal(lenOld / len); 842*63a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject) dm), "Function convergence rate at coarsening %" PetscInt_FMT ": %.2f\n", c, (double)p)); 843c4762a1bSJed Brown rel = errorDer/errorDerOld; 844c4762a1bSJed Brown p = PetscLogReal(rel) / PetscLogReal(lenOld / len); 845*63a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject) dm), "Derivative convergence rate at coarsening %" PetscInt_FMT ": %.2f\n", c, (double)p)); 8469566063dSJacob Faibussowitsch PetscCall(DMDestroy(&odm)); 847c4762a1bSJed Brown odm = cdm; 848c4762a1bSJed Brown errorOld = error; 849c4762a1bSJed Brown errorDerOld = errorDer; 850c4762a1bSJed Brown lenOld = len; 851c4762a1bSJed Brown } 852c4762a1bSJed Brown } 8539566063dSJacob Faibussowitsch PetscCall(DMDestroy(&odm)); 854c4762a1bSJed Brown PetscFunctionReturn(0); 855c4762a1bSJed Brown } 856c4762a1bSJed Brown 857c4762a1bSJed Brown int main(int argc, char **argv) 858c4762a1bSJed Brown { 859c4762a1bSJed Brown DM dm; 860c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 86130602db0SMatthew G. Knepley PetscInt dim = 2; 86230602db0SMatthew G. Knepley PetscBool simplex = PETSC_FALSE; 863c4762a1bSJed Brown 8649566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 8659566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 8669566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 86730602db0SMatthew G. Knepley if (!user.useDA) { 8689566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 8699566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 87030602db0SMatthew G. Knepley } 8719566063dSJacob Faibussowitsch PetscCall(DMPlexMetricSetFromOptions(dm)); 87230602db0SMatthew G. Knepley user.numComponents = user.numComponents < 0 ? dim : user.numComponents; 8739566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_WORLD, dim, user.numComponents, simplex, NULL, user.qorder, &user.fe)); 8749566063dSJacob Faibussowitsch PetscCall(SetupSection(dm, &user)); 8759566063dSJacob Faibussowitsch if (user.testFEjacobian) PetscCall(TestFEJacobian(dm, &user)); 8769566063dSJacob Faibussowitsch if (user.testFVgrad) PetscCall(TestFVGrad(dm, &user)); 8779566063dSJacob Faibussowitsch if (user.testInjector) PetscCall(TestInjector(dm, &user)); 8789566063dSJacob Faibussowitsch PetscCall(CheckFunctions(dm, user.porder, &user)); 879c4762a1bSJed Brown { 880c4762a1bSJed Brown PetscDualSpace dsp; 881c4762a1bSJed Brown PetscInt k; 882c4762a1bSJed Brown 8839566063dSJacob Faibussowitsch PetscCall(PetscFEGetDualSpace(user.fe, &dsp)); 8849566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDeRahm(dsp, &k)); 88530602db0SMatthew G. Knepley if (dim == 2 && simplex == PETSC_TRUE && user.tree == PETSC_FALSE && k == 0) { 8869566063dSJacob Faibussowitsch PetscCall(CheckInterpolation(dm, PETSC_FALSE, user.porder, &user)); 8879566063dSJacob Faibussowitsch PetscCall(CheckInterpolation(dm, PETSC_TRUE, user.porder, &user)); 888c4762a1bSJed Brown } 889c4762a1bSJed Brown } 8909566063dSJacob Faibussowitsch PetscCall(CheckConvergence(dm, 3, &user)); 8919566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&user.fe)); 8929566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 8939566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 894b122ec5aSJacob Faibussowitsch return 0; 895c4762a1bSJed Brown } 896c4762a1bSJed Brown 897c4762a1bSJed Brown /*TEST 898c4762a1bSJed Brown 899c4762a1bSJed Brown test: 900c4762a1bSJed Brown suffix: 1 901c4762a1bSJed Brown requires: triangle 902c4762a1bSJed Brown 903c4762a1bSJed Brown # 2D P_1 on a triangle 904c4762a1bSJed Brown test: 905c4762a1bSJed Brown suffix: p1_2d_0 906c4762a1bSJed Brown requires: triangle 907c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -convergence 908c4762a1bSJed Brown test: 909c4762a1bSJed Brown suffix: p1_2d_1 910c4762a1bSJed Brown requires: triangle 911c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -porder 1 912c4762a1bSJed Brown test: 913c4762a1bSJed Brown suffix: p1_2d_2 914c4762a1bSJed Brown requires: triangle 915c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -porder 2 916c4762a1bSJed Brown test: 917c4762a1bSJed Brown suffix: p1_2d_3 918e788613bSJoe Wallwork requires: triangle mmg 919c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0 920c4762a1bSJed Brown test: 921c4762a1bSJed Brown suffix: p1_2d_4 922e788613bSJoe Wallwork requires: triangle mmg 923c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0 924c4762a1bSJed Brown test: 925c4762a1bSJed Brown suffix: p1_2d_5 926e788613bSJoe Wallwork requires: triangle mmg 927c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0 928c4762a1bSJed Brown 929c4762a1bSJed Brown # 3D P_1 on a tetrahedron 930c4762a1bSJed Brown test: 931c4762a1bSJed Brown suffix: p1_3d_0 932c4762a1bSJed Brown requires: ctetgen 93330602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -convergence 934c4762a1bSJed Brown test: 935c4762a1bSJed Brown suffix: p1_3d_1 936c4762a1bSJed Brown requires: ctetgen 93730602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -porder 1 938c4762a1bSJed Brown test: 939c4762a1bSJed Brown suffix: p1_3d_2 940c4762a1bSJed Brown requires: ctetgen 94130602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -porder 2 942c4762a1bSJed Brown test: 943c4762a1bSJed Brown suffix: p1_3d_3 944e788613bSJoe Wallwork requires: ctetgen mmg 94530602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0 946c4762a1bSJed Brown test: 947c4762a1bSJed Brown suffix: p1_3d_4 948e788613bSJoe Wallwork requires: ctetgen mmg 94930602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0 950c4762a1bSJed Brown test: 951c4762a1bSJed Brown suffix: p1_3d_5 952e788613bSJoe Wallwork requires: ctetgen mmg 95330602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0 954c4762a1bSJed Brown 955c4762a1bSJed Brown # 2D P_2 on a triangle 956c4762a1bSJed Brown test: 957c4762a1bSJed Brown suffix: p2_2d_0 958c4762a1bSJed Brown requires: triangle 959c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 2 -convergence 960c4762a1bSJed Brown test: 961c4762a1bSJed Brown suffix: p2_2d_1 962c4762a1bSJed Brown requires: triangle 963c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 2 -porder 1 964c4762a1bSJed Brown test: 965c4762a1bSJed Brown suffix: p2_2d_2 966c4762a1bSJed Brown requires: triangle 967c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 2 -porder 2 968c4762a1bSJed Brown test: 969c4762a1bSJed Brown suffix: p2_2d_3 970e788613bSJoe Wallwork requires: triangle mmg 971c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -convergence -conv_refine 0 972c4762a1bSJed Brown test: 973c4762a1bSJed Brown suffix: p2_2d_4 974e788613bSJoe Wallwork requires: triangle mmg 975c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 1 -conv_refine 0 976c4762a1bSJed Brown test: 977c4762a1bSJed Brown suffix: p2_2d_5 978e788613bSJoe Wallwork requires: triangle mmg 979c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 2 -conv_refine 0 980c4762a1bSJed Brown 981c4762a1bSJed Brown # 3D P_2 on a tetrahedron 982c4762a1bSJed Brown test: 983c4762a1bSJed Brown suffix: p2_3d_0 984c4762a1bSJed Brown requires: ctetgen 98530602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -convergence 986c4762a1bSJed Brown test: 987c4762a1bSJed Brown suffix: p2_3d_1 988c4762a1bSJed Brown requires: ctetgen 98930602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -porder 1 990c4762a1bSJed Brown test: 991c4762a1bSJed Brown suffix: p2_3d_2 992c4762a1bSJed Brown requires: ctetgen 99330602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -porder 2 994c4762a1bSJed Brown test: 995c4762a1bSJed Brown suffix: p2_3d_3 996e788613bSJoe Wallwork requires: ctetgen mmg 99730602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -convergence -conv_refine 0 998c4762a1bSJed Brown test: 999c4762a1bSJed Brown suffix: p2_3d_4 1000e788613bSJoe Wallwork requires: ctetgen mmg 100130602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 1 -conv_refine 0 1002c4762a1bSJed Brown test: 1003c4762a1bSJed Brown suffix: p2_3d_5 1004e788613bSJoe Wallwork requires: ctetgen mmg 100530602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 2 -conv_refine 0 1006c4762a1bSJed Brown 1007c4762a1bSJed Brown # 2D Q_1 on a quadrilaterial DA 1008c4762a1bSJed Brown test: 1009c4762a1bSJed Brown suffix: q1_2d_da_0 101099a192c5SJunchao Zhang requires: broken 101130602db0SMatthew G. Knepley args: -use_da 1 -petscspace_degree 1 -qorder 1 -convergence 1012c4762a1bSJed Brown test: 1013c4762a1bSJed Brown suffix: q1_2d_da_1 101499a192c5SJunchao Zhang requires: broken 101530602db0SMatthew G. Knepley args: -use_da 1 -petscspace_degree 1 -qorder 1 -porder 1 1016c4762a1bSJed Brown test: 1017c4762a1bSJed Brown suffix: q1_2d_da_2 101899a192c5SJunchao Zhang requires: broken 101930602db0SMatthew G. Knepley args: -use_da 1 -petscspace_degree 1 -qorder 1 -porder 2 1020c4762a1bSJed Brown 1021c4762a1bSJed Brown # 2D Q_1 on a quadrilaterial Plex 1022c4762a1bSJed Brown test: 1023c4762a1bSJed Brown suffix: q1_2d_plex_0 102430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -convergence 1025c4762a1bSJed Brown test: 1026c4762a1bSJed Brown suffix: q1_2d_plex_1 102730602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 1 1028c4762a1bSJed Brown test: 1029c4762a1bSJed Brown suffix: q1_2d_plex_2 103030602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 2 1031c4762a1bSJed Brown test: 1032c4762a1bSJed Brown suffix: q1_2d_plex_3 103330602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 1 -shear_coords 1034c4762a1bSJed Brown test: 1035c4762a1bSJed Brown suffix: q1_2d_plex_4 103630602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 2 -shear_coords 1037c4762a1bSJed Brown test: 1038c4762a1bSJed Brown suffix: q1_2d_plex_5 103930602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 0 -non_affine_coords -convergence 1040c4762a1bSJed Brown test: 1041c4762a1bSJed Brown suffix: q1_2d_plex_6 104230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 1 -non_affine_coords -convergence 1043c4762a1bSJed Brown test: 1044c4762a1bSJed Brown suffix: q1_2d_plex_7 104530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 2 -non_affine_coords -convergence 1046c4762a1bSJed Brown 1047c4762a1bSJed Brown # 2D Q_2 on a quadrilaterial 1048c4762a1bSJed Brown test: 1049c4762a1bSJed Brown suffix: q2_2d_plex_0 105030602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -convergence 1051c4762a1bSJed Brown test: 1052c4762a1bSJed Brown suffix: q2_2d_plex_1 105330602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 1054c4762a1bSJed Brown test: 1055c4762a1bSJed Brown suffix: q2_2d_plex_2 105630602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 1057c4762a1bSJed Brown test: 1058c4762a1bSJed Brown suffix: q2_2d_plex_3 105930602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 -shear_coords 1060c4762a1bSJed Brown test: 1061c4762a1bSJed Brown suffix: q2_2d_plex_4 106230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 -shear_coords 1063c4762a1bSJed Brown test: 1064c4762a1bSJed Brown suffix: q2_2d_plex_5 106530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 0 -non_affine_coords -convergence 1066c4762a1bSJed Brown test: 1067c4762a1bSJed Brown suffix: q2_2d_plex_6 106830602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 1 -non_affine_coords -convergence 1069c4762a1bSJed Brown test: 1070c4762a1bSJed Brown suffix: q2_2d_plex_7 107130602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 2 -non_affine_coords -convergence 1072c4762a1bSJed Brown 1073c4762a1bSJed Brown # 2D P_3 on a triangle 1074c4762a1bSJed Brown test: 1075c4762a1bSJed Brown suffix: p3_2d_0 1076c4762a1bSJed Brown requires: triangle !single 1077c4762a1bSJed Brown args: -petscspace_degree 3 -qorder 3 -convergence 1078c4762a1bSJed Brown test: 1079c4762a1bSJed Brown suffix: p3_2d_1 1080c4762a1bSJed Brown requires: triangle !single 1081c4762a1bSJed Brown args: -petscspace_degree 3 -qorder 3 -porder 1 1082c4762a1bSJed Brown test: 1083c4762a1bSJed Brown suffix: p3_2d_2 1084c4762a1bSJed Brown requires: triangle !single 1085c4762a1bSJed Brown args: -petscspace_degree 3 -qorder 3 -porder 2 1086c4762a1bSJed Brown test: 1087c4762a1bSJed Brown suffix: p3_2d_3 1088c4762a1bSJed Brown requires: triangle !single 1089c4762a1bSJed Brown args: -petscspace_degree 3 -qorder 3 -porder 3 1090c4762a1bSJed Brown test: 1091c4762a1bSJed Brown suffix: p3_2d_4 1092e788613bSJoe Wallwork requires: triangle mmg 1093c4762a1bSJed Brown args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -convergence -conv_refine 0 1094c4762a1bSJed Brown test: 1095c4762a1bSJed Brown suffix: p3_2d_5 1096e788613bSJoe Wallwork requires: triangle mmg 1097c4762a1bSJed Brown args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -porder 1 -conv_refine 0 1098c4762a1bSJed Brown test: 1099c4762a1bSJed Brown suffix: p3_2d_6 1100e788613bSJoe Wallwork requires: triangle mmg 1101c4762a1bSJed Brown args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -porder 3 -conv_refine 0 1102c4762a1bSJed Brown 1103c4762a1bSJed Brown # 2D Q_3 on a quadrilaterial 1104c4762a1bSJed Brown test: 1105c4762a1bSJed Brown suffix: q3_2d_0 110699a192c5SJunchao Zhang requires: !single 110730602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -convergence 1108c4762a1bSJed Brown test: 1109c4762a1bSJed Brown suffix: q3_2d_1 111099a192c5SJunchao Zhang requires: !single 111130602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 1 1112c4762a1bSJed Brown test: 1113c4762a1bSJed Brown suffix: q3_2d_2 111499a192c5SJunchao Zhang requires: !single 111530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 2 1116c4762a1bSJed Brown test: 1117c4762a1bSJed Brown suffix: q3_2d_3 111899a192c5SJunchao Zhang requires: !single 111930602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 3 1120c4762a1bSJed Brown 1121c4762a1bSJed Brown # 2D P_1disc on a triangle/quadrilateral 1122c4762a1bSJed Brown test: 1123c4762a1bSJed Brown suffix: p1d_2d_0 1124c4762a1bSJed Brown requires: triangle 1125c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -convergence 1126c4762a1bSJed Brown test: 1127c4762a1bSJed Brown suffix: p1d_2d_1 1128c4762a1bSJed Brown requires: triangle 1129c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 1 1130c4762a1bSJed Brown test: 1131c4762a1bSJed Brown suffix: p1d_2d_2 1132c4762a1bSJed Brown requires: triangle 1133c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 2 1134c4762a1bSJed Brown test: 1135c4762a1bSJed Brown suffix: p1d_2d_3 1136c4762a1bSJed Brown requires: triangle 113730602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -convergence 1138c4762a1bSJed Brown filter: sed -e "s/convergence rate at refinement 0: 2/convergence rate at refinement 0: 1.9/g" 1139c4762a1bSJed Brown test: 1140c4762a1bSJed Brown suffix: p1d_2d_4 1141c4762a1bSJed Brown requires: triangle 114230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 1 1143c4762a1bSJed Brown test: 1144c4762a1bSJed Brown suffix: p1d_2d_5 1145c4762a1bSJed Brown requires: triangle 114630602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 2 1147c4762a1bSJed Brown 1148c4762a1bSJed Brown # 2D BDM_1 on a triangle 1149c4762a1bSJed Brown test: 1150c4762a1bSJed Brown suffix: bdm1_2d_0 1151c4762a1bSJed Brown requires: triangle 1152c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_type bdm \ 1153c4762a1bSJed Brown -num_comp 2 -qorder 1 -convergence 1154c4762a1bSJed Brown test: 1155c4762a1bSJed Brown suffix: bdm1_2d_1 1156c4762a1bSJed Brown requires: triangle 1157c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_type bdm \ 1158c4762a1bSJed Brown -num_comp 2 -qorder 1 -porder 1 1159c4762a1bSJed Brown test: 1160c4762a1bSJed Brown suffix: bdm1_2d_2 1161c4762a1bSJed Brown requires: triangle 1162c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_type bdm \ 1163c4762a1bSJed Brown -num_comp 2 -qorder 1 -porder 2 1164c4762a1bSJed Brown 1165c4762a1bSJed Brown # 2D BDM_1 on a quadrilateral 1166c4762a1bSJed Brown test: 1167c4762a1bSJed Brown suffix: bdm1q_2d_0 1168c4762a1bSJed Brown requires: triangle 1169c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_type bdm \ 11703f27d899SToby Isaac -petscdualspace_lagrange_tensor 1 \ 117130602db0SMatthew G. Knepley -dm_plex_simplex 0 -num_comp 2 -qorder 1 -convergence 1172c4762a1bSJed Brown test: 1173c4762a1bSJed Brown suffix: bdm1q_2d_1 1174c4762a1bSJed Brown requires: triangle 1175c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_type bdm \ 11763f27d899SToby Isaac -petscdualspace_lagrange_tensor 1 \ 117730602db0SMatthew G. Knepley -dm_plex_simplex 0 -num_comp 2 -qorder 1 -porder 1 1178c4762a1bSJed Brown test: 1179c4762a1bSJed Brown suffix: bdm1q_2d_2 1180c4762a1bSJed Brown requires: triangle 1181c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_type bdm \ 11823f27d899SToby Isaac -petscdualspace_lagrange_tensor 1 \ 118330602db0SMatthew G. Knepley -dm_plex_simplex 0 -num_comp 2 -qorder 1 -porder 2 1184c4762a1bSJed Brown 1185c4762a1bSJed Brown # Test high order quadrature 1186c4762a1bSJed Brown test: 1187c4762a1bSJed Brown suffix: p1_quad_2 1188c4762a1bSJed Brown requires: triangle 1189c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 2 -porder 1 1190c4762a1bSJed Brown test: 1191c4762a1bSJed Brown suffix: p1_quad_5 1192c4762a1bSJed Brown requires: triangle 1193c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 5 -porder 1 1194c4762a1bSJed Brown test: 1195c4762a1bSJed Brown suffix: p2_quad_3 1196c4762a1bSJed Brown requires: triangle 1197c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 3 -porder 2 1198c4762a1bSJed Brown test: 1199c4762a1bSJed Brown suffix: p2_quad_5 1200c4762a1bSJed Brown requires: triangle 1201c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 5 -porder 2 1202c4762a1bSJed Brown test: 1203c4762a1bSJed Brown suffix: q1_quad_2 120430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 2 -porder 1 1205c4762a1bSJed Brown test: 1206c4762a1bSJed Brown suffix: q1_quad_5 120730602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 5 -porder 1 1208c4762a1bSJed Brown test: 1209c4762a1bSJed Brown suffix: q2_quad_3 121030602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 3 -porder 1 1211c4762a1bSJed Brown test: 1212c4762a1bSJed Brown suffix: q2_quad_5 121330602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 5 -porder 1 1214c4762a1bSJed Brown 1215c4762a1bSJed Brown # Nonconforming tests 1216c4762a1bSJed Brown test: 1217c4762a1bSJed Brown suffix: constraints 121830602db0SMatthew G. Knepley args: -dm_coord_space 0 -dm_plex_simplex 0 -petscspace_type tensor -petscspace_degree 1 -qorder 0 -constraints 1219c4762a1bSJed Brown test: 1220c4762a1bSJed Brown suffix: nonconforming_tensor_2 1221c4762a1bSJed Brown nsize: 4 122230602db0SMatthew G. Knepley 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 1223c4762a1bSJed Brown test: 1224c4762a1bSJed Brown suffix: nonconforming_tensor_3 1225c4762a1bSJed Brown nsize: 4 122630602db0SMatthew G. Knepley 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 1227c4762a1bSJed Brown test: 1228c4762a1bSJed Brown suffix: nonconforming_tensor_2_fv 1229c4762a1bSJed Brown nsize: 4 123030602db0SMatthew G. Knepley args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -dm_plex_simplex 0 -num_comp 2 1231c4762a1bSJed Brown test: 1232c4762a1bSJed Brown suffix: nonconforming_tensor_3_fv 1233c4762a1bSJed Brown nsize: 4 123430602db0SMatthew G. Knepley 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 1235c4762a1bSJed Brown test: 1236c4762a1bSJed Brown suffix: nonconforming_tensor_2_hi 1237c4762a1bSJed Brown requires: !single 1238c4762a1bSJed Brown nsize: 4 123930602db0SMatthew G. Knepley 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 1240c4762a1bSJed Brown test: 1241c4762a1bSJed Brown suffix: nonconforming_tensor_3_hi 1242c4762a1bSJed Brown requires: !single skip 1243c4762a1bSJed Brown nsize: 4 124430602db0SMatthew G. Knepley 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 1245c4762a1bSJed Brown test: 1246c4762a1bSJed Brown suffix: nonconforming_simplex_2 1247c4762a1bSJed Brown requires: triangle 1248c4762a1bSJed Brown nsize: 4 124930602db0SMatthew G. Knepley 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 1250c4762a1bSJed Brown test: 1251c4762a1bSJed Brown suffix: nonconforming_simplex_2_hi 1252c4762a1bSJed Brown requires: triangle !single 1253c4762a1bSJed Brown nsize: 4 125430602db0SMatthew G. Knepley args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_max_projection_height 1 -petscspace_degree 4 -qorder 4 1255c4762a1bSJed Brown test: 1256c4762a1bSJed Brown suffix: nonconforming_simplex_2_fv 1257c4762a1bSJed Brown requires: triangle 1258c4762a1bSJed Brown nsize: 4 125930602db0SMatthew G. Knepley args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -num_comp 2 1260c4762a1bSJed Brown test: 1261c4762a1bSJed Brown suffix: nonconforming_simplex_3 1262c4762a1bSJed Brown requires: ctetgen 1263c4762a1bSJed Brown nsize: 4 126430602db0SMatthew G. Knepley 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 1265c4762a1bSJed Brown test: 1266c4762a1bSJed Brown suffix: nonconforming_simplex_3_hi 1267c4762a1bSJed Brown requires: ctetgen skip 1268c4762a1bSJed Brown nsize: 4 126930602db0SMatthew G. Knepley 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 1270c4762a1bSJed Brown test: 1271c4762a1bSJed Brown suffix: nonconforming_simplex_3_fv 1272c4762a1bSJed Brown requires: ctetgen 1273c4762a1bSJed Brown nsize: 4 127430602db0SMatthew G. Knepley args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -dm_plex_dim 3 -num_comp 3 1275c4762a1bSJed Brown 1276d21efd2eSMatthew G. Knepley # 3D WXY on a triangular prism 1277d21efd2eSMatthew G. Knepley test: 1278d21efd2eSMatthew G. Knepley suffix: wxy_0 1279d21efd2eSMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -qorder 2 -porder 0 \ 1280e239af90SMatthew G. Knepley -petscspace_type sum \ 1281e239af90SMatthew G. Knepley -petscspace_variables 3 \ 1282e239af90SMatthew G. Knepley -petscspace_components 3 \ 1283e239af90SMatthew G. Knepley -petscspace_sum_spaces 2 \ 1284e239af90SMatthew G. Knepley -petscspace_sum_concatenate false \ 1285e239af90SMatthew G. Knepley -sumcomp_0_petscspace_variables 3 \ 1286e239af90SMatthew G. Knepley -sumcomp_0_petscspace_components 3 \ 1287e239af90SMatthew G. Knepley -sumcomp_0_petscspace_degree 1 \ 1288e239af90SMatthew G. Knepley -sumcomp_1_petscspace_variables 3 \ 1289e239af90SMatthew G. Knepley -sumcomp_1_petscspace_components 3 \ 1290e239af90SMatthew G. Knepley -sumcomp_1_petscspace_type wxy \ 1291e239af90SMatthew G. Knepley -petscdualspace_refcell triangular_prism \ 1292e239af90SMatthew G. Knepley -petscdualspace_form_degree 0 \ 1293e239af90SMatthew G. Knepley -petscdualspace_order 1 \ 1294e239af90SMatthew G. Knepley -petscdualspace_components 3 1295d21efd2eSMatthew G. Knepley 1296c4762a1bSJed Brown TEST*/ 1297c4762a1bSJed Brown 1298c4762a1bSJed Brown /* 1299c4762a1bSJed Brown # 2D Q_2 on a quadrilaterial Plex 1300c4762a1bSJed Brown test: 1301c4762a1bSJed Brown suffix: q2_2d_plex_0 130230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -convergence 1303c4762a1bSJed Brown test: 1304c4762a1bSJed Brown suffix: q2_2d_plex_1 130530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 1306c4762a1bSJed Brown test: 1307c4762a1bSJed Brown suffix: q2_2d_plex_2 130830602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 1309c4762a1bSJed Brown test: 1310c4762a1bSJed Brown suffix: q2_2d_plex_3 131130602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 -shear_coords 1312c4762a1bSJed Brown test: 1313c4762a1bSJed Brown suffix: q2_2d_plex_4 131430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 -shear_coords 1315c4762a1bSJed Brown test: 1316c4762a1bSJed Brown suffix: q2_2d_plex_5 131730602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 0 -non_affine_coords 1318c4762a1bSJed Brown test: 1319c4762a1bSJed Brown suffix: q2_2d_plex_6 132030602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 1 -non_affine_coords 1321c4762a1bSJed Brown test: 1322c4762a1bSJed Brown suffix: q2_2d_plex_7 132330602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 2 -non_affine_coords 1324c4762a1bSJed Brown 1325c4762a1bSJed Brown test: 1326c4762a1bSJed Brown suffix: p1d_2d_6 1327e788613bSJoe Wallwork requires: mmg 1328c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0 1329c4762a1bSJed Brown test: 1330c4762a1bSJed Brown suffix: p1d_2d_7 1331e788613bSJoe Wallwork requires: mmg 1332c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0 1333c4762a1bSJed Brown test: 1334c4762a1bSJed Brown suffix: p1d_2d_8 1335e788613bSJoe Wallwork requires: mmg 1336c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0 1337c4762a1bSJed Brown */ 1338