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; 89566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 99566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 109566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 119566063dSJacob Faibussowitsch PetscCall(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; 219566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 229566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "cells")); 239566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "cells", &label)); 249566063dSJacob Faibussowitsch PetscCall(DMLabelClearStratum(label, 1)); 25c4762a1bSJed Brown if (start) {cStartSub = cStart; cEndSub = cEnd/2;} 26c4762a1bSJed Brown else {cStartSub = cEnd/2; cEndSub = cEnd;} 279566063dSJacob Faibussowitsch for (c = cStartSub; c < cEndSub; ++c) PetscCall(DMLabelSetValue(label, c, 1)); 289566063dSJacob Faibussowitsch PetscCall(DMPlexFilter(dm, label, 1, subdm)); 299566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) *subdm, "Submesh")); 309566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*subdm, NULL, "-dm_view")); 319566063dSJacob Faibussowitsch PetscCall(DMPlexGetSubpointMap(*subdm, &map)); 329566063dSJacob Faibussowitsch PetscCall(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 40*327415f7SBarry Smith PetscFunctionBeginUser; 419566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 429566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm)); 439566063dSJacob Faibussowitsch PetscCall(CreateSubmesh(dm, PETSC_TRUE, &subdm)); 449566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(subdm)); 459566063dSJacob Faibussowitsch PetscCall(DMDestroy(&subdm)); 469566063dSJacob Faibussowitsch PetscCall(CreateSubmesh(dm, PETSC_FALSE, &subdm)); 479566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(subdm)); 489566063dSJacob Faibussowitsch PetscCall(DMDestroy(&subdm)); 499566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 509566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 51b122ec5aSJacob Faibussowitsch return 0; 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