1 static char help[] = "Tests for creation of submeshes\n\n"; 2 3 #include <petscdmplex.h> 4 5 PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm) 6 { 7 PetscFunctionBegin; 8 PetscCall(DMCreate(comm, dm)); 9 PetscCall(DMSetType(*dm, DMPLEX)); 10 PetscCall(DMSetFromOptions(*dm)); 11 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 12 PetscFunctionReturn(0); 13 } 14 15 PetscErrorCode CreateSubmesh(DM dm, PetscBool start, DM *subdm) 16 { 17 DMLabel label, map; 18 PetscInt cStart, cEnd, cStartSub, cEndSub, c; 19 20 PetscFunctionBegin; 21 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 22 PetscCall(DMCreateLabel(dm, "cells")); 23 PetscCall(DMGetLabel(dm, "cells", &label)); 24 PetscCall(DMLabelClearStratum(label, 1)); 25 if (start) {cStartSub = cStart; cEndSub = cEnd/2;} 26 else {cStartSub = cEnd/2; cEndSub = cEnd;} 27 for (c = cStartSub; c < cEndSub; ++c) PetscCall(DMLabelSetValue(label, c, 1)); 28 PetscCall(DMPlexFilter(dm, label, 1, subdm)); 29 PetscCall(PetscObjectSetName((PetscObject) *subdm, "Submesh")); 30 PetscCall(DMViewFromOptions(*subdm, NULL, "-dm_view")); 31 PetscCall(DMPlexGetSubpointMap(*subdm, &map)); 32 PetscCall(DMLabelView(map, PETSC_VIEWER_STDOUT_WORLD)); 33 PetscFunctionReturn(0); 34 } 35 36 int main(int argc, char **argv) 37 { 38 DM dm, subdm; 39 40 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 41 PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm)); 42 PetscCall(CreateSubmesh(dm, PETSC_TRUE, &subdm)); 43 PetscCall(DMSetFromOptions(subdm)); 44 PetscCall(DMDestroy(&subdm)); 45 PetscCall(CreateSubmesh(dm, PETSC_FALSE, &subdm)); 46 PetscCall(DMSetFromOptions(subdm)); 47 PetscCall(DMDestroy(&subdm)); 48 PetscCall(DMDestroy(&dm)); 49 PetscCall(PetscFinalize()); 50 return 0; 51 } 52 53 /*TEST 54 55 test: 56 suffix: 0 57 requires: triangle 58 args: -dm_coord_space 0 -dm_view ascii::ascii_info_detail -dm_plex_check_all 59 60 TEST*/ 61