1c4762a1bSJed Brown static char help[] = "Tests for creation of submeshes\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 4c4762a1bSJed Brown 5d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm) 6d71ae5a4SJacob Faibussowitsch { 7d3ef4daaSMatthew G. Knepley PetscFunctionBeginUser; 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")); 123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13c4762a1bSJed Brown } 14c4762a1bSJed Brown 15d3ef4daaSMatthew G. Knepley // Label half of the cells 16d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateHalfCellsLabel(DM dm, PetscBool lower, DMLabel *label) 17d71ae5a4SJacob Faibussowitsch { 18d3ef4daaSMatthew G. Knepley PetscInt cStart, cEnd, cStartSub, cEndSub; 19d3ef4daaSMatthew G. Knepley 20d3ef4daaSMatthew G. Knepley PetscFunctionBeginUser; 21d3ef4daaSMatthew G. Knepley PetscCall(DMCreateLabel(dm, "cells")); 22d3ef4daaSMatthew G. Knepley PetscCall(DMGetLabel(dm, "cells", label)); 23d3ef4daaSMatthew G. Knepley PetscCall(DMLabelClearStratum(*label, 1)); 24d3ef4daaSMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 259371c9d4SSatish Balay if (lower) { 269371c9d4SSatish Balay cStartSub = cStart; 279371c9d4SSatish Balay cEndSub = cEnd / 2; 289371c9d4SSatish Balay } else { 299371c9d4SSatish Balay cStartSub = cEnd / 2; 309371c9d4SSatish Balay cEndSub = cEnd; 319371c9d4SSatish Balay } 32d3ef4daaSMatthew G. Knepley for (PetscInt c = cStartSub; c < cEndSub; ++c) PetscCall(DMLabelSetValue(*label, c, 1)); 33d3ef4daaSMatthew G. Knepley PetscCall(DMPlexLabelComplete(dm, *label)); 343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 35d3ef4daaSMatthew G. Knepley } 36d3ef4daaSMatthew G. Knepley 37d3ef4daaSMatthew G. Knepley // Label everything on the right half of the domain 38d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateHalfDomainLabel(DM dm, PetscBool lower, PetscReal height, DMLabel *label) 39d71ae5a4SJacob Faibussowitsch { 40d3ef4daaSMatthew G. Knepley PetscReal centroid[3]; 41d3ef4daaSMatthew G. Knepley PetscInt cStart, cEnd, cdim; 42d3ef4daaSMatthew G. Knepley 43d3ef4daaSMatthew G. Knepley PetscFunctionBeginUser; 44a44be8dcSMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(dm)); 45d3ef4daaSMatthew G. Knepley PetscCall(DMCreateLabel(dm, "cells")); 46d3ef4daaSMatthew G. Knepley PetscCall(DMGetLabel(dm, "cells", label)); 47d3ef4daaSMatthew G. Knepley PetscCall(DMLabelClearStratum(*label, 1)); 48d3ef4daaSMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 49d3ef4daaSMatthew G. Knepley PetscCall(DMGetCoordinateDim(dm, &cdim)); 50d3ef4daaSMatthew G. Knepley for (PetscInt c = cStart; c < cEnd; ++c) { 51d3ef4daaSMatthew G. Knepley PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL)); 52a44be8dcSMatthew G. Knepley if (height > 0.0 && PetscAbsReal(centroid[1] - height) > PETSC_SMALL) continue; 539371c9d4SSatish Balay if (lower) { 549371c9d4SSatish Balay if (centroid[0] < 0.5) PetscCall(DMLabelSetValue(*label, c, 1)); 559371c9d4SSatish Balay } else { 569371c9d4SSatish Balay if (centroid[0] > 0.5) PetscCall(DMLabelSetValue(*label, c, 1)); 579371c9d4SSatish Balay } 58d3ef4daaSMatthew G. Knepley } 59d3ef4daaSMatthew G. Knepley PetscCall(DMPlexLabelComplete(dm, *label)); 603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 61d3ef4daaSMatthew G. Knepley } 62d3ef4daaSMatthew G. Knepley 63836a3779SMatthew G. Knepley // Create a line of faces at a given x value 64d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateLineLabel(DM dm, PetscReal x, DMLabel *label) 65d71ae5a4SJacob Faibussowitsch { 66836a3779SMatthew G. Knepley PetscReal centroid[3]; 67836a3779SMatthew G. Knepley PetscInt fStart, fEnd; 68836a3779SMatthew G. Knepley 69836a3779SMatthew G. Knepley PetscFunctionBeginUser; 70836a3779SMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(dm)); 71836a3779SMatthew G. Knepley PetscCall(DMCreateLabel(dm, "faces")); 72836a3779SMatthew G. Knepley PetscCall(DMGetLabel(dm, "faces", label)); 73836a3779SMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 74836a3779SMatthew G. Knepley for (PetscInt f = fStart; f < fEnd; ++f) { 75836a3779SMatthew G. Knepley PetscCall(DMPlexComputeCellGeometryFVM(dm, f, NULL, centroid, NULL)); 76836a3779SMatthew G. Knepley if (PetscAbsReal(centroid[0] - x) < PETSC_SMALL) PetscCall(DMLabelSetValue(*label, f, 1)); 77836a3779SMatthew G. Knepley } 78836a3779SMatthew G. Knepley PetscCall(DMPlexLabelComplete(dm, *label)); 793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 80836a3779SMatthew G. Knepley } 81836a3779SMatthew G. Knepley 82d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateVolumeSubmesh(DM dm, PetscBool domain, PetscBool lower, PetscReal height, DM *subdm) 83d71ae5a4SJacob Faibussowitsch { 84c4762a1bSJed Brown DMLabel label, map; 85c4762a1bSJed Brown 86c4762a1bSJed Brown PetscFunctionBegin; 87a44be8dcSMatthew G. Knepley if (domain) PetscCall(CreateHalfDomainLabel(dm, lower, height, &label)); 88d3ef4daaSMatthew G. Knepley else PetscCall(CreateHalfCellsLabel(dm, lower, &label)); 899566063dSJacob Faibussowitsch PetscCall(DMPlexFilter(dm, label, 1, subdm)); 909566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*subdm, "Submesh")); 91d3ef4daaSMatthew G. Knepley PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*subdm, "sub_")); 929566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*subdm, NULL, "-dm_view")); 939566063dSJacob Faibussowitsch PetscCall(DMPlexGetSubpointMap(*subdm, &map)); 94d3ef4daaSMatthew G. Knepley PetscCall(PetscObjectViewFromOptions((PetscObject)map, NULL, "-map_view")); 953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown 98d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestBoundaryField(DM dm) 99d71ae5a4SJacob Faibussowitsch { 100836a3779SMatthew G. Knepley DM subdm; 101836a3779SMatthew G. Knepley DMLabel label, map; 102836a3779SMatthew G. Knepley PetscFV fvm; 103836a3779SMatthew G. Knepley Vec gv; 104836a3779SMatthew G. Knepley PetscInt cdim; 105836a3779SMatthew G. Knepley 106836a3779SMatthew G. Knepley PetscFunctionBeginUser; 107836a3779SMatthew G. Knepley PetscCall(CreateLineLabel(dm, 0.5, &label)); 108836a3779SMatthew G. Knepley PetscCall(DMPlexFilter(dm, label, 1, &subdm)); 109836a3779SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)subdm, "Submesh")); 110836a3779SMatthew G. Knepley PetscCall(PetscObjectSetOptionsPrefix((PetscObject)subdm, "sub_")); 111836a3779SMatthew G. Knepley PetscCall(DMViewFromOptions(subdm, NULL, "-dm_view")); 112836a3779SMatthew G. Knepley PetscCall(DMPlexGetSubpointMap(subdm, &map)); 113836a3779SMatthew G. Knepley PetscCall(PetscObjectViewFromOptions((PetscObject)map, NULL, "-map_view")); 114836a3779SMatthew G. Knepley 115836a3779SMatthew G. Knepley PetscCall(PetscFVCreate(PetscObjectComm((PetscObject)subdm), &fvm)); 116836a3779SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)fvm, "testField")); 117836a3779SMatthew G. Knepley PetscCall(PetscFVSetNumComponents(fvm, 1)); 118836a3779SMatthew G. Knepley PetscCall(DMGetCoordinateDim(subdm, &cdim)); 119836a3779SMatthew G. Knepley PetscCall(PetscFVSetSpatialDimension(fvm, cdim)); 120836a3779SMatthew G. Knepley PetscCall(PetscFVSetFromOptions(fvm)); 121836a3779SMatthew G. Knepley 122836a3779SMatthew G. Knepley PetscCall(DMAddField(subdm, NULL, (PetscObject)fvm)); 123836a3779SMatthew G. Knepley PetscCall(PetscFVDestroy(&fvm)); 124836a3779SMatthew G. Knepley PetscCall(DMCreateDS(subdm)); 125836a3779SMatthew G. Knepley 126836a3779SMatthew G. Knepley PetscCall(DMCreateGlobalVector(subdm, &gv)); 127836a3779SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)gv, "potential")); 128836a3779SMatthew G. Knepley PetscCall(VecSet(gv, 1.)); 129836a3779SMatthew G. Knepley PetscCall(VecViewFromOptions(gv, NULL, "-vec_view")); 130836a3779SMatthew G. Knepley PetscCall(VecDestroy(&gv)); 131836a3779SMatthew G. Knepley PetscCall(DMDestroy(&subdm)); 1323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 133836a3779SMatthew G. Knepley } 134836a3779SMatthew G. Knepley 135d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 136d71ae5a4SJacob Faibussowitsch { 137c4762a1bSJed Brown DM dm, subdm; 138a44be8dcSMatthew G. Knepley PetscReal height = -1.0; 139836a3779SMatthew G. Knepley PetscBool volume = PETSC_TRUE, domain = PETSC_FALSE; 140c4762a1bSJed Brown 141327415f7SBarry Smith PetscFunctionBeginUser; 1429566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 143a44be8dcSMatthew G. Knepley PetscCall(PetscOptionsGetReal(NULL, NULL, "-height", &height, NULL)); 144836a3779SMatthew G. Knepley PetscCall(PetscOptionsGetBool(NULL, NULL, "-volume", &volume, NULL)); 145d3ef4daaSMatthew G. Knepley PetscCall(PetscOptionsGetBool(NULL, NULL, "-domain", &domain, NULL)); 146836a3779SMatthew G. Knepley 1479566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm)); 148836a3779SMatthew G. Knepley if (volume) { 149a44be8dcSMatthew G. Knepley PetscCall(CreateVolumeSubmesh(dm, domain, PETSC_TRUE, height, &subdm)); 1509566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(subdm)); 1519566063dSJacob Faibussowitsch PetscCall(DMDestroy(&subdm)); 152a44be8dcSMatthew G. Knepley PetscCall(CreateVolumeSubmesh(dm, domain, PETSC_FALSE, height, &subdm)); 1539566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(subdm)); 1549566063dSJacob Faibussowitsch PetscCall(DMDestroy(&subdm)); 155836a3779SMatthew G. Knepley } else { 156836a3779SMatthew G. Knepley PetscCall(TestBoundaryField(dm)); 157836a3779SMatthew G. Knepley } 1589566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 1599566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 160b122ec5aSJacob Faibussowitsch return 0; 161c4762a1bSJed Brown } 162c4762a1bSJed Brown 163c4762a1bSJed Brown /*TEST 164c4762a1bSJed Brown 165c4762a1bSJed Brown test: 166c4762a1bSJed Brown suffix: 0 167c4762a1bSJed Brown requires: triangle 168d3ef4daaSMatthew G. Knepley args: -dm_coord_space 0 -sub_dm_plex_check_all \ 169d3ef4daaSMatthew G. Knepley -dm_view ascii::ascii_info_detail -sub_dm_view ascii::ascii_info_detail -map_view 170d3ef4daaSMatthew G. Knepley 171*d0de8544SStefano Zampini test: 172*d0de8544SStefano Zampini suffix: 0_vtk 173*d0de8544SStefano Zampini requires: triangle 174*d0de8544SStefano Zampini args: -dm_coord_space 0 -sub_dm_plex_check_all \ 175*d0de8544SStefano Zampini -dm_view vtk: -sub_dm_view vtk: -map_view 176*d0de8544SStefano Zampini 177d3ef4daaSMatthew G. Knepley # These tests check that filtering is stable when boundary point ownership could change, so it needs 3 processes 178d3ef4daaSMatthew G. Knepley testset: 179d3ef4daaSMatthew G. Knepley nsize: 3 180d3ef4daaSMatthew G. Knepley requires: parmetis 181d3ef4daaSMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 20,20 -petscpartitioner_type parmetis -dm_distribute_overlap 1 -sub_dm_distribute 0 \ 182d3ef4daaSMatthew G. Knepley -sub_dm_plex_check_all -dm_view -sub_dm_view 183d3ef4daaSMatthew G. Knepley 184d3ef4daaSMatthew G. Knepley test: 185d3ef4daaSMatthew G. Knepley suffix: 1 186d3ef4daaSMatthew G. Knepley 187d3ef4daaSMatthew G. Knepley test: 188d3ef4daaSMatthew G. Knepley suffix: 2 189d3ef4daaSMatthew G. Knepley args: -domain 190c4762a1bSJed Brown 191a44be8dcSMatthew G. Knepley # This set tests that global numberings can be made when some strata are missing on a process 192a44be8dcSMatthew G. Knepley testset: 193a44be8dcSMatthew G. Knepley nsize: 3 194a44be8dcSMatthew G. Knepley requires: hdf5 195a44be8dcSMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 4,4 -petscpartitioner_type simple -sub_dm_distribute 0 \ 196a44be8dcSMatthew G. Knepley -sub_dm_plex_check_all -sub_dm_view hdf5:subdm.h5 197a44be8dcSMatthew G. Knepley 198a44be8dcSMatthew G. Knepley test: 199a44be8dcSMatthew G. Knepley suffix: 3 200a44be8dcSMatthew G. Knepley args: -domain -height 0.625 201a44be8dcSMatthew G. Knepley 202*d0de8544SStefano Zampini # This set tests that global numberings can be made when some strata are missing on a process 203*d0de8544SStefano Zampini testset: 204*d0de8544SStefano Zampini nsize: 3 205*d0de8544SStefano Zampini args: -dm_plex_simplex 0 -dm_plex_box_faces 4,4 -petscpartitioner_type simple -sub_dm_distribute 0 \ 206*d0de8544SStefano Zampini -sub_dm_plex_check_all -sub_dm_view {{vtk:subdm.vtk: vtk:subdm.vtu :subdm.txt :subdm_d.txt:ascii_info_detail}} 207*d0de8544SStefano Zampini 208*d0de8544SStefano Zampini test: 209*d0de8544SStefano Zampini suffix: 3_vtk 210*d0de8544SStefano Zampini args: -domain -height 0.625 211*d0de8544SStefano Zampini 212836a3779SMatthew G. Knepley # This test checks whether filter can extract a lower-dimensional manifold and output a field on it 213836a3779SMatthew G. Knepley testset: 214836a3779SMatthew G. Knepley requires: hdf5 215*d0de8544SStefano Zampini args: -volume 0 -dm_plex_simplex 0 -sub_dm_view hdf5:subdm.h5 -vec_view hdf5:subdm.h5::append 216836a3779SMatthew G. Knepley 217836a3779SMatthew G. Knepley test: 218836a3779SMatthew G. Knepley suffix: surface_2d 219836a3779SMatthew G. Knepley args: -dm_plex_box_faces 5,5 220836a3779SMatthew G. Knepley 221836a3779SMatthew G. Knepley test: 222836a3779SMatthew G. Knepley suffix: surface_3d 223836a3779SMatthew G. Knepley args: -dm_plex_box_faces 5,5,5 224836a3779SMatthew G. Knepley 22533905a2dSMatthew G. Knepley # This test checks that if cells are present in the SF, the dm is marked as having an overlap 22633905a2dSMatthew G. Knepley test: 22733905a2dSMatthew G. Knepley suffix: surface_2d_2 22833905a2dSMatthew G. Knepley nsize: 3 22933905a2dSMatthew G. Knepley args: -dm_plex_box_faces 5,5 -domain -height 0.625 23033905a2dSMatthew G. Knepley 231c4762a1bSJed Brown TEST*/ 232