1c4762a1bSJed Brown static char help[] = "Tests for creation of submeshes\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown 5*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm) { 6d3ef4daaSMatthew G. Knepley PetscFunctionBeginUser; 79566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 89566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 99566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 109566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 11c4762a1bSJed Brown PetscFunctionReturn(0); 12c4762a1bSJed Brown } 13c4762a1bSJed Brown 14d3ef4daaSMatthew G. Knepley // Label half of the cells 15*9371c9d4SSatish Balay static PetscErrorCode CreateHalfCellsLabel(DM dm, PetscBool lower, DMLabel *label) { 16d3ef4daaSMatthew G. Knepley PetscInt cStart, cEnd, cStartSub, cEndSub; 17d3ef4daaSMatthew G. Knepley 18d3ef4daaSMatthew G. Knepley PetscFunctionBeginUser; 19d3ef4daaSMatthew G. Knepley PetscCall(DMCreateLabel(dm, "cells")); 20d3ef4daaSMatthew G. Knepley PetscCall(DMGetLabel(dm, "cells", label)); 21d3ef4daaSMatthew G. Knepley PetscCall(DMLabelClearStratum(*label, 1)); 22d3ef4daaSMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 23*9371c9d4SSatish Balay if (lower) { 24*9371c9d4SSatish Balay cStartSub = cStart; 25*9371c9d4SSatish Balay cEndSub = cEnd / 2; 26*9371c9d4SSatish Balay } else { 27*9371c9d4SSatish Balay cStartSub = cEnd / 2; 28*9371c9d4SSatish Balay cEndSub = cEnd; 29*9371c9d4SSatish Balay } 30d3ef4daaSMatthew G. Knepley for (PetscInt c = cStartSub; c < cEndSub; ++c) PetscCall(DMLabelSetValue(*label, c, 1)); 31d3ef4daaSMatthew G. Knepley PetscCall(DMPlexLabelComplete(dm, *label)); 32d3ef4daaSMatthew G. Knepley PetscFunctionReturn(0); 33d3ef4daaSMatthew G. Knepley } 34d3ef4daaSMatthew G. Knepley 35d3ef4daaSMatthew G. Knepley // Label everything on the right half of the domain 36*9371c9d4SSatish Balay static PetscErrorCode CreateHalfDomainLabel(DM dm, PetscBool lower, DMLabel *label) { 37d3ef4daaSMatthew G. Knepley PetscReal centroid[3]; 38d3ef4daaSMatthew G. Knepley PetscInt cStart, cEnd, cdim; 39d3ef4daaSMatthew G. Knepley 40d3ef4daaSMatthew G. Knepley PetscFunctionBeginUser; 41d3ef4daaSMatthew G. Knepley PetscCall(DMCreateLabel(dm, "cells")); 42d3ef4daaSMatthew G. Knepley PetscCall(DMGetLabel(dm, "cells", label)); 43d3ef4daaSMatthew G. Knepley PetscCall(DMLabelClearStratum(*label, 1)); 44d3ef4daaSMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 45d3ef4daaSMatthew G. Knepley PetscCall(DMGetCoordinateDim(dm, &cdim)); 46d3ef4daaSMatthew G. Knepley for (PetscInt c = cStart; c < cEnd; ++c) { 47d3ef4daaSMatthew G. Knepley PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL)); 48*9371c9d4SSatish Balay if (lower) { 49*9371c9d4SSatish Balay if (centroid[0] < 0.5) PetscCall(DMLabelSetValue(*label, c, 1)); 50*9371c9d4SSatish Balay } else { 51*9371c9d4SSatish Balay if (centroid[0] > 0.5) PetscCall(DMLabelSetValue(*label, c, 1)); 52*9371c9d4SSatish Balay } 53d3ef4daaSMatthew G. Knepley } 54d3ef4daaSMatthew G. Knepley PetscCall(DMPlexLabelComplete(dm, *label)); 55d3ef4daaSMatthew G. Knepley PetscFunctionReturn(0); 56d3ef4daaSMatthew G. Knepley } 57d3ef4daaSMatthew G. Knepley 58*9371c9d4SSatish Balay PetscErrorCode CreateSubmesh(DM dm, PetscBool domain, PetscBool lower, DM *subdm) { 59c4762a1bSJed Brown DMLabel label, map; 60c4762a1bSJed Brown 61c4762a1bSJed Brown PetscFunctionBegin; 62d3ef4daaSMatthew G. Knepley if (domain) PetscCall(CreateHalfDomainLabel(dm, lower, &label)); 63d3ef4daaSMatthew G. Knepley else PetscCall(CreateHalfCellsLabel(dm, lower, &label)); 649566063dSJacob Faibussowitsch PetscCall(DMPlexFilter(dm, label, 1, subdm)); 659566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*subdm, "Submesh")); 66d3ef4daaSMatthew G. Knepley PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*subdm, "sub_")); 679566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*subdm, NULL, "-dm_view")); 689566063dSJacob Faibussowitsch PetscCall(DMPlexGetSubpointMap(*subdm, &map)); 69d3ef4daaSMatthew G. Knepley PetscCall(PetscObjectViewFromOptions((PetscObject)map, NULL, "-map_view")); 70c4762a1bSJed Brown PetscFunctionReturn(0); 71c4762a1bSJed Brown } 72c4762a1bSJed Brown 73*9371c9d4SSatish Balay int main(int argc, char **argv) { 74c4762a1bSJed Brown DM dm, subdm; 75d3ef4daaSMatthew G. Knepley PetscBool domain = PETSC_FALSE; 76c4762a1bSJed Brown 77327415f7SBarry Smith PetscFunctionBeginUser; 789566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 79d3ef4daaSMatthew G. Knepley PetscCall(PetscOptionsGetBool(NULL, NULL, "-domain", &domain, NULL)); 809566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm)); 81d3ef4daaSMatthew G. Knepley PetscCall(CreateSubmesh(dm, domain, PETSC_TRUE, &subdm)); 829566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(subdm)); 839566063dSJacob Faibussowitsch PetscCall(DMDestroy(&subdm)); 84d3ef4daaSMatthew G. Knepley PetscCall(CreateSubmesh(dm, domain, PETSC_FALSE, &subdm)); 859566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(subdm)); 869566063dSJacob Faibussowitsch PetscCall(DMDestroy(&subdm)); 879566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 889566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 89b122ec5aSJacob Faibussowitsch return 0; 90c4762a1bSJed Brown } 91c4762a1bSJed Brown 92c4762a1bSJed Brown /*TEST 93c4762a1bSJed Brown 94c4762a1bSJed Brown test: 95c4762a1bSJed Brown suffix: 0 96c4762a1bSJed Brown requires: triangle 97d3ef4daaSMatthew G. Knepley args: -dm_coord_space 0 -sub_dm_plex_check_all \ 98d3ef4daaSMatthew G. Knepley -dm_view ascii::ascii_info_detail -sub_dm_view ascii::ascii_info_detail -map_view 99d3ef4daaSMatthew G. Knepley 100d3ef4daaSMatthew G. Knepley # These tests check that filtering is stable when boundary point ownership could change, so it needs 3 processes 101d3ef4daaSMatthew G. Knepley testset: 102d3ef4daaSMatthew G. Knepley nsize: 3 103d3ef4daaSMatthew G. Knepley requires: parmetis 104d3ef4daaSMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 20,20 -petscpartitioner_type parmetis -dm_distribute_overlap 1 -sub_dm_distribute 0 \ 105d3ef4daaSMatthew G. Knepley -sub_dm_plex_check_all -dm_view -sub_dm_view 106d3ef4daaSMatthew G. Knepley 107d3ef4daaSMatthew G. Knepley test: 108d3ef4daaSMatthew G. Knepley suffix: 1 109d3ef4daaSMatthew G. Knepley 110d3ef4daaSMatthew G. Knepley test: 111d3ef4daaSMatthew G. Knepley suffix: 2 112d3ef4daaSMatthew G. Knepley args: -domain 113c4762a1bSJed Brown 114c4762a1bSJed Brown TEST*/ 115