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