1*df2803a8Sksagiyam static char help[] = "Tests for submesh creation for periodic meshes\n\n"; 2*df2803a8Sksagiyam 3*df2803a8Sksagiyam #include <petscdmplex.h> 4*df2803a8Sksagiyam #include <petsc/private/dmimpl.h> 5*df2803a8Sksagiyam #include <petscsf.h> 6*df2803a8Sksagiyam 7*df2803a8Sksagiyam /* 3 x 2 2D mesh periodic in x-direction using nproc = 1. 8*df2803a8Sksagiyam 9*df2803a8Sksagiyam Cell submesh (subdim = 2): 10*df2803a8Sksagiyam 11*df2803a8Sksagiyam 12--21--13--22--14--23---~ 6--12---7--13---8--14---~ 12*df2803a8Sksagiyam | | | ~ | | | ~ 13*df2803a8Sksagiyam 25 3 27 4 29 5 ~ 15 0 16 1 17 2 ~ 14*df2803a8Sksagiyam | | | ~ | | | ~ 15*df2803a8Sksagiyam 9--18--10--19--11--20---~ -> 3---9---4--10---5--11---~ 16*df2803a8Sksagiyam | | | ~ 17*df2803a8Sksagiyam 24 0 26 1 28 2 ~ 18*df2803a8Sksagiyam | | | ~ 19*df2803a8Sksagiyam 6--15---7--16---8--17---~ 20*df2803a8Sksagiyam 21*df2803a8Sksagiyam Facet submesh (subdim = 1): 22*df2803a8Sksagiyam 23*df2803a8Sksagiyam 12--21--13--22--14--23---~ 3---0---4---1---5---2---~ 24*df2803a8Sksagiyam | | | ~ 25*df2803a8Sksagiyam 25 3 27 4 29 5 ~ 26*df2803a8Sksagiyam | | | ~ 27*df2803a8Sksagiyam 9--18--10--19--11--20---~ -> 28*df2803a8Sksagiyam | | | ~ 29*df2803a8Sksagiyam 24 0 26 1 28 2 ~ 30*df2803a8Sksagiyam | | | ~ 31*df2803a8Sksagiyam 6--15---7--16---8--17---~ 32*df2803a8Sksagiyam 33*df2803a8Sksagiyam */ 34*df2803a8Sksagiyam 35*df2803a8Sksagiyam typedef struct { 36*df2803a8Sksagiyam PetscInt subdim; /* The topological dimension of the submesh */ 37*df2803a8Sksagiyam } AppCtx; 38*df2803a8Sksagiyam 39*df2803a8Sksagiyam PetscErrorCode ProcessOptions(AppCtx *options) 40*df2803a8Sksagiyam { 41*df2803a8Sksagiyam PetscFunctionBegin; 42*df2803a8Sksagiyam options->subdim = 2; 43*df2803a8Sksagiyam 44*df2803a8Sksagiyam PetscOptionsBegin(PETSC_COMM_SELF, "", "Filtering Problem Options", "DMPLEX"); 45*df2803a8Sksagiyam PetscCall(PetscOptionsBoundedInt("-subdim", "The topological dimension of the submesh", "ex74.c", options->subdim, &options->subdim, NULL, 0)); 46*df2803a8Sksagiyam PetscOptionsEnd(); 47*df2803a8Sksagiyam PetscFunctionReturn(PETSC_SUCCESS); 48*df2803a8Sksagiyam } 49*df2803a8Sksagiyam 50*df2803a8Sksagiyam int main(int argc, char **argv) 51*df2803a8Sksagiyam { 52*df2803a8Sksagiyam PetscInt dim = 2; 53*df2803a8Sksagiyam DM dm, subdm, coorddm; 54*df2803a8Sksagiyam PetscFE coordfe; 55*df2803a8Sksagiyam const PetscInt faces[2] = {3, 2}; 56*df2803a8Sksagiyam const PetscReal lower[2] = {0., 0.}, upper[2] = {3., 2.}; 57*df2803a8Sksagiyam const DMBoundaryType periodicity[2] = {DM_BOUNDARY_PERIODIC, DM_BOUNDARY_NONE}; 58*df2803a8Sksagiyam DMLabel filter; 59*df2803a8Sksagiyam const PetscInt filterValue = 1; 60*df2803a8Sksagiyam MPI_Comm comm; 61*df2803a8Sksagiyam PetscMPIInt size; 62*df2803a8Sksagiyam AppCtx user; 63*df2803a8Sksagiyam 64*df2803a8Sksagiyam PetscFunctionBeginUser; 65*df2803a8Sksagiyam PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 66*df2803a8Sksagiyam PetscCall(ProcessOptions(&user)); 67*df2803a8Sksagiyam comm = PETSC_COMM_WORLD; 68*df2803a8Sksagiyam PetscCallMPI(MPI_Comm_size(comm, &size)); 69*df2803a8Sksagiyam if (size != 1) { 70*df2803a8Sksagiyam PetscCall(PetscPrintf(comm, "This example is specifically designed for comm size == 1.\n")); 71*df2803a8Sksagiyam PetscCall(PetscFinalize()); 72*df2803a8Sksagiyam return 0; 73*df2803a8Sksagiyam } 74*df2803a8Sksagiyam PetscCall(DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, faces, lower, upper, periodicity, PETSC_TRUE, 0, PETSC_TRUE, &dm)); 75*df2803a8Sksagiyam /* Reset DG coordinates so that DMLocalizeCoordinates() will run again */ 76*df2803a8Sksagiyam /* in DMSetFromOptions() after CG coordinate FE is set. */ 77*df2803a8Sksagiyam dm->coordinates[1].dim = -1; 78*df2803a8Sksagiyam /* Localize coordinates on facets, too, if we are to make a facet submesh. */ 79*df2803a8Sksagiyam /* Otherwise, DG coordinate values will not be copied from the parent mesh. */ 80*df2803a8Sksagiyam PetscCall(DMSetFromOptions(dm)); 81*df2803a8Sksagiyam PetscCall(DMLabelCreate(PETSC_COMM_SELF, "filter", &filter)); 82*df2803a8Sksagiyam switch (user.subdim) { 83*df2803a8Sksagiyam case 2: 84*df2803a8Sksagiyam PetscCall(DMLabelSetValue(filter, 3, filterValue)); 85*df2803a8Sksagiyam PetscCall(DMLabelSetValue(filter, 4, filterValue)); 86*df2803a8Sksagiyam PetscCall(DMLabelSetValue(filter, 5, filterValue)); 87*df2803a8Sksagiyam break; 88*df2803a8Sksagiyam case 1: 89*df2803a8Sksagiyam PetscCall(DMLabelSetValue(filter, 21, filterValue)); 90*df2803a8Sksagiyam PetscCall(DMLabelSetValue(filter, 22, filterValue)); 91*df2803a8Sksagiyam PetscCall(DMLabelSetValue(filter, 23, filterValue)); 92*df2803a8Sksagiyam break; 93*df2803a8Sksagiyam default: 94*df2803a8Sksagiyam PetscCall(PetscPrintf(comm, "This example is only for subdim == {1, 2}.\n")); 95*df2803a8Sksagiyam PetscCall(PetscFinalize()); 96*df2803a8Sksagiyam return 0; 97*df2803a8Sksagiyam } 98*df2803a8Sksagiyam PetscCall(DMPlexFilter(dm, filter, filterValue, PETSC_FALSE, PETSC_FALSE, NULL, &subdm)); 99*df2803a8Sksagiyam PetscCall(DMLabelDestroy(&filter)); 100*df2803a8Sksagiyam PetscCall(DMDestroy(&dm)); 101*df2803a8Sksagiyam PetscCall(PetscObjectSetName((PetscObject)subdm, "Example_SubDM")); 102*df2803a8Sksagiyam PetscCall(DMViewFromOptions(subdm, NULL, "-dm_view")); 103*df2803a8Sksagiyam PetscCall(DMGetCellCoordinateDM(subdm, &coorddm)); 104*df2803a8Sksagiyam PetscCall(DMGetField(coorddm, 0, NULL, (PetscObject *)&coordfe)); 105*df2803a8Sksagiyam PetscCall(PetscFEView(coordfe, PETSC_VIEWER_STDOUT_WORLD)); 106*df2803a8Sksagiyam PetscCall(DMDestroy(&subdm)); 107*df2803a8Sksagiyam PetscCall(PetscFinalize()); 108*df2803a8Sksagiyam return 0; 109*df2803a8Sksagiyam } 110*df2803a8Sksagiyam 111*df2803a8Sksagiyam /*TEST 112*df2803a8Sksagiyam 113*df2803a8Sksagiyam testset: 114*df2803a8Sksagiyam args: -dm_coord_space 1 -dm_coord_petscspace_degree 1 -dm_localize 1 -dm_view ascii::ascii_info_detail 115*df2803a8Sksagiyam 116*df2803a8Sksagiyam test: 117*df2803a8Sksagiyam suffix: 0 118*df2803a8Sksagiyam args: -dm_sparse_localize 0 -dm_localize_height 0 -subdim 2 119*df2803a8Sksagiyam 120*df2803a8Sksagiyam test: 121*df2803a8Sksagiyam suffix: 1 122*df2803a8Sksagiyam args: -dm_sparse_localize 1 -dm_localize_height 0 -subdim 2 123*df2803a8Sksagiyam 124*df2803a8Sksagiyam test: 125*df2803a8Sksagiyam suffix: 2 126*df2803a8Sksagiyam args: -dm_sparse_localize 0 -dm_localize_height 1 -subdim 1 127*df2803a8Sksagiyam 128*df2803a8Sksagiyam test: 129*df2803a8Sksagiyam suffix: 3 130*df2803a8Sksagiyam args: -dm_sparse_localize 1 -dm_localize_height 1 -subdim 1 131*df2803a8Sksagiyam 132*df2803a8Sksagiyam TEST*/ 133