1*607e733fSJames Wright static char help[] = "Test CGNS parallel load-save-reload cycle, including data and DMLabels\n\n";
2*607e733fSJames Wright // This is a modification of src/dm/impls/plex/tutorials/ex15.c, but with additional tests that don't make sense for a tutorial problem (such as verify FaceLabels)
3*607e733fSJames Wright
4*607e733fSJames Wright #include <petscdmlabel.h>
5*607e733fSJames Wright #include <petscdmplex.h>
6*607e733fSJames Wright #include <petscsf.h>
7*607e733fSJames Wright #include <petscviewerhdf5.h>
8*607e733fSJames Wright #define EX "ex103.c"
9*607e733fSJames Wright
10*607e733fSJames Wright typedef struct {
11*607e733fSJames Wright char infile[PETSC_MAX_PATH_LEN]; /* Input mesh filename */
12*607e733fSJames Wright char outfile[PETSC_MAX_PATH_LEN]; /* Dump/reload mesh filename */
13*607e733fSJames Wright PetscBool heterogeneous; /* Test save on N / load on M */
14*607e733fSJames Wright PetscInt ntimes; /* How many times do the cycle */
15*607e733fSJames Wright } AppCtx;
16*607e733fSJames Wright
ProcessOptions(MPI_Comm comm,AppCtx * options)17*607e733fSJames Wright static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
18*607e733fSJames Wright {
19*607e733fSJames Wright PetscBool flg;
20*607e733fSJames Wright
21*607e733fSJames Wright PetscFunctionBeginUser;
22*607e733fSJames Wright options->infile[0] = '\0';
23*607e733fSJames Wright options->outfile[0] = '\0';
24*607e733fSJames Wright options->heterogeneous = PETSC_FALSE;
25*607e733fSJames Wright options->ntimes = 2;
26*607e733fSJames Wright PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
27*607e733fSJames Wright PetscCall(PetscOptionsString("-infile", "The input CGNS file", EX, options->infile, options->infile, sizeof(options->infile), &flg));
28*607e733fSJames Wright PetscCall(PetscOptionsString("-outfile", "The output CGNS file", EX, options->outfile, options->outfile, sizeof(options->outfile), &flg));
29*607e733fSJames Wright PetscCall(PetscOptionsBool("-heterogeneous", "Test save on N / load on M", EX, options->heterogeneous, &options->heterogeneous, NULL));
30*607e733fSJames Wright PetscCall(PetscOptionsInt("-ntimes", "How many times do the cycle", EX, options->ntimes, &options->ntimes, NULL));
31*607e733fSJames Wright PetscOptionsEnd();
32*607e733fSJames Wright PetscCheck(flg, comm, PETSC_ERR_USER_INPUT, "-infile needs to be specified");
33*607e733fSJames Wright PetscCheck(flg, comm, PETSC_ERR_USER_INPUT, "-outfile needs to be specified");
34*607e733fSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
35*607e733fSJames Wright }
36*607e733fSJames Wright
37*607e733fSJames Wright // @brief Create DM from CGNS file and setup PetscFE to VecLoad solution from that file
ReadCGNSDM(MPI_Comm comm,const char filename[],DM * dm)38*607e733fSJames Wright PetscErrorCode ReadCGNSDM(MPI_Comm comm, const char filename[], DM *dm)
39*607e733fSJames Wright {
40*607e733fSJames Wright PetscInt degree;
41*607e733fSJames Wright
42*607e733fSJames Wright PetscFunctionBeginUser;
43*607e733fSJames Wright PetscCall(DMPlexCreateFromFile(comm, filename, "ex15_plex", PETSC_TRUE, dm));
44*607e733fSJames Wright PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
45*607e733fSJames Wright PetscCall(DMSetFromOptions(*dm));
46*607e733fSJames Wright PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
47*607e733fSJames Wright
48*607e733fSJames Wright /* Redistribute */
49*607e733fSJames Wright PetscCall(DMSetOptionsPrefix(*dm, "redistributed_"));
50*607e733fSJames Wright PetscCall(DMSetFromOptions(*dm));
51*607e733fSJames Wright PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
52*607e733fSJames Wright
53*607e733fSJames Wright { // Get degree of the natural section
54*607e733fSJames Wright PetscFE fe_natural;
55*607e733fSJames Wright PetscDualSpace dual_space_natural;
56*607e733fSJames Wright
57*607e733fSJames Wright PetscCall(DMGetField(*dm, 0, NULL, (PetscObject *)&fe_natural));
58*607e733fSJames Wright PetscCall(PetscFEGetDualSpace(fe_natural, &dual_space_natural));
59*607e733fSJames Wright PetscCall(PetscDualSpaceGetOrder(dual_space_natural, °ree));
60*607e733fSJames Wright PetscCall(DMClearFields(*dm));
61*607e733fSJames Wright PetscCall(DMSetLocalSection(*dm, NULL));
62*607e733fSJames Wright }
63*607e733fSJames Wright
64*607e733fSJames Wright { // Setup fe to load in the initial condition data
65*607e733fSJames Wright PetscFE fe;
66*607e733fSJames Wright PetscInt dim;
67*607e733fSJames Wright
68*607e733fSJames Wright PetscCall(DMGetDimension(*dm, &dim));
69*607e733fSJames Wright PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, 5, PETSC_FALSE, degree, PETSC_DETERMINE, &fe));
70*607e733fSJames Wright PetscCall(PetscObjectSetName((PetscObject)fe, "FE for VecLoad"));
71*607e733fSJames Wright PetscCall(DMAddField(*dm, NULL, (PetscObject)fe));
72*607e733fSJames Wright PetscCall(DMCreateDS(*dm));
73*607e733fSJames Wright PetscCall(PetscFEDestroy(&fe));
74*607e733fSJames Wright }
75*607e733fSJames Wright
76*607e733fSJames Wright // Set section component names, used when writing out CGNS files
77*607e733fSJames Wright PetscSection section;
78*607e733fSJames Wright PetscCall(DMGetLocalSection(*dm, §ion));
79*607e733fSJames Wright PetscCall(PetscSectionSetFieldName(section, 0, ""));
80*607e733fSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "Pressure"));
81*607e733fSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "VelocityX"));
82*607e733fSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "VelocityY"));
83*607e733fSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "VelocityZ"));
84*607e733fSJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "Temperature"));
85*607e733fSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
86*607e733fSJames Wright }
87*607e733fSJames Wright
88*607e733fSJames Wright // Verify that V_load is equivalent to V_serial, even if V_load is distributed
VerifyLoadedSolution(DM dm_serial,Vec V_serial,DM dm_load,Vec V_load,PetscScalar tol)89*607e733fSJames Wright PetscErrorCode VerifyLoadedSolution(DM dm_serial, Vec V_serial, DM dm_load, Vec V_load, PetscScalar tol)
90*607e733fSJames Wright {
91*607e733fSJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)dm_load);
92*607e733fSJames Wright PetscSF load_to_serial_sf_;
93*607e733fSJames Wright PetscScalar *array_load_bcast = NULL;
94*607e733fSJames Wright PetscInt num_comps = 5;
95*607e733fSJames Wright
96*607e733fSJames Wright PetscFunctionBeginUser;
97*607e733fSJames Wright { // Create SF to broadcast loaded vec nodes to serial vec nodes
98*607e733fSJames Wright PetscInt dim, num_local_serial = 0, num_local_load;
99*607e733fSJames Wright Vec coord_Vec_serial, coord_Vec_load;
100*607e733fSJames Wright const PetscScalar *coord_serial = NULL, *coord_load;
101*607e733fSJames Wright
102*607e733fSJames Wright PetscCall(DMGetCoordinateDim(dm_load, &dim));
103*607e733fSJames Wright PetscCall(DMGetCoordinates(dm_load, &coord_Vec_load));
104*607e733fSJames Wright PetscCall(VecGetLocalSize(coord_Vec_load, &num_local_load));
105*607e733fSJames Wright num_local_load /= dim;
106*607e733fSJames Wright
107*607e733fSJames Wright PetscCall(VecGetArrayRead(coord_Vec_load, &coord_load));
108*607e733fSJames Wright
109*607e733fSJames Wright if (dm_serial) {
110*607e733fSJames Wright PetscCall(DMGetCoordinates(dm_serial, &coord_Vec_serial));
111*607e733fSJames Wright PetscCall(VecGetLocalSize(coord_Vec_serial, &num_local_serial));
112*607e733fSJames Wright num_local_serial /= dim;
113*607e733fSJames Wright PetscCall(VecGetArrayRead(coord_Vec_serial, &coord_serial));
114*607e733fSJames Wright }
115*607e733fSJames Wright
116*607e733fSJames Wright PetscCall(PetscSFCreate(comm, &load_to_serial_sf_));
117*607e733fSJames Wright PetscCall(PetscSFSetGraphFromCoordinates(load_to_serial_sf_, num_local_load, num_local_serial, dim, 100 * PETSC_MACHINE_EPSILON, coord_load, coord_serial));
118*607e733fSJames Wright PetscCall(PetscSFViewFromOptions(load_to_serial_sf_, NULL, "-verify_sf_view"));
119*607e733fSJames Wright
120*607e733fSJames Wright PetscCall(VecRestoreArrayRead(coord_Vec_load, &coord_load));
121*607e733fSJames Wright if (dm_serial) PetscCall(VecRestoreArrayRead(coord_Vec_serial, &coord_serial));
122*607e733fSJames Wright }
123*607e733fSJames Wright
124*607e733fSJames Wright { // Broadcast the loaded vector data into a serial-sized array
125*607e733fSJames Wright PetscInt size_local_serial = 0;
126*607e733fSJames Wright const PetscScalar *array_load;
127*607e733fSJames Wright MPI_Datatype unit;
128*607e733fSJames Wright
129*607e733fSJames Wright PetscCall(VecGetArrayRead(V_load, &array_load));
130*607e733fSJames Wright if (V_serial) {
131*607e733fSJames Wright PetscCall(VecGetLocalSize(V_serial, &size_local_serial));
132*607e733fSJames Wright PetscCall(PetscMalloc1(size_local_serial, &array_load_bcast));
133*607e733fSJames Wright }
134*607e733fSJames Wright
135*607e733fSJames Wright PetscCallMPI(MPI_Type_contiguous(num_comps, MPIU_REAL, &unit));
136*607e733fSJames Wright PetscCallMPI(MPI_Type_commit(&unit));
137*607e733fSJames Wright PetscCall(PetscSFBcastBegin(load_to_serial_sf_, unit, array_load, array_load_bcast, MPI_REPLACE));
138*607e733fSJames Wright PetscCall(PetscSFBcastEnd(load_to_serial_sf_, unit, array_load, array_load_bcast, MPI_REPLACE));
139*607e733fSJames Wright PetscCallMPI(MPI_Type_free(&unit));
140*607e733fSJames Wright PetscCall(VecRestoreArrayRead(V_load, &array_load));
141*607e733fSJames Wright }
142*607e733fSJames Wright
143*607e733fSJames Wright if (V_serial) {
144*607e733fSJames Wright const PetscScalar *array_serial;
145*607e733fSJames Wright PetscInt size_local_serial;
146*607e733fSJames Wright
147*607e733fSJames Wright PetscCall(VecGetArrayRead(V_serial, &array_serial));
148*607e733fSJames Wright PetscCall(VecGetLocalSize(V_serial, &size_local_serial));
149*607e733fSJames Wright for (PetscInt i = 0; i < size_local_serial; i++) {
150*607e733fSJames Wright if (PetscAbs(array_serial[i] - array_load_bcast[i]) > tol) PetscCall(PetscPrintf(comm, "DoF %" PetscInt_FMT " is inconsistent. Found %g, expected %g\n", i, (double)array_load_bcast[i], (double)array_serial[i]));
151*607e733fSJames Wright }
152*607e733fSJames Wright PetscCall(VecRestoreArrayRead(V_serial, &array_serial));
153*607e733fSJames Wright }
154*607e733fSJames Wright
155*607e733fSJames Wright PetscCall(PetscFree(array_load_bcast));
156*607e733fSJames Wright PetscCall(PetscSFDestroy(&load_to_serial_sf_));
157*607e733fSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
158*607e733fSJames Wright }
159*607e733fSJames Wright
160*607e733fSJames Wright // Get centroids associated with every Plex point
DMPlexGetPointsCentroids(DM dm,PetscReal * points_centroids[])161*607e733fSJames Wright PetscErrorCode DMPlexGetPointsCentroids(DM dm, PetscReal *points_centroids[])
162*607e733fSJames Wright {
163*607e733fSJames Wright PetscInt coords_dim, pStart, pEnd, depth = 0;
164*607e733fSJames Wright PetscSection coord_section;
165*607e733fSJames Wright Vec coords_vec;
166*607e733fSJames Wright PetscReal *points_centroids_;
167*607e733fSJames Wright
168*607e733fSJames Wright PetscFunctionBeginUser;
169*607e733fSJames Wright PetscCall(DMGetCoordinateDim(dm, &coords_dim));
170*607e733fSJames Wright PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
171*607e733fSJames Wright PetscCall(DMGetCoordinateSection(dm, &coord_section));
172*607e733fSJames Wright PetscCall(DMGetCoordinatesLocal(dm, &coords_vec));
173*607e733fSJames Wright
174*607e733fSJames Wright PetscCall(PetscCalloc1((pEnd - pStart) * coords_dim, &points_centroids_));
175*607e733fSJames Wright for (PetscInt p = pStart; p < pEnd; p++) {
176*607e733fSJames Wright PetscScalar *coords = NULL;
177*607e733fSJames Wright PetscInt coords_size, num_coords;
178*607e733fSJames Wright
179*607e733fSJames Wright PetscCall(DMPlexVecGetClosureAtDepth(dm, coord_section, coords_vec, p, depth, &coords_size, &coords));
180*607e733fSJames Wright num_coords = coords_size / coords_dim;
181*607e733fSJames Wright for (PetscInt c = 0; c < num_coords; c++) {
182*607e733fSJames Wright for (PetscInt d = 0; d < coords_dim; d++) points_centroids_[p * coords_dim + d] += PetscRealPart(coords[c * coords_dim + d]);
183*607e733fSJames Wright }
184*607e733fSJames Wright for (PetscInt d = 0; d < coords_dim; d++) points_centroids_[p * coords_dim + d] /= num_coords;
185*607e733fSJames Wright PetscCall(DMPlexVecRestoreClosure(dm, coord_section, coords_vec, p, &coords_size, &coords));
186*607e733fSJames Wright }
187*607e733fSJames Wright *points_centroids = points_centroids_;
188*607e733fSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
189*607e733fSJames Wright }
190*607e733fSJames Wright
191*607e733fSJames Wright // Verify that the DMLabel in dm_serial is identical to that in dm_load
VerifyDMLabels(DM dm_serial,DM dm_load,const char label_name[],PetscSF * serial2loadPointSF)192*607e733fSJames Wright PetscErrorCode VerifyDMLabels(DM dm_serial, DM dm_load, const char label_name[], PetscSF *serial2loadPointSF)
193*607e733fSJames Wright {
194*607e733fSJames Wright PetscMPIInt rank;
195*607e733fSJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)dm_load);
196*607e733fSJames Wright PetscInt num_values_serial = 0, dim;
197*607e733fSJames Wright PetscInt *values_serial = NULL;
198*607e733fSJames Wright DMLabel label_serial = NULL, label_load;
199*607e733fSJames Wright
200*607e733fSJames Wright PetscFunctionBeginUser;
201*607e733fSJames Wright PetscCall(DMGetCoordinateDim(dm_load, &dim));
202*607e733fSJames Wright PetscCallMPI(MPI_Comm_rank(comm, &rank));
203*607e733fSJames Wright if (dm_serial) { // Communicate valid label values to all ranks
204*607e733fSJames Wright IS serialValuesIS;
205*607e733fSJames Wright const PetscInt *values_serial_is;
206*607e733fSJames Wright
207*607e733fSJames Wright PetscCheck(rank == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Rank with serial DM not rank 0");
208*607e733fSJames Wright PetscCall(DMGetLabel(dm_serial, label_name, &label_serial));
209*607e733fSJames Wright PetscCall(DMLabelGetNonEmptyStratumValuesIS(label_serial, &serialValuesIS));
210*607e733fSJames Wright PetscCall(ISGetLocalSize(serialValuesIS, &num_values_serial));
211*607e733fSJames Wright
212*607e733fSJames Wright PetscCall(PetscMalloc1(num_values_serial, &values_serial));
213*607e733fSJames Wright PetscCall(ISGetIndices(serialValuesIS, &values_serial_is));
214*607e733fSJames Wright PetscCall(PetscArraycpy(values_serial, values_serial_is, num_values_serial));
215*607e733fSJames Wright PetscCall(PetscSortInt(num_values_serial, values_serial));
216*607e733fSJames Wright PetscCall(ISRestoreIndices(serialValuesIS, &values_serial_is));
217*607e733fSJames Wright PetscCall(ISDestroy(&serialValuesIS));
218*607e733fSJames Wright }
219*607e733fSJames Wright PetscCallMPI(MPI_Bcast(&num_values_serial, 1, MPIU_INT, 0, comm));
220*607e733fSJames Wright if (values_serial == NULL) PetscCall(PetscMalloc1(num_values_serial, &values_serial));
221*607e733fSJames Wright PetscCallMPI(MPI_Bcast(values_serial, num_values_serial, MPIU_INT, 0, comm));
222*607e733fSJames Wright
223*607e733fSJames Wright IS loadValuesIS;
224*607e733fSJames Wright PetscInt num_values_global;
225*607e733fSJames Wright const PetscInt *values_global;
226*607e733fSJames Wright PetscBool are_values_same = PETSC_TRUE;
227*607e733fSJames Wright
228*607e733fSJames Wright PetscCall(DMGetLabel(dm_load, label_name, &label_load));
229*607e733fSJames Wright PetscCall(DMLabelGetValueISGlobal(comm, label_load, PETSC_TRUE, &loadValuesIS));
230*607e733fSJames Wright PetscCall(ISSort(loadValuesIS));
231*607e733fSJames Wright PetscCall(ISGetLocalSize(loadValuesIS, &num_values_global));
232*607e733fSJames Wright PetscCall(ISGetIndices(loadValuesIS, &values_global));
233*607e733fSJames Wright if (num_values_serial != num_values_global) {
234*607e733fSJames Wright PetscCall(PetscPrintf(comm, "DMLabel '%s': Number of nonempty values in serial DM (%" PetscInt_FMT ") not equal to number of values in global DM (%" PetscInt_FMT ")\n", label_name, num_values_serial, num_values_global));
235*607e733fSJames Wright are_values_same = PETSC_FALSE;
236*607e733fSJames Wright }
237*607e733fSJames Wright PetscCall(PetscPrintf(comm, "DMLabel '%s': serial values:\n", label_name));
238*607e733fSJames Wright PetscCall(PetscIntView(num_values_serial, values_serial, PETSC_VIEWER_STDOUT_(comm)));
239*607e733fSJames Wright PetscCall(PetscPrintf(comm, "DMLabel '%s': global values:\n", label_name));
240*607e733fSJames Wright PetscCall(PetscIntView(num_values_global, values_global, PETSC_VIEWER_STDOUT_(comm)));
241*607e733fSJames Wright for (PetscInt i = 0; i < num_values_serial; i++) {
242*607e733fSJames Wright PetscInt loc;
243*607e733fSJames Wright PetscCall(PetscFindInt(values_serial[i], num_values_global, values_global, &loc));
244*607e733fSJames Wright if (loc < 0) {
245*607e733fSJames Wright PetscCall(PetscPrintf(comm, "DMLabel '%s': Label value %" PetscInt_FMT " in serial DM not found in global DM\n", label_name, values_serial[i]));
246*607e733fSJames Wright are_values_same = PETSC_FALSE;
247*607e733fSJames Wright }
248*607e733fSJames Wright }
249*607e733fSJames Wright PetscCheck(are_values_same, comm, PETSC_ERR_PLIB, "The values in DMLabel are not the same");
250*607e733fSJames Wright PetscCall(PetscFree(values_serial));
251*607e733fSJames Wright
252*607e733fSJames Wright PetscSF serial2loadPointSF_;
253*607e733fSJames Wright PetscInt pStart, pEnd, pStartSerial = -1, pEndSerial = -2;
254*607e733fSJames Wright PetscInt num_points_load, num_points_serial = 0;
255*607e733fSJames Wright { // Create SF mapping serialDM points to loadDM points
256*607e733fSJames Wright PetscReal *points_centroid_load, *points_centroid_serial = NULL;
257*607e733fSJames Wright
258*607e733fSJames Wright if (rank == 0) {
259*607e733fSJames Wright PetscCall(DMPlexGetChart(dm_serial, &pStartSerial, &pEndSerial));
260*607e733fSJames Wright num_points_serial = pEndSerial - pStartSerial;
261*607e733fSJames Wright PetscCall(DMPlexGetPointsCentroids(dm_serial, &points_centroid_serial));
262*607e733fSJames Wright }
263*607e733fSJames Wright PetscCall(DMPlexGetChart(dm_load, &pStart, &pEnd));
264*607e733fSJames Wright num_points_load = pEnd - pStart;
265*607e733fSJames Wright PetscCall(DMPlexGetPointsCentroids(dm_load, &points_centroid_load));
266*607e733fSJames Wright
267*607e733fSJames Wright PetscCall(PetscSFCreate(comm, &serial2loadPointSF_));
268*607e733fSJames Wright PetscCall(PetscSFSetGraphFromCoordinates(serial2loadPointSF_, num_points_serial, num_points_load, dim, 100 * PETSC_MACHINE_EPSILON, points_centroid_serial, points_centroid_load));
269*607e733fSJames Wright PetscCall(PetscObjectSetName((PetscObject)serial2loadPointSF_, "Serial To Loaded DM Points SF"));
270*607e733fSJames Wright PetscCall(PetscSFViewFromOptions(serial2loadPointSF_, NULL, "-verify_points_sf_view"));
271*607e733fSJames Wright PetscCall(PetscFree(points_centroid_load));
272*607e733fSJames Wright PetscCall(PetscFree(points_centroid_serial));
273*607e733fSJames Wright }
274*607e733fSJames Wright
275*607e733fSJames Wright PetscSection pointSerialSection = NULL;
276*607e733fSJames Wright PetscInt npointMaskSerial = 0;
277*607e733fSJames Wright PetscBool *pointMask, *pointMaskSerial = NULL;
278*607e733fSJames Wright
279*607e733fSJames Wright if (rank == 0) {
280*607e733fSJames Wright const PetscInt *root_degree;
281*607e733fSJames Wright PetscCall(PetscSectionCreate(PETSC_COMM_SELF, &pointSerialSection));
282*607e733fSJames Wright PetscCall(PetscSectionSetChart(pointSerialSection, pStartSerial, pEndSerial));
283*607e733fSJames Wright PetscCall(PetscSFComputeDegreeBegin(serial2loadPointSF_, &root_degree));
284*607e733fSJames Wright PetscCall(PetscSFComputeDegreeEnd(serial2loadPointSF_, &root_degree));
285*607e733fSJames Wright for (PetscInt p = 0; p < num_points_serial; p++) PetscCall(PetscSectionSetDof(pointSerialSection, p, root_degree[p]));
286*607e733fSJames Wright PetscCall(PetscSectionSetUp(pointSerialSection));
287*607e733fSJames Wright PetscCall(PetscSectionGetStorageSize(pointSerialSection, &npointMaskSerial));
288*607e733fSJames Wright
289*607e733fSJames Wright PetscCall(PetscMalloc1(npointMaskSerial, &pointMaskSerial));
290*607e733fSJames Wright }
291*607e733fSJames Wright PetscCall(PetscMalloc1(num_points_load, &pointMask));
292*607e733fSJames Wright
293*607e733fSJames Wright for (PetscInt v = 0; v < num_values_global; v++) {
294*607e733fSJames Wright PetscInt value = values_global[v];
295*607e733fSJames Wright IS stratumIS = NULL;
296*607e733fSJames Wright
297*607e733fSJames Wright if (pointMaskSerial) PetscCall(PetscArrayzero(pointMaskSerial, npointMaskSerial));
298*607e733fSJames Wright PetscCall(PetscArrayzero(pointMask, num_points_load));
299*607e733fSJames Wright PetscCall(DMLabelGetStratumIS(label_load, value, &stratumIS));
300*607e733fSJames Wright if (stratumIS) {
301*607e733fSJames Wright PetscInt num_points;
302*607e733fSJames Wright const PetscInt *points;
303*607e733fSJames Wright
304*607e733fSJames Wright PetscCall(ISGetLocalSize(stratumIS, &num_points));
305*607e733fSJames Wright PetscCall(ISGetIndices(stratumIS, &points));
306*607e733fSJames Wright for (PetscInt p = 0; p < num_points; p++) pointMask[points[p]] = PETSC_TRUE;
307*607e733fSJames Wright PetscCall(ISRestoreIndices(stratumIS, &points));
308*607e733fSJames Wright PetscCall(ISDestroy(&stratumIS));
309*607e733fSJames Wright }
310*607e733fSJames Wright PetscCall(PetscSFGatherBegin(serial2loadPointSF_, MPI_C_BOOL, pointMask, pointMaskSerial));
311*607e733fSJames Wright PetscCall(PetscSFGatherEnd(serial2loadPointSF_, MPI_C_BOOL, pointMask, pointMaskSerial));
312*607e733fSJames Wright
313*607e733fSJames Wright if (rank == 0) {
314*607e733fSJames Wright IS stratumIS = NULL;
315*607e733fSJames Wright
316*607e733fSJames Wright PetscCall(DMLabelGetStratumIS(label_serial, value, &stratumIS));
317*607e733fSJames Wright if (stratumIS) {
318*607e733fSJames Wright PetscInt num_points;
319*607e733fSJames Wright const PetscInt *points;
320*607e733fSJames Wright
321*607e733fSJames Wright PetscCall(ISSort(stratumIS));
322*607e733fSJames Wright PetscCall(ISGetLocalSize(stratumIS, &num_points));
323*607e733fSJames Wright PetscCall(ISGetIndices(stratumIS, &points));
324*607e733fSJames Wright for (PetscInt p = 0; p < num_points_serial; p++) {
325*607e733fSJames Wright PetscInt ndof, offset, loc;
326*607e733fSJames Wright
327*607e733fSJames Wright PetscCall(PetscSectionGetDof(pointSerialSection, p, &ndof));
328*607e733fSJames Wright PetscCall(PetscSectionGetOffset(pointSerialSection, p, &offset));
329*607e733fSJames Wright PetscCall(PetscFindInt(p, num_points, points, &loc));
330*607e733fSJames Wright PetscBool serial_has_point = loc >= 0;
331*607e733fSJames Wright
332*607e733fSJames Wright for (PetscInt d = 0; d < ndof; d++) {
333*607e733fSJames Wright if (serial_has_point != pointMaskSerial[offset + d]) PetscCall(PetscPrintf(comm, "DMLabel '%s': Serial and global DM disagree on point %" PetscInt_FMT " valid for label value %" PetscInt_FMT "\n", label_name, p, value));
334*607e733fSJames Wright }
335*607e733fSJames Wright }
336*607e733fSJames Wright PetscCall(ISRestoreIndices(stratumIS, &points));
337*607e733fSJames Wright PetscCall(ISDestroy(&stratumIS));
338*607e733fSJames Wright }
339*607e733fSJames Wright }
340*607e733fSJames Wright }
341*607e733fSJames Wright if (serial2loadPointSF && !*serial2loadPointSF) *serial2loadPointSF = serial2loadPointSF_;
342*607e733fSJames Wright else PetscCall(PetscSFDestroy(&serial2loadPointSF_));
343*607e733fSJames Wright
344*607e733fSJames Wright PetscCall(ISRestoreIndices(loadValuesIS, &values_global));
345*607e733fSJames Wright PetscCall(ISDestroy(&loadValuesIS));
346*607e733fSJames Wright PetscCall(PetscSectionDestroy(&pointSerialSection));
347*607e733fSJames Wright PetscCall(PetscFree(pointMaskSerial));
348*607e733fSJames Wright PetscCall(PetscFree(pointMask));
349*607e733fSJames Wright PetscFunctionReturn(PETSC_SUCCESS);
350*607e733fSJames Wright }
351*607e733fSJames Wright
main(int argc,char ** argv)352*607e733fSJames Wright int main(int argc, char **argv)
353*607e733fSJames Wright {
354*607e733fSJames Wright AppCtx user;
355*607e733fSJames Wright MPI_Comm comm;
356*607e733fSJames Wright PetscMPIInt gsize, grank, mycolor;
357*607e733fSJames Wright PetscBool flg;
358*607e733fSJames Wright const char *infilename;
359*607e733fSJames Wright DM dm_serial = NULL;
360*607e733fSJames Wright Vec V_serial = NULL;
361*607e733fSJames Wright
362*607e733fSJames Wright PetscFunctionBeginUser;
363*607e733fSJames Wright PetscCall(PetscInitialize(&argc, &argv, NULL, help));
364*607e733fSJames Wright PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
365*607e733fSJames Wright PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &gsize));
366*607e733fSJames Wright PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &grank));
367*607e733fSJames Wright
368*607e733fSJames Wright { // Read infile in serial
369*607e733fSJames Wright PetscViewer viewer;
370*607e733fSJames Wright PetscMPIInt gsize_serial;
371*607e733fSJames Wright
372*607e733fSJames Wright mycolor = grank == 0 ? 0 : 1;
373*607e733fSJames Wright PetscCallMPI(MPI_Comm_split(PETSC_COMM_WORLD, mycolor, grank, &comm));
374*607e733fSJames Wright
375*607e733fSJames Wright if (grank == 0) {
376*607e733fSJames Wright PetscCallMPI(MPI_Comm_size(comm, &gsize_serial));
377*607e733fSJames Wright
378*607e733fSJames Wright PetscCall(ReadCGNSDM(comm, user.infile, &dm_serial));
379*607e733fSJames Wright PetscCall(DMSetOptionsPrefix(dm_serial, "serial_"));
380*607e733fSJames Wright PetscCall(DMSetFromOptions(dm_serial));
381*607e733fSJames Wright
382*607e733fSJames Wright /* We just test/demonstrate DM is indeed distributed - unneeded in the application code */
383*607e733fSJames Wright PetscCall(DMPlexIsDistributed(dm_serial, &flg));
384*607e733fSJames Wright PetscCall(PetscPrintf(comm, "Loaded mesh distributed? %s\n", PetscBools[flg]));
385*607e733fSJames Wright
386*607e733fSJames Wright PetscCall(DMViewFromOptions(dm_serial, NULL, "-dm_view"));
387*607e733fSJames Wright PetscCall(PetscViewerCGNSOpen(comm, user.infile, FILE_MODE_READ, &viewer));
388*607e733fSJames Wright PetscCall(DMGetGlobalVector(dm_serial, &V_serial));
389*607e733fSJames Wright PetscCall(VecLoad(V_serial, viewer));
390*607e733fSJames Wright PetscCall(PetscViewerDestroy(&viewer));
391*607e733fSJames Wright
392*607e733fSJames Wright PetscCallMPI(MPI_Comm_free(&comm));
393*607e733fSJames Wright }
394*607e733fSJames Wright PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
395*607e733fSJames Wright }
396*607e733fSJames Wright
397*607e733fSJames Wright for (PetscInt i = 0; i < user.ntimes; i++) {
398*607e733fSJames Wright if (i == 0) {
399*607e733fSJames Wright /* Use infile for the initial load */
400*607e733fSJames Wright infilename = user.infile;
401*607e733fSJames Wright } else {
402*607e733fSJames Wright /* Use outfile for all I/O except the very initial load */
403*607e733fSJames Wright infilename = user.outfile;
404*607e733fSJames Wright }
405*607e733fSJames Wright
406*607e733fSJames Wright if (user.heterogeneous) {
407*607e733fSJames Wright mycolor = (PetscMPIInt)(grank < (gsize - i) ? 0 : 1);
408*607e733fSJames Wright } else {
409*607e733fSJames Wright mycolor = (PetscMPIInt)0;
410*607e733fSJames Wright }
411*607e733fSJames Wright PetscCallMPI(MPI_Comm_split(PETSC_COMM_WORLD, mycolor, grank, &comm));
412*607e733fSJames Wright
413*607e733fSJames Wright if (mycolor == 0) {
414*607e733fSJames Wright /* Load/Save only on processes with mycolor == 0 */
415*607e733fSJames Wright DM dm;
416*607e733fSJames Wright Vec V;
417*607e733fSJames Wright PetscViewer viewer;
418*607e733fSJames Wright PetscMPIInt comm_size;
419*607e733fSJames Wright const char *name;
420*607e733fSJames Wright PetscReal time;
421*607e733fSJames Wright PetscBool set;
422*607e733fSJames Wright
423*607e733fSJames Wright PetscCallMPI(MPI_Comm_size(comm, &comm_size));
424*607e733fSJames Wright PetscCall(PetscPrintf(comm, "Begin cycle %" PetscInt_FMT ", comm size %d\n", i, comm_size));
425*607e733fSJames Wright
426*607e733fSJames Wright // Load DM from CGNS file
427*607e733fSJames Wright PetscCall(ReadCGNSDM(comm, infilename, &dm));
428*607e733fSJames Wright PetscCall(DMSetOptionsPrefix(dm, "loaded_"));
429*607e733fSJames Wright PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
430*607e733fSJames Wright
431*607e733fSJames Wright // Load solution from CGNS file
432*607e733fSJames Wright PetscCall(PetscViewerCGNSOpen(comm, infilename, FILE_MODE_READ, &viewer));
433*607e733fSJames Wright PetscCall(DMGetGlobalVector(dm, &V));
434*607e733fSJames Wright PetscCall(PetscViewerCGNSSetSolutionIndex(viewer, 1));
435*607e733fSJames Wright { // Test GetSolutionIndex, not needed in application code
436*607e733fSJames Wright PetscInt solution_index;
437*607e733fSJames Wright PetscCall(PetscViewerCGNSGetSolutionIndex(viewer, &solution_index));
438*607e733fSJames Wright PetscCheck(solution_index == 1, comm, PETSC_ERR_ARG_INCOMP, "Returned solution index wrong.");
439*607e733fSJames Wright }
440*607e733fSJames Wright PetscCall(PetscViewerCGNSGetSolutionName(viewer, &name));
441*607e733fSJames Wright PetscCall(PetscViewerCGNSGetSolutionTime(viewer, &time, &set));
442*607e733fSJames Wright PetscCheck(set, comm, PETSC_ERR_RETURN, "Time data wasn't set!");
443*607e733fSJames Wright PetscCall(PetscPrintf(comm, "Solution Name: %s, and time %g\n", name, time));
444*607e733fSJames Wright PetscCall(VecLoad(V, viewer));
445*607e733fSJames Wright PetscCall(PetscViewerDestroy(&viewer));
446*607e733fSJames Wright
447*607e733fSJames Wright // Verify loaded solution against serial solution
448*607e733fSJames Wright PetscCall(VerifyLoadedSolution(dm_serial, V_serial, dm, V, 100 * PETSC_MACHINE_EPSILON));
449*607e733fSJames Wright
450*607e733fSJames Wright // Verify DMLabel values against the serial DM
451*607e733fSJames Wright PetscCall(VerifyDMLabels(dm_serial, dm, "Face Sets", NULL));
452*607e733fSJames Wright
453*607e733fSJames Wright { // Complete the label so that the writer must sort through non-face points
454*607e733fSJames Wright DMLabel label;
455*607e733fSJames Wright PetscCall(DMGetLabel(dm, "Face Sets", &label));
456*607e733fSJames Wright PetscCall(DMPlexLabelComplete(dm, label));
457*607e733fSJames Wright }
458*607e733fSJames Wright
459*607e733fSJames Wright // Write loaded solution to CGNS file
460*607e733fSJames Wright PetscCall(PetscViewerCGNSOpen(comm, user.outfile, FILE_MODE_WRITE, &viewer));
461*607e733fSJames Wright PetscCall(VecView(V, viewer));
462*607e733fSJames Wright PetscCall(PetscViewerDestroy(&viewer));
463*607e733fSJames Wright
464*607e733fSJames Wright PetscCall(DMRestoreGlobalVector(dm, &V));
465*607e733fSJames Wright PetscCall(DMDestroy(&dm));
466*607e733fSJames Wright PetscCall(PetscPrintf(comm, "End cycle %" PetscInt_FMT "\n--------\n", i));
467*607e733fSJames Wright }
468*607e733fSJames Wright PetscCallMPI(MPI_Comm_free(&comm));
469*607e733fSJames Wright PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
470*607e733fSJames Wright }
471*607e733fSJames Wright
472*607e733fSJames Wright if (V_serial && dm_serial) PetscCall(DMRestoreGlobalVector(dm_serial, &V_serial));
473*607e733fSJames Wright if (dm_serial) PetscCall(DMDestroy(&dm_serial));
474*607e733fSJames Wright PetscCall(PetscFinalize());
475*607e733fSJames Wright return 0;
476*607e733fSJames Wright }
477*607e733fSJames Wright
478*607e733fSJames Wright /*TEST
479*607e733fSJames Wright build:
480*607e733fSJames Wright requires: cgns
481*607e733fSJames Wright
482*607e733fSJames Wright testset:
483*607e733fSJames Wright suffix: cgns_3x3
484*607e733fSJames Wright requires: !complex
485*607e733fSJames Wright nsize: 4
486*607e733fSJames Wright args: -infile ${DATAFILESPATH}/meshes/3x3_Q1.cgns -outfile 3x3_Q1_output.cgns -dm_plex_cgns_parallel -ntimes 3 -heterogeneous -serial_dm_view -loaded_dm_view -redistributed_dm_distribute
487*607e733fSJames Wright test:
488*607e733fSJames Wright # this partitioner should not shuffle anything, it should yield the same partitioning as the XDMF reader - added just for testing
489*607e733fSJames Wright suffix: simple
490*607e733fSJames Wright args: -petscpartitioner_type simple
491*607e733fSJames Wright test:
492*607e733fSJames Wright suffix: parmetis
493*607e733fSJames Wright requires: parmetis
494*607e733fSJames Wright args: -petscpartitioner_type parmetis
495*607e733fSJames Wright test:
496*607e733fSJames Wright suffix: ptscotch
497*607e733fSJames Wright requires: ptscotch
498*607e733fSJames Wright args: -petscpartitioner_type ptscotch
499*607e733fSJames Wright
500*607e733fSJames Wright testset:
501*607e733fSJames Wright suffix: cgns_2x2x2
502*607e733fSJames Wright requires: !complex
503*607e733fSJames Wright nsize: 4
504*607e733fSJames Wright args: -infile ${DATAFILESPATH}/meshes/2x2x2_Q3_wave.cgns -outfile 2x2x2_Q3_wave_output.cgns -dm_plex_cgns_parallel -ntimes 3 -heterogeneous -serial_dm_view -loaded_dm_view -redistributed_dm_distribute
505*607e733fSJames Wright test:
506*607e733fSJames Wright # this partitioner should not shuffle anything, it should yield the same partitioning as the XDMF reader - added just for testing
507*607e733fSJames Wright suffix: simple
508*607e733fSJames Wright args: -petscpartitioner_type simple
509*607e733fSJames Wright test:
510*607e733fSJames Wright suffix: parmetis
511*607e733fSJames Wright requires: parmetis
512*607e733fSJames Wright args: -petscpartitioner_type parmetis
513*607e733fSJames Wright test:
514*607e733fSJames Wright suffix: ptscotch
515*607e733fSJames Wright requires: ptscotch
516*607e733fSJames Wright args: -petscpartitioner_type ptscotch
517*607e733fSJames Wright
518*607e733fSJames Wright test:
519*607e733fSJames Wright # This file is meant to explicitly have a special case where a partition completely surrounds a boundary element, but does not own it
520*607e733fSJames Wright requires: !complex
521*607e733fSJames Wright suffix: cgns_3x3_2
522*607e733fSJames Wright args: -infile ${DATAFILESPATH}/meshes/3x3_Q1.cgns -outfile 3x3_Q1_output.cgns -dm_plex_cgns_parallel -serial_dm_view -loaded_dm_view -redistributed_dm_distribute -petscpartitioner_type simple
523*607e733fSJames Wright TEST*/
524