1*c4762a1bSJed Brown static char help[] = "Test section ordering for FEM discretizations\n\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscdmplex.h> 4*c4762a1bSJed Brown #include <petscds.h> 5*c4762a1bSJed Brown 6*c4762a1bSJed Brown typedef struct { 7*c4762a1bSJed Brown PetscInt dim; /* The topological mesh dimension */ 8*c4762a1bSJed Brown PetscInt faces[3]; /* Number of faces per dimension */ 9*c4762a1bSJed Brown PetscBool simplex; /* Use simplices or hexes */ 10*c4762a1bSJed Brown } AppCtx; 11*c4762a1bSJed Brown 12*c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 13*c4762a1bSJed Brown { 14*c4762a1bSJed Brown PetscInt dim; 15*c4762a1bSJed Brown PetscErrorCode ierr; 16*c4762a1bSJed Brown 17*c4762a1bSJed Brown PetscFunctionBegin; 18*c4762a1bSJed Brown options->dim = 1; 19*c4762a1bSJed Brown options->faces[0] = 1; 20*c4762a1bSJed Brown options->faces[1] = 1; 21*c4762a1bSJed Brown options->faces[2] = 1; 22*c4762a1bSJed Brown options->simplex = PETSC_TRUE; 23*c4762a1bSJed Brown 24*c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Plex ex27", "DMPLEX");CHKERRQ(ierr); 25*c4762a1bSJed Brown ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex27.c", options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr); 26*c4762a1bSJed Brown ierr = PetscOptionsBool("-simplex", "Use simplices if true, otherwise tensor product cells", "ex27.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr); 27*c4762a1bSJed Brown dim = options->dim; 28*c4762a1bSJed Brown ierr = PetscOptionsIntArray("-faces", "Number of faces per dimension", "ex27.c", options->faces, &dim, NULL);CHKERRQ(ierr); 29*c4762a1bSJed Brown if (dim) options->dim = dim; 30*c4762a1bSJed Brown ierr = PetscOptionsEnd(); 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown if (options->dim > 3) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %d must in [0, 4)", options->dim); 33*c4762a1bSJed Brown PetscFunctionReturn(0); 34*c4762a1bSJed Brown } 35*c4762a1bSJed Brown 36*c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *options, DM *dm) 37*c4762a1bSJed Brown { 38*c4762a1bSJed Brown PetscInt dim = options->dim; 39*c4762a1bSJed Brown PetscInt *faces = options->faces; 40*c4762a1bSJed Brown PetscBool simplex = options->simplex; 41*c4762a1bSJed Brown PetscErrorCode ierr; 42*c4762a1bSJed Brown 43*c4762a1bSJed Brown PetscFunctionBegin; 44*c4762a1bSJed Brown ierr = DMPlexCreateBoxMesh(comm, dim, simplex, faces, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr); 45*c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) *dm, "Serial Mesh");CHKERRQ(ierr); 46*c4762a1bSJed Brown ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 47*c4762a1bSJed Brown PetscFunctionReturn(0); 48*c4762a1bSJed Brown } 49*c4762a1bSJed Brown 50*c4762a1bSJed Brown static PetscErrorCode TestLocalDofOrder(DM dm, AppCtx *ctx) 51*c4762a1bSJed Brown { 52*c4762a1bSJed Brown PetscFE fe[3]; 53*c4762a1bSJed Brown PetscSection s; 54*c4762a1bSJed Brown PetscInt dim, Nf, f; 55*c4762a1bSJed Brown PetscErrorCode ierr; 56*c4762a1bSJed Brown 57*c4762a1bSJed Brown PetscFunctionBegin; 58*c4762a1bSJed Brown ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 59*c4762a1bSJed Brown ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, dim, ctx->simplex, "field0_", -1, &fe[0]);CHKERRQ(ierr); 60*c4762a1bSJed Brown ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, ctx->simplex, "field1_", -1, &fe[1]);CHKERRQ(ierr); 61*c4762a1bSJed Brown ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, ctx->simplex, "field2_", -1, &fe[2]);CHKERRQ(ierr); 62*c4762a1bSJed Brown 63*c4762a1bSJed Brown ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = DMSetField(dm, 2, NULL, (PetscObject) fe[2]);CHKERRQ(ierr); 66*c4762a1bSJed Brown ierr = DMCreateDS(dm);CHKERRQ(ierr); 67*c4762a1bSJed Brown ierr = DMGetLocalSection(dm, &s);CHKERRQ(ierr); 68*c4762a1bSJed Brown ierr = PetscObjectViewFromOptions((PetscObject) s, NULL, "-dof_view");CHKERRQ(ierr); 69*c4762a1bSJed Brown 70*c4762a1bSJed Brown ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr); 71*c4762a1bSJed Brown for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&fe[f]);CHKERRQ(ierr);} 72*c4762a1bSJed Brown PetscFunctionReturn(0); 73*c4762a1bSJed Brown } 74*c4762a1bSJed Brown 75*c4762a1bSJed Brown int main(int argc, char **argv) 76*c4762a1bSJed Brown { 77*c4762a1bSJed Brown MPI_Comm comm; 78*c4762a1bSJed Brown DM dm; 79*c4762a1bSJed Brown AppCtx ctx; 80*c4762a1bSJed Brown PetscErrorCode ierr; 81*c4762a1bSJed Brown 82*c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 83*c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 84*c4762a1bSJed Brown ierr = ProcessOptions(comm, &ctx);CHKERRQ(ierr); 85*c4762a1bSJed Brown ierr = CreateMesh(comm, &ctx, &dm);CHKERRQ(ierr); 86*c4762a1bSJed Brown ierr = TestLocalDofOrder(dm, &ctx);CHKERRQ(ierr); 87*c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 88*c4762a1bSJed Brown ierr = PetscFinalize(); 89*c4762a1bSJed Brown return ierr; 90*c4762a1bSJed Brown } 91*c4762a1bSJed Brown 92*c4762a1bSJed Brown /*TEST 93*c4762a1bSJed Brown 94*c4762a1bSJed Brown test: 95*c4762a1bSJed Brown suffix: tri_pm 96*c4762a1bSJed Brown requires: triangle 97*c4762a1bSJed Brown args: -dim 2 -field0_petscspace_degree 2 -field1_petscspace_degree 1 -field2_petscspace_degree 1 -dm_view -dof_view 98*c4762a1bSJed Brown 99*c4762a1bSJed Brown test: 100*c4762a1bSJed Brown suffix: quad_pm 101*c4762a1bSJed Brown requires: 102*c4762a1bSJed Brown args: -dim 2 -simplex 0 -field0_petscspace_degree 2 -field1_petscspace_degree 1 -field2_petscspace_degree 1 -dm_view -dof_view 103*c4762a1bSJed Brown 104*c4762a1bSJed Brown test: 105*c4762a1bSJed Brown suffix: tri_fm 106*c4762a1bSJed Brown requires: triangle 107*c4762a1bSJed Brown args: -dim 2 -field0_petscspace_degree 2 -field1_petscspace_degree 1 -field2_petscspace_degree 1 -petscsection_point_major 0 -dm_view -dof_view 108*c4762a1bSJed Brown 109*c4762a1bSJed Brown test: 110*c4762a1bSJed Brown suffix: quad_fm 111*c4762a1bSJed Brown requires: 112*c4762a1bSJed Brown args: -dim 2 -simplex 0 -field0_petscspace_degree 2 -field1_petscspace_degree 1 -field2_petscspace_degree 1 -petscsection_point_major 0 -dm_view -dof_view 113*c4762a1bSJed Brown 114*c4762a1bSJed Brown TEST*/ 115