1*90d1c1a4SMatthew G. Knepley static char help[] = "Tests dof numberings for external integrators such as LibCEED.\n\n"; 2*90d1c1a4SMatthew G. Knepley 3*90d1c1a4SMatthew G. Knepley #include <petscdmplex.h> 4*90d1c1a4SMatthew G. Knepley #include <petscds.h> 5*90d1c1a4SMatthew G. Knepley 6*90d1c1a4SMatthew G. Knepley typedef struct { 7*90d1c1a4SMatthew G. Knepley PetscInt dummy; 8*90d1c1a4SMatthew G. Knepley } AppCtx; 9*90d1c1a4SMatthew G. Knepley 10*90d1c1a4SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 11*90d1c1a4SMatthew G. Knepley { 12*90d1c1a4SMatthew G. Knepley PetscFunctionBeginUser; 13*90d1c1a4SMatthew G. Knepley options->dummy = 1; 14*90d1c1a4SMatthew G. Knepley PetscOptionsBegin(comm, "", "Dof Ordering Options", "DMPLEX"); 15*90d1c1a4SMatthew G. Knepley //PetscCall(PetscOptionsInt("-num_refine", "Refine cycles", "ex46.c", options->Nr, &options->Nr, NULL)); 16*90d1c1a4SMatthew G. Knepley PetscOptionsEnd(); 17*90d1c1a4SMatthew G. Knepley PetscFunctionReturn(0); 18*90d1c1a4SMatthew G. Knepley } 19*90d1c1a4SMatthew G. Knepley 20*90d1c1a4SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 21*90d1c1a4SMatthew G. Knepley { 22*90d1c1a4SMatthew G. Knepley PetscFunctionBeginUser; 23*90d1c1a4SMatthew G. Knepley PetscCall(DMCreate(comm, dm)); 24*90d1c1a4SMatthew G. Knepley PetscCall(DMSetType(*dm, DMPLEX)); 25*90d1c1a4SMatthew G. Knepley PetscCall(DMSetFromOptions(*dm)); 26*90d1c1a4SMatthew G. Knepley PetscCall(DMSetApplicationContext(*dm, user)); 27*90d1c1a4SMatthew G. Knepley PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 28*90d1c1a4SMatthew G. Knepley PetscFunctionReturn(0); 29*90d1c1a4SMatthew G. Knepley } 30*90d1c1a4SMatthew G. Knepley 31*90d1c1a4SMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 32*90d1c1a4SMatthew G. Knepley { 33*90d1c1a4SMatthew G. Knepley DM cdm = dm; 34*90d1c1a4SMatthew G. Knepley PetscFE fe; 35*90d1c1a4SMatthew G. Knepley PetscInt dim; 36*90d1c1a4SMatthew G. Knepley 37*90d1c1a4SMatthew G. Knepley PetscFunctionBeginUser; 38*90d1c1a4SMatthew G. Knepley PetscCall(DMGetDimension(dm, &dim)); 39*90d1c1a4SMatthew G. Knepley PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, NULL, -1, &fe)); 40*90d1c1a4SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject) fe, "scalar")); 41*90d1c1a4SMatthew G. Knepley PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe)); 42*90d1c1a4SMatthew G. Knepley PetscCall(DMSetField(dm, 1, NULL, (PetscObject) fe)); 43*90d1c1a4SMatthew G. Knepley PetscCall(PetscFEDestroy(&fe)); 44*90d1c1a4SMatthew G. Knepley PetscCall(DMCreateDS(dm)); 45*90d1c1a4SMatthew G. Knepley while (cdm) { 46*90d1c1a4SMatthew G. Knepley PetscCall(DMCopyDisc(dm,cdm)); 47*90d1c1a4SMatthew G. Knepley PetscCall(DMGetCoarseDM(cdm, &cdm)); 48*90d1c1a4SMatthew G. Knepley } 49*90d1c1a4SMatthew G. Knepley PetscFunctionReturn(0); 50*90d1c1a4SMatthew G. Knepley } 51*90d1c1a4SMatthew G. Knepley 52*90d1c1a4SMatthew G. Knepley static PetscErrorCode CheckOffsets(DM dm, const char *domain_name, PetscInt label_value, PetscInt height) 53*90d1c1a4SMatthew G. Knepley { 54*90d1c1a4SMatthew G. Knepley const char *height_name[] = {"cells", "faces"}; 55*90d1c1a4SMatthew G. Knepley DMLabel domain_label = NULL; 56*90d1c1a4SMatthew G. Knepley DM cdm; 57*90d1c1a4SMatthew G. Knepley IS offIS; 58*90d1c1a4SMatthew G. Knepley PetscInt *offsets, Ncell, Ncl, Nc, n; 59*90d1c1a4SMatthew G. Knepley PetscInt Nf, f; 60*90d1c1a4SMatthew G. Knepley 61*90d1c1a4SMatthew G. Knepley PetscFunctionBeginUser; 62*90d1c1a4SMatthew G. Knepley if (domain_name) PetscCall(DMGetLabel(dm, domain_name, &domain_label)); 63*90d1c1a4SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_WORLD, "## %s: '%s' {%" PetscInt_FMT "}\n", height_name[height], domain_name?domain_name:"default", label_value)); 64*90d1c1a4SMatthew G. Knepley // Offsets for cell closures 65*90d1c1a4SMatthew G. Knepley PetscCall(DMGetNumFields(dm, &Nf)); 66*90d1c1a4SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 67*90d1c1a4SMatthew G. Knepley char name[PETSC_MAX_PATH_LEN]; 68*90d1c1a4SMatthew G. Knepley 69*90d1c1a4SMatthew G. Knepley PetscCall(DMPlexGetLocalOffsets(dm, domain_label, label_value, height, f, &Ncell, &Ncl, &Nc, &n, &offsets)); 70*90d1c1a4SMatthew G. Knepley PetscCall(ISCreateGeneral(PETSC_COMM_SELF, Ncell*Ncl, offsets, PETSC_OWN_POINTER, &offIS)); 71*90d1c1a4SMatthew G. Knepley PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "Field %" PetscInt_FMT " Offsets", f)); 72*90d1c1a4SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject) offIS, name)); 73*90d1c1a4SMatthew G. Knepley PetscCall(ISViewFromOptions(offIS, NULL, "-offsets_view")); 74*90d1c1a4SMatthew G. Knepley PetscCall(ISDestroy(&offIS)); 75*90d1c1a4SMatthew G. Knepley } 76*90d1c1a4SMatthew G. Knepley // Offsets for coordinates 77*90d1c1a4SMatthew G. Knepley { 78*90d1c1a4SMatthew G. Knepley Vec X; 79*90d1c1a4SMatthew G. Knepley PetscSection s; 80*90d1c1a4SMatthew G. Knepley const PetscScalar *x; 81*90d1c1a4SMatthew G. Knepley const char *cname; 82*90d1c1a4SMatthew G. Knepley PetscInt cdim; 83*90d1c1a4SMatthew G. Knepley PetscBool isDG = PETSC_FALSE; 84*90d1c1a4SMatthew G. Knepley 85*90d1c1a4SMatthew G. Knepley PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 86*90d1c1a4SMatthew G. Knepley if (!cdm) { 87*90d1c1a4SMatthew G. Knepley PetscCall(DMGetCoordinateDM(dm, &cdm)); cname = "Coordinates"; 88*90d1c1a4SMatthew G. Knepley PetscCall(DMGetCoordinatesLocal(dm, &X)); 89*90d1c1a4SMatthew G. Knepley } else { 90*90d1c1a4SMatthew G. Knepley isDG = PETSC_TRUE; 91*90d1c1a4SMatthew G. Knepley cname = "DG Coordinates"; 92*90d1c1a4SMatthew G. Knepley PetscCall(DMGetCellCoordinatesLocal(dm, &X)); 93*90d1c1a4SMatthew G. Knepley } 94*90d1c1a4SMatthew G. Knepley if (isDG && height) PetscFunctionReturn(0); 95*90d1c1a4SMatthew G. Knepley if (domain_name) PetscCall(DMGetLabel(cdm, domain_name, &domain_label)); 96*90d1c1a4SMatthew G. Knepley PetscCall(DMPlexGetLocalOffsets(cdm, domain_label, label_value, height, 0, &Ncell, &Ncl, &Nc, &n, &offsets)); 97*90d1c1a4SMatthew G. Knepley PetscCall(DMGetCoordinateDim(dm, &cdim)); 98*90d1c1a4SMatthew G. Knepley PetscCheck(Nc == cdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Geometric dimension %" PetscInt_FMT " should be %" PetscInt_FMT, Nc, cdim); 99*90d1c1a4SMatthew G. Knepley PetscCall(DMGetLocalSection(cdm, &s)); 100*90d1c1a4SMatthew G. Knepley PetscCall(VecGetArrayRead(X, &x)); 101*90d1c1a4SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s by element\n", cname)); 102*90d1c1a4SMatthew G. Knepley for (PetscInt c = 0; c < Ncell; ++c) { 103*90d1c1a4SMatthew G. Knepley for (PetscInt v = 0; v < Ncl; ++v) { 104*90d1c1a4SMatthew G. Knepley PetscInt off = offsets[c*Ncl+v], dgdof; 105*90d1c1a4SMatthew G. Knepley const PetscScalar *vx = &x[off]; 106*90d1c1a4SMatthew G. Knepley 107*90d1c1a4SMatthew G. Knepley if (isDG) { 108*90d1c1a4SMatthew G. Knepley PetscCall(PetscSectionGetDof(s, c, &dgdof)); 109*90d1c1a4SMatthew G. Knepley PetscCheck(Ncl*Nc == dgdof, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Offset size %" PetscInt_FMT " should be %" PetscInt_FMT, Ncl*Nc, dgdof); 110*90d1c1a4SMatthew G. Knepley } 111*90d1c1a4SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[%" PetscInt_FMT "] %" PetscInt_FMT " <-- %2" PetscInt_FMT " (% 4.2f, % 4.2f)\n", c, v, off, (double) PetscRealPart(vx[0]), (double) PetscRealPart(vx[1]))); 112*90d1c1a4SMatthew G. Knepley } 113*90d1c1a4SMatthew G. Knepley } 114*90d1c1a4SMatthew G. Knepley PetscCall(VecRestoreArrayRead(X, &x)); 115*90d1c1a4SMatthew G. Knepley PetscCall(PetscFree(offsets)); 116*90d1c1a4SMatthew G. Knepley } 117*90d1c1a4SMatthew G. Knepley PetscFunctionReturn(0); 118*90d1c1a4SMatthew G. Knepley } 119*90d1c1a4SMatthew G. Knepley 120*90d1c1a4SMatthew G. Knepley int main(int argc, char **argv) 121*90d1c1a4SMatthew G. Knepley { 122*90d1c1a4SMatthew G. Knepley DM dm; 123*90d1c1a4SMatthew G. Knepley AppCtx user; 124*90d1c1a4SMatthew G. Knepley 125*90d1c1a4SMatthew G. Knepley PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 126*90d1c1a4SMatthew G. Knepley PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 127*90d1c1a4SMatthew G. Knepley PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 128*90d1c1a4SMatthew G. Knepley PetscCall(SetupDiscretization(dm, &user)); 129*90d1c1a4SMatthew G. Knepley PetscCall(CheckOffsets(dm, NULL, 0, 0)); 130*90d1c1a4SMatthew G. Knepley PetscCall(CheckOffsets(dm, "Face Sets", 1, 1)); 131*90d1c1a4SMatthew G. Knepley PetscCall(DMDestroy(&dm)); 132*90d1c1a4SMatthew G. Knepley PetscCall(PetscFinalize()); 133*90d1c1a4SMatthew G. Knepley return 0; 134*90d1c1a4SMatthew G. Knepley } 135*90d1c1a4SMatthew G. Knepley 136*90d1c1a4SMatthew G. Knepley /*TEST 137*90d1c1a4SMatthew G. Knepley 138*90d1c1a4SMatthew G. Knepley test: 139*90d1c1a4SMatthew G. Knepley suffix: 0 140*90d1c1a4SMatthew G. Knepley requires: triangle 141*90d1c1a4SMatthew G. Knepley args: -dm_refine 1 -petscspace_degree 1 -dm_view -offsets_view 142*90d1c1a4SMatthew G. Knepley 143*90d1c1a4SMatthew G. Knepley test: 144*90d1c1a4SMatthew G. Knepley suffix: 1 145*90d1c1a4SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_bd periodic,none -dm_plex_box_faces 3,3 -dm_sparse_localize 0 -petscspace_degree 1 \ 146*90d1c1a4SMatthew G. Knepley -dm_view -offsets_view 147*90d1c1a4SMatthew G. Knepley 148*90d1c1a4SMatthew G. Knepley test: 149*90d1c1a4SMatthew G. Knepley suffix: cg_2d 150*90d1c1a4SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_bd none,none -dm_plex_box_faces 3,3 -petscspace_degree 1 \ 151*90d1c1a4SMatthew G. Knepley -dm_view -offsets_view 152*90d1c1a4SMatthew G. Knepley TEST*/ 153