190a64276Sksagiyam static char help[] = "Tests for submesh creation\n\n"; 290a64276Sksagiyam 390a64276Sksagiyam #include <petscdmplex.h> 490a64276Sksagiyam #include <petscsf.h> 590a64276Sksagiyam 690a64276Sksagiyam /* Submesh of a 2 x 2 mesh using 3 processes. 790a64276Sksagiyam 890a64276Sksagiyam Local numbering on each rank: 990a64276Sksagiyam 1090a64276Sksagiyam (6)(16)-(7)(17)-(8) 5--14---6--15---7 (10)(20)(11)(21)(12) 1190a64276Sksagiyam | | | | | | | | | 1290a64276Sksagiyam (18) (1)(19) (2)(20) 16 0 17 1 18 (22) (2)(23) (3)(24) 1390a64276Sksagiyam | | | | | | | | | 1490a64276Sksagiyam (5)(15)(11)(22)(12) 4--13-(11)(22)(12) (9)(19)--6--14---7 1590a64276Sksagiyam | | | | | | | | | 1690a64276Sksagiyam 14 0 (23) (3)(24) (20) (2)(23) (3)(24) (18) (1) 15 0 16 1790a64276Sksagiyam | | | | | | | | | 1890a64276Sksagiyam 4--13--(9)(21)(10) (8)(19)-(9)(21)(10) (8)(17)--4--13---5 1990a64276Sksagiyam 2090a64276Sksagiyam mesh_0 mesh_1 mesh_2 2190a64276Sksagiyam 2290a64276Sksagiyam where () represents ghost points. We extract the left 2 cells. 2390a64276Sksagiyam With sanitize_submesh = PETSC_FALSE, we get: 2490a64276Sksagiyam 2590a64276Sksagiyam (4)(11)-(5) 3---9---4 2690a64276Sksagiyam | | | | 2790a64276Sksagiyam (12) (1)(13) 10 0 11 2890a64276Sksagiyam | | | | 2990a64276Sksagiyam (3)(10) (7) 2---8---7 3090a64276Sksagiyam | | | | 3190a64276Sksagiyam 9 0 (14) (13) (1) 14 3290a64276Sksagiyam | | | | 3390a64276Sksagiyam 2---8--(6) (5)(12)--6 3490a64276Sksagiyam 3590a64276Sksagiyam On the other hand, with sanitize_submesh = PETSC_TRUE, we get: 3690a64276Sksagiyam 3790a64276Sksagiyam (4)(11)-(5) 3---9---4 3890a64276Sksagiyam | | | | 3990a64276Sksagiyam (12) (1)(13) 10 0 11 4090a64276Sksagiyam | | | | 4190a64276Sksagiyam (3)(10) (7) 2---8---7 4290a64276Sksagiyam | | | | 4390a64276Sksagiyam 9 0 14 (13) (1)(14) 4490a64276Sksagiyam | | | | 4590a64276Sksagiyam 2---8---6 (5)(12)-(6) 4690a64276Sksagiyam 4790a64276Sksagiyam submesh_0 submesh_1 submesh_2 4890a64276Sksagiyam 4990a64276Sksagiyam as points 15 and 4 of mesh_2 are in the closure of a submesh cell owned by rank 0 (point 0 of submesh_0), 5090a64276Sksagiyam and not in the closure of any submesh cell owned by rank 1. 5190a64276Sksagiyam 5290a64276Sksagiyam */ 5390a64276Sksagiyam 5490a64276Sksagiyam typedef struct { 5590a64276Sksagiyam PetscBool ignoreLabelHalo; /* Ignore filter values in the halo. */ 5690a64276Sksagiyam PetscBool sanitizeSubmesh; /* Sanitize submesh. */ 5790a64276Sksagiyam } AppCtx; 5890a64276Sksagiyam 5990a64276Sksagiyam PetscErrorCode ProcessOptions(AppCtx *options) 6090a64276Sksagiyam { 6190a64276Sksagiyam PetscFunctionBegin; 6290a64276Sksagiyam options->ignoreLabelHalo = PETSC_FALSE; 6390a64276Sksagiyam options->sanitizeSubmesh = PETSC_FALSE; 6490a64276Sksagiyam 6590a64276Sksagiyam PetscOptionsBegin(PETSC_COMM_SELF, "", "Filtering Problem Options", "DMPLEX"); 6690a64276Sksagiyam PetscCall(PetscOptionsBool("-ignore_label_halo", "Ignore filter values in the halo", "ex80.c", options->ignoreLabelHalo, &options->ignoreLabelHalo, NULL)); 6790a64276Sksagiyam PetscCall(PetscOptionsBool("-sanitize_submesh", "Sanitize submesh", "ex80.c", options->sanitizeSubmesh, &options->sanitizeSubmesh, NULL)); 6890a64276Sksagiyam PetscOptionsEnd(); 6990a64276Sksagiyam PetscFunctionReturn(PETSC_SUCCESS); 7090a64276Sksagiyam } 7190a64276Sksagiyam int main(int argc, char **argv) 7290a64276Sksagiyam { 7390a64276Sksagiyam DM dm, subdm; 7490a64276Sksagiyam PetscSF ownershipTransferSF; 7590a64276Sksagiyam DMLabel filter; 7690a64276Sksagiyam const PetscInt filterValue = 1; 7790a64276Sksagiyam MPI_Comm comm; 7890a64276Sksagiyam PetscMPIInt size, rank; 7990a64276Sksagiyam AppCtx user; 8090a64276Sksagiyam 8190a64276Sksagiyam PetscFunctionBeginUser; 8290a64276Sksagiyam PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 8390a64276Sksagiyam PetscCall(ProcessOptions(&user)); 8490a64276Sksagiyam comm = PETSC_COMM_WORLD; 8590a64276Sksagiyam PetscCallMPI(MPI_Comm_size(comm, &size)); 8690a64276Sksagiyam if (size != 3) { 8790a64276Sksagiyam PetscCall(PetscPrintf(comm, "This example is specifically designed for comm size == 3.\n")); 8890a64276Sksagiyam PetscCall(PetscFinalize()); 8990a64276Sksagiyam return 0; 9090a64276Sksagiyam } 9190a64276Sksagiyam PetscCallMPI(MPI_Comm_rank(comm, &rank)); 9290a64276Sksagiyam { 9390a64276Sksagiyam DM pdm; 9490a64276Sksagiyam const PetscInt faces[2] = {2, 2}; 9590a64276Sksagiyam PetscInt overlap = 1; 9690a64276Sksagiyam 97*42108689Sksagiyam PetscCall(DMPlexCreateBoxMesh(comm, 2, PETSC_FALSE, faces, NULL, NULL, NULL, PETSC_TRUE, 0, PETSC_TRUE, &dm)); 9890a64276Sksagiyam { 9990a64276Sksagiyam PetscPartitioner part; 10090a64276Sksagiyam PetscInt *sizes = NULL; 10190a64276Sksagiyam PetscInt *points = NULL; 10290a64276Sksagiyam 10390a64276Sksagiyam if (rank == 0) { 10490a64276Sksagiyam PetscInt sizes1[3] = {1, 2, 1}; 10590a64276Sksagiyam PetscInt points1[4] = {0, 2, 3, 1}; 10690a64276Sksagiyam 10790a64276Sksagiyam PetscCall(PetscMalloc2(3, &sizes, 4, &points)); 10890a64276Sksagiyam PetscCall(PetscArraycpy(sizes, sizes1, 3)); 10990a64276Sksagiyam PetscCall(PetscArraycpy(points, points1, 4)); 11090a64276Sksagiyam } 11190a64276Sksagiyam PetscCall(DMPlexGetPartitioner(dm, &part)); 11290a64276Sksagiyam PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL)); 11390a64276Sksagiyam PetscCall(PetscPartitionerShellSetPartition(part, size, sizes, points)); 11490a64276Sksagiyam PetscCall(PetscFree2(sizes, points)); 11590a64276Sksagiyam } 11690a64276Sksagiyam PetscCall(DMSetAdjacency(dm, -1, PETSC_FALSE, PETSC_TRUE)); 11790a64276Sksagiyam PetscCall(DMPlexDistribute(dm, overlap, NULL, &pdm)); 11890a64276Sksagiyam if (pdm) { 11990a64276Sksagiyam PetscCall(DMDestroy(&dm)); 12090a64276Sksagiyam dm = pdm; 12190a64276Sksagiyam } 12290a64276Sksagiyam } 12390a64276Sksagiyam PetscCall(DMLabelCreate(PETSC_COMM_SELF, "filter", &filter)); 12490a64276Sksagiyam switch (rank) { 12590a64276Sksagiyam case 0: 12690a64276Sksagiyam PetscCall(DMLabelSetValue(filter, 0, filterValue)); 12790a64276Sksagiyam PetscCall(DMLabelSetValue(filter, 1, filterValue)); 12890a64276Sksagiyam break; 12990a64276Sksagiyam case 1: 13090a64276Sksagiyam PetscCall(DMLabelSetValue(filter, 0, filterValue)); 13190a64276Sksagiyam PetscCall(DMLabelSetValue(filter, 2, filterValue)); 13290a64276Sksagiyam break; 13390a64276Sksagiyam case 2: 13490a64276Sksagiyam break; 13590a64276Sksagiyam } 13690a64276Sksagiyam PetscCall(PetscObjectSetName((PetscObject)dm, "Example_DM")); 13790a64276Sksagiyam PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 13890a64276Sksagiyam PetscCall(DMPlexFilter(dm, filter, filterValue, user.ignoreLabelHalo, user.sanitizeSubmesh, &ownershipTransferSF, &subdm)); 13990a64276Sksagiyam PetscCall(DMLabelDestroy(&filter)); 14090a64276Sksagiyam PetscCall(DMDestroy(&dm)); 14190a64276Sksagiyam PetscCall(PetscObjectSetName((PetscObject)subdm, "Example_SubDM")); 14290a64276Sksagiyam PetscCall(DMViewFromOptions(subdm, NULL, "-dm_view")); 14390a64276Sksagiyam PetscCall(DMDestroy(&subdm)); 14490a64276Sksagiyam PetscCall(PetscObjectSetName((PetscObject)ownershipTransferSF, "Example_Ownership_Transfer_SF")); 14590a64276Sksagiyam PetscCall(PetscSFView(ownershipTransferSF, PETSC_VIEWER_STDOUT_WORLD)); 14690a64276Sksagiyam PetscCall(PetscSFDestroy(&ownershipTransferSF)); 14790a64276Sksagiyam PetscCall(PetscFinalize()); 14890a64276Sksagiyam return 0; 14990a64276Sksagiyam } 15090a64276Sksagiyam 15190a64276Sksagiyam /*TEST 15290a64276Sksagiyam 15390a64276Sksagiyam testset: 15490a64276Sksagiyam nsize: 3 15590a64276Sksagiyam args: -dm_view ascii::ascii_info_detail 15690a64276Sksagiyam 15790a64276Sksagiyam test: 15890a64276Sksagiyam suffix: 0 15990a64276Sksagiyam args: 16090a64276Sksagiyam 16190a64276Sksagiyam test: 16290a64276Sksagiyam suffix: 1 16390a64276Sksagiyam args: -sanitize_submesh 16490a64276Sksagiyam 16590a64276Sksagiyam test: 16690a64276Sksagiyam suffix: 2 16790a64276Sksagiyam args: -ignore_label_halo 16890a64276Sksagiyam 16990a64276Sksagiyam TEST*/ 170