xref: /petsc/src/vec/is/tests/ex5.c (revision 1e1ea65d8de51fde77ce8a787efbef25e407badc)
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   PetscErrorCode  ierr;
993ceb54acSksagiyam 
1003ceb54acSksagiyam   PetscFunctionBegin;
1013ceb54acSksagiyam   options->fname[0] = '\0';
1023ceb54acSksagiyam   options->includes_constraints = PETSC_TRUE;
1033ceb54acSksagiyam   ierr = PetscOptionsBegin(comm, "", "PetscSectionView()/Load() in HDF5 Test Options", "DMPLEX");CHKERRQ(ierr);
1043ceb54acSksagiyam   ierr = PetscOptionsString("-fname", "The output file", "ex5.c", options->fname, options->fname, sizeof(options->fname), NULL);CHKERRQ(ierr);
1053ceb54acSksagiyam   ierr = PetscOptionsBool("-includes_constraints", "Flag for if global section is to include constrained DoFs or not", "ex5.c", options->includes_constraints, &options->includes_constraints, NULL);CHKERRQ(ierr);
106*1e1ea65dSPierre Jolivet   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1073ceb54acSksagiyam   PetscFunctionReturn(0);
1083ceb54acSksagiyam }
1093ceb54acSksagiyam 
1103ceb54acSksagiyam int main(int argc, char **argv)
1113ceb54acSksagiyam {
1123ceb54acSksagiyam   MPI_Comm        comm;
1133ceb54acSksagiyam   PetscMPIInt     size, rank, mycolor;
1143ceb54acSksagiyam   AppCtx          user;
1153ceb54acSksagiyam   PetscErrorCode  ierr;
1163ceb54acSksagiyam 
1173ceb54acSksagiyam   ierr = PetscInitialize(&argc, &argv, NULL, help); if (ierr) return ierr;
1183ceb54acSksagiyam   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
1193ceb54acSksagiyam   ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRMPI(ierr);
1203ceb54acSksagiyam   ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRMPI(ierr);
1213ceb54acSksagiyam   if (size < 3) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Example only works with three or more processes");
1223ceb54acSksagiyam 
1233ceb54acSksagiyam   /* Save */
1243ceb54acSksagiyam   mycolor = (PetscMPIInt)(rank >= 2);
1253ceb54acSksagiyam   ierr = MPI_Comm_split(PETSC_COMM_WORLD, mycolor, rank, &comm);CHKERRMPI(ierr);
1263ceb54acSksagiyam   if (mycolor == 0) {
1273ceb54acSksagiyam     PetscSection  section, gsection;
1283ceb54acSksagiyam     PetscSF       sf;
1293ceb54acSksagiyam     PetscInt      nroots = -1, nleaves = -1, *ilocal;
1303ceb54acSksagiyam     PetscSFNode  *iremote;
1313ceb54acSksagiyam     PetscViewer   viewer;
1323ceb54acSksagiyam 
1333ceb54acSksagiyam     /* Create section */
1343ceb54acSksagiyam     ierr = PetscSectionCreate(comm, &section);CHKERRQ(ierr);
1353ceb54acSksagiyam     ierr = PetscSectionSetNumFields(section, 2);CHKERRQ(ierr);
1363ceb54acSksagiyam     switch (rank) {
1373ceb54acSksagiyam     case 0:
1383ceb54acSksagiyam       ierr = PetscSectionSetChart(section, 0, 4);CHKERRQ(ierr);
1393ceb54acSksagiyam       ierr = PetscSectionSetDof(section, 0, 3);CHKERRQ(ierr);
1403ceb54acSksagiyam       ierr = PetscSectionSetDof(section, 1, 3);CHKERRQ(ierr);
1413ceb54acSksagiyam       ierr = PetscSectionSetDof(section, 2, 5);CHKERRQ(ierr);
1423ceb54acSksagiyam       ierr = PetscSectionSetDof(section, 3, 7);CHKERRQ(ierr);
1433ceb54acSksagiyam       ierr = PetscSectionSetFieldDof(section, 0, 0, 2);CHKERRQ(ierr);
1443ceb54acSksagiyam       ierr = PetscSectionSetFieldDof(section, 1, 0, 3);CHKERRQ(ierr);
1453ceb54acSksagiyam       ierr = PetscSectionSetFieldDof(section, 2, 0, 5);CHKERRQ(ierr);
1463ceb54acSksagiyam       ierr = PetscSectionSetFieldDof(section, 3, 0, 7);CHKERRQ(ierr);
1473ceb54acSksagiyam       ierr = PetscSectionSetFieldDof(section, 0, 1, 1);CHKERRQ(ierr);
1483ceb54acSksagiyam       break;
1493ceb54acSksagiyam     case 1:
1503ceb54acSksagiyam       ierr = PetscSectionSetChart(section, 0, 3);CHKERRQ(ierr);
1513ceb54acSksagiyam       ierr = PetscSectionSetDof(section, 0, 7);CHKERRQ(ierr);
1523ceb54acSksagiyam       ierr = PetscSectionSetDof(section, 1, 5);CHKERRQ(ierr);
1533ceb54acSksagiyam       ierr = PetscSectionSetDof(section, 2, 13);CHKERRQ(ierr);
1543ceb54acSksagiyam       ierr = PetscSectionSetConstraintDof(section, 2, 1);CHKERRQ(ierr);
1553ceb54acSksagiyam       ierr = PetscSectionSetFieldDof(section, 0, 0, 7);CHKERRQ(ierr);
1563ceb54acSksagiyam       ierr = PetscSectionSetFieldDof(section, 1, 0, 5);CHKERRQ(ierr);
1573ceb54acSksagiyam       ierr = PetscSectionSetFieldDof(section, 2, 0, 11);CHKERRQ(ierr);
1583ceb54acSksagiyam       ierr = PetscSectionSetFieldDof(section, 2, 1, 2);CHKERRQ(ierr);
1593ceb54acSksagiyam       ierr = PetscSectionSetFieldConstraintDof(section, 2, 0, 1);CHKERRQ(ierr);
1603ceb54acSksagiyam       break;
1613ceb54acSksagiyam     }
1623ceb54acSksagiyam     ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
1633ceb54acSksagiyam     if (rank == 1)
1643ceb54acSksagiyam     {
1653ceb54acSksagiyam       const PetscInt indices[] = {7};
1663ceb54acSksagiyam       const PetscInt indices0[] = {7};
1673ceb54acSksagiyam 
1683ceb54acSksagiyam       ierr = PetscSectionSetConstraintIndices(section, 2, indices);CHKERRQ(ierr);
1693ceb54acSksagiyam       ierr = PetscSectionSetFieldConstraintIndices(section, 2, 0, indices0);CHKERRQ(ierr);
1703ceb54acSksagiyam     }
1713ceb54acSksagiyam     /* Create sf */
1723ceb54acSksagiyam     switch (rank) {
1733ceb54acSksagiyam     case 0:
1743ceb54acSksagiyam       nroots = 4;
1753ceb54acSksagiyam       nleaves = 1;
1763ceb54acSksagiyam       ierr = PetscMalloc1(nleaves, &ilocal);CHKERRQ(ierr);
1773ceb54acSksagiyam       ierr = PetscMalloc1(nleaves, &iremote);CHKERRQ(ierr);
1783ceb54acSksagiyam       ilocal[0] = 3;
1793ceb54acSksagiyam       iremote[0].rank = 1;
1803ceb54acSksagiyam       iremote[0].index = 0;
1813ceb54acSksagiyam       break;
1823ceb54acSksagiyam     case 1:
1833ceb54acSksagiyam       nroots = 3;
1843ceb54acSksagiyam       nleaves = 1;
1853ceb54acSksagiyam       ierr = PetscMalloc1(nleaves, &ilocal);CHKERRQ(ierr);
1863ceb54acSksagiyam       ierr = PetscMalloc1(nleaves, &iremote);CHKERRQ(ierr);
1873ceb54acSksagiyam       ilocal[0] = 1;
1883ceb54acSksagiyam       iremote[0].rank = 0;
1893ceb54acSksagiyam       iremote[0].index = 2;
1903ceb54acSksagiyam       break;
1913ceb54acSksagiyam     }
1923ceb54acSksagiyam     ierr = PetscSFCreate(comm, &sf);CHKERRQ(ierr);
1933ceb54acSksagiyam     ierr = PetscSFSetGraph(sf, nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr);
1943ceb54acSksagiyam     /* Create global section*/
1953ceb54acSksagiyam     ierr = PetscSectionCreateGlobalSection(section, sf, user.includes_constraints, PETSC_FALSE, &gsection);CHKERRQ(ierr);
1963ceb54acSksagiyam     ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
1973ceb54acSksagiyam     /* View */
1983ceb54acSksagiyam     ierr = PetscViewerHDF5Open(comm, user.fname, FILE_MODE_WRITE, &viewer);CHKERRQ(ierr);
1993ceb54acSksagiyam     ierr = PetscSectionView(gsection, viewer);CHKERRQ(ierr);
2003ceb54acSksagiyam     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
2013ceb54acSksagiyam     ierr = PetscObjectSetName((PetscObject)section, "Save: local section");CHKERRQ(ierr);
2023ceb54acSksagiyam     ierr = PetscSectionView(section, PETSC_VIEWER_STDOUT_(comm));CHKERRQ(ierr);
2033ceb54acSksagiyam     ierr = PetscObjectSetName((PetscObject)gsection, "Save: global section");CHKERRQ(ierr);
2043ceb54acSksagiyam     ierr = PetscSectionView(gsection, PETSC_VIEWER_STDOUT_(comm));CHKERRQ(ierr);
2053ceb54acSksagiyam     ierr = PetscSectionDestroy(&gsection);CHKERRQ(ierr);
2063ceb54acSksagiyam     ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
2073ceb54acSksagiyam   }
2083ceb54acSksagiyam   ierr = MPI_Comm_free(&comm);CHKERRMPI(ierr);
2093ceb54acSksagiyam 
2103ceb54acSksagiyam   /* Load */
2113ceb54acSksagiyam   mycolor = (PetscMPIInt)(rank >= 3);
2123ceb54acSksagiyam   ierr = MPI_Comm_split(PETSC_COMM_WORLD, mycolor, rank, &comm);CHKERRMPI(ierr);
2133ceb54acSksagiyam   if (mycolor == 0) {
2143ceb54acSksagiyam     PetscSection  section;
2153ceb54acSksagiyam     PetscInt      chartSize = -1;
2163ceb54acSksagiyam     PetscViewer   viewer;
2173ceb54acSksagiyam 
2183ceb54acSksagiyam     ierr = PetscSectionCreate(comm, &section);CHKERRQ(ierr);
2193ceb54acSksagiyam     switch (rank) {
2203ceb54acSksagiyam     case 0:
2213ceb54acSksagiyam       chartSize = 4;
2223ceb54acSksagiyam       break;
2233ceb54acSksagiyam     case 1:
2243ceb54acSksagiyam       chartSize = 0;
2253ceb54acSksagiyam       break;
2263ceb54acSksagiyam     case 2:
2273ceb54acSksagiyam       chartSize = 1;
2283ceb54acSksagiyam       break;
2293ceb54acSksagiyam     }
2303ceb54acSksagiyam     ierr = PetscSectionSetChart(section, 0, chartSize);CHKERRQ(ierr);
2313ceb54acSksagiyam     ierr = PetscViewerHDF5Open(comm, user.fname, FILE_MODE_READ, &viewer);CHKERRQ(ierr);
2323ceb54acSksagiyam     ierr = PetscSectionLoad(section, viewer);CHKERRQ(ierr);
2333ceb54acSksagiyam     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
2343ceb54acSksagiyam     ierr = PetscObjectSetName((PetscObject)section, "Load: section");CHKERRQ(ierr);
2353ceb54acSksagiyam     ierr = PetscSectionView(section, PETSC_VIEWER_STDOUT_(comm));CHKERRQ(ierr);
2363ceb54acSksagiyam     ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
2373ceb54acSksagiyam   }
2383ceb54acSksagiyam   ierr = MPI_Comm_free(&comm);CHKERRMPI(ierr);
2393ceb54acSksagiyam 
2403ceb54acSksagiyam   /* Finalize */
2413ceb54acSksagiyam   ierr = PetscFinalize();
2423ceb54acSksagiyam   return ierr;
2433ceb54acSksagiyam }
2443ceb54acSksagiyam 
2453ceb54acSksagiyam /*TEST
2463ceb54acSksagiyam 
2473ceb54acSksagiyam   build:
2483ceb54acSksagiyam     requires: hdf5
2493ceb54acSksagiyam     requires: !complex
2503ceb54acSksagiyam   testset:
2513ceb54acSksagiyam     nsize: 4
2523ceb54acSksagiyam     test:
2533ceb54acSksagiyam       suffix: 0
2543ceb54acSksagiyam       args: -fname ex5_dump.h5 -includes_constraints 0
2553ceb54acSksagiyam     test:
2563ceb54acSksagiyam       suffix: 1
2573ceb54acSksagiyam       args: -fname ex5_dump.h5 -includes_constraints 1
2583ceb54acSksagiyam 
2593ceb54acSksagiyam TEST*/
260