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