1*df2db7a0SMatthew G. Knepley static char help[] = "Tests save/load of plex/section/vec in HDF5.\n\n"; 2*df2db7a0SMatthew G. Knepley 3*df2db7a0SMatthew G. Knepley #include <petscdmshell.h> 4*df2db7a0SMatthew G. Knepley #include <petscdmplex.h> 5*df2db7a0SMatthew G. Knepley #include <petscsection.h> 6*df2db7a0SMatthew G. Knepley #include <petscsf.h> 7*df2db7a0SMatthew G. Knepley #include <petsclayouthdf5.h> 8*df2db7a0SMatthew G. Knepley 9*df2db7a0SMatthew G. Knepley /* A six-element mesh 10*df2db7a0SMatthew G. Knepley 11*df2db7a0SMatthew G. Knepley ===================== 12*df2db7a0SMatthew G. Knepley Save on 2 processes 13*df2db7a0SMatthew G. Knepley ===================== 14*df2db7a0SMatthew G. Knepley 15*df2db7a0SMatthew G. Knepley exampleDMPlex: Local numbering: 16*df2db7a0SMatthew G. Knepley 17*df2db7a0SMatthew G. Knepley 7---17--8---18--9--19--(12)(24)(13) 18*df2db7a0SMatthew G. Knepley | | | | | 19*df2db7a0SMatthew G. Knepley rank 0: 20 0 21 1 22 2 (25) (3)(26) 20*df2db7a0SMatthew G. Knepley | | | | | 21*df2db7a0SMatthew G. Knepley 4---14--5---15--6--16--(10)(23)(11) 22*df2db7a0SMatthew G. Knepley 23*df2db7a0SMatthew G. Knepley (13)(25)--8--17---9--18--10--19--11 24*df2db7a0SMatthew G. Knepley | | | | | 25*df2db7a0SMatthew G. Knepley rank 1: (26) (3) 20 0 21 1 22 2 23 26*df2db7a0SMatthew G. Knepley | | | | | 27*df2db7a0SMatthew G. Knepley (12)(24)--4--14---5--15---6--16---7 28*df2db7a0SMatthew G. Knepley 29*df2db7a0SMatthew G. Knepley exampleDMPlex: globalPointNumbering: 30*df2db7a0SMatthew G. Knepley 31*df2db7a0SMatthew G. Knepley 9--23--10--24--11--25--16--32--17--33--18--34--19 32*df2db7a0SMatthew G. Knepley | | | | | | | 33*df2db7a0SMatthew G. Knepley 26 0 27 1 28 2 35 3 36 4 37 5 38 34*df2db7a0SMatthew G. Knepley | | | | | | | 35*df2db7a0SMatthew G. Knepley 6--20---7--21---8--22--12--29--13--30--14--31--15 36*df2db7a0SMatthew G. Knepley 37*df2db7a0SMatthew G. Knepley exampleSectionDM: 38*df2db7a0SMatthew G. Knepley - includesConstraints = TRUE for local section (default) 39*df2db7a0SMatthew G. Knepley - includesConstraints = FALSE for global section (default) 40*df2db7a0SMatthew G. Knepley 41*df2db7a0SMatthew G. Knepley exampleSectionDM: Dofs (Field 0): 42*df2db7a0SMatthew G. Knepley 43*df2db7a0SMatthew G. Knepley 0---0---0---0---0---0---2---0---0---0---0---0---0 44*df2db7a0SMatthew G. Knepley | | | | | | | 45*df2db7a0SMatthew G. Knepley 0 0 0 0 0 0 0 2 0 0 0 0 0 46*df2db7a0SMatthew G. Knepley | | | | | | | 47*df2db7a0SMatthew G. Knepley 0---0---0---0---0---0---0---0---0---0---0---0---0 48*df2db7a0SMatthew G. Knepley 49*df2db7a0SMatthew G. Knepley exampleSectionDM: Dofs (Field 1): constrained 50*df2db7a0SMatthew G. Knepley / 51*df2db7a0SMatthew G. Knepley 0---0---0---0---0---0---1---0---0---0---0---0---0 52*df2db7a0SMatthew G. Knepley | | | | | | | 53*df2db7a0SMatthew G. Knepley 0 0 0 0 0 0 2 0 0 1 0 0 0 54*df2db7a0SMatthew G. Knepley | | | | | | | 55*df2db7a0SMatthew G. Knepley 0---0---0---0---0---0---0---0---0---0---0---0---0 56*df2db7a0SMatthew G. Knepley 57*df2db7a0SMatthew G. Knepley exampleSectionDM: Offsets (total) in global section: 58*df2db7a0SMatthew G. Knepley 59*df2db7a0SMatthew G. Knepley 0---0---0---0---0---0---3---5---5---5---5---5---5 60*df2db7a0SMatthew G. Knepley | | | | | | | 61*df2db7a0SMatthew G. Knepley 0 0 0 0 0 0 5 0 7 2 7 3 7 62*df2db7a0SMatthew G. Knepley | | | | | | | 63*df2db7a0SMatthew G. Knepley 0---0---0---0---0---0---3---5---3---5---3---5---3 64*df2db7a0SMatthew G. Knepley 65*df2db7a0SMatthew G. Knepley exampleVec: Values (Field 0): (1.3, 1.4) 66*df2db7a0SMatthew G. Knepley / 67*df2db7a0SMatthew G. Knepley +-------+-------+-------*-------+-------+-------+ 68*df2db7a0SMatthew G. Knepley | | | | | | | 69*df2db7a0SMatthew G. Knepley | | | | * (1.0, 1.1)| | 70*df2db7a0SMatthew G. Knepley | | | | | | | 71*df2db7a0SMatthew G. Knepley +-------+-------+-------+-------+-------+-------+ 72*df2db7a0SMatthew G. Knepley 73*df2db7a0SMatthew G. Knepley exampleVec: Values (Field 1): (1.5,) constrained 74*df2db7a0SMatthew G. Knepley / 75*df2db7a0SMatthew G. Knepley +-------+-------+-------*-------+-------+-------+ 76*df2db7a0SMatthew G. Knepley | | | | | | | 77*df2db7a0SMatthew G. Knepley | | (1.6, 1.7) * | * (1.2,) | 78*df2db7a0SMatthew G. Knepley | | | | | | | 79*df2db7a0SMatthew G. Knepley +-------+-------+-------+-------+-------+-------+ 80*df2db7a0SMatthew G. Knepley 81*df2db7a0SMatthew G. Knepley exampleVec: as global vector 82*df2db7a0SMatthew G. Knepley 83*df2db7a0SMatthew G. Knepley rank 0: [] 84*df2db7a0SMatthew G. Knepley rakn 1: [1.0, 1.1, 1.2, 1.3, 1.4, 1.6, 1.7] 85*df2db7a0SMatthew G. Knepley 86*df2db7a0SMatthew G. Knepley ===================== 87*df2db7a0SMatthew G. Knepley Load on 3 Processes 88*df2db7a0SMatthew G. Knepley ===================== 89*df2db7a0SMatthew G. Knepley 90*df2db7a0SMatthew G. Knepley exampleDMPlex: Loaded/Distributed: 91*df2db7a0SMatthew G. Knepley 92*df2db7a0SMatthew G. Knepley 5--13---6--14--(8)(18)(10) 93*df2db7a0SMatthew G. Knepley | | | | 94*df2db7a0SMatthew G. Knepley rank 0: 15 0 16 1 (19)(2)(20) 95*df2db7a0SMatthew G. Knepley | | | | 96*df2db7a0SMatthew G. Knepley 3--11---4--12--(7)(17)-(9) 97*df2db7a0SMatthew G. Knepley 98*df2db7a0SMatthew G. Knepley (9)(21)--5--15---7--18-(12)(24)(13) 99*df2db7a0SMatthew G. Knepley | | | | | 100*df2db7a0SMatthew G. Knepley rank 1: (22) (2) 16 0 19 1 (25) (3)(26) 101*df2db7a0SMatthew G. Knepley | | | | | 102*df2db7a0SMatthew G. Knepley (8)(20)--4--14---6--17-(10)(23)(11) 103*df2db7a0SMatthew G. Knepley 104*df2db7a0SMatthew G. Knepley +-> (10)(19)--6--13---7--14---8 105*df2db7a0SMatthew G. Knepley permute | | | | | 106*df2db7a0SMatthew G. Knepley rank 2: +-> (20) (2) 15 0 16 1 17 107*df2db7a0SMatthew G. Knepley | | | | 108*df2db7a0SMatthew G. Knepley (9)(18)--3--11---4--12---5 109*df2db7a0SMatthew G. Knepley 110*df2db7a0SMatthew G. Knepley exampleSectionDM: 111*df2db7a0SMatthew G. Knepley - includesConstraints = TRUE for local section (default) 112*df2db7a0SMatthew G. Knepley - includesConstraints = FALSE for global section (default) 113*df2db7a0SMatthew G. Knepley 114*df2db7a0SMatthew G. Knepley exampleVec: as local vector: 115*df2db7a0SMatthew G. Knepley 116*df2db7a0SMatthew G. Knepley rank 0: [1.3, 1.4, 1.5, 1.6, 1.7] 117*df2db7a0SMatthew G. Knepley rank 1: [1.0, 1.1, 1.2, 1.3, 1.4, 1.5, 1.6, 1.7] 118*df2db7a0SMatthew G. Knepley rank 2: [1.2, 1.0, 1.1, 1.6, 1.7, 1.3, 1.4, 1.5] 119*df2db7a0SMatthew G. Knepley 120*df2db7a0SMatthew G. Knepley exampleVec: as global vector: 121*df2db7a0SMatthew G. Knepley 122*df2db7a0SMatthew G. Knepley rank 0: [] 123*df2db7a0SMatthew G. Knepley rank 1: [1.0, 1.1, 1.3, 1.4, 1.6, 1.7] 124*df2db7a0SMatthew G. Knepley rank 2: [1.2] 125*df2db7a0SMatthew G. Knepley 126*df2db7a0SMatthew G. Knepley */ 127*df2db7a0SMatthew G. Knepley 128*df2db7a0SMatthew G. Knepley typedef struct { 129*df2db7a0SMatthew G. Knepley char fname[PETSC_MAX_PATH_LEN]; /* Output mesh filename */ 130*df2db7a0SMatthew G. Knepley PetscBool shell; /* Use DMShell to wrap sections */ 131*df2db7a0SMatthew G. Knepley } AppCtx; 132*df2db7a0SMatthew G. Knepley 133*df2db7a0SMatthew G. Knepley PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 134*df2db7a0SMatthew G. Knepley { 135*df2db7a0SMatthew G. Knepley PetscBool flg; 136*df2db7a0SMatthew G. Knepley PetscErrorCode ierr; 137*df2db7a0SMatthew G. Knepley 138*df2db7a0SMatthew G. Knepley PetscFunctionBegin; 139*df2db7a0SMatthew G. Knepley options->fname[0] = '\0'; 140*df2db7a0SMatthew G. Knepley ierr = PetscOptionsBegin(comm, "", "DMPlex View/Load Test Options", "DMPLEX");PetscCall(ierr); 141*df2db7a0SMatthew G. Knepley PetscCall(PetscOptionsString("-fname", "The output mesh file", "ex12.c", options->fname, options->fname, sizeof(options->fname), &flg)); 142*df2db7a0SMatthew G. Knepley PetscCall(PetscOptionsBool("-shell", "Use DMShell to wrap sections", "ex12.c", options->shell, &options->shell, NULL)); 143*df2db7a0SMatthew G. Knepley ierr = PetscOptionsEnd();PetscCall(ierr); 144*df2db7a0SMatthew G. Knepley PetscFunctionReturn(0); 145*df2db7a0SMatthew G. Knepley } 146*df2db7a0SMatthew G. Knepley 147*df2db7a0SMatthew G. Knepley int main(int argc, char **argv) 148*df2db7a0SMatthew G. Knepley { 149*df2db7a0SMatthew G. Knepley MPI_Comm comm; 150*df2db7a0SMatthew G. Knepley PetscMPIInt size, rank, mycolor; 151*df2db7a0SMatthew G. Knepley const char exampleDMPlexName[] = "exampleDMPlex"; 152*df2db7a0SMatthew G. Knepley const char exampleSectionDMName[] = "exampleSectionDM"; 153*df2db7a0SMatthew G. Knepley const char exampleVecName[] = "exampleVec"; 154*df2db7a0SMatthew G. Knepley PetscScalar constraintValue = 1.5; 155*df2db7a0SMatthew G. Knepley PetscViewerFormat format = PETSC_VIEWER_HDF5_PETSC; 156*df2db7a0SMatthew G. Knepley AppCtx user; 157*df2db7a0SMatthew G. Knepley 158*df2db7a0SMatthew G. Knepley PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 159*df2db7a0SMatthew G. Knepley PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 160*df2db7a0SMatthew G. Knepley PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 161*df2db7a0SMatthew G. Knepley PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 162*df2db7a0SMatthew G. Knepley PetscCheckFalse(size < 3,PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Example only works with three or more processes"); 163*df2db7a0SMatthew G. Knepley 164*df2db7a0SMatthew G. Knepley /* Save */ 165*df2db7a0SMatthew G. Knepley mycolor = (PetscMPIInt)(rank >= 2); 166*df2db7a0SMatthew G. Knepley PetscCallMPI(MPI_Comm_split(PETSC_COMM_WORLD, mycolor, rank, &comm)); 167*df2db7a0SMatthew G. Knepley if (mycolor == 0) { 168*df2db7a0SMatthew G. Knepley DM dm; 169*df2db7a0SMatthew G. Knepley PetscViewer viewer; 170*df2db7a0SMatthew G. Knepley 171*df2db7a0SMatthew G. Knepley PetscCall(PetscViewerHDF5Open(comm, user.fname, FILE_MODE_WRITE, &viewer)); 172*df2db7a0SMatthew G. Knepley /* Save exampleDMPlex */ 173*df2db7a0SMatthew G. Knepley { 174*df2db7a0SMatthew G. Knepley DM pdm; 175*df2db7a0SMatthew G. Knepley const PetscInt faces[2] = {6, 1}; 176*df2db7a0SMatthew G. Knepley PetscSF sf; 177*df2db7a0SMatthew G. Knepley PetscInt overlap = 1; 178*df2db7a0SMatthew G. Knepley 179*df2db7a0SMatthew G. Knepley PetscCall(DMPlexCreateBoxMesh(comm, 2, PETSC_FALSE, faces, NULL, NULL, NULL, PETSC_TRUE, &dm)); 180*df2db7a0SMatthew G. Knepley PetscCall(DMPlexDistribute(dm, overlap, &sf, &pdm)); 181*df2db7a0SMatthew G. Knepley if (pdm) { 182*df2db7a0SMatthew G. Knepley PetscCall(DMDestroy(&dm)); 183*df2db7a0SMatthew G. Knepley dm = pdm; 184*df2db7a0SMatthew G. Knepley } 185*df2db7a0SMatthew G. Knepley PetscCall(PetscSFDestroy(&sf)); 186*df2db7a0SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)dm, exampleDMPlexName)); 187*df2db7a0SMatthew G. Knepley PetscCall(PetscViewerPushFormat(viewer, format)); 188*df2db7a0SMatthew G. Knepley PetscCall(DMPlexTopologyView(dm, viewer)); 189*df2db7a0SMatthew G. Knepley PetscCall(DMPlexLabelsView(dm, viewer)); 190*df2db7a0SMatthew G. Knepley PetscCall(PetscViewerPopFormat(viewer)); 191*df2db7a0SMatthew G. Knepley } 192*df2db7a0SMatthew G. Knepley /* Save coordinates */ 193*df2db7a0SMatthew G. Knepley PetscCall(PetscViewerPushFormat(viewer, format)); 194*df2db7a0SMatthew G. Knepley PetscCall(DMPlexCoordinatesView(dm, viewer)); 195*df2db7a0SMatthew G. Knepley PetscCall(PetscViewerPopFormat(viewer)); 196*df2db7a0SMatthew G. Knepley /* Save exampleVec */ 197*df2db7a0SMatthew G. Knepley { 198*df2db7a0SMatthew G. Knepley PetscInt pStart = -1, pEnd = -1; 199*df2db7a0SMatthew G. Knepley DM sdm; 200*df2db7a0SMatthew G. Knepley PetscSection section, gsection; 201*df2db7a0SMatthew G. Knepley PetscBool includesConstraints = PETSC_FALSE; 202*df2db7a0SMatthew G. Knepley Vec vec; 203*df2db7a0SMatthew G. Knepley PetscScalar *array = NULL; 204*df2db7a0SMatthew G. Knepley 205*df2db7a0SMatthew G. Knepley /* Create section */ 206*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionCreate(comm, §ion)); 207*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetNumFields(section, 2)); 208*df2db7a0SMatthew G. Knepley PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 209*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetChart(section, pStart, pEnd)); 210*df2db7a0SMatthew G. Knepley switch (rank) { 211*df2db7a0SMatthew G. Knepley case 0: 212*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetDof(section, 3, 2)); 213*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetDof(section, 12, 3)); 214*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetDof(section, 25, 2)); 215*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetConstraintDof(section, 12, 1)); 216*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetFieldDof(section, 3, 0, 2)); 217*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetFieldDof(section, 12, 0, 2)); 218*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetFieldDof(section, 12, 1, 1)); 219*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetFieldDof(section, 25, 1, 2)); 220*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetFieldConstraintDof(section, 12, 1, 1)); 221*df2db7a0SMatthew G. Knepley break; 222*df2db7a0SMatthew G. Knepley case 1: 223*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetDof(section, 0, 2)); 224*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetDof(section, 1, 1)); 225*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetDof(section, 8, 3)); 226*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetDof(section, 20, 2)); 227*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetConstraintDof(section, 8, 1)); 228*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetFieldDof(section, 0, 0, 2)); 229*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetFieldDof(section, 8, 0, 2)); 230*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetFieldDof(section, 1, 1, 1)); 231*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetFieldDof(section, 8, 1, 1)); 232*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetFieldDof(section, 20, 1, 2)); 233*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetFieldConstraintDof(section, 8, 1, 1)); 234*df2db7a0SMatthew G. Knepley break; 235*df2db7a0SMatthew G. Knepley } 236*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetUp(section)); 237*df2db7a0SMatthew G. Knepley { 238*df2db7a0SMatthew G. Knepley const PetscInt indices[] = {2}; 239*df2db7a0SMatthew G. Knepley const PetscInt indices1[] = {0}; 240*df2db7a0SMatthew G. Knepley 241*df2db7a0SMatthew G. Knepley switch (rank) { 242*df2db7a0SMatthew G. Knepley case 0: 243*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetConstraintIndices(section, 12, indices)); 244*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetFieldConstraintIndices(section, 12, 1, indices1)); 245*df2db7a0SMatthew G. Knepley break; 246*df2db7a0SMatthew G. Knepley case 1: 247*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetConstraintIndices(section, 8, indices)); 248*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetFieldConstraintIndices(section, 8, 1, indices1)); 249*df2db7a0SMatthew G. Knepley break; 250*df2db7a0SMatthew G. Knepley } 251*df2db7a0SMatthew G. Knepley } 252*df2db7a0SMatthew G. Knepley if (user.shell) { 253*df2db7a0SMatthew G. Knepley PetscSF sf; 254*df2db7a0SMatthew G. Knepley 255*df2db7a0SMatthew G. Knepley PetscCall(DMShellCreate(comm, &sdm)); 256*df2db7a0SMatthew G. Knepley PetscCall(DMGetPointSF(dm, &sf)); 257*df2db7a0SMatthew G. Knepley PetscCall(DMSetPointSF(sdm, sf)); 258*df2db7a0SMatthew G. Knepley } 259*df2db7a0SMatthew G. Knepley else { 260*df2db7a0SMatthew G. Knepley PetscCall(DMClone(dm, &sdm)); 261*df2db7a0SMatthew G. Knepley } 262*df2db7a0SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)sdm, exampleSectionDMName)); 263*df2db7a0SMatthew G. Knepley PetscCall(DMSetLocalSection(sdm, section)); 264*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionDestroy(§ion)); 265*df2db7a0SMatthew G. Knepley PetscCall(DMPlexSectionView(dm, viewer, sdm)); 266*df2db7a0SMatthew G. Knepley /* Create global vector */ 267*df2db7a0SMatthew G. Knepley PetscCall(DMGetGlobalSection(sdm, &gsection)); 268*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionGetIncludesConstraints(gsection, &includesConstraints)); 269*df2db7a0SMatthew G. Knepley if (user.shell) { 270*df2db7a0SMatthew G. Knepley PetscInt n = -1; 271*df2db7a0SMatthew G. Knepley 272*df2db7a0SMatthew G. Knepley PetscCall(VecCreate(comm, &vec)); 273*df2db7a0SMatthew G. Knepley if (includesConstraints) PetscCall(PetscSectionGetStorageSize(gsection, &n)); 274*df2db7a0SMatthew G. Knepley else PetscCall(PetscSectionGetConstrainedStorageSize(gsection, &n)); 275*df2db7a0SMatthew G. Knepley PetscCall(VecSetSizes(vec, n, PETSC_DECIDE)); 276*df2db7a0SMatthew G. Knepley PetscCall(VecSetUp(vec)); 277*df2db7a0SMatthew G. Knepley } else { 278*df2db7a0SMatthew G. Knepley PetscCall(DMGetGlobalVector(sdm, &vec)); 279*df2db7a0SMatthew G. Knepley } 280*df2db7a0SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)vec, exampleVecName)); 281*df2db7a0SMatthew G. Knepley PetscCall(VecGetArrayWrite(vec, &array)); 282*df2db7a0SMatthew G. Knepley if (includesConstraints) { 283*df2db7a0SMatthew G. Knepley switch (rank) { 284*df2db7a0SMatthew G. Knepley case 0: 285*df2db7a0SMatthew G. Knepley break; 286*df2db7a0SMatthew G. Knepley case 1: 287*df2db7a0SMatthew G. Knepley array[0] = 1.0; 288*df2db7a0SMatthew G. Knepley array[1] = 1.1; 289*df2db7a0SMatthew G. Knepley array[2] = 1.2; 290*df2db7a0SMatthew G. Knepley array[3] = 1.3; 291*df2db7a0SMatthew G. Knepley array[4] = 1.4; 292*df2db7a0SMatthew G. Knepley array[5] = 1.5; 293*df2db7a0SMatthew G. Knepley array[6] = 1.6; 294*df2db7a0SMatthew G. Knepley array[7] = 1.7; 295*df2db7a0SMatthew G. Knepley break; 296*df2db7a0SMatthew G. Knepley } 297*df2db7a0SMatthew G. Knepley } else { 298*df2db7a0SMatthew G. Knepley switch (rank) { 299*df2db7a0SMatthew G. Knepley case 0: 300*df2db7a0SMatthew G. Knepley break; 301*df2db7a0SMatthew G. Knepley case 1: 302*df2db7a0SMatthew G. Knepley array[0] = 1.0; 303*df2db7a0SMatthew G. Knepley array[1] = 1.1; 304*df2db7a0SMatthew G. Knepley array[2] = 1.2; 305*df2db7a0SMatthew G. Knepley array[3] = 1.3; 306*df2db7a0SMatthew G. Knepley array[4] = 1.4; 307*df2db7a0SMatthew G. Knepley array[5] = 1.6; 308*df2db7a0SMatthew G. Knepley array[6] = 1.7; 309*df2db7a0SMatthew G. Knepley break; 310*df2db7a0SMatthew G. Knepley } 311*df2db7a0SMatthew G. Knepley } 312*df2db7a0SMatthew G. Knepley PetscCall(VecRestoreArrayWrite(vec, &array)); 313*df2db7a0SMatthew G. Knepley PetscCall(DMPlexGlobalVectorView(dm, viewer, sdm, vec)); 314*df2db7a0SMatthew G. Knepley if (user.shell) { 315*df2db7a0SMatthew G. Knepley PetscCall(VecDestroy(&vec)); 316*df2db7a0SMatthew G. Knepley } else { 317*df2db7a0SMatthew G. Knepley PetscCall(DMRestoreGlobalVector(sdm, &vec)); 318*df2db7a0SMatthew G. Knepley } 319*df2db7a0SMatthew G. Knepley PetscCall(DMDestroy(&sdm)); 320*df2db7a0SMatthew G. Knepley } 321*df2db7a0SMatthew G. Knepley PetscCall(PetscViewerDestroy(&viewer)); 322*df2db7a0SMatthew G. Knepley PetscCall(DMDestroy(&dm)); 323*df2db7a0SMatthew G. Knepley } 324*df2db7a0SMatthew G. Knepley PetscCallMPI(MPI_Comm_free(&comm)); 325*df2db7a0SMatthew G. Knepley /* Load */ 326*df2db7a0SMatthew G. Knepley mycolor = (PetscMPIInt)(rank >= 3); 327*df2db7a0SMatthew G. Knepley PetscCallMPI(MPI_Comm_split(PETSC_COMM_WORLD, mycolor, rank, &comm)); 328*df2db7a0SMatthew G. Knepley if (mycolor == 0) { 329*df2db7a0SMatthew G. Knepley DM dm; 330*df2db7a0SMatthew G. Knepley PetscSF sfXC; 331*df2db7a0SMatthew G. Knepley PetscViewer viewer; 332*df2db7a0SMatthew G. Knepley 333*df2db7a0SMatthew G. Knepley PetscCall(PetscViewerHDF5Open(comm, user.fname, FILE_MODE_READ, &viewer)); 334*df2db7a0SMatthew G. Knepley /* Load exampleDMPlex */ 335*df2db7a0SMatthew G. Knepley { 336*df2db7a0SMatthew G. Knepley PetscSF sfXB, sfBC; 337*df2db7a0SMatthew G. Knepley 338*df2db7a0SMatthew G. Knepley PetscCall(DMCreate(comm, &dm)); 339*df2db7a0SMatthew G. Knepley PetscCall(DMSetType(dm, DMPLEX)); 340*df2db7a0SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)dm, exampleDMPlexName)); 341*df2db7a0SMatthew G. Knepley /* sfXB: X -> B */ 342*df2db7a0SMatthew G. Knepley /* X: set of globalPointNumbers, [0, N) */ 343*df2db7a0SMatthew G. Knepley /* B: loaded naive in-memory plex */ 344*df2db7a0SMatthew G. Knepley PetscCall(PetscViewerPushFormat(viewer, format)); 345*df2db7a0SMatthew G. Knepley PetscCall(DMPlexTopologyLoad(dm, viewer, &sfXB)); 346*df2db7a0SMatthew G. Knepley PetscCall(PetscViewerPopFormat(viewer)); 347*df2db7a0SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)dm, exampleDMPlexName)); 348*df2db7a0SMatthew G. Knepley { 349*df2db7a0SMatthew G. Knepley DM distributedDM; 350*df2db7a0SMatthew G. Knepley PetscInt overlap = 1; 351*df2db7a0SMatthew G. Knepley PetscPartitioner part; 352*df2db7a0SMatthew G. Knepley 353*df2db7a0SMatthew G. Knepley PetscCall(DMPlexGetPartitioner(dm, &part)); 354*df2db7a0SMatthew G. Knepley PetscCall(PetscPartitionerSetFromOptions(part)); 355*df2db7a0SMatthew G. Knepley /* sfBC: B -> C */ 356*df2db7a0SMatthew G. Knepley /* B: loaded naive in-memory plex */ 357*df2db7a0SMatthew G. Knepley /* C: redistributed good in-memory */ 358*df2db7a0SMatthew G. Knepley PetscCall(DMPlexDistribute(dm, overlap, &sfBC, &distributedDM)); 359*df2db7a0SMatthew G. Knepley if (distributedDM) { 360*df2db7a0SMatthew G. Knepley PetscCall(DMDestroy(&dm)); 361*df2db7a0SMatthew G. Knepley dm = distributedDM; 362*df2db7a0SMatthew G. Knepley } 363*df2db7a0SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)dm, exampleDMPlexName)); 364*df2db7a0SMatthew G. Knepley } 365*df2db7a0SMatthew G. Knepley /* sfXC: X -> C */ 366*df2db7a0SMatthew G. Knepley PetscCall(PetscSFCompose(sfXB, sfBC, &sfXC)); 367*df2db7a0SMatthew G. Knepley PetscCall(PetscSFDestroy(&sfXB)); 368*df2db7a0SMatthew G. Knepley PetscCall(PetscSFDestroy(&sfBC)); 369*df2db7a0SMatthew G. Knepley } 370*df2db7a0SMatthew G. Knepley /* Load labels */ 371*df2db7a0SMatthew G. Knepley PetscCall(PetscViewerPushFormat(viewer, format)); 372*df2db7a0SMatthew G. Knepley PetscCall(DMPlexLabelsLoad(dm, viewer, sfXC)); 373*df2db7a0SMatthew G. Knepley PetscCall(PetscViewerPopFormat(viewer)); 374*df2db7a0SMatthew G. Knepley /* Load coordinates */ 375*df2db7a0SMatthew G. Knepley PetscCall(PetscViewerPushFormat(viewer, format)); 376*df2db7a0SMatthew G. Knepley PetscCall(DMPlexCoordinatesLoad(dm, viewer, sfXC)); 377*df2db7a0SMatthew G. Knepley PetscCall(PetscViewerPopFormat(viewer)); 378*df2db7a0SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)dm, "Load: DM (with coordinates)")); 379*df2db7a0SMatthew G. Knepley PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 380*df2db7a0SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)dm, exampleDMPlexName)); 381*df2db7a0SMatthew G. Knepley /* Load exampleVec */ 382*df2db7a0SMatthew G. Knepley { 383*df2db7a0SMatthew G. Knepley DM sdm; 384*df2db7a0SMatthew G. Knepley PetscSection section, gsection; 385*df2db7a0SMatthew G. Knepley IS perm; 386*df2db7a0SMatthew G. Knepley PetscBool includesConstraints = PETSC_FALSE; 387*df2db7a0SMatthew G. Knepley Vec vec; 388*df2db7a0SMatthew G. Knepley PetscSF lsf, gsf; 389*df2db7a0SMatthew G. Knepley 390*df2db7a0SMatthew G. Knepley if (user.shell) { 391*df2db7a0SMatthew G. Knepley PetscSF sf; 392*df2db7a0SMatthew G. Knepley 393*df2db7a0SMatthew G. Knepley PetscCall(DMShellCreate(comm, &sdm)); 394*df2db7a0SMatthew G. Knepley PetscCall(DMGetPointSF(dm, &sf)); 395*df2db7a0SMatthew G. Knepley PetscCall(DMSetPointSF(sdm, sf)); 396*df2db7a0SMatthew G. Knepley } else { 397*df2db7a0SMatthew G. Knepley PetscCall(DMClone(dm, &sdm)); 398*df2db7a0SMatthew G. Knepley } 399*df2db7a0SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)sdm, exampleSectionDMName)); 400*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionCreate(comm, §ion)); 401*df2db7a0SMatthew G. Knepley { 402*df2db7a0SMatthew G. Knepley PetscInt pStart = -1, pEnd = -1, p = -1; 403*df2db7a0SMatthew G. Knepley PetscInt *pinds = NULL; 404*df2db7a0SMatthew G. Knepley 405*df2db7a0SMatthew G. Knepley PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 406*df2db7a0SMatthew G. Knepley PetscCall(PetscMalloc1(pEnd - pStart, &pinds)); 407*df2db7a0SMatthew G. Knepley for (p = 0; p < pEnd - pStart; ++p) pinds[p] = p; 408*df2db7a0SMatthew G. Knepley if (rank == 2) {pinds[10] = 20; pinds[20] = 10;} 409*df2db7a0SMatthew G. Knepley PetscCall(ISCreateGeneral(comm, pEnd - pStart, pinds, PETSC_OWN_POINTER, &perm)); 410*df2db7a0SMatthew G. Knepley } 411*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionSetPermutation(section, perm)); 412*df2db7a0SMatthew G. Knepley PetscCall(ISDestroy(&perm)); 413*df2db7a0SMatthew G. Knepley PetscCall(DMSetLocalSection(sdm, section)); 414*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionDestroy(§ion)); 415*df2db7a0SMatthew G. Knepley PetscCall(DMPlexSectionLoad(dm, viewer, sdm, sfXC, &gsf, &lsf)); 416*df2db7a0SMatthew G. Knepley /* Load as local vector */ 417*df2db7a0SMatthew G. Knepley PetscCall(DMGetLocalSection(sdm, §ion)); 418*df2db7a0SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)section, "Load: local section")); 419*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionView(section, PETSC_VIEWER_STDOUT_(comm))); 420*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionGetIncludesConstraints(section, &includesConstraints)); 421*df2db7a0SMatthew G. Knepley if (user.shell) { 422*df2db7a0SMatthew G. Knepley PetscInt m = -1; 423*df2db7a0SMatthew G. Knepley 424*df2db7a0SMatthew G. Knepley PetscCall(VecCreate(comm, &vec)); 425*df2db7a0SMatthew G. Knepley if (includesConstraints) PetscCall(PetscSectionGetStorageSize(section, &m)); 426*df2db7a0SMatthew G. Knepley else PetscCall(PetscSectionGetConstrainedStorageSize(section, &m)); 427*df2db7a0SMatthew G. Knepley PetscCall(VecSetSizes(vec, m, PETSC_DECIDE)); 428*df2db7a0SMatthew G. Knepley PetscCall(VecSetUp(vec)); 429*df2db7a0SMatthew G. Knepley } else { 430*df2db7a0SMatthew G. Knepley PetscCall(DMGetLocalVector(sdm, &vec)); 431*df2db7a0SMatthew G. Knepley } 432*df2db7a0SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)vec, exampleVecName)); 433*df2db7a0SMatthew G. Knepley PetscCall(VecSet(vec, constraintValue)); 434*df2db7a0SMatthew G. Knepley PetscCall(DMPlexLocalVectorLoad(dm, viewer, sdm, lsf, vec)); 435*df2db7a0SMatthew G. Knepley PetscCall(PetscSFDestroy(&lsf)); 436*df2db7a0SMatthew G. Knepley if (user.shell) { 437*df2db7a0SMatthew G. Knepley PetscCall(VecView(vec, PETSC_VIEWER_STDOUT_(comm))); 438*df2db7a0SMatthew G. Knepley PetscCall(VecDestroy(&vec)); 439*df2db7a0SMatthew G. Knepley } else { 440*df2db7a0SMatthew G. Knepley PetscCall(DMRestoreLocalVector(sdm, &vec)); 441*df2db7a0SMatthew G. Knepley } 442*df2db7a0SMatthew G. Knepley /* Load as global vector */ 443*df2db7a0SMatthew G. Knepley PetscCall(DMGetGlobalSection(sdm, &gsection)); 444*df2db7a0SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)gsection, "Load: global section")); 445*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionView(gsection, PETSC_VIEWER_STDOUT_(comm))); 446*df2db7a0SMatthew G. Knepley PetscCall(PetscSectionGetIncludesConstraints(gsection, &includesConstraints)); 447*df2db7a0SMatthew G. Knepley if (user.shell) { 448*df2db7a0SMatthew G. Knepley PetscInt m = -1; 449*df2db7a0SMatthew G. Knepley 450*df2db7a0SMatthew G. Knepley PetscCall(VecCreate(comm, &vec)); 451*df2db7a0SMatthew G. Knepley if (includesConstraints) PetscCall(PetscSectionGetStorageSize(gsection, &m)); 452*df2db7a0SMatthew G. Knepley else PetscCall(PetscSectionGetConstrainedStorageSize(gsection, &m)); 453*df2db7a0SMatthew G. Knepley PetscCall(VecSetSizes(vec, m, PETSC_DECIDE)); 454*df2db7a0SMatthew G. Knepley PetscCall(VecSetUp(vec)); 455*df2db7a0SMatthew G. Knepley } else { 456*df2db7a0SMatthew G. Knepley PetscCall(DMGetGlobalVector(sdm, &vec)); 457*df2db7a0SMatthew G. Knepley } 458*df2db7a0SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)vec, exampleVecName)); 459*df2db7a0SMatthew G. Knepley PetscCall(DMPlexGlobalVectorLoad(dm, viewer, sdm, gsf, vec)); 460*df2db7a0SMatthew G. Knepley PetscCall(PetscSFDestroy(&gsf)); 461*df2db7a0SMatthew G. Knepley PetscCall(VecView(vec, PETSC_VIEWER_STDOUT_(comm))); 462*df2db7a0SMatthew G. Knepley if (user.shell) { 463*df2db7a0SMatthew G. Knepley PetscCall(VecDestroy(&vec)); 464*df2db7a0SMatthew G. Knepley } else { 465*df2db7a0SMatthew G. Knepley PetscCall(DMRestoreGlobalVector(sdm, &vec)); 466*df2db7a0SMatthew G. Knepley } 467*df2db7a0SMatthew G. Knepley PetscCall(DMDestroy(&sdm)); 468*df2db7a0SMatthew G. Knepley } 469*df2db7a0SMatthew G. Knepley PetscCall(PetscViewerDestroy(&viewer)); 470*df2db7a0SMatthew G. Knepley PetscCall(PetscSFDestroy(&sfXC)); 471*df2db7a0SMatthew G. Knepley PetscCall(DMDestroy(&dm)); 472*df2db7a0SMatthew G. Knepley } 473*df2db7a0SMatthew G. Knepley PetscCallMPI(MPI_Comm_free(&comm)); 474*df2db7a0SMatthew G. Knepley 475*df2db7a0SMatthew G. Knepley /* Finalize */ 476*df2db7a0SMatthew G. Knepley PetscCall(PetscFinalize()); 477*df2db7a0SMatthew G. Knepley return 0; 478*df2db7a0SMatthew G. Knepley } 479*df2db7a0SMatthew G. Knepley 480*df2db7a0SMatthew G. Knepley /*TEST 481*df2db7a0SMatthew G. Knepley 482*df2db7a0SMatthew G. Knepley build: 483*df2db7a0SMatthew G. Knepley requires: hdf5 484*df2db7a0SMatthew G. Knepley testset: 485*df2db7a0SMatthew G. Knepley suffix: 0 486*df2db7a0SMatthew G. Knepley requires: !complex 487*df2db7a0SMatthew G. Knepley nsize: 4 488*df2db7a0SMatthew G. Knepley args: -fname ex12_dump.h5 -shell {{True False}separate output} -dm_view ascii::ascii_info_detail 489*df2db7a0SMatthew G. Knepley args: -dm_plex_view_hdf5_storage_version 2.0.0 490*df2db7a0SMatthew G. Knepley test: 491*df2db7a0SMatthew G. Knepley suffix: parmetis 492*df2db7a0SMatthew G. Knepley requires: parmetis 493*df2db7a0SMatthew G. Knepley args: -petscpartitioner_type parmetis 494*df2db7a0SMatthew G. Knepley test: 495*df2db7a0SMatthew G. Knepley suffix: ptscotch 496*df2db7a0SMatthew G. Knepley requires: ptscotch 497*df2db7a0SMatthew G. Knepley args: -petscpartitioner_type ptscotch 498*df2db7a0SMatthew G. Knepley 499*df2db7a0SMatthew G. Knepley TEST*/ 500