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 CHKERRQ(DMCreate(comm, dm)); 9 CHKERRQ(DMSetType(*dm, DMPLEX)); 10 CHKERRQ(DMSetFromOptions(*dm)); 11 CHKERRQ(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 CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 22 CHKERRQ(DMCreateLabel(dm, "cells")); 23 CHKERRQ(DMGetLabel(dm, "cells", &label)); 24 CHKERRQ(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) CHKERRQ(DMLabelSetValue(label, c, 1)); 28 CHKERRQ(DMPlexFilter(dm, label, 1, subdm)); 29 CHKERRQ(PetscObjectSetName((PetscObject) *subdm, "Submesh")); 30 CHKERRQ(DMViewFromOptions(*subdm, NULL, "-dm_view")); 31 CHKERRQ(DMPlexGetSubpointMap(*subdm, &map)); 32 CHKERRQ(DMLabelView(map, PETSC_VIEWER_STDOUT_WORLD)); 33 PetscFunctionReturn(0); 34 } 35 36 int main(int argc, char **argv) 37 { 38 DM dm, subdm; 39 PetscErrorCode ierr; 40 41 ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 42 CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &dm)); 43 CHKERRQ(CreateSubmesh(dm, PETSC_TRUE, &subdm)); 44 CHKERRQ(DMSetFromOptions(subdm)); 45 CHKERRQ(DMDestroy(&subdm)); 46 CHKERRQ(CreateSubmesh(dm, PETSC_FALSE, &subdm)); 47 CHKERRQ(DMSetFromOptions(subdm)); 48 CHKERRQ(DMDestroy(&subdm)); 49 CHKERRQ(DMDestroy(&dm)); 50 ierr = PetscFinalize(); 51 return ierr; 52 } 53 54 /*TEST 55 56 test: 57 suffix: 0 58 requires: triangle 59 args: -dm_coord_space 0 -dm_view ascii::ascii_info_detail -dm_plex_check_all 60 61 TEST*/ 62