190d1c1a4SMatthew G. Knepley static char help[] = "Tests dof numberings for external integrators such as LibCEED.\n\n"; 290d1c1a4SMatthew G. Knepley 390d1c1a4SMatthew G. Knepley #include <petscdmplex.h> 490d1c1a4SMatthew G. Knepley #include <petscds.h> 590d1c1a4SMatthew G. Knepley 690d1c1a4SMatthew G. Knepley typedef struct { 790d1c1a4SMatthew G. Knepley PetscInt dummy; 890d1c1a4SMatthew G. Knepley } AppCtx; 990d1c1a4SMatthew G. Knepley 10*d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 11*d71ae5a4SJacob Faibussowitsch { 1290d1c1a4SMatthew G. Knepley PetscFunctionBeginUser; 1390d1c1a4SMatthew G. Knepley options->dummy = 1; 1490d1c1a4SMatthew G. Knepley PetscOptionsBegin(comm, "", "Dof Ordering Options", "DMPLEX"); 1590d1c1a4SMatthew G. Knepley //PetscCall(PetscOptionsInt("-num_refine", "Refine cycles", "ex46.c", options->Nr, &options->Nr, NULL)); 1690d1c1a4SMatthew G. Knepley PetscOptionsEnd(); 1790d1c1a4SMatthew G. Knepley PetscFunctionReturn(0); 1890d1c1a4SMatthew G. Knepley } 1990d1c1a4SMatthew G. Knepley 20*d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 21*d71ae5a4SJacob Faibussowitsch { 2290d1c1a4SMatthew G. Knepley PetscFunctionBeginUser; 2390d1c1a4SMatthew G. Knepley PetscCall(DMCreate(comm, dm)); 2490d1c1a4SMatthew G. Knepley PetscCall(DMSetType(*dm, DMPLEX)); 2590d1c1a4SMatthew G. Knepley PetscCall(DMSetFromOptions(*dm)); 2690d1c1a4SMatthew G. Knepley PetscCall(DMSetApplicationContext(*dm, user)); 2790d1c1a4SMatthew G. Knepley PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 2890d1c1a4SMatthew G. Knepley PetscFunctionReturn(0); 2990d1c1a4SMatthew G. Knepley } 3090d1c1a4SMatthew G. Knepley 31*d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) 32*d71ae5a4SJacob Faibussowitsch { 3390d1c1a4SMatthew G. Knepley DM cdm = dm; 3490d1c1a4SMatthew G. Knepley PetscFE fe; 3590d1c1a4SMatthew G. Knepley PetscInt dim; 3690d1c1a4SMatthew G. Knepley 3790d1c1a4SMatthew G. Knepley PetscFunctionBeginUser; 3890d1c1a4SMatthew G. Knepley PetscCall(DMGetDimension(dm, &dim)); 3990d1c1a4SMatthew G. Knepley PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, NULL, -1, &fe)); 4090d1c1a4SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)fe, "scalar")); 4190d1c1a4SMatthew G. Knepley PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 4290d1c1a4SMatthew G. Knepley PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe)); 4390d1c1a4SMatthew G. Knepley PetscCall(PetscFEDestroy(&fe)); 4490d1c1a4SMatthew G. Knepley PetscCall(DMCreateDS(dm)); 4590d1c1a4SMatthew G. Knepley while (cdm) { 4690d1c1a4SMatthew G. Knepley PetscCall(DMCopyDisc(dm, cdm)); 4790d1c1a4SMatthew G. Knepley PetscCall(DMGetCoarseDM(cdm, &cdm)); 4890d1c1a4SMatthew G. Knepley } 4990d1c1a4SMatthew G. Knepley PetscFunctionReturn(0); 5090d1c1a4SMatthew G. Knepley } 5190d1c1a4SMatthew G. Knepley 52*d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckOffsets(DM dm, const char *domain_name, PetscInt label_value, PetscInt height) 53*d71ae5a4SJacob Faibussowitsch { 5490d1c1a4SMatthew G. Knepley const char *height_name[] = {"cells", "faces"}; 5590d1c1a4SMatthew G. Knepley DMLabel domain_label = NULL; 5690d1c1a4SMatthew G. Knepley DM cdm; 5790d1c1a4SMatthew G. Knepley IS offIS; 5890d1c1a4SMatthew G. Knepley PetscInt *offsets, Ncell, Ncl, Nc, n; 5990d1c1a4SMatthew G. Knepley PetscInt Nf, f; 6090d1c1a4SMatthew G. Knepley 6190d1c1a4SMatthew G. Knepley PetscFunctionBeginUser; 6290d1c1a4SMatthew G. Knepley if (domain_name) PetscCall(DMGetLabel(dm, domain_name, &domain_label)); 6390d1c1a4SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_WORLD, "## %s: '%s' {%" PetscInt_FMT "}\n", height_name[height], domain_name ? domain_name : "default", label_value)); 6490d1c1a4SMatthew G. Knepley // Offsets for cell closures 6590d1c1a4SMatthew G. Knepley PetscCall(DMGetNumFields(dm, &Nf)); 6690d1c1a4SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 6790d1c1a4SMatthew G. Knepley char name[PETSC_MAX_PATH_LEN]; 6890d1c1a4SMatthew G. Knepley 6990d1c1a4SMatthew G. Knepley PetscCall(DMPlexGetLocalOffsets(dm, domain_label, label_value, height, f, &Ncell, &Ncl, &Nc, &n, &offsets)); 7090d1c1a4SMatthew G. Knepley PetscCall(ISCreateGeneral(PETSC_COMM_SELF, Ncell * Ncl, offsets, PETSC_OWN_POINTER, &offIS)); 7190d1c1a4SMatthew G. Knepley PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN, "Field %" PetscInt_FMT " Offsets", f)); 7290d1c1a4SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)offIS, name)); 7390d1c1a4SMatthew G. Knepley PetscCall(ISViewFromOptions(offIS, NULL, "-offsets_view")); 7490d1c1a4SMatthew G. Knepley PetscCall(ISDestroy(&offIS)); 7590d1c1a4SMatthew G. Knepley } 7690d1c1a4SMatthew G. Knepley // Offsets for coordinates 7790d1c1a4SMatthew G. Knepley { 7890d1c1a4SMatthew G. Knepley Vec X; 7990d1c1a4SMatthew G. Knepley PetscSection s; 8090d1c1a4SMatthew G. Knepley const PetscScalar *x; 8190d1c1a4SMatthew G. Knepley const char *cname; 8290d1c1a4SMatthew G. Knepley PetscInt cdim; 8390d1c1a4SMatthew G. Knepley PetscBool isDG = PETSC_FALSE; 8490d1c1a4SMatthew G. Knepley 8590d1c1a4SMatthew G. Knepley PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 8690d1c1a4SMatthew G. Knepley if (!cdm) { 879371c9d4SSatish Balay PetscCall(DMGetCoordinateDM(dm, &cdm)); 889371c9d4SSatish Balay cname = "Coordinates"; 8990d1c1a4SMatthew G. Knepley PetscCall(DMGetCoordinatesLocal(dm, &X)); 9090d1c1a4SMatthew G. Knepley } else { 9190d1c1a4SMatthew G. Knepley isDG = PETSC_TRUE; 9290d1c1a4SMatthew G. Knepley cname = "DG Coordinates"; 9390d1c1a4SMatthew G. Knepley PetscCall(DMGetCellCoordinatesLocal(dm, &X)); 9490d1c1a4SMatthew G. Knepley } 9590d1c1a4SMatthew G. Knepley if (isDG && height) PetscFunctionReturn(0); 9690d1c1a4SMatthew G. Knepley if (domain_name) PetscCall(DMGetLabel(cdm, domain_name, &domain_label)); 9790d1c1a4SMatthew G. Knepley PetscCall(DMPlexGetLocalOffsets(cdm, domain_label, label_value, height, 0, &Ncell, &Ncl, &Nc, &n, &offsets)); 9890d1c1a4SMatthew G. Knepley PetscCall(DMGetCoordinateDim(dm, &cdim)); 9990d1c1a4SMatthew G. Knepley PetscCheck(Nc == cdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Geometric dimension %" PetscInt_FMT " should be %" PetscInt_FMT, Nc, cdim); 10090d1c1a4SMatthew G. Knepley PetscCall(DMGetLocalSection(cdm, &s)); 10190d1c1a4SMatthew G. Knepley PetscCall(VecGetArrayRead(X, &x)); 10290d1c1a4SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s by element\n", cname)); 10390d1c1a4SMatthew G. Knepley for (PetscInt c = 0; c < Ncell; ++c) { 10490d1c1a4SMatthew G. Knepley for (PetscInt v = 0; v < Ncl; ++v) { 10590d1c1a4SMatthew G. Knepley PetscInt off = offsets[c * Ncl + v], dgdof; 10690d1c1a4SMatthew G. Knepley const PetscScalar *vx = &x[off]; 10790d1c1a4SMatthew G. Knepley 10890d1c1a4SMatthew G. Knepley if (isDG) { 10990d1c1a4SMatthew G. Knepley PetscCall(PetscSectionGetDof(s, c, &dgdof)); 11090d1c1a4SMatthew G. Knepley PetscCheck(Ncl * Nc == dgdof, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Offset size %" PetscInt_FMT " should be %" PetscInt_FMT, Ncl * Nc, dgdof); 11190d1c1a4SMatthew G. Knepley } 112eefd97a2SJed Brown switch (cdim) { 113*d71ae5a4SJacob Faibussowitsch case 2: 114*d71ae5a4SJacob Faibussowitsch 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]))); 115*d71ae5a4SJacob Faibussowitsch break; 116eefd97a2SJed Brown case 3: 117eefd97a2SJed Brown PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[%" PetscInt_FMT "] %" PetscInt_FMT " <-- %2" PetscInt_FMT " (% 4.2f, % 4.2f, % 4.2f)\n", c, v, off, (double)PetscRealPart(vx[0]), (double)PetscRealPart(vx[1]), (double)PetscRealPart(vx[2]))); 118eefd97a2SJed Brown break; 119eefd97a2SJed Brown } 12090d1c1a4SMatthew G. Knepley } 12190d1c1a4SMatthew G. Knepley } 12290d1c1a4SMatthew G. Knepley PetscCall(VecRestoreArrayRead(X, &x)); 12390d1c1a4SMatthew G. Knepley PetscCall(PetscFree(offsets)); 12490d1c1a4SMatthew G. Knepley } 12590d1c1a4SMatthew G. Knepley PetscFunctionReturn(0); 12690d1c1a4SMatthew G. Knepley } 12790d1c1a4SMatthew G. Knepley 128*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 129*d71ae5a4SJacob Faibussowitsch { 13090d1c1a4SMatthew G. Knepley DM dm; 13190d1c1a4SMatthew G. Knepley AppCtx user; 13290d1c1a4SMatthew G. Knepley 133327415f7SBarry Smith PetscFunctionBeginUser; 13490d1c1a4SMatthew G. Knepley PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 13590d1c1a4SMatthew G. Knepley PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 13690d1c1a4SMatthew G. Knepley PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 13790d1c1a4SMatthew G. Knepley PetscCall(SetupDiscretization(dm, &user)); 13890d1c1a4SMatthew G. Knepley PetscCall(CheckOffsets(dm, NULL, 0, 0)); 13990d1c1a4SMatthew G. Knepley PetscCall(CheckOffsets(dm, "Face Sets", 1, 1)); 14090d1c1a4SMatthew G. Knepley PetscCall(DMDestroy(&dm)); 14190d1c1a4SMatthew G. Knepley PetscCall(PetscFinalize()); 14290d1c1a4SMatthew G. Knepley return 0; 14390d1c1a4SMatthew G. Knepley } 14490d1c1a4SMatthew G. Knepley 14590d1c1a4SMatthew G. Knepley /*TEST 14690d1c1a4SMatthew G. Knepley 14790d1c1a4SMatthew G. Knepley test: 14890d1c1a4SMatthew G. Knepley suffix: 0 14990d1c1a4SMatthew G. Knepley requires: triangle 15090d1c1a4SMatthew G. Knepley args: -dm_refine 1 -petscspace_degree 1 -dm_view -offsets_view 15190d1c1a4SMatthew G. Knepley 15290d1c1a4SMatthew G. Knepley test: 15390d1c1a4SMatthew G. Knepley suffix: 1 15490d1c1a4SMatthew 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 \ 15590d1c1a4SMatthew G. Knepley -dm_view -offsets_view 15690d1c1a4SMatthew G. Knepley 15790d1c1a4SMatthew G. Knepley test: 15890d1c1a4SMatthew G. Knepley suffix: cg_2d 15990d1c1a4SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_bd none,none -dm_plex_box_faces 3,3 -petscspace_degree 1 \ 16090d1c1a4SMatthew G. Knepley -dm_view -offsets_view 16190d1c1a4SMatthew G. Knepley TEST*/ 162