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 PetscErrorCode ierr; 115c4762a1bSJed Brown 116c4762a1bSJed Brown PetscFunctionBeginUser; 11730602db0SMatthew G. Knepley options->useDA = PETSC_FALSE; 118c4762a1bSJed Brown options->shearCoords = PETSC_FALSE; 119c4762a1bSJed Brown options->nonaffineCoords = PETSC_FALSE; 120c4762a1bSJed Brown options->qorder = 0; 121c4762a1bSJed Brown options->numComponents = PETSC_DEFAULT; 122c4762a1bSJed Brown options->porder = 0; 123c4762a1bSJed Brown options->convergence = PETSC_FALSE; 124c4762a1bSJed Brown options->convRefine = PETSC_TRUE; 125c4762a1bSJed Brown options->constraints = PETSC_FALSE; 126c4762a1bSJed Brown options->tree = PETSC_FALSE; 127c4762a1bSJed Brown options->treeCell = 0; 128c4762a1bSJed Brown options->testFEjacobian = PETSC_FALSE; 129c4762a1bSJed Brown options->testFVgrad = PETSC_FALSE; 130c4762a1bSJed Brown options->testInjector = PETSC_FALSE; 131c4762a1bSJed Brown options->constants[0] = 1.0; 132c4762a1bSJed Brown options->constants[1] = 2.0; 133c4762a1bSJed Brown options->constants[2] = 3.0; 134c4762a1bSJed Brown 135c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Projection Test Options", "DMPlex");CHKERRQ(ierr); 1365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-use_da", "Flag for DMDA mesh", "ex3.c", options->useDA, &options->useDA, NULL)); 1375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-shear_coords", "Transform coordinates with a shear", "ex3.c", options->shearCoords, &options->shearCoords, NULL)); 1385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-non_affine_coords", "Transform coordinates with a non-affine transform", "ex3.c", options->nonaffineCoords, &options->nonaffineCoords, NULL)); 1395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBoundedInt("-qorder", "The quadrature order", "ex3.c", options->qorder, &options->qorder, NULL,0)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBoundedInt("-num_comp", "The number of field components", "ex3.c", options->numComponents, &options->numComponents, NULL,PETSC_DEFAULT)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBoundedInt("-porder", "The order of polynomials to test", "ex3.c", options->porder, &options->porder, NULL,0)); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-convergence", "Check the convergence rate", "ex3.c", options->convergence, &options->convergence, NULL)); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-conv_refine", "Use refinement for the convergence rate", "ex3.c", options->convRefine, &options->convRefine, NULL)); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-constraints", "Test local constraints (serial only)", "ex3.c", options->constraints, &options->constraints, NULL)); 1455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-tree", "Test tree routines", "ex3.c", options->tree, &options->tree, NULL)); 1465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBoundedInt("-tree_cell", "cell to refine in tree test", "ex3.c", options->treeCell, &options->treeCell, NULL,0)); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-test_fe_jacobian", "Test finite element Jacobian assembly", "ex3.c", options->testFEjacobian, &options->testFEjacobian, NULL)); 1485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-test_fv_grad", "Test finite volume gradient reconstruction", "ex3.c", options->testFVgrad, &options->testFVgrad, NULL)); 1495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-test_injector","Test finite element injection", "ex3.c", options->testInjector, &options->testInjector,NULL)); 1505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsRealArray("-constants","Set the constant values", "ex3.c", options->constants, &n,NULL)); 151c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 152c4762a1bSJed Brown 153c4762a1bSJed Brown PetscFunctionReturn(0); 154c4762a1bSJed Brown } 155c4762a1bSJed Brown 156c4762a1bSJed Brown static PetscErrorCode TransformCoordinates(DM dm, AppCtx *user) 157c4762a1bSJed Brown { 158c4762a1bSJed Brown PetscSection coordSection; 159c4762a1bSJed Brown Vec coordinates; 160c4762a1bSJed Brown PetscScalar *coords; 161c4762a1bSJed Brown PetscInt vStart, vEnd, v; 162c4762a1bSJed Brown 163c4762a1bSJed Brown PetscFunctionBeginUser; 164c4762a1bSJed Brown if (user->nonaffineCoords) { 165c4762a1bSJed Brown /* x' = r^(1/p) (x/r), y' = r^(1/p) (y/r), z' = z */ 1665f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateSection(dm, &coordSection)); 1675f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocal(dm, &coordinates)); 1685f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 1695f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(coordinates, &coords)); 170c4762a1bSJed Brown for (v = vStart; v < vEnd; ++v) { 171c4762a1bSJed Brown PetscInt dof, off; 172c4762a1bSJed Brown PetscReal p = 4.0, r; 173c4762a1bSJed Brown 1745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetDof(coordSection, v, &dof)); 1755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetOffset(coordSection, v, &off)); 176c4762a1bSJed Brown switch (dof) { 177c4762a1bSJed Brown case 2: 178c4762a1bSJed Brown r = PetscSqr(PetscRealPart(coords[off+0])) + PetscSqr(PetscRealPart(coords[off+1])); 179c4762a1bSJed Brown coords[off+0] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p)/(2*p))*coords[off+0]; 180c4762a1bSJed Brown coords[off+1] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p)/(2*p))*coords[off+1]; 181c4762a1bSJed Brown break; 182c4762a1bSJed Brown case 3: 183c4762a1bSJed Brown r = PetscSqr(PetscRealPart(coords[off+0])) + PetscSqr(PetscRealPart(coords[off+1])); 184c4762a1bSJed Brown coords[off+0] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p)/(2*p))*coords[off+0]; 185c4762a1bSJed Brown coords[off+1] = r == 0.0 ? 0.0 : PetscPowReal(r, (1 - p)/(2*p))*coords[off+1]; 186c4762a1bSJed Brown coords[off+2] = coords[off+2]; 187c4762a1bSJed Brown break; 188c4762a1bSJed Brown } 189c4762a1bSJed Brown } 1905f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(coordinates, &coords)); 191c4762a1bSJed Brown } 192c4762a1bSJed Brown if (user->shearCoords) { 193c4762a1bSJed Brown /* x' = x + m y + m z, y' = y + m z, z' = z */ 1945f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateSection(dm, &coordSection)); 1955f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocal(dm, &coordinates)); 1965f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 1975f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(coordinates, &coords)); 198c4762a1bSJed Brown for (v = vStart; v < vEnd; ++v) { 199c4762a1bSJed Brown PetscInt dof, off; 200c4762a1bSJed Brown PetscReal m = 1.0; 201c4762a1bSJed Brown 2025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetDof(coordSection, v, &dof)); 2035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetOffset(coordSection, v, &off)); 204c4762a1bSJed Brown switch (dof) { 205c4762a1bSJed Brown case 2: 206c4762a1bSJed Brown coords[off+0] = coords[off+0] + m*coords[off+1]; 207c4762a1bSJed Brown coords[off+1] = coords[off+1]; 208c4762a1bSJed Brown break; 209c4762a1bSJed Brown case 3: 210c4762a1bSJed Brown coords[off+0] = coords[off+0] + m*coords[off+1] + m*coords[off+2]; 211c4762a1bSJed Brown coords[off+1] = coords[off+1] + m*coords[off+2]; 212c4762a1bSJed Brown coords[off+2] = coords[off+2]; 213c4762a1bSJed Brown break; 214c4762a1bSJed Brown } 215c4762a1bSJed Brown } 2165f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(coordinates, &coords)); 217c4762a1bSJed Brown } 218c4762a1bSJed Brown PetscFunctionReturn(0); 219c4762a1bSJed Brown } 220c4762a1bSJed Brown 221c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 222c4762a1bSJed Brown { 22330602db0SMatthew G. Knepley PetscInt dim = 2; 22430602db0SMatthew G. Knepley PetscBool simplex; 225c4762a1bSJed Brown 226c4762a1bSJed Brown PetscFunctionBeginUser; 22730602db0SMatthew G. Knepley if (user->useDA) { 22830602db0SMatthew G. Knepley switch (dim) { 229c4762a1bSJed Brown case 2: 2305f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 2, 2, PETSC_DETERMINE, PETSC_DETERMINE, 1, 1, NULL, NULL, dm)); 2315f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 2325f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(*dm)); 2335f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetVertexCoordinates(*dm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 234c4762a1bSJed Brown break; 235c4762a1bSJed Brown default: 23698921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot create structured mesh of dimension %d", dim); 237c4762a1bSJed Brown } 2385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) *dm, "Hexahedral Mesh")); 23930602db0SMatthew G. Knepley } else { 2405f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 2415f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 2425f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 2435f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 244c4762a1bSJed Brown 2455f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(*dm, &dim)); 2465f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexIsSimplex(*dm, &simplex)); 2475f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Bcast(&simplex, 1, MPIU_BOOL, 0, comm)); 248c4762a1bSJed Brown if (user->tree) { 24930602db0SMatthew G. Knepley DM refTree, ncdm = NULL; 250c4762a1bSJed Brown 2515f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateDefaultReferenceTree(comm,dim,simplex,&refTree)); 2525f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(refTree,NULL,"-reftree_dm_view")); 2535f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetReferenceTree(*dm,refTree)); 2545f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&refTree)); 2555f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexTreeRefineCell(*dm,user->treeCell,&ncdm)); 256c4762a1bSJed Brown if (ncdm) { 2575f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 258c4762a1bSJed Brown *dm = ncdm; 2595f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetRefinementUniform(*dm, PETSC_FALSE)); 260c4762a1bSJed Brown } 2615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) *dm, "tree_")); 2625f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 2635f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 2645f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm,NULL,"-dm_view")); 265c4762a1bSJed Brown } else { 2665f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetRefinementUniform(*dm, PETSC_TRUE)); 267c4762a1bSJed Brown } 2685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) *dm, "dist_")); 2695f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistributeSetDefault(*dm, PETSC_FALSE)); 2705f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 2715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) *dm, NULL)); 2725f80ce2aSJacob Faibussowitsch if (simplex) CHKERRQ(PetscObjectSetName((PetscObject) *dm, "Simplicial Mesh")); 2735f80ce2aSJacob Faibussowitsch else CHKERRQ(PetscObjectSetName((PetscObject) *dm, "Hexahedral Mesh")); 274c4762a1bSJed Brown } 2755f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 2765f80ce2aSJacob Faibussowitsch CHKERRQ(TransformCoordinates(*dm, user)); 2775f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm,NULL,"-dm_view")); 278c4762a1bSJed Brown PetscFunctionReturn(0); 279c4762a1bSJed Brown } 280c4762a1bSJed Brown 281c4762a1bSJed Brown static void simple_mass(PetscInt dim, PetscInt Nf, PetscInt NfAux, 282c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 283c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 284c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 285c4762a1bSJed Brown { 286c4762a1bSJed Brown PetscInt d, e; 287c4762a1bSJed Brown for (d = 0, e = 0; d < dim; d++, e+=dim+1) { 288c4762a1bSJed Brown g0[e] = 1.; 289c4762a1bSJed Brown } 290c4762a1bSJed Brown } 291c4762a1bSJed Brown 292c4762a1bSJed Brown /* < \nabla v, 1/2(\nabla u + {\nabla u}^T) > */ 293c4762a1bSJed Brown static void symmetric_gradient_inner_product(PetscInt dim, PetscInt Nf, PetscInt NfAux, 294c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 295c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 296c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar C[]) 297c4762a1bSJed Brown { 298c4762a1bSJed Brown PetscInt compI, compJ, d, e; 299c4762a1bSJed Brown 300c4762a1bSJed Brown for (compI = 0; compI < dim; ++compI) { 301c4762a1bSJed Brown for (compJ = 0; compJ < dim; ++compJ) { 302c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 303c4762a1bSJed Brown for (e = 0; e < dim; e++) { 304c4762a1bSJed Brown if (d == e && d == compI && d == compJ) { 305c4762a1bSJed Brown C[((compI*dim+compJ)*dim+d)*dim+e] = 1.0; 306c4762a1bSJed Brown } else if ((d == compJ && e == compI) || (d == e && compI == compJ)) { 307c4762a1bSJed Brown C[((compI*dim+compJ)*dim+d)*dim+e] = 0.5; 308c4762a1bSJed Brown } else { 309c4762a1bSJed Brown C[((compI*dim+compJ)*dim+d)*dim+e] = 0.0; 310c4762a1bSJed Brown } 311c4762a1bSJed Brown } 312c4762a1bSJed Brown } 313c4762a1bSJed Brown } 314c4762a1bSJed Brown } 315c4762a1bSJed Brown } 316c4762a1bSJed Brown 317c4762a1bSJed Brown static PetscErrorCode SetupSection(DM dm, AppCtx *user) 318c4762a1bSJed Brown { 319c4762a1bSJed Brown PetscFunctionBeginUser; 32030602db0SMatthew G. Knepley if (user->constraints) { 321c4762a1bSJed Brown /* test local constraints */ 322c4762a1bSJed Brown DM coordDM; 323c4762a1bSJed Brown PetscInt fStart, fEnd, f, vStart, vEnd, v; 324c4762a1bSJed Brown PetscInt edgesx = 2, vertsx; 325c4762a1bSJed Brown PetscInt edgesy = 2, vertsy; 326c4762a1bSJed Brown PetscMPIInt size; 327c4762a1bSJed Brown PetscInt numConst; 328c4762a1bSJed Brown PetscSection aSec; 329c4762a1bSJed Brown PetscInt *anchors; 330c4762a1bSJed Brown PetscInt offset; 331c4762a1bSJed Brown IS aIS; 332c4762a1bSJed Brown MPI_Comm comm = PetscObjectComm((PetscObject)dm); 333c4762a1bSJed Brown 3345f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(comm,&size)); 3352c71b3e2SJacob Faibussowitsch PetscCheckFalse(size > 1,comm,PETSC_ERR_SUP,"Local constraint test can only be performed in serial"); 336c4762a1bSJed Brown 337c4762a1bSJed Brown /* we are going to test constraints by using them to enforce periodicity 338c4762a1bSJed Brown * in one direction, and comparing to the existing method of enforcing 339c4762a1bSJed Brown * periodicity */ 340c4762a1bSJed Brown 341c4762a1bSJed Brown /* first create the coordinate section so that it does not clone the 342c4762a1bSJed Brown * constraints */ 3435f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(dm,&coordDM)); 344c4762a1bSJed Brown 345c4762a1bSJed Brown /* create the constrained-to-anchor section */ 3465f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetDepthStratum(dm,0,&vStart,&vEnd)); 3475f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetDepthStratum(dm,1,&fStart,&fEnd)); 3485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionCreate(PETSC_COMM_SELF,&aSec)); 3495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetChart(aSec,PetscMin(fStart,vStart),PetscMax(fEnd,vEnd))); 350c4762a1bSJed Brown 351c4762a1bSJed Brown /* define the constraints */ 3525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-da_grid_x",&edgesx,NULL)); 3535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-da_grid_y",&edgesy,NULL)); 354c4762a1bSJed Brown vertsx = edgesx + 1; 355c4762a1bSJed Brown vertsy = edgesy + 1; 356c4762a1bSJed Brown numConst = vertsy + edgesy; 3575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(numConst,&anchors)); 358c4762a1bSJed Brown offset = 0; 359c4762a1bSJed Brown for (v = vStart + edgesx; v < vEnd; v+= vertsx) { 3605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetDof(aSec,v,1)); 361c4762a1bSJed Brown anchors[offset++] = v - edgesx; 362c4762a1bSJed Brown } 363c4762a1bSJed Brown for (f = fStart + edgesx * vertsy + edgesx * edgesy; f < fEnd; f++) { 3645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetDof(aSec,f,1)); 365c4762a1bSJed Brown anchors[offset++] = f - edgesx * edgesy; 366c4762a1bSJed Brown } 3675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetUp(aSec)); 3685f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,numConst,anchors,PETSC_OWN_POINTER,&aIS)); 369c4762a1bSJed Brown 3705f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetAnchors(dm,aSec,aIS)); 3715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionDestroy(&aSec)); 3725f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&aIS)); 373c4762a1bSJed Brown } 3745f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetNumFields(dm, 1)); 3755f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) user->fe)); 3765f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(dm)); 37730602db0SMatthew G. Knepley if (user->constraints) { 378c4762a1bSJed Brown /* test getting local constraint matrix that matches section */ 379c4762a1bSJed Brown PetscSection aSec; 380c4762a1bSJed Brown IS aIS; 381c4762a1bSJed Brown 3825f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetAnchors(dm,&aSec,&aIS)); 383c4762a1bSJed Brown if (aSec) { 384c4762a1bSJed Brown PetscDS ds; 385c4762a1bSJed Brown PetscSection cSec, section; 386c4762a1bSJed Brown PetscInt cStart, cEnd, c, numComp; 387c4762a1bSJed Brown Mat cMat, mass; 388c4762a1bSJed Brown Vec local; 389c4762a1bSJed Brown const PetscInt *anchors; 390c4762a1bSJed Brown 3915f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalSection(dm,§ion)); 392c4762a1bSJed Brown /* this creates the matrix and preallocates the matrix structure: we 393c4762a1bSJed Brown * just have to fill in the values */ 3945f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDefaultConstraints(dm,&cSec,&cMat,NULL)); 3955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetChart(cSec,&cStart,&cEnd)); 3965f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(aIS,&anchors)); 3975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEGetNumComponents(user->fe, &numComp)); 398c4762a1bSJed Brown for (c = cStart; c < cEnd; c++) { 399c4762a1bSJed Brown PetscInt cDof; 400c4762a1bSJed Brown 401c4762a1bSJed Brown /* is this point constrained? (does it have an anchor?) */ 4025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetDof(aSec,c,&cDof)); 403c4762a1bSJed Brown if (cDof) { 404c4762a1bSJed Brown PetscInt cOff, a, aDof, aOff, j; 4052c71b3e2SJacob Faibussowitsch PetscCheckFalse(cDof != 1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Found %d anchor points: should be just one",cDof); 406c4762a1bSJed Brown 407c4762a1bSJed Brown /* find the anchor point */ 4085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetOffset(aSec,c,&cOff)); 409c4762a1bSJed Brown a = anchors[cOff]; 410c4762a1bSJed Brown 411c4762a1bSJed Brown /* find the constrained dofs (row in constraint matrix) */ 4125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetDof(cSec,c,&cDof)); 4135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetOffset(cSec,c,&cOff)); 414c4762a1bSJed Brown 415c4762a1bSJed Brown /* find the anchor dofs (column in constraint matrix) */ 4165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetDof(section,a,&aDof)); 4175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetOffset(section,a,&aOff)); 418c4762a1bSJed Brown 4192c71b3e2SJacob Faibussowitsch PetscCheckFalse(cDof != aDof,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Point and anchor have different number of dofs: %d, %d",cDof,aDof); 4202c71b3e2SJacob Faibussowitsch PetscCheckFalse(cDof % numComp,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Point dofs not divisible by field components: %d, %d",cDof,numComp); 421c4762a1bSJed Brown 422c4762a1bSJed Brown /* put in a simple equality constraint */ 423c4762a1bSJed Brown for (j = 0; j < cDof; j++) { 4245f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(cMat,cOff+j,aOff+j,1.,INSERT_VALUES)); 425c4762a1bSJed Brown } 426c4762a1bSJed Brown } 427c4762a1bSJed Brown } 4285f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(cMat,MAT_FINAL_ASSEMBLY)); 4295f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(cMat,MAT_FINAL_ASSEMBLY)); 4305f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(aIS,&anchors)); 431c4762a1bSJed Brown 432c4762a1bSJed Brown /* Now that we have constructed the constraint matrix, any FE matrix 433c4762a1bSJed Brown * that we construct will apply the constraints during construction */ 434c4762a1bSJed Brown 4355f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(dm,&mass)); 436c4762a1bSJed Brown /* get a dummy local variable to serve as the solution */ 4375f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dm,&local)); 4385f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm,&ds)); 439c4762a1bSJed Brown /* set the jacobian to be the mass matrix */ 4405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds, 0, 0, simple_mass, NULL, NULL, NULL)); 441c4762a1bSJed Brown /* build the mass matrix */ 4425f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSNESComputeJacobianFEM(dm,local,mass,mass,NULL)); 4435f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(mass,PETSC_VIEWER_STDOUT_WORLD)); 4445f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mass)); 4455f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dm,&local)); 446c4762a1bSJed Brown } 447c4762a1bSJed Brown } 448c4762a1bSJed Brown PetscFunctionReturn(0); 449c4762a1bSJed Brown } 450c4762a1bSJed Brown 451c4762a1bSJed Brown static PetscErrorCode TestFEJacobian(DM dm, AppCtx *user) 452c4762a1bSJed Brown { 453c4762a1bSJed Brown PetscFunctionBeginUser; 45430602db0SMatthew G. Knepley if (!user->useDA) { 455c4762a1bSJed Brown Vec local; 456c4762a1bSJed Brown const Vec *vecs; 457c4762a1bSJed Brown Mat E; 458c4762a1bSJed Brown MatNullSpace sp; 459c4762a1bSJed Brown PetscBool isNullSpace, hasConst; 46030602db0SMatthew G. Knepley PetscInt dim, n, i; 461c4762a1bSJed Brown Vec res = NULL, localX, localRes; 462c4762a1bSJed Brown PetscDS ds; 463c4762a1bSJed Brown 4645f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 4652c71b3e2SJacob Faibussowitsch PetscCheckFalse(user->numComponents != dim,PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "The number of components %d must be equal to the dimension %d for this test", user->numComponents, dim); 4665f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm,&ds)); 4675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetJacobian(ds,0,0,NULL,NULL,NULL,symmetric_gradient_inner_product)); 4685f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(dm,&E)); 4695f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dm,&local)); 4705f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSNESComputeJacobianFEM(dm,local,E,E,NULL)); 4715f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateRigidBody(dm,0,&sp)); 4725f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceGetVecs(sp,&hasConst,&n,&vecs)); 4735f80ce2aSJacob Faibussowitsch if (n) CHKERRQ(VecDuplicate(vecs[0],&res)); 4745f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLocalVector(dm,&localX)); 4755f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLocalVector(dm,&localRes)); 476c4762a1bSJed Brown for (i = 0; i < n; i++) { /* also test via matrix-free Jacobian application */ 477c4762a1bSJed Brown PetscReal resNorm; 478c4762a1bSJed Brown 4795f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(localRes,0.)); 4805f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(localX,0.)); 4815f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(local,0.)); 4825f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(res,0.)); 4835f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(dm,vecs[i],INSERT_VALUES,localX)); 4845f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(dm,vecs[i],INSERT_VALUES,localX)); 4855f80ce2aSJacob Faibussowitsch CHKERRQ(DMSNESComputeJacobianAction(dm,local,localX,localRes,NULL)); 4865f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalBegin(dm,localRes,ADD_VALUES,res)); 4875f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalEnd(dm,localRes,ADD_VALUES,res)); 4885f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(res,NORM_2,&resNorm)); 489c4762a1bSJed Brown if (resNorm > PETSC_SMALL) { 4905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)dm),"Symmetric gradient action null space vector %D residual: %E\n",i,resNorm)); 491c4762a1bSJed Brown } 492c4762a1bSJed Brown } 4935f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&localRes)); 4945f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&localX)); 4955f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&res)); 4965f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceTest(sp,E,&isNullSpace)); 497c4762a1bSJed Brown if (isNullSpace) { 4985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)dm),"Symmetric gradient null space: PASS\n")); 499c4762a1bSJed Brown } else { 5005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)dm),"Symmetric gradient null space: FAIL\n")); 501c4762a1bSJed Brown } 5025f80ce2aSJacob Faibussowitsch CHKERRQ(MatNullSpaceDestroy(&sp)); 5035f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&E)); 5045f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dm,&local)); 505c4762a1bSJed Brown } 506c4762a1bSJed Brown PetscFunctionReturn(0); 507c4762a1bSJed Brown } 508c4762a1bSJed Brown 509c4762a1bSJed Brown static PetscErrorCode TestInjector(DM dm, AppCtx *user) 510c4762a1bSJed Brown { 511c4762a1bSJed Brown DM refTree; 512c4762a1bSJed Brown PetscMPIInt rank; 513c4762a1bSJed Brown 514c4762a1bSJed Brown PetscFunctionBegin; 5155f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetReferenceTree(dm,&refTree)); 516c4762a1bSJed Brown if (refTree) { 517c4762a1bSJed Brown Mat inj; 518c4762a1bSJed Brown 5195f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexComputeInjectorReferenceTree(refTree,&inj)); 5205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)inj,"Reference Tree Injector")); 5215f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 522dd400576SPatrick Sanan if (rank == 0) { 5235f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(inj,PETSC_VIEWER_STDOUT_SELF)); 524c4762a1bSJed Brown } 5255f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&inj)); 526c4762a1bSJed Brown } 527c4762a1bSJed Brown PetscFunctionReturn(0); 528c4762a1bSJed Brown } 529c4762a1bSJed Brown 530c4762a1bSJed Brown static PetscErrorCode TestFVGrad(DM dm, AppCtx *user) 531c4762a1bSJed Brown { 532c4762a1bSJed Brown MPI_Comm comm; 533c4762a1bSJed Brown DM dmRedist, dmfv, dmgrad, dmCell, refTree; 534c4762a1bSJed Brown PetscFV fv; 53530602db0SMatthew G. Knepley PetscInt dim, nvecs, v, cStart, cEnd, cEndInterior; 536c4762a1bSJed Brown PetscMPIInt size; 537c4762a1bSJed Brown Vec cellgeom, grad, locGrad; 538c4762a1bSJed Brown const PetscScalar *cgeom; 539c4762a1bSJed Brown PetscReal allVecMaxDiff = 0., fvTol = 100. * PETSC_MACHINE_EPSILON; 540c4762a1bSJed Brown 541c4762a1bSJed Brown PetscFunctionBeginUser; 542c4762a1bSJed Brown comm = PetscObjectComm((PetscObject)dm); 5435f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 544c4762a1bSJed Brown /* duplicate DM, give dup. a FV discretization */ 5455f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetBasicAdjacency(dm,PETSC_TRUE,PETSC_FALSE)); 5465f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(comm,&size)); 547c4762a1bSJed Brown dmRedist = NULL; 548c4762a1bSJed Brown if (size > 1) { 5495f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistributeOverlap(dm,1,NULL,&dmRedist)); 550c4762a1bSJed Brown } 551c4762a1bSJed Brown if (!dmRedist) { 552c4762a1bSJed Brown dmRedist = dm; 5535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)dmRedist)); 554c4762a1bSJed Brown } 5555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFVCreate(comm,&fv)); 5565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFVSetType(fv,PETSCFVLEASTSQUARES)); 5575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFVSetNumComponents(fv,user->numComponents)); 5585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFVSetSpatialDimension(fv,dim)); 5595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFVSetFromOptions(fv)); 5605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFVSetUp(fv)); 5615f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexConstructGhostCells(dmRedist,NULL,NULL,&dmfv)); 5625f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dmRedist)); 5635f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetNumFields(dmfv,1)); 5645f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetField(dmfv, 0, NULL, (PetscObject) fv)); 5655f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(dmfv)); 5665f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetReferenceTree(dm,&refTree)); 5675f80ce2aSJacob Faibussowitsch if (refTree) CHKERRQ(DMCopyDisc(dmfv,refTree)); 5685f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetGradientDM(dmfv, fv, &dmgrad)); 5695f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dmfv,0,&cStart,&cEnd)); 57030602db0SMatthew G. Knepley nvecs = dim * (dim+1) / 2; 5715f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetGeometryFVM(dmfv,NULL,&cellgeom,NULL)); 5725f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetDM(cellgeom,&dmCell)); 5735f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(cellgeom,&cgeom)); 5745f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(dmgrad,&grad)); 5755f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dmgrad,&locGrad)); 5765f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetGhostCellStratum(dmgrad,&cEndInterior,NULL)); 577c4762a1bSJed Brown cEndInterior = (cEndInterior < 0) ? cEnd: cEndInterior; 578c4762a1bSJed Brown for (v = 0; v < nvecs; v++) { 579c4762a1bSJed Brown Vec locX; 580c4762a1bSJed Brown PetscInt c; 581c4762a1bSJed Brown PetscScalar trueGrad[3][3] = {{0.}}; 582c4762a1bSJed Brown const PetscScalar *gradArray; 583c4762a1bSJed Brown PetscReal maxDiff, maxDiffGlob; 584c4762a1bSJed Brown 5855f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dmfv,&locX)); 586c4762a1bSJed Brown /* get the local projection of the rigid body mode */ 587c4762a1bSJed Brown for (c = cStart; c < cEnd; c++) { 588c4762a1bSJed Brown PetscFVCellGeom *cg; 589c4762a1bSJed Brown PetscScalar cx[3] = {0.,0.,0.}; 590c4762a1bSJed Brown 5915f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexPointLocalRead(dmCell, c, cgeom, &cg)); 59230602db0SMatthew G. Knepley if (v < dim) { 593c4762a1bSJed Brown cx[v] = 1.; 594c4762a1bSJed Brown } else { 59530602db0SMatthew G. Knepley PetscInt w = v - dim; 596c4762a1bSJed Brown 59730602db0SMatthew G. Knepley cx[(w + 1) % dim] = cg->centroid[(w + 2) % dim]; 59830602db0SMatthew G. Knepley cx[(w + 2) % dim] = -cg->centroid[(w + 1) % dim]; 599c4762a1bSJed Brown } 6005f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecSetClosure(dmfv,NULL,locX,c,cx,INSERT_ALL_VALUES)); 601c4762a1bSJed Brown } 602c4762a1bSJed Brown /* TODO: this isn't in any header */ 6035f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexReconstructGradientsFVM(dmfv,locX,grad)); 6045f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(dmgrad,grad,INSERT_VALUES,locGrad)); 6055f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(dmgrad,grad,INSERT_VALUES,locGrad)); 6065f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(locGrad,&gradArray)); 607c4762a1bSJed Brown /* compare computed gradient to exact gradient */ 60830602db0SMatthew G. Knepley if (v >= dim) { 60930602db0SMatthew G. Knepley PetscInt w = v - dim; 610c4762a1bSJed Brown 61130602db0SMatthew G. Knepley trueGrad[(w + 1) % dim][(w + 2) % dim] = 1.; 61230602db0SMatthew G. Knepley trueGrad[(w + 2) % dim][(w + 1) % dim] = -1.; 613c4762a1bSJed Brown } 614c4762a1bSJed Brown maxDiff = 0.; 615c4762a1bSJed Brown for (c = cStart; c < cEndInterior; c++) { 616c4762a1bSJed Brown PetscScalar *compGrad; 617c4762a1bSJed Brown PetscInt i, j, k; 618c4762a1bSJed Brown PetscReal FrobDiff = 0.; 619c4762a1bSJed Brown 6205f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexPointLocalRead(dmgrad, c, gradArray, &compGrad)); 621c4762a1bSJed Brown 62230602db0SMatthew G. Knepley for (i = 0, k = 0; i < dim; i++) { 62330602db0SMatthew G. Knepley for (j = 0; j < dim; j++, k++) { 624c4762a1bSJed Brown PetscScalar diff = compGrad[k] - trueGrad[i][j]; 625c4762a1bSJed Brown FrobDiff += PetscRealPart(diff * PetscConj(diff)); 626c4762a1bSJed Brown } 627c4762a1bSJed Brown } 628c4762a1bSJed Brown FrobDiff = PetscSqrtReal(FrobDiff); 629c4762a1bSJed Brown maxDiff = PetscMax(maxDiff,FrobDiff); 630c4762a1bSJed Brown } 6315f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Allreduce(&maxDiff,&maxDiffGlob,1,MPIU_REAL,MPIU_MAX,comm)); 632c4762a1bSJed Brown allVecMaxDiff = PetscMax(allVecMaxDiff,maxDiffGlob); 6335f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(locGrad,&gradArray)); 6345f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dmfv,&locX)); 635c4762a1bSJed Brown } 636c4762a1bSJed Brown if (allVecMaxDiff < fvTol) { 6375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)dm),"Finite volume gradient reconstruction: PASS\n")); 638c4762a1bSJed Brown } else { 6395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject)dm),"Finite volume gradient reconstruction: FAIL at tolerance %g with max difference %g\n",fvTol,allVecMaxDiff)); 640c4762a1bSJed Brown } 6415f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dmgrad,&locGrad)); 6425f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(dmgrad,&grad)); 6435f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(cellgeom,&cgeom)); 6445f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dmfv)); 6455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFVDestroy(&fv)); 646c4762a1bSJed Brown PetscFunctionReturn(0); 647c4762a1bSJed Brown } 648c4762a1bSJed Brown 649c4762a1bSJed Brown static PetscErrorCode ComputeError(DM dm, PetscErrorCode (**exactFuncs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), 650c4762a1bSJed Brown PetscErrorCode (**exactFuncDers)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), 651c4762a1bSJed Brown void **exactCtxs, PetscReal *error, PetscReal *errorDer, AppCtx *user) 652c4762a1bSJed Brown { 653c4762a1bSJed Brown Vec u; 654c4762a1bSJed Brown PetscReal n[3] = {1.0, 1.0, 1.0}; 655c4762a1bSJed Brown 656c4762a1bSJed Brown PetscFunctionBeginUser; 6575f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(dm, &u)); 658c4762a1bSJed Brown /* Project function into FE function space */ 6595f80ce2aSJacob Faibussowitsch CHKERRQ(DMProjectFunction(dm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, u)); 6605f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(u, NULL, "-projection_view")); 661c4762a1bSJed Brown /* Compare approximation to exact in L_2 */ 6625f80ce2aSJacob Faibussowitsch CHKERRQ(DMComputeL2Diff(dm, 0.0, exactFuncs, exactCtxs, u, error)); 6635f80ce2aSJacob Faibussowitsch CHKERRQ(DMComputeL2GradientDiff(dm, 0.0, exactFuncDers, exactCtxs, u, n, errorDer)); 6645f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(dm, &u)); 665c4762a1bSJed Brown PetscFunctionReturn(0); 666c4762a1bSJed Brown } 667c4762a1bSJed Brown 668c4762a1bSJed Brown static PetscErrorCode CheckFunctions(DM dm, PetscInt order, AppCtx *user) 669c4762a1bSJed Brown { 670c4762a1bSJed Brown PetscErrorCode (*exactFuncs[1]) (PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 671c4762a1bSJed Brown PetscErrorCode (*exactFuncDers[1]) (PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx); 672c4762a1bSJed Brown void *exactCtxs[3]; 673c4762a1bSJed Brown MPI_Comm comm; 674c4762a1bSJed Brown PetscReal error, errorDer, tol = PETSC_SMALL; 675c4762a1bSJed Brown 676c4762a1bSJed Brown PetscFunctionBeginUser; 677c4762a1bSJed Brown exactCtxs[0] = user; 678c4762a1bSJed Brown exactCtxs[1] = user; 679c4762a1bSJed Brown exactCtxs[2] = user; 6805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)dm, &comm)); 681c4762a1bSJed Brown /* Setup functions to approximate */ 682c4762a1bSJed Brown switch (order) { 683c4762a1bSJed Brown case 0: 684c4762a1bSJed Brown exactFuncs[0] = constant; 685c4762a1bSJed Brown exactFuncDers[0] = constantDer; 686c4762a1bSJed Brown break; 687c4762a1bSJed Brown case 1: 688c4762a1bSJed Brown exactFuncs[0] = linear; 689c4762a1bSJed Brown exactFuncDers[0] = linearDer; 690c4762a1bSJed Brown break; 691c4762a1bSJed Brown case 2: 692c4762a1bSJed Brown exactFuncs[0] = quadratic; 693c4762a1bSJed Brown exactFuncDers[0] = quadraticDer; 694c4762a1bSJed Brown break; 695c4762a1bSJed Brown case 3: 696c4762a1bSJed Brown exactFuncs[0] = cubic; 697c4762a1bSJed Brown exactFuncDers[0] = cubicDer; 698c4762a1bSJed Brown break; 699c4762a1bSJed Brown default: 70098921bdaSJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Could not determine functions to test for order %d", order); 701c4762a1bSJed Brown } 7025f80ce2aSJacob Faibussowitsch CHKERRQ(ComputeError(dm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user)); 703c4762a1bSJed Brown /* Report result */ 7045f80ce2aSJacob Faibussowitsch if (error > tol) CHKERRQ(PetscPrintf(comm, "Function tests FAIL for order %D at tolerance %g error %g\n", order, (double)tol,(double) error)); 7055f80ce2aSJacob Faibussowitsch else CHKERRQ(PetscPrintf(comm, "Function tests pass for order %D at tolerance %g\n", order, (double)tol)); 7065f80ce2aSJacob Faibussowitsch if (errorDer > tol) CHKERRQ(PetscPrintf(comm, "Function tests FAIL for order %D derivatives at tolerance %g error %g\n", order, (double)tol, (double)errorDer)); 7075f80ce2aSJacob Faibussowitsch else CHKERRQ(PetscPrintf(comm, "Function tests pass for order %D derivatives at tolerance %g\n", order, (double)tol)); 708c4762a1bSJed Brown PetscFunctionReturn(0); 709c4762a1bSJed Brown } 710c4762a1bSJed Brown 711c4762a1bSJed Brown static PetscErrorCode CheckInterpolation(DM dm, PetscBool checkRestrict, PetscInt order, AppCtx *user) 712c4762a1bSJed Brown { 713c4762a1bSJed Brown PetscErrorCode (*exactFuncs[1]) (PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx); 714c4762a1bSJed Brown PetscErrorCode (*exactFuncDers[1]) (PetscInt, PetscReal, const PetscReal x[], const PetscReal n[], PetscInt, PetscScalar *u, void *ctx); 715c4762a1bSJed Brown PetscReal n[3] = {1.0, 1.0, 1.0}; 716c4762a1bSJed Brown void *exactCtxs[3]; 717c4762a1bSJed Brown DM rdm, idm, fdm; 718c4762a1bSJed Brown Mat Interp; 719c4762a1bSJed Brown Vec iu, fu, scaling; 720c4762a1bSJed Brown MPI_Comm comm; 72130602db0SMatthew G. Knepley PetscInt dim; 722c4762a1bSJed Brown PetscReal error, errorDer, tol = PETSC_SMALL; 723c4762a1bSJed Brown 724c4762a1bSJed Brown PetscFunctionBeginUser; 725c4762a1bSJed Brown exactCtxs[0] = user; 726c4762a1bSJed Brown exactCtxs[1] = user; 727c4762a1bSJed Brown exactCtxs[2] = user; 7285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)dm,&comm)); 7295f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 7305f80ce2aSJacob Faibussowitsch CHKERRQ(DMRefine(dm, comm, &rdm)); 7315f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetCoarseDM(rdm, dm)); 7325f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetRegularRefinement(rdm, user->convRefine)); 73330602db0SMatthew G. Knepley if (user->tree) { 734c4762a1bSJed Brown DM refTree; 7355f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetReferenceTree(dm,&refTree)); 7365f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetReferenceTree(rdm,refTree)); 737c4762a1bSJed Brown } 7385f80ce2aSJacob Faibussowitsch if (user->useDA) CHKERRQ(DMDASetVertexCoordinates(rdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 7395f80ce2aSJacob Faibussowitsch CHKERRQ(SetupSection(rdm, user)); 740c4762a1bSJed Brown /* Setup functions to approximate */ 741c4762a1bSJed Brown switch (order) { 742c4762a1bSJed Brown case 0: 743c4762a1bSJed Brown exactFuncs[0] = constant; 744c4762a1bSJed Brown exactFuncDers[0] = constantDer; 745c4762a1bSJed Brown break; 746c4762a1bSJed Brown case 1: 747c4762a1bSJed Brown exactFuncs[0] = linear; 748c4762a1bSJed Brown exactFuncDers[0] = linearDer; 749c4762a1bSJed Brown break; 750c4762a1bSJed Brown case 2: 751c4762a1bSJed Brown exactFuncs[0] = quadratic; 752c4762a1bSJed Brown exactFuncDers[0] = quadraticDer; 753c4762a1bSJed Brown break; 754c4762a1bSJed Brown case 3: 755c4762a1bSJed Brown exactFuncs[0] = cubic; 756c4762a1bSJed Brown exactFuncDers[0] = cubicDer; 757c4762a1bSJed Brown break; 758c4762a1bSJed Brown default: 75998921bdaSJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Could not determine functions to test for dimension %D order %D", dim, order); 760c4762a1bSJed Brown } 761c4762a1bSJed Brown idm = checkRestrict ? rdm : dm; 762c4762a1bSJed Brown fdm = checkRestrict ? dm : rdm; 7635f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(idm, &iu)); 7645f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(fdm, &fu)); 7655f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetApplicationContext(dm, user)); 7665f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetApplicationContext(rdm, user)); 7675f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateInterpolation(dm, rdm, &Interp, &scaling)); 768c4762a1bSJed Brown /* Project function into initial FE function space */ 7695f80ce2aSJacob Faibussowitsch CHKERRQ(DMProjectFunction(idm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, iu)); 770c4762a1bSJed Brown /* Interpolate function into final FE function space */ 7715f80ce2aSJacob Faibussowitsch if (checkRestrict) {CHKERRQ(MatRestrict(Interp, iu, fu));CHKERRQ(VecPointwiseMult(fu, scaling, fu));} 7725f80ce2aSJacob Faibussowitsch else CHKERRQ(MatInterpolate(Interp, iu, fu)); 773c4762a1bSJed Brown /* Compare approximation to exact in L_2 */ 7745f80ce2aSJacob Faibussowitsch CHKERRQ(DMComputeL2Diff(fdm, 0.0, exactFuncs, exactCtxs, fu, &error)); 7755f80ce2aSJacob Faibussowitsch CHKERRQ(DMComputeL2GradientDiff(fdm, 0.0, exactFuncDers, exactCtxs, fu, n, &errorDer)); 776c4762a1bSJed Brown /* Report result */ 7775f80ce2aSJacob Faibussowitsch if (error > tol) CHKERRQ(PetscPrintf(comm, "Interpolation tests FAIL for order %D at tolerance %g error %g\n", order, (double)tol, (double)error)); 7785f80ce2aSJacob Faibussowitsch else CHKERRQ(PetscPrintf(comm, "Interpolation tests pass for order %D at tolerance %g\n", order, (double)tol)); 7795f80ce2aSJacob Faibussowitsch if (errorDer > tol) CHKERRQ(PetscPrintf(comm, "Interpolation tests FAIL for order %D derivatives at tolerance %g error %g\n", order, (double)tol, (double)errorDer)); 7805f80ce2aSJacob Faibussowitsch else CHKERRQ(PetscPrintf(comm, "Interpolation tests pass for order %D derivatives at tolerance %g\n", order, (double)tol)); 7815f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(idm, &iu)); 7825f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(fdm, &fu)); 7835f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Interp)); 7845f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&scaling)); 7855f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&rdm)); 786c4762a1bSJed Brown PetscFunctionReturn(0); 787c4762a1bSJed Brown } 788c4762a1bSJed Brown 789c4762a1bSJed Brown static PetscErrorCode CheckConvergence(DM dm, PetscInt Nr, AppCtx *user) 790c4762a1bSJed Brown { 791c4762a1bSJed Brown DM odm = dm, rdm = NULL, cdm = NULL; 792c4762a1bSJed Brown PetscErrorCode (*exactFuncs[1]) (PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) = {trig}; 793c4762a1bSJed Brown PetscErrorCode (*exactFuncDers[1]) (PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx) = {trigDer}; 794c4762a1bSJed Brown void *exactCtxs[3]; 795c4762a1bSJed Brown PetscInt r, c, cStart, cEnd; 796c4762a1bSJed Brown PetscReal errorOld, errorDerOld, error, errorDer, rel, len, lenOld; 797c4762a1bSJed Brown double p; 798c4762a1bSJed Brown 799c4762a1bSJed Brown PetscFunctionBeginUser; 800c4762a1bSJed Brown if (!user->convergence) PetscFunctionReturn(0); 801c4762a1bSJed Brown exactCtxs[0] = user; 802c4762a1bSJed Brown exactCtxs[1] = user; 803c4762a1bSJed Brown exactCtxs[2] = user; 8045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject) odm)); 805c4762a1bSJed Brown if (!user->convRefine) { 806c4762a1bSJed Brown for (r = 0; r < Nr; ++r) { 8075f80ce2aSJacob Faibussowitsch CHKERRQ(DMRefine(odm, PetscObjectComm((PetscObject) dm), &rdm)); 8085f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&odm)); 809c4762a1bSJed Brown odm = rdm; 810c4762a1bSJed Brown } 8115f80ce2aSJacob Faibussowitsch CHKERRQ(SetupSection(odm, user)); 812c4762a1bSJed Brown } 8135f80ce2aSJacob Faibussowitsch CHKERRQ(ComputeError(odm, exactFuncs, exactFuncDers, exactCtxs, &errorOld, &errorDerOld, user)); 814c4762a1bSJed Brown if (user->convRefine) { 815c4762a1bSJed Brown for (r = 0; r < Nr; ++r) { 8165f80ce2aSJacob Faibussowitsch CHKERRQ(DMRefine(odm, PetscObjectComm((PetscObject) dm), &rdm)); 8175f80ce2aSJacob Faibussowitsch if (user->useDA) CHKERRQ(DMDASetVertexCoordinates(rdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 8185f80ce2aSJacob Faibussowitsch CHKERRQ(SetupSection(rdm, user)); 8195f80ce2aSJacob Faibussowitsch CHKERRQ(ComputeError(rdm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user)); 820c4762a1bSJed Brown p = PetscLog2Real(errorOld/error); 8215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) dm), "Function convergence rate at refinement %D: %.2f\n", r, (double)p)); 822c4762a1bSJed Brown p = PetscLog2Real(errorDerOld/errorDer); 8235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) dm), "Derivative convergence rate at refinement %D: %.2f\n", r, (double)p)); 8245f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&odm)); 825c4762a1bSJed Brown odm = rdm; 826c4762a1bSJed Brown errorOld = error; 827c4762a1bSJed Brown errorDerOld = errorDer; 828c4762a1bSJed Brown } 829c4762a1bSJed Brown } else { 8305f80ce2aSJacob Faibussowitsch /* CHKERRQ(ComputeLongestEdge(dm, &lenOld)); */ 8315f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(odm, 0, &cStart, &cEnd)); 832c4762a1bSJed Brown lenOld = cEnd - cStart; 833c4762a1bSJed Brown for (c = 0; c < Nr; ++c) { 8345f80ce2aSJacob Faibussowitsch CHKERRQ(DMCoarsen(odm, PetscObjectComm((PetscObject) dm), &cdm)); 8355f80ce2aSJacob Faibussowitsch if (user->useDA) CHKERRQ(DMDASetVertexCoordinates(cdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 8365f80ce2aSJacob Faibussowitsch CHKERRQ(SetupSection(cdm, user)); 8375f80ce2aSJacob Faibussowitsch CHKERRQ(ComputeError(cdm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user)); 8385f80ce2aSJacob Faibussowitsch /* CHKERRQ(ComputeLongestEdge(cdm, &len)); */ 8395f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(cdm, 0, &cStart, &cEnd)); 840c4762a1bSJed Brown len = cEnd - cStart; 841c4762a1bSJed Brown rel = error/errorOld; 842c4762a1bSJed Brown p = PetscLogReal(rel) / PetscLogReal(lenOld / len); 8435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) dm), "Function convergence rate at coarsening %D: %.2f\n", c, (double)p)); 844c4762a1bSJed Brown rel = errorDer/errorDerOld; 845c4762a1bSJed Brown p = PetscLogReal(rel) / PetscLogReal(lenOld / len); 8465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) dm), "Derivative convergence rate at coarsening %D: %.2f\n", c, (double)p)); 8475f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&odm)); 848c4762a1bSJed Brown odm = cdm; 849c4762a1bSJed Brown errorOld = error; 850c4762a1bSJed Brown errorDerOld = errorDer; 851c4762a1bSJed Brown lenOld = len; 852c4762a1bSJed Brown } 853c4762a1bSJed Brown } 8545f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&odm)); 855c4762a1bSJed Brown PetscFunctionReturn(0); 856c4762a1bSJed Brown } 857c4762a1bSJed Brown 858c4762a1bSJed Brown int main(int argc, char **argv) 859c4762a1bSJed Brown { 860c4762a1bSJed Brown DM dm; 861c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 86230602db0SMatthew G. Knepley PetscInt dim = 2; 86330602db0SMatthew G. Knepley PetscBool simplex = PETSC_FALSE; 864c4762a1bSJed Brown 865*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, NULL, help)); 8665f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user)); 8675f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 86830602db0SMatthew G. Knepley if (!user.useDA) { 8695f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &dim)); 8705f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexIsSimplex(dm, &simplex)); 87130602db0SMatthew G. Knepley } 8725f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexMetricSetFromOptions(dm)); 87330602db0SMatthew G. Knepley user.numComponents = user.numComponents < 0 ? dim : user.numComponents; 8745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateDefault(PETSC_COMM_WORLD, dim, user.numComponents, simplex, NULL, user.qorder, &user.fe)); 8755f80ce2aSJacob Faibussowitsch CHKERRQ(SetupSection(dm, &user)); 8765f80ce2aSJacob Faibussowitsch if (user.testFEjacobian) CHKERRQ(TestFEJacobian(dm, &user)); 8775f80ce2aSJacob Faibussowitsch if (user.testFVgrad) CHKERRQ(TestFVGrad(dm, &user)); 8785f80ce2aSJacob Faibussowitsch if (user.testInjector) CHKERRQ(TestInjector(dm, &user)); 8795f80ce2aSJacob Faibussowitsch CHKERRQ(CheckFunctions(dm, user.porder, &user)); 880c4762a1bSJed Brown { 881c4762a1bSJed Brown PetscDualSpace dsp; 882c4762a1bSJed Brown PetscInt k; 883c4762a1bSJed Brown 8845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEGetDualSpace(user.fe, &dsp)); 8855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDualSpaceGetDeRahm(dsp, &k)); 88630602db0SMatthew G. Knepley if (dim == 2 && simplex == PETSC_TRUE && user.tree == PETSC_FALSE && k == 0) { 8875f80ce2aSJacob Faibussowitsch CHKERRQ(CheckInterpolation(dm, PETSC_FALSE, user.porder, &user)); 8885f80ce2aSJacob Faibussowitsch CHKERRQ(CheckInterpolation(dm, PETSC_TRUE, user.porder, &user)); 889c4762a1bSJed Brown } 890c4762a1bSJed Brown } 8915f80ce2aSJacob Faibussowitsch CHKERRQ(CheckConvergence(dm, 3, &user)); 8925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&user.fe)); 8935f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 894*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 895*b122ec5aSJacob Faibussowitsch return 0; 896c4762a1bSJed Brown } 897c4762a1bSJed Brown 898c4762a1bSJed Brown /*TEST 899c4762a1bSJed Brown 900c4762a1bSJed Brown test: 901c4762a1bSJed Brown suffix: 1 902c4762a1bSJed Brown requires: triangle 903c4762a1bSJed Brown 904c4762a1bSJed Brown # 2D P_1 on a triangle 905c4762a1bSJed Brown test: 906c4762a1bSJed Brown suffix: p1_2d_0 907c4762a1bSJed Brown requires: triangle 908c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -convergence 909c4762a1bSJed Brown test: 910c4762a1bSJed Brown suffix: p1_2d_1 911c4762a1bSJed Brown requires: triangle 912c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -porder 1 913c4762a1bSJed Brown test: 914c4762a1bSJed Brown suffix: p1_2d_2 915c4762a1bSJed Brown requires: triangle 916c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -porder 2 917c4762a1bSJed Brown test: 918c4762a1bSJed Brown suffix: p1_2d_3 919e788613bSJoe Wallwork requires: triangle mmg 920c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0 921c4762a1bSJed Brown test: 922c4762a1bSJed Brown suffix: p1_2d_4 923e788613bSJoe Wallwork requires: triangle mmg 924c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0 925c4762a1bSJed Brown test: 926c4762a1bSJed Brown suffix: p1_2d_5 927e788613bSJoe Wallwork requires: triangle mmg 928c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0 929c4762a1bSJed Brown 930c4762a1bSJed Brown # 3D P_1 on a tetrahedron 931c4762a1bSJed Brown test: 932c4762a1bSJed Brown suffix: p1_3d_0 933c4762a1bSJed Brown requires: ctetgen 93430602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -convergence 935c4762a1bSJed Brown test: 936c4762a1bSJed Brown suffix: p1_3d_1 937c4762a1bSJed Brown requires: ctetgen 93830602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -porder 1 939c4762a1bSJed Brown test: 940c4762a1bSJed Brown suffix: p1_3d_2 941c4762a1bSJed Brown requires: ctetgen 94230602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -porder 2 943c4762a1bSJed Brown test: 944c4762a1bSJed Brown suffix: p1_3d_3 945e788613bSJoe Wallwork requires: ctetgen mmg 94630602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0 947c4762a1bSJed Brown test: 948c4762a1bSJed Brown suffix: p1_3d_4 949e788613bSJoe Wallwork requires: ctetgen mmg 95030602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0 951c4762a1bSJed Brown test: 952c4762a1bSJed Brown suffix: p1_3d_5 953e788613bSJoe Wallwork requires: ctetgen mmg 95430602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0 955c4762a1bSJed Brown 956c4762a1bSJed Brown # 2D P_2 on a triangle 957c4762a1bSJed Brown test: 958c4762a1bSJed Brown suffix: p2_2d_0 959c4762a1bSJed Brown requires: triangle 960c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 2 -convergence 961c4762a1bSJed Brown test: 962c4762a1bSJed Brown suffix: p2_2d_1 963c4762a1bSJed Brown requires: triangle 964c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 2 -porder 1 965c4762a1bSJed Brown test: 966c4762a1bSJed Brown suffix: p2_2d_2 967c4762a1bSJed Brown requires: triangle 968c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 2 -porder 2 969c4762a1bSJed Brown test: 970c4762a1bSJed Brown suffix: p2_2d_3 971e788613bSJoe Wallwork requires: triangle mmg 972c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -convergence -conv_refine 0 973c4762a1bSJed Brown test: 974c4762a1bSJed Brown suffix: p2_2d_4 975e788613bSJoe Wallwork requires: triangle mmg 976c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 1 -conv_refine 0 977c4762a1bSJed Brown test: 978c4762a1bSJed Brown suffix: p2_2d_5 979e788613bSJoe Wallwork requires: triangle mmg 980c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 2 -conv_refine 0 981c4762a1bSJed Brown 982c4762a1bSJed Brown # 3D P_2 on a tetrahedron 983c4762a1bSJed Brown test: 984c4762a1bSJed Brown suffix: p2_3d_0 985c4762a1bSJed Brown requires: ctetgen 98630602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -convergence 987c4762a1bSJed Brown test: 988c4762a1bSJed Brown suffix: p2_3d_1 989c4762a1bSJed Brown requires: ctetgen 99030602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -porder 1 991c4762a1bSJed Brown test: 992c4762a1bSJed Brown suffix: p2_3d_2 993c4762a1bSJed Brown requires: ctetgen 99430602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -porder 2 995c4762a1bSJed Brown test: 996c4762a1bSJed Brown suffix: p2_3d_3 997e788613bSJoe Wallwork requires: ctetgen mmg 99830602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -convergence -conv_refine 0 999c4762a1bSJed Brown test: 1000c4762a1bSJed Brown suffix: p2_3d_4 1001e788613bSJoe Wallwork requires: ctetgen mmg 100230602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 1 -conv_refine 0 1003c4762a1bSJed Brown test: 1004c4762a1bSJed Brown suffix: p2_3d_5 1005e788613bSJoe Wallwork requires: ctetgen mmg 100630602db0SMatthew G. Knepley args: -dm_plex_dim 3 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder 2 -conv_refine 0 1007c4762a1bSJed Brown 1008c4762a1bSJed Brown # 2D Q_1 on a quadrilaterial DA 1009c4762a1bSJed Brown test: 1010c4762a1bSJed Brown suffix: q1_2d_da_0 101199a192c5SJunchao Zhang requires: broken 101230602db0SMatthew G. Knepley args: -use_da 1 -petscspace_degree 1 -qorder 1 -convergence 1013c4762a1bSJed Brown test: 1014c4762a1bSJed Brown suffix: q1_2d_da_1 101599a192c5SJunchao Zhang requires: broken 101630602db0SMatthew G. Knepley args: -use_da 1 -petscspace_degree 1 -qorder 1 -porder 1 1017c4762a1bSJed Brown test: 1018c4762a1bSJed Brown suffix: q1_2d_da_2 101999a192c5SJunchao Zhang requires: broken 102030602db0SMatthew G. Knepley args: -use_da 1 -petscspace_degree 1 -qorder 1 -porder 2 1021c4762a1bSJed Brown 1022c4762a1bSJed Brown # 2D Q_1 on a quadrilaterial Plex 1023c4762a1bSJed Brown test: 1024c4762a1bSJed Brown suffix: q1_2d_plex_0 102530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -convergence 1026c4762a1bSJed Brown test: 1027c4762a1bSJed Brown suffix: q1_2d_plex_1 102830602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 1 1029c4762a1bSJed Brown test: 1030c4762a1bSJed Brown suffix: q1_2d_plex_2 103130602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 2 1032c4762a1bSJed Brown test: 1033c4762a1bSJed Brown suffix: q1_2d_plex_3 103430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 1 -shear_coords 1035c4762a1bSJed Brown test: 1036c4762a1bSJed Brown suffix: q1_2d_plex_4 103730602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder 2 -shear_coords 1038c4762a1bSJed Brown test: 1039c4762a1bSJed Brown suffix: q1_2d_plex_5 104030602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 0 -non_affine_coords -convergence 1041c4762a1bSJed Brown test: 1042c4762a1bSJed Brown suffix: q1_2d_plex_6 104330602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 1 -non_affine_coords -convergence 1044c4762a1bSJed Brown test: 1045c4762a1bSJed Brown suffix: q1_2d_plex_7 104630602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -petscspace_type tensor -qorder 1 -porder 2 -non_affine_coords -convergence 1047c4762a1bSJed Brown 1048c4762a1bSJed Brown # 2D Q_2 on a quadrilaterial 1049c4762a1bSJed Brown test: 1050c4762a1bSJed Brown suffix: q2_2d_plex_0 105130602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -convergence 1052c4762a1bSJed Brown test: 1053c4762a1bSJed Brown suffix: q2_2d_plex_1 105430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 1055c4762a1bSJed Brown test: 1056c4762a1bSJed Brown suffix: q2_2d_plex_2 105730602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 1058c4762a1bSJed Brown test: 1059c4762a1bSJed Brown suffix: q2_2d_plex_3 106030602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 -shear_coords 1061c4762a1bSJed Brown test: 1062c4762a1bSJed Brown suffix: q2_2d_plex_4 106330602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 -shear_coords 1064c4762a1bSJed Brown test: 1065c4762a1bSJed Brown suffix: q2_2d_plex_5 106630602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 0 -non_affine_coords -convergence 1067c4762a1bSJed Brown test: 1068c4762a1bSJed Brown suffix: q2_2d_plex_6 106930602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 1 -non_affine_coords -convergence 1070c4762a1bSJed Brown test: 1071c4762a1bSJed Brown suffix: q2_2d_plex_7 107230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_type tensor -qorder 2 -porder 2 -non_affine_coords -convergence 1073c4762a1bSJed Brown 1074c4762a1bSJed Brown # 2D P_3 on a triangle 1075c4762a1bSJed Brown test: 1076c4762a1bSJed Brown suffix: p3_2d_0 1077c4762a1bSJed Brown requires: triangle !single 1078c4762a1bSJed Brown args: -petscspace_degree 3 -qorder 3 -convergence 1079c4762a1bSJed Brown test: 1080c4762a1bSJed Brown suffix: p3_2d_1 1081c4762a1bSJed Brown requires: triangle !single 1082c4762a1bSJed Brown args: -petscspace_degree 3 -qorder 3 -porder 1 1083c4762a1bSJed Brown test: 1084c4762a1bSJed Brown suffix: p3_2d_2 1085c4762a1bSJed Brown requires: triangle !single 1086c4762a1bSJed Brown args: -petscspace_degree 3 -qorder 3 -porder 2 1087c4762a1bSJed Brown test: 1088c4762a1bSJed Brown suffix: p3_2d_3 1089c4762a1bSJed Brown requires: triangle !single 1090c4762a1bSJed Brown args: -petscspace_degree 3 -qorder 3 -porder 3 1091c4762a1bSJed Brown test: 1092c4762a1bSJed Brown suffix: p3_2d_4 1093e788613bSJoe Wallwork requires: triangle mmg 1094c4762a1bSJed Brown args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -convergence -conv_refine 0 1095c4762a1bSJed Brown test: 1096c4762a1bSJed Brown suffix: p3_2d_5 1097e788613bSJoe Wallwork requires: triangle mmg 1098c4762a1bSJed Brown args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -porder 1 -conv_refine 0 1099c4762a1bSJed Brown test: 1100c4762a1bSJed Brown suffix: p3_2d_6 1101e788613bSJoe Wallwork requires: triangle mmg 1102c4762a1bSJed Brown args: -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -porder 3 -conv_refine 0 1103c4762a1bSJed Brown 1104c4762a1bSJed Brown # 2D Q_3 on a quadrilaterial 1105c4762a1bSJed Brown test: 1106c4762a1bSJed Brown suffix: q3_2d_0 110799a192c5SJunchao Zhang requires: !single 110830602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -convergence 1109c4762a1bSJed Brown test: 1110c4762a1bSJed Brown suffix: q3_2d_1 111199a192c5SJunchao Zhang requires: !single 111230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 1 1113c4762a1bSJed Brown test: 1114c4762a1bSJed Brown suffix: q3_2d_2 111599a192c5SJunchao Zhang requires: !single 111630602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 2 1117c4762a1bSJed Brown test: 1118c4762a1bSJed Brown suffix: q3_2d_3 111999a192c5SJunchao Zhang requires: !single 112030602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder 3 1121c4762a1bSJed Brown 1122c4762a1bSJed Brown # 2D P_1disc on a triangle/quadrilateral 1123c4762a1bSJed Brown test: 1124c4762a1bSJed Brown suffix: p1d_2d_0 1125c4762a1bSJed Brown requires: triangle 1126c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -convergence 1127c4762a1bSJed Brown test: 1128c4762a1bSJed Brown suffix: p1d_2d_1 1129c4762a1bSJed Brown requires: triangle 1130c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 1 1131c4762a1bSJed Brown test: 1132c4762a1bSJed Brown suffix: p1d_2d_2 1133c4762a1bSJed Brown requires: triangle 1134c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 2 1135c4762a1bSJed Brown test: 1136c4762a1bSJed Brown suffix: p1d_2d_3 1137c4762a1bSJed Brown requires: triangle 113830602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -convergence 1139c4762a1bSJed Brown filter: sed -e "s/convergence rate at refinement 0: 2/convergence rate at refinement 0: 1.9/g" 1140c4762a1bSJed Brown test: 1141c4762a1bSJed Brown suffix: p1d_2d_4 1142c4762a1bSJed Brown requires: triangle 114330602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 1 1144c4762a1bSJed Brown test: 1145c4762a1bSJed Brown suffix: p1d_2d_5 1146c4762a1bSJed Brown requires: triangle 114730602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder 2 1148c4762a1bSJed Brown 1149c4762a1bSJed Brown # 2D BDM_1 on a triangle 1150c4762a1bSJed Brown test: 1151c4762a1bSJed Brown suffix: bdm1_2d_0 1152c4762a1bSJed Brown requires: triangle 1153c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_type bdm \ 1154c4762a1bSJed Brown -num_comp 2 -qorder 1 -convergence 1155c4762a1bSJed Brown test: 1156c4762a1bSJed Brown suffix: bdm1_2d_1 1157c4762a1bSJed Brown requires: triangle 1158c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_type bdm \ 1159c4762a1bSJed Brown -num_comp 2 -qorder 1 -porder 1 1160c4762a1bSJed Brown test: 1161c4762a1bSJed Brown suffix: bdm1_2d_2 1162c4762a1bSJed Brown requires: triangle 1163c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_type bdm \ 1164c4762a1bSJed Brown -num_comp 2 -qorder 1 -porder 2 1165c4762a1bSJed Brown 1166c4762a1bSJed Brown # 2D BDM_1 on a quadrilateral 1167c4762a1bSJed Brown test: 1168c4762a1bSJed Brown suffix: bdm1q_2d_0 1169c4762a1bSJed Brown requires: triangle 1170c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_type bdm \ 11713f27d899SToby Isaac -petscdualspace_lagrange_tensor 1 \ 117230602db0SMatthew G. Knepley -dm_plex_simplex 0 -num_comp 2 -qorder 1 -convergence 1173c4762a1bSJed Brown test: 1174c4762a1bSJed Brown suffix: bdm1q_2d_1 1175c4762a1bSJed Brown requires: triangle 1176c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_type bdm \ 11773f27d899SToby Isaac -petscdualspace_lagrange_tensor 1 \ 117830602db0SMatthew G. Knepley -dm_plex_simplex 0 -num_comp 2 -qorder 1 -porder 1 1179c4762a1bSJed Brown test: 1180c4762a1bSJed Brown suffix: bdm1q_2d_2 1181c4762a1bSJed Brown requires: triangle 1182c4762a1bSJed Brown args: -petscspace_degree 1 -petscdualspace_type bdm \ 11833f27d899SToby Isaac -petscdualspace_lagrange_tensor 1 \ 118430602db0SMatthew G. Knepley -dm_plex_simplex 0 -num_comp 2 -qorder 1 -porder 2 1185c4762a1bSJed Brown 1186c4762a1bSJed Brown # Test high order quadrature 1187c4762a1bSJed Brown test: 1188c4762a1bSJed Brown suffix: p1_quad_2 1189c4762a1bSJed Brown requires: triangle 1190c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 2 -porder 1 1191c4762a1bSJed Brown test: 1192c4762a1bSJed Brown suffix: p1_quad_5 1193c4762a1bSJed Brown requires: triangle 1194c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 5 -porder 1 1195c4762a1bSJed Brown test: 1196c4762a1bSJed Brown suffix: p2_quad_3 1197c4762a1bSJed Brown requires: triangle 1198c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 3 -porder 2 1199c4762a1bSJed Brown test: 1200c4762a1bSJed Brown suffix: p2_quad_5 1201c4762a1bSJed Brown requires: triangle 1202c4762a1bSJed Brown args: -petscspace_degree 2 -qorder 5 -porder 2 1203c4762a1bSJed Brown test: 1204c4762a1bSJed Brown suffix: q1_quad_2 120530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 2 -porder 1 1206c4762a1bSJed Brown test: 1207c4762a1bSJed Brown suffix: q1_quad_5 120830602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 1 -qorder 5 -porder 1 1209c4762a1bSJed Brown test: 1210c4762a1bSJed Brown suffix: q2_quad_3 121130602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 3 -porder 1 1212c4762a1bSJed Brown test: 1213c4762a1bSJed Brown suffix: q2_quad_5 121430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 5 -porder 1 1215c4762a1bSJed Brown 1216c4762a1bSJed Brown # Nonconforming tests 1217c4762a1bSJed Brown test: 1218c4762a1bSJed Brown suffix: constraints 121930602db0SMatthew G. Knepley args: -dm_coord_space 0 -dm_plex_simplex 0 -petscspace_type tensor -petscspace_degree 1 -qorder 0 -constraints 1220c4762a1bSJed Brown test: 1221c4762a1bSJed Brown suffix: nonconforming_tensor_2 1222c4762a1bSJed Brown nsize: 4 122330602db0SMatthew 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 1224c4762a1bSJed Brown test: 1225c4762a1bSJed Brown suffix: nonconforming_tensor_3 1226c4762a1bSJed Brown nsize: 4 122730602db0SMatthew 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 1228c4762a1bSJed Brown test: 1229c4762a1bSJed Brown suffix: nonconforming_tensor_2_fv 1230c4762a1bSJed Brown nsize: 4 123130602db0SMatthew 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 1232c4762a1bSJed Brown test: 1233c4762a1bSJed Brown suffix: nonconforming_tensor_3_fv 1234c4762a1bSJed Brown nsize: 4 123530602db0SMatthew 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 1236c4762a1bSJed Brown test: 1237c4762a1bSJed Brown suffix: nonconforming_tensor_2_hi 1238c4762a1bSJed Brown requires: !single 1239c4762a1bSJed Brown nsize: 4 124030602db0SMatthew 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 1241c4762a1bSJed Brown test: 1242c4762a1bSJed Brown suffix: nonconforming_tensor_3_hi 1243c4762a1bSJed Brown requires: !single skip 1244c4762a1bSJed Brown nsize: 4 124530602db0SMatthew 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 1246c4762a1bSJed Brown test: 1247c4762a1bSJed Brown suffix: nonconforming_simplex_2 1248c4762a1bSJed Brown requires: triangle 1249c4762a1bSJed Brown nsize: 4 125030602db0SMatthew 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 1251c4762a1bSJed Brown test: 1252c4762a1bSJed Brown suffix: nonconforming_simplex_2_hi 1253c4762a1bSJed Brown requires: triangle !single 1254c4762a1bSJed Brown nsize: 4 125530602db0SMatthew G. Knepley args: -dist_dm_distribute -test_fe_jacobian -petscpartitioner_type simple -tree -dm_plex_max_projection_height 1 -petscspace_degree 4 -qorder 4 1256c4762a1bSJed Brown test: 1257c4762a1bSJed Brown suffix: nonconforming_simplex_2_fv 1258c4762a1bSJed Brown requires: triangle 1259c4762a1bSJed Brown nsize: 4 126030602db0SMatthew G. Knepley args: -dist_dm_distribute -test_fv_grad -test_injector -petsclimiter_type none -petscpartitioner_type simple -tree -num_comp 2 1261c4762a1bSJed Brown test: 1262c4762a1bSJed Brown suffix: nonconforming_simplex_3 1263c4762a1bSJed Brown requires: ctetgen 1264c4762a1bSJed Brown nsize: 4 126530602db0SMatthew 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 1266c4762a1bSJed Brown test: 1267c4762a1bSJed Brown suffix: nonconforming_simplex_3_hi 1268c4762a1bSJed Brown requires: ctetgen skip 1269c4762a1bSJed Brown nsize: 4 127030602db0SMatthew 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 1271c4762a1bSJed Brown test: 1272c4762a1bSJed Brown suffix: nonconforming_simplex_3_fv 1273c4762a1bSJed Brown requires: ctetgen 1274c4762a1bSJed Brown nsize: 4 127530602db0SMatthew 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 1276c4762a1bSJed Brown 1277d21efd2eSMatthew G. Knepley # 3D WXY on a triangular prism 1278d21efd2eSMatthew G. Knepley test: 1279d21efd2eSMatthew G. Knepley suffix: wxy_0 1280d21efd2eSMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -qorder 2 -porder 0 \ 1281e239af90SMatthew G. Knepley -petscspace_type sum \ 1282e239af90SMatthew G. Knepley -petscspace_variables 3 \ 1283e239af90SMatthew G. Knepley -petscspace_components 3 \ 1284e239af90SMatthew G. Knepley -petscspace_sum_spaces 2 \ 1285e239af90SMatthew G. Knepley -petscspace_sum_concatenate false \ 1286e239af90SMatthew G. Knepley -sumcomp_0_petscspace_variables 3 \ 1287e239af90SMatthew G. Knepley -sumcomp_0_petscspace_components 3 \ 1288e239af90SMatthew G. Knepley -sumcomp_0_petscspace_degree 1 \ 1289e239af90SMatthew G. Knepley -sumcomp_1_petscspace_variables 3 \ 1290e239af90SMatthew G. Knepley -sumcomp_1_petscspace_components 3 \ 1291e239af90SMatthew G. Knepley -sumcomp_1_petscspace_type wxy \ 1292e239af90SMatthew G. Knepley -petscdualspace_refcell triangular_prism \ 1293e239af90SMatthew G. Knepley -petscdualspace_form_degree 0 \ 1294e239af90SMatthew G. Knepley -petscdualspace_order 1 \ 1295e239af90SMatthew G. Knepley -petscdualspace_components 3 1296d21efd2eSMatthew G. Knepley 1297c4762a1bSJed Brown TEST*/ 1298c4762a1bSJed Brown 1299c4762a1bSJed Brown /* 1300c4762a1bSJed Brown # 2D Q_2 on a quadrilaterial Plex 1301c4762a1bSJed Brown test: 1302c4762a1bSJed Brown suffix: q2_2d_plex_0 130330602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -convergence 1304c4762a1bSJed Brown test: 1305c4762a1bSJed Brown suffix: q2_2d_plex_1 130630602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 1307c4762a1bSJed Brown test: 1308c4762a1bSJed Brown suffix: q2_2d_plex_2 130930602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 1310c4762a1bSJed Brown test: 1311c4762a1bSJed Brown suffix: q2_2d_plex_3 131230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 1 -shear_coords 1313c4762a1bSJed Brown test: 1314c4762a1bSJed Brown suffix: q2_2d_plex_4 131530602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder 2 -shear_coords 1316c4762a1bSJed Brown test: 1317c4762a1bSJed Brown suffix: q2_2d_plex_5 131830602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 0 -non_affine_coords 1319c4762a1bSJed Brown test: 1320c4762a1bSJed Brown suffix: q2_2d_plex_6 132130602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 1 -non_affine_coords 1322c4762a1bSJed Brown test: 1323c4762a1bSJed Brown suffix: q2_2d_plex_7 132430602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -petscspace_degree 2 -petscspace_poly_tensor 1 -qorder 2 -porder 2 -non_affine_coords 1325c4762a1bSJed Brown 1326c4762a1bSJed Brown test: 1327c4762a1bSJed Brown suffix: p1d_2d_6 1328e788613bSJoe Wallwork requires: mmg 1329c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -convergence -conv_refine 0 1330c4762a1bSJed Brown test: 1331c4762a1bSJed Brown suffix: p1d_2d_7 1332e788613bSJoe Wallwork requires: mmg 1333c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 1 -conv_refine 0 1334c4762a1bSJed Brown test: 1335c4762a1bSJed Brown suffix: p1d_2d_8 1336e788613bSJoe Wallwork requires: mmg 1337c4762a1bSJed Brown args: -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder 2 -conv_refine 0 1338c4762a1bSJed Brown */ 1339