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