1 static char help[] = "Tests PetscSectionView()/Load() with HDF5.\n\n"; 2 3 #include <petscdmshell.h> 4 #include <petscdmplex.h> 5 #include <petscsection.h> 6 #include <petscsf.h> 7 #include <petsclayouthdf5.h> 8 9 /* Save/Load abstract sections 10 11 ===================== 12 Save on 2 processes 13 ===================== 14 15 section: 16 0 1 2 3 17 rank 0: Dof (Field 0) 2 3 5 7 18 Dof (Field 1) 1 0 0 0 19 20 0 1 2 21 rank 1: Dof (Field 0) 7 5 11 <- DoF 7 is constrained 22 Dof (Field 1) 0 0 2 23 24 sf: 25 [0] 3 <- (1, 0) 26 [1] 1 <- (0, 2) 27 28 global section (includesConstraints = PETSC_FALSE): 29 0 1 2 3 30 rank 0: Dof (Field 0) 2 3 5 -8 31 Off (Field 0) 0 3 6 -12 32 Dof (Field 1) 1 0 0 -1 33 Off (Field 1) 2 6 11 -19 34 35 0 1 2 36 rank 1: Dof (Field 0) 7 -6 11 37 Off (Field 0) 11 -7 18 38 Dof (Field 1) 0 -1 2 39 Off (Field 1) 18 -12 28 40 41 global section (includesConstraints = PETSC_TRUE): 42 0 1 2 3 43 rank 0: Dof (Field 0) 2 3 5 -8 44 Off (Field 0) 0 3 6 -12 45 Dof (Field 1) 1 0 0 -1 46 Off (Field 1) 2 6 11 -19 47 48 0 1 2 49 rank 1: Dof (Field 0) 7 -6 11 50 Off (Field 0) 11 -7 18 51 Dof (Field 1) 0 -1 2 52 Off (Field 1) 18 -12 29 53 54 ===================== 55 Load on 3 Processes 56 ===================== 57 58 (Set chartSize = 4, 0, 1 for rank 0, 1, 2, respectively) 59 60 global section (includesConstraints = PETSC_FALSE): 61 62 rank 0: Dof (Field 0) 2 3 5 7 63 Off (Field 0) 0 3 6 11 64 Dof (Field 1) 1 0 0 0 65 Off (Field 1) 2 6 11 18 66 67 rank 1: Dof (Field 0) 68 Dof (Field 1) 69 70 rank 2: Dof (Field 0) 11 71 Off (Field 0) 18 72 Dof (Field 1) 2 73 Off (Field 1) 28 74 75 global section (includesConstraints = PETSC_TRUE): 76 77 rank 0: Dof (Field 0) 2 3 5 7 78 Off (Field 0) 0 3 6 11 79 Dof (Field 1) 1 0 0 0 80 Off (Field 1) 2 6 11 18 81 82 rank 1: Dof (Field 0) 83 Dof (Field 1) 84 85 rank 2: Dof (Field 0) 11 86 Off (Field 0) 18 87 Dof (Field 1) 2 88 Off (Field 1) 29 89 */ 90 91 typedef struct { 92 char fname[PETSC_MAX_PATH_LEN]; /* Output mesh filename */ 93 PetscBool includes_constraints; /* Flag for if global section is to include constrained DoFs or not */ 94 } AppCtx; 95 96 PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 97 { 98 PetscErrorCode ierr; 99 100 PetscFunctionBegin; 101 options->fname[0] = '\0'; 102 options->includes_constraints = PETSC_TRUE; 103 ierr = PetscOptionsBegin(comm, "", "PetscSectionView()/Load() in HDF5 Test Options", "DMPLEX");CHKERRQ(ierr); 104 CHKERRQ(PetscOptionsString("-fname", "The output file", "ex5.c", options->fname, options->fname, sizeof(options->fname), NULL)); 105 CHKERRQ(PetscOptionsBool("-includes_constraints", "Flag for if global section is to include constrained DoFs or not", "ex5.c", options->includes_constraints, &options->includes_constraints, NULL)); 106 ierr = PetscOptionsEnd();CHKERRQ(ierr); 107 PetscFunctionReturn(0); 108 } 109 110 int main(int argc, char **argv) 111 { 112 MPI_Comm comm; 113 PetscMPIInt size, rank, mycolor; 114 AppCtx user; 115 PetscErrorCode ierr; 116 117 ierr = PetscInitialize(&argc, &argv, NULL, help); if (ierr) return ierr; 118 CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user)); 119 CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 120 CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 121 PetscCheckFalse(size < 3,PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Example only works with three or more processes"); 122 123 /* Save */ 124 mycolor = (PetscMPIInt)(rank >= 2); 125 CHKERRMPI(MPI_Comm_split(PETSC_COMM_WORLD, mycolor, rank, &comm)); 126 if (mycolor == 0) { 127 PetscSection section, gsection; 128 PetscSF sf; 129 PetscInt nroots = -1, nleaves = -1, *ilocal; 130 PetscSFNode *iremote; 131 PetscViewer viewer; 132 133 /* Create section */ 134 CHKERRQ(PetscSectionCreate(comm, §ion)); 135 CHKERRQ(PetscSectionSetNumFields(section, 2)); 136 switch (rank) { 137 case 0: 138 CHKERRQ(PetscSectionSetChart(section, 0, 4)); 139 CHKERRQ(PetscSectionSetDof(section, 0, 3)); 140 CHKERRQ(PetscSectionSetDof(section, 1, 3)); 141 CHKERRQ(PetscSectionSetDof(section, 2, 5)); 142 CHKERRQ(PetscSectionSetDof(section, 3, 7)); 143 CHKERRQ(PetscSectionSetFieldDof(section, 0, 0, 2)); 144 CHKERRQ(PetscSectionSetFieldDof(section, 1, 0, 3)); 145 CHKERRQ(PetscSectionSetFieldDof(section, 2, 0, 5)); 146 CHKERRQ(PetscSectionSetFieldDof(section, 3, 0, 7)); 147 CHKERRQ(PetscSectionSetFieldDof(section, 0, 1, 1)); 148 break; 149 case 1: 150 CHKERRQ(PetscSectionSetChart(section, 0, 3)); 151 CHKERRQ(PetscSectionSetDof(section, 0, 7)); 152 CHKERRQ(PetscSectionSetDof(section, 1, 5)); 153 CHKERRQ(PetscSectionSetDof(section, 2, 13)); 154 CHKERRQ(PetscSectionSetConstraintDof(section, 2, 1)); 155 CHKERRQ(PetscSectionSetFieldDof(section, 0, 0, 7)); 156 CHKERRQ(PetscSectionSetFieldDof(section, 1, 0, 5)); 157 CHKERRQ(PetscSectionSetFieldDof(section, 2, 0, 11)); 158 CHKERRQ(PetscSectionSetFieldDof(section, 2, 1, 2)); 159 CHKERRQ(PetscSectionSetFieldConstraintDof(section, 2, 0, 1)); 160 break; 161 } 162 CHKERRQ(PetscSectionSetUp(section)); 163 if (rank == 1) 164 { 165 const PetscInt indices[] = {7}; 166 const PetscInt indices0[] = {7}; 167 168 CHKERRQ(PetscSectionSetConstraintIndices(section, 2, indices)); 169 CHKERRQ(PetscSectionSetFieldConstraintIndices(section, 2, 0, indices0)); 170 } 171 /* Create sf */ 172 switch (rank) { 173 case 0: 174 nroots = 4; 175 nleaves = 1; 176 CHKERRQ(PetscMalloc1(nleaves, &ilocal)); 177 CHKERRQ(PetscMalloc1(nleaves, &iremote)); 178 ilocal[0] = 3; 179 iremote[0].rank = 1; 180 iremote[0].index = 0; 181 break; 182 case 1: 183 nroots = 3; 184 nleaves = 1; 185 CHKERRQ(PetscMalloc1(nleaves, &ilocal)); 186 CHKERRQ(PetscMalloc1(nleaves, &iremote)); 187 ilocal[0] = 1; 188 iremote[0].rank = 0; 189 iremote[0].index = 2; 190 break; 191 } 192 CHKERRQ(PetscSFCreate(comm, &sf)); 193 CHKERRQ(PetscSFSetGraph(sf, nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER)); 194 /* Create global section*/ 195 CHKERRQ(PetscSectionCreateGlobalSection(section, sf, user.includes_constraints, PETSC_FALSE, &gsection)); 196 CHKERRQ(PetscSFDestroy(&sf)); 197 /* View */ 198 CHKERRQ(PetscViewerHDF5Open(comm, user.fname, FILE_MODE_WRITE, &viewer)); 199 CHKERRQ(PetscSectionView(gsection, viewer)); 200 CHKERRQ(PetscViewerDestroy(&viewer)); 201 CHKERRQ(PetscObjectSetName((PetscObject)section, "Save: local section")); 202 CHKERRQ(PetscSectionView(section, PETSC_VIEWER_STDOUT_(comm))); 203 CHKERRQ(PetscObjectSetName((PetscObject)gsection, "Save: global section")); 204 CHKERRQ(PetscSectionView(gsection, PETSC_VIEWER_STDOUT_(comm))); 205 CHKERRQ(PetscSectionDestroy(&gsection)); 206 CHKERRQ(PetscSectionDestroy(§ion)); 207 } 208 CHKERRMPI(MPI_Comm_free(&comm)); 209 210 /* Load */ 211 mycolor = (PetscMPIInt)(rank >= 3); 212 CHKERRMPI(MPI_Comm_split(PETSC_COMM_WORLD, mycolor, rank, &comm)); 213 if (mycolor == 0) { 214 PetscSection section; 215 PetscInt chartSize = -1; 216 PetscViewer viewer; 217 218 CHKERRQ(PetscSectionCreate(comm, §ion)); 219 switch (rank) { 220 case 0: 221 chartSize = 4; 222 break; 223 case 1: 224 chartSize = 0; 225 break; 226 case 2: 227 chartSize = 1; 228 break; 229 } 230 CHKERRQ(PetscSectionSetChart(section, 0, chartSize)); 231 CHKERRQ(PetscViewerHDF5Open(comm, user.fname, FILE_MODE_READ, &viewer)); 232 CHKERRQ(PetscSectionLoad(section, viewer)); 233 CHKERRQ(PetscViewerDestroy(&viewer)); 234 CHKERRQ(PetscObjectSetName((PetscObject)section, "Load: section")); 235 CHKERRQ(PetscSectionView(section, PETSC_VIEWER_STDOUT_(comm))); 236 CHKERRQ(PetscSectionDestroy(§ion)); 237 } 238 CHKERRMPI(MPI_Comm_free(&comm)); 239 240 /* Finalize */ 241 ierr = PetscFinalize(); 242 return ierr; 243 } 244 245 /*TEST 246 247 build: 248 requires: hdf5 249 requires: !complex 250 testset: 251 nsize: 4 252 test: 253 suffix: 0 254 args: -fname ex5_dump.h5 -includes_constraints 0 255 test: 256 suffix: 1 257 args: -fname ex5_dump.h5 -includes_constraints 1 258 259 TEST*/ 260