13ead50e0SBlaise Bourdin static char help[] = "Test DMPlexGetCellType\n\n"; 23ead50e0SBlaise Bourdin 33ead50e0SBlaise Bourdin #include <petsc.h> 43ead50e0SBlaise Bourdin 5*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 6*d71ae5a4SJacob Faibussowitsch { 73ead50e0SBlaise Bourdin DM dm, pdm; 83ead50e0SBlaise Bourdin char ifilename[PETSC_MAX_PATH_LEN]; 93ead50e0SBlaise Bourdin PetscInt pStart, pEnd, p; 103ead50e0SBlaise Bourdin DMPolytopeType cellType; 113ead50e0SBlaise Bourdin DMLabel label; 123ead50e0SBlaise Bourdin 133ead50e0SBlaise Bourdin PetscFunctionBeginUser; 143ead50e0SBlaise Bourdin PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 153ead50e0SBlaise Bourdin PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FEM Layout Options", "ex97"); 163ead50e0SBlaise Bourdin PetscCall(PetscOptionsString("-i", "Filename to read", "ex97", ifilename, ifilename, sizeof(ifilename), NULL)); 173ead50e0SBlaise Bourdin PetscOptionsEnd(); 183ead50e0SBlaise Bourdin 193ead50e0SBlaise Bourdin PetscCall(DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename, NULL, PETSC_TRUE, &dm)); 203ead50e0SBlaise Bourdin PetscCall(DMPlexDistributeSetDefault(dm, PETSC_FALSE)); 213ead50e0SBlaise Bourdin PetscCall(DMSetFromOptions(dm)); 223ead50e0SBlaise Bourdin 233ead50e0SBlaise Bourdin PetscCall(DMPlexDistribute(dm, 0, NULL, &pdm)); 243ead50e0SBlaise Bourdin if (pdm) { 253ead50e0SBlaise Bourdin PetscCall(DMDestroy(&dm)); 263ead50e0SBlaise Bourdin dm = pdm; 273ead50e0SBlaise Bourdin } 283ead50e0SBlaise Bourdin PetscCall(PetscObjectSetName((PetscObject)dm, "ex97")); 293ead50e0SBlaise Bourdin PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 303ead50e0SBlaise Bourdin 313ead50e0SBlaise Bourdin PetscCall(DMGetLabel(dm, "celltype", &label)); 323ead50e0SBlaise Bourdin PetscCall(DMLabelView(label, PETSC_VIEWER_STDOUT_WORLD)); 333ead50e0SBlaise Bourdin PetscCall(DMPlexGetHeightStratum(dm, 0, &pStart, &pEnd)); 343ead50e0SBlaise Bourdin for (p = pStart; p < pEnd; ++p) { 353ead50e0SBlaise Bourdin PetscCall(DMPlexGetCellType(dm, p, &cellType)); 363ead50e0SBlaise Bourdin PetscCall(PetscPrintf(PETSC_COMM_SELF, "cell: %" PetscInt_FMT " type: %d\n", p, cellType)); 373ead50e0SBlaise Bourdin } 383ead50e0SBlaise Bourdin PetscCall(DMDestroy(&dm)); 393ead50e0SBlaise Bourdin 403ead50e0SBlaise Bourdin PetscCall(PetscFinalize()); 413ead50e0SBlaise Bourdin return 0; 423ead50e0SBlaise Bourdin } 433ead50e0SBlaise Bourdin 443ead50e0SBlaise Bourdin /*TEST 453ead50e0SBlaise Bourdin build: 463ead50e0SBlaise Bourdin requires: !complex 473ead50e0SBlaise Bourdin testset: 483ead50e0SBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/doublet-tet.msh -dm_view 493ead50e0SBlaise Bourdin nsize: 1 503ead50e0SBlaise Bourdin test: 513ead50e0SBlaise Bourdin suffix: 0 523ead50e0SBlaise Bourdin args: 533ead50e0SBlaise Bourdin TEST*/ 54