1c4762a1bSJed Brown static char help[] = "Tests for creation of submeshes\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown 59371c9d4SSatish 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 159371c9d4SSatish 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)); 239371c9d4SSatish Balay if (lower) { 249371c9d4SSatish Balay cStartSub = cStart; 259371c9d4SSatish Balay cEndSub = cEnd / 2; 269371c9d4SSatish Balay } else { 279371c9d4SSatish Balay cStartSub = cEnd / 2; 289371c9d4SSatish Balay cEndSub = cEnd; 299371c9d4SSatish 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 369371c9d4SSatish 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)); 489371c9d4SSatish Balay if (lower) { 499371c9d4SSatish Balay if (centroid[0] < 0.5) PetscCall(DMLabelSetValue(*label, c, 1)); 509371c9d4SSatish Balay } else { 519371c9d4SSatish Balay if (centroid[0] > 0.5) PetscCall(DMLabelSetValue(*label, c, 1)); 529371c9d4SSatish Balay } 53d3ef4daaSMatthew G. Knepley } 54d3ef4daaSMatthew G. Knepley PetscCall(DMPlexLabelComplete(dm, *label)); 55d3ef4daaSMatthew G. Knepley PetscFunctionReturn(0); 56d3ef4daaSMatthew G. Knepley } 57d3ef4daaSMatthew G. Knepley 58*836a3779SMatthew G. Knepley // Create a line of faces at a given x value 59*836a3779SMatthew G. Knepley static PetscErrorCode CreateLineLabel(DM dm, PetscReal x, DMLabel *label) { 60*836a3779SMatthew G. Knepley PetscReal centroid[3]; 61*836a3779SMatthew G. Knepley PetscInt fStart, fEnd; 62*836a3779SMatthew G. Knepley 63*836a3779SMatthew G. Knepley PetscFunctionBeginUser; 64*836a3779SMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(dm)); 65*836a3779SMatthew G. Knepley PetscCall(DMCreateLabel(dm, "faces")); 66*836a3779SMatthew G. Knepley PetscCall(DMGetLabel(dm, "faces", label)); 67*836a3779SMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 68*836a3779SMatthew G. Knepley for (PetscInt f = fStart; f < fEnd; ++f) { 69*836a3779SMatthew G. Knepley PetscCall(DMPlexComputeCellGeometryFVM(dm, f, NULL, centroid, NULL)); 70*836a3779SMatthew G. Knepley if (PetscAbsReal(centroid[0] - x) < PETSC_SMALL) PetscCall(DMLabelSetValue(*label, f, 1)); 71*836a3779SMatthew G. Knepley } 72*836a3779SMatthew G. Knepley PetscCall(DMPlexLabelComplete(dm, *label)); 73*836a3779SMatthew G. Knepley PetscFunctionReturn(0); 74*836a3779SMatthew G. Knepley } 75*836a3779SMatthew G. Knepley 76*836a3779SMatthew G. Knepley static PetscErrorCode CreateVolumeSubmesh(DM dm, PetscBool domain, PetscBool lower, DM *subdm) { 77c4762a1bSJed Brown DMLabel label, map; 78c4762a1bSJed Brown 79c4762a1bSJed Brown PetscFunctionBegin; 80d3ef4daaSMatthew G. Knepley if (domain) PetscCall(CreateHalfDomainLabel(dm, lower, &label)); 81d3ef4daaSMatthew G. Knepley else PetscCall(CreateHalfCellsLabel(dm, lower, &label)); 829566063dSJacob Faibussowitsch PetscCall(DMPlexFilter(dm, label, 1, subdm)); 839566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*subdm, "Submesh")); 84d3ef4daaSMatthew G. Knepley PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*subdm, "sub_")); 859566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*subdm, NULL, "-dm_view")); 869566063dSJacob Faibussowitsch PetscCall(DMPlexGetSubpointMap(*subdm, &map)); 87d3ef4daaSMatthew G. Knepley PetscCall(PetscObjectViewFromOptions((PetscObject)map, NULL, "-map_view")); 88c4762a1bSJed Brown PetscFunctionReturn(0); 89c4762a1bSJed Brown } 90c4762a1bSJed Brown 91*836a3779SMatthew G. Knepley static PetscErrorCode TestBoundaryField(DM dm) { 92*836a3779SMatthew G. Knepley DM subdm; 93*836a3779SMatthew G. Knepley DMLabel label, map; 94*836a3779SMatthew G. Knepley PetscFV fvm; 95*836a3779SMatthew G. Knepley Vec gv; 96*836a3779SMatthew G. Knepley PetscInt cdim; 97*836a3779SMatthew G. Knepley 98*836a3779SMatthew G. Knepley PetscFunctionBeginUser; 99*836a3779SMatthew G. Knepley PetscCall(CreateLineLabel(dm, 0.5, &label)); 100*836a3779SMatthew G. Knepley PetscCall(DMPlexFilter(dm, label, 1, &subdm)); 101*836a3779SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)subdm, "Submesh")); 102*836a3779SMatthew G. Knepley PetscCall(PetscObjectSetOptionsPrefix((PetscObject)subdm, "sub_")); 103*836a3779SMatthew G. Knepley PetscCall(DMViewFromOptions(subdm, NULL, "-dm_view")); 104*836a3779SMatthew G. Knepley PetscCall(DMPlexGetSubpointMap(subdm, &map)); 105*836a3779SMatthew G. Knepley PetscCall(PetscObjectViewFromOptions((PetscObject)map, NULL, "-map_view")); 106*836a3779SMatthew G. Knepley 107*836a3779SMatthew G. Knepley PetscCall(PetscFVCreate(PetscObjectComm((PetscObject)subdm), &fvm)); 108*836a3779SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)fvm, "testField")); 109*836a3779SMatthew G. Knepley PetscCall(PetscFVSetNumComponents(fvm, 1)); 110*836a3779SMatthew G. Knepley PetscCall(DMGetCoordinateDim(subdm, &cdim)); 111*836a3779SMatthew G. Knepley PetscCall(PetscFVSetSpatialDimension(fvm, cdim)); 112*836a3779SMatthew G. Knepley PetscCall(PetscFVSetFromOptions(fvm)); 113*836a3779SMatthew G. Knepley 114*836a3779SMatthew G. Knepley PetscCall(DMAddField(subdm, NULL, (PetscObject)fvm)); 115*836a3779SMatthew G. Knepley PetscCall(PetscFVDestroy(&fvm)); 116*836a3779SMatthew G. Knepley PetscCall(DMCreateDS(subdm)); 117*836a3779SMatthew G. Knepley 118*836a3779SMatthew G. Knepley PetscCall(DMCreateGlobalVector(subdm, &gv)); 119*836a3779SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)gv, "potential")); 120*836a3779SMatthew G. Knepley PetscCall(VecSet(gv, 1.)); 121*836a3779SMatthew G. Knepley PetscCall(VecViewFromOptions(gv, NULL, "-vec_view")); 122*836a3779SMatthew G. Knepley PetscCall(VecDestroy(&gv)); 123*836a3779SMatthew G. Knepley PetscCall(DMDestroy(&subdm)); 124*836a3779SMatthew G. Knepley PetscFunctionReturn(0); 125*836a3779SMatthew G. Knepley } 126*836a3779SMatthew G. Knepley 1279371c9d4SSatish Balay int main(int argc, char **argv) { 128c4762a1bSJed Brown DM dm, subdm; 129*836a3779SMatthew G. Knepley PetscBool volume = PETSC_TRUE, domain = PETSC_FALSE; 130c4762a1bSJed Brown 131327415f7SBarry Smith PetscFunctionBeginUser; 1329566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 133*836a3779SMatthew G. Knepley PetscCall(PetscOptionsGetBool(NULL, NULL, "-volume", &volume, NULL)); 134d3ef4daaSMatthew G. Knepley PetscCall(PetscOptionsGetBool(NULL, NULL, "-domain", &domain, NULL)); 135*836a3779SMatthew G. Knepley 1369566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm)); 137*836a3779SMatthew G. Knepley if (volume) { 138*836a3779SMatthew G. Knepley PetscCall(CreateVolumeSubmesh(dm, domain, PETSC_TRUE, &subdm)); 1399566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(subdm)); 1409566063dSJacob Faibussowitsch PetscCall(DMDestroy(&subdm)); 141*836a3779SMatthew G. Knepley PetscCall(CreateVolumeSubmesh(dm, domain, PETSC_FALSE, &subdm)); 1429566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(subdm)); 1439566063dSJacob Faibussowitsch PetscCall(DMDestroy(&subdm)); 144*836a3779SMatthew G. Knepley } else { 145*836a3779SMatthew G. Knepley PetscCall(TestBoundaryField(dm)); 146*836a3779SMatthew G. Knepley } 1479566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 1489566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 149b122ec5aSJacob Faibussowitsch return 0; 150c4762a1bSJed Brown } 151c4762a1bSJed Brown 152c4762a1bSJed Brown /*TEST 153c4762a1bSJed Brown 154c4762a1bSJed Brown test: 155c4762a1bSJed Brown suffix: 0 156c4762a1bSJed Brown requires: triangle 157d3ef4daaSMatthew G. Knepley args: -dm_coord_space 0 -sub_dm_plex_check_all \ 158d3ef4daaSMatthew G. Knepley -dm_view ascii::ascii_info_detail -sub_dm_view ascii::ascii_info_detail -map_view 159d3ef4daaSMatthew G. Knepley 160d3ef4daaSMatthew G. Knepley # These tests check that filtering is stable when boundary point ownership could change, so it needs 3 processes 161d3ef4daaSMatthew G. Knepley testset: 162d3ef4daaSMatthew G. Knepley nsize: 3 163d3ef4daaSMatthew G. Knepley requires: parmetis 164d3ef4daaSMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 20,20 -petscpartitioner_type parmetis -dm_distribute_overlap 1 -sub_dm_distribute 0 \ 165d3ef4daaSMatthew G. Knepley -sub_dm_plex_check_all -dm_view -sub_dm_view 166d3ef4daaSMatthew G. Knepley 167d3ef4daaSMatthew G. Knepley test: 168d3ef4daaSMatthew G. Knepley suffix: 1 169d3ef4daaSMatthew G. Knepley 170d3ef4daaSMatthew G. Knepley test: 171d3ef4daaSMatthew G. Knepley suffix: 2 172d3ef4daaSMatthew G. Knepley args: -domain 173c4762a1bSJed Brown 174*836a3779SMatthew G. Knepley # This test checks whether filter can extract a lower-dimensional manifold and output a field on it 175*836a3779SMatthew G. Knepley testset: 176*836a3779SMatthew G. Knepley args: -volume 0 -dm_plex_simplex 0 -sub_dm_view hdf5:subdm.h5 -vec_view hdf5:subdm.h5::append 177*836a3779SMatthew G. Knepley requires: hdf5 178*836a3779SMatthew G. Knepley 179*836a3779SMatthew G. Knepley test: 180*836a3779SMatthew G. Knepley suffix: surface_2d 181*836a3779SMatthew G. Knepley args: -dm_plex_box_faces 5,5 182*836a3779SMatthew G. Knepley 183*836a3779SMatthew G. Knepley test: 184*836a3779SMatthew G. Knepley suffix: surface_3d 185*836a3779SMatthew G. Knepley args: -dm_plex_box_faces 5,5,5 186*836a3779SMatthew G. Knepley 187c4762a1bSJed Brown TEST*/ 188