13ceb54acSksagiyam static char help[] = "Tests PetscSectionView()/Load() with HDF5.\n\n"; 23ceb54acSksagiyam 33ceb54acSksagiyam #include <petscdmshell.h> 43ceb54acSksagiyam #include <petscdmplex.h> 53ceb54acSksagiyam #include <petscsection.h> 63ceb54acSksagiyam #include <petscsf.h> 73ceb54acSksagiyam #include <petsclayouthdf5.h> 83ceb54acSksagiyam 93ceb54acSksagiyam /* Save/Load abstract sections 103ceb54acSksagiyam 113ceb54acSksagiyam ===================== 123ceb54acSksagiyam Save on 2 processes 133ceb54acSksagiyam ===================== 143ceb54acSksagiyam 153ceb54acSksagiyam section: 163ceb54acSksagiyam 0 1 2 3 173ceb54acSksagiyam rank 0: Dof (Field 0) 2 3 5 7 183ceb54acSksagiyam Dof (Field 1) 1 0 0 0 193ceb54acSksagiyam 203ceb54acSksagiyam 0 1 2 213ceb54acSksagiyam rank 1: Dof (Field 0) 7 5 11 <- DoF 7 is constrained 223ceb54acSksagiyam Dof (Field 1) 0 0 2 233ceb54acSksagiyam 243ceb54acSksagiyam sf: 253ceb54acSksagiyam [0] 3 <- (1, 0) 263ceb54acSksagiyam [1] 1 <- (0, 2) 273ceb54acSksagiyam 283ceb54acSksagiyam global section (includesConstraints = PETSC_FALSE): 293ceb54acSksagiyam 0 1 2 3 303ceb54acSksagiyam rank 0: Dof (Field 0) 2 3 5 -8 313ceb54acSksagiyam Off (Field 0) 0 3 6 -12 323ceb54acSksagiyam Dof (Field 1) 1 0 0 -1 333ceb54acSksagiyam Off (Field 1) 2 6 11 -19 343ceb54acSksagiyam 353ceb54acSksagiyam 0 1 2 363ceb54acSksagiyam rank 1: Dof (Field 0) 7 -6 11 373ceb54acSksagiyam Off (Field 0) 11 -7 18 383ceb54acSksagiyam Dof (Field 1) 0 -1 2 393ceb54acSksagiyam Off (Field 1) 18 -12 28 403ceb54acSksagiyam 413ceb54acSksagiyam global section (includesConstraints = PETSC_TRUE): 423ceb54acSksagiyam 0 1 2 3 433ceb54acSksagiyam rank 0: Dof (Field 0) 2 3 5 -8 443ceb54acSksagiyam Off (Field 0) 0 3 6 -12 453ceb54acSksagiyam Dof (Field 1) 1 0 0 -1 463ceb54acSksagiyam Off (Field 1) 2 6 11 -19 473ceb54acSksagiyam 483ceb54acSksagiyam 0 1 2 493ceb54acSksagiyam rank 1: Dof (Field 0) 7 -6 11 503ceb54acSksagiyam Off (Field 0) 11 -7 18 513ceb54acSksagiyam Dof (Field 1) 0 -1 2 523ceb54acSksagiyam Off (Field 1) 18 -12 29 533ceb54acSksagiyam 543ceb54acSksagiyam ===================== 553ceb54acSksagiyam Load on 3 Processes 563ceb54acSksagiyam ===================== 573ceb54acSksagiyam 583ceb54acSksagiyam (Set chartSize = 4, 0, 1 for rank 0, 1, 2, respectively) 593ceb54acSksagiyam 603ceb54acSksagiyam global section (includesConstraints = PETSC_FALSE): 613ceb54acSksagiyam 623ceb54acSksagiyam rank 0: Dof (Field 0) 2 3 5 7 633ceb54acSksagiyam Off (Field 0) 0 3 6 11 643ceb54acSksagiyam Dof (Field 1) 1 0 0 0 653ceb54acSksagiyam Off (Field 1) 2 6 11 18 663ceb54acSksagiyam 673ceb54acSksagiyam rank 1: Dof (Field 0) 683ceb54acSksagiyam Dof (Field 1) 693ceb54acSksagiyam 703ceb54acSksagiyam rank 2: Dof (Field 0) 11 713ceb54acSksagiyam Off (Field 0) 18 723ceb54acSksagiyam Dof (Field 1) 2 733ceb54acSksagiyam Off (Field 1) 28 743ceb54acSksagiyam 753ceb54acSksagiyam global section (includesConstraints = PETSC_TRUE): 763ceb54acSksagiyam 773ceb54acSksagiyam rank 0: Dof (Field 0) 2 3 5 7 783ceb54acSksagiyam Off (Field 0) 0 3 6 11 793ceb54acSksagiyam Dof (Field 1) 1 0 0 0 803ceb54acSksagiyam Off (Field 1) 2 6 11 18 813ceb54acSksagiyam 823ceb54acSksagiyam rank 1: Dof (Field 0) 833ceb54acSksagiyam Dof (Field 1) 843ceb54acSksagiyam 853ceb54acSksagiyam rank 2: Dof (Field 0) 11 863ceb54acSksagiyam Off (Field 0) 18 873ceb54acSksagiyam Dof (Field 1) 2 883ceb54acSksagiyam Off (Field 1) 29 893ceb54acSksagiyam */ 903ceb54acSksagiyam 913ceb54acSksagiyam typedef struct { 923ceb54acSksagiyam char fname[PETSC_MAX_PATH_LEN]; /* Output mesh filename */ 933ceb54acSksagiyam PetscBool includes_constraints; /* Flag for if global section is to include constrained DoFs or not */ 943ceb54acSksagiyam } AppCtx; 953ceb54acSksagiyam 963ceb54acSksagiyam PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 973ceb54acSksagiyam { 983ceb54acSksagiyam PetscFunctionBegin; 993ceb54acSksagiyam options->fname[0] = '\0'; 1003ceb54acSksagiyam options->includes_constraints = PETSC_TRUE; 101d0609cedSBarry Smith PetscOptionsBegin(comm, "", "PetscSectionView()/Load() in HDF5 Test Options", "DMPLEX"); 1029566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-fname", "The output file", "ex5.c", options->fname, options->fname, sizeof(options->fname), NULL)); 1039566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-includes_constraints", "Flag for if global section is to include constrained DoFs or not", "ex5.c", options->includes_constraints, &options->includes_constraints, NULL)); 104d0609cedSBarry Smith PetscOptionsEnd(); 1053ceb54acSksagiyam PetscFunctionReturn(0); 1063ceb54acSksagiyam } 1073ceb54acSksagiyam 1083ceb54acSksagiyam int main(int argc, char **argv) 1093ceb54acSksagiyam { 1103ceb54acSksagiyam MPI_Comm comm; 1113ceb54acSksagiyam PetscMPIInt size, rank, mycolor; 1123ceb54acSksagiyam AppCtx user; 1133ceb54acSksagiyam 114*327415f7SBarry Smith PetscFunctionBeginUser; 1159566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1169566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 1179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 1189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 119be096a46SBarry Smith PetscCheck(size >= 3,PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Example only works with three or more processes"); 1203ceb54acSksagiyam 1213ceb54acSksagiyam /* Save */ 1223ceb54acSksagiyam mycolor = (PetscMPIInt)(rank >= 2); 1239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_split(PETSC_COMM_WORLD, mycolor, rank, &comm)); 1243ceb54acSksagiyam if (mycolor == 0) { 1253ceb54acSksagiyam PetscSection section, gsection; 1263ceb54acSksagiyam PetscSF sf; 1273ceb54acSksagiyam PetscInt nroots = -1, nleaves = -1, *ilocal; 1283ceb54acSksagiyam PetscSFNode *iremote; 1293ceb54acSksagiyam PetscViewer viewer; 1303ceb54acSksagiyam 1313ceb54acSksagiyam /* Create section */ 1329566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(comm, §ion)); 1339566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(section, 2)); 1343ceb54acSksagiyam switch (rank) { 1353ceb54acSksagiyam case 0: 1369566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(section, 0, 4)); 1379566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(section, 0, 3)); 1389566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(section, 1, 3)); 1399566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(section, 2, 5)); 1409566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(section, 3, 7)); 1419566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(section, 0, 0, 2)); 1429566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(section, 1, 0, 3)); 1439566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(section, 2, 0, 5)); 1449566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(section, 3, 0, 7)); 1459566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(section, 0, 1, 1)); 1463ceb54acSksagiyam break; 1473ceb54acSksagiyam case 1: 1489566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(section, 0, 3)); 1499566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(section, 0, 7)); 1509566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(section, 1, 5)); 1519566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(section, 2, 13)); 1529566063dSJacob Faibussowitsch PetscCall(PetscSectionSetConstraintDof(section, 2, 1)); 1539566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(section, 0, 0, 7)); 1549566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(section, 1, 0, 5)); 1559566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(section, 2, 0, 11)); 1569566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(section, 2, 1, 2)); 1579566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldConstraintDof(section, 2, 0, 1)); 1583ceb54acSksagiyam break; 1593ceb54acSksagiyam } 1609566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(section)); 1613ceb54acSksagiyam if (rank == 1) 1623ceb54acSksagiyam { 1633ceb54acSksagiyam const PetscInt indices[] = {7}; 1643ceb54acSksagiyam const PetscInt indices0[] = {7}; 1653ceb54acSksagiyam 1669566063dSJacob Faibussowitsch PetscCall(PetscSectionSetConstraintIndices(section, 2, indices)); 1679566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldConstraintIndices(section, 2, 0, indices0)); 1683ceb54acSksagiyam } 1693ceb54acSksagiyam /* Create sf */ 1703ceb54acSksagiyam switch (rank) { 1713ceb54acSksagiyam case 0: 1723ceb54acSksagiyam nroots = 4; 1733ceb54acSksagiyam nleaves = 1; 1749566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &ilocal)); 1759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &iremote)); 1763ceb54acSksagiyam ilocal[0] = 3; 1773ceb54acSksagiyam iremote[0].rank = 1; 1783ceb54acSksagiyam iremote[0].index = 0; 1793ceb54acSksagiyam break; 1803ceb54acSksagiyam case 1: 1813ceb54acSksagiyam nroots = 3; 1823ceb54acSksagiyam nleaves = 1; 1839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &ilocal)); 1849566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleaves, &iremote)); 1853ceb54acSksagiyam ilocal[0] = 1; 1863ceb54acSksagiyam iremote[0].rank = 0; 1873ceb54acSksagiyam iremote[0].index = 2; 1883ceb54acSksagiyam break; 1893ceb54acSksagiyam } 1909566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(comm, &sf)); 1919566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER)); 1923ceb54acSksagiyam /* Create global section*/ 1939566063dSJacob Faibussowitsch PetscCall(PetscSectionCreateGlobalSection(section, sf, user.includes_constraints, PETSC_FALSE, &gsection)); 1949566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf)); 1953ceb54acSksagiyam /* View */ 1969566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Open(comm, user.fname, FILE_MODE_WRITE, &viewer)); 1979566063dSJacob Faibussowitsch PetscCall(PetscSectionView(gsection, viewer)); 1989566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 1999566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)section, "Save: local section")); 2009566063dSJacob Faibussowitsch PetscCall(PetscSectionView(section, PETSC_VIEWER_STDOUT_(comm))); 2019566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)gsection, "Save: global section")); 2029566063dSJacob Faibussowitsch PetscCall(PetscSectionView(gsection, PETSC_VIEWER_STDOUT_(comm))); 2039566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&gsection)); 2049566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(§ion)); 2053ceb54acSksagiyam } 2069566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&comm)); 2073ceb54acSksagiyam 2083ceb54acSksagiyam /* Load */ 2093ceb54acSksagiyam mycolor = (PetscMPIInt)(rank >= 3); 2109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_split(PETSC_COMM_WORLD, mycolor, rank, &comm)); 2113ceb54acSksagiyam if (mycolor == 0) { 2123ceb54acSksagiyam PetscSection section; 2133ceb54acSksagiyam PetscInt chartSize = -1; 2143ceb54acSksagiyam PetscViewer viewer; 2153ceb54acSksagiyam 2169566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(comm, §ion)); 2173ceb54acSksagiyam switch (rank) { 2183ceb54acSksagiyam case 0: 2193ceb54acSksagiyam chartSize = 4; 2203ceb54acSksagiyam break; 2213ceb54acSksagiyam case 1: 2223ceb54acSksagiyam chartSize = 0; 2233ceb54acSksagiyam break; 2243ceb54acSksagiyam case 2: 2253ceb54acSksagiyam chartSize = 1; 2263ceb54acSksagiyam break; 2273ceb54acSksagiyam } 2289566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(section, 0, chartSize)); 2299566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Open(comm, user.fname, FILE_MODE_READ, &viewer)); 2309566063dSJacob Faibussowitsch PetscCall(PetscSectionLoad(section, viewer)); 2319566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 2329566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)section, "Load: section")); 2339566063dSJacob Faibussowitsch PetscCall(PetscSectionView(section, PETSC_VIEWER_STDOUT_(comm))); 2349566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(§ion)); 2353ceb54acSksagiyam } 2369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&comm)); 2373ceb54acSksagiyam 2383ceb54acSksagiyam /* Finalize */ 2399566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 240b122ec5aSJacob Faibussowitsch return 0; 2413ceb54acSksagiyam } 2423ceb54acSksagiyam 2433ceb54acSksagiyam /*TEST 2443ceb54acSksagiyam 2453ceb54acSksagiyam build: 2463ceb54acSksagiyam requires: hdf5 2473ceb54acSksagiyam requires: !complex 2483ceb54acSksagiyam testset: 2493ceb54acSksagiyam nsize: 4 2503ceb54acSksagiyam test: 2513ceb54acSksagiyam suffix: 0 2523ceb54acSksagiyam args: -fname ex5_dump.h5 -includes_constraints 0 2533ceb54acSksagiyam test: 2543ceb54acSksagiyam suffix: 1 2553ceb54acSksagiyam args: -fname ex5_dump.h5 -includes_constraints 1 2563ceb54acSksagiyam 2573ceb54acSksagiyam TEST*/ 258