xref: /petsc/src/dm/impls/plex/tests/ex26.c (revision 2fa40bb9206b96114faa7cb222621ec184d31cd2)
1 static char help[] = "Test FEM layout with DM and ExodusII storage\n\n";
2 
3 /*
4   In order to see the vectors which are being tested, use
5 
6      -ua_vec_view -s_vec_view
7 */
8 
9 #include <petsc.h>
10 #include <exodusII.h>
11 
12 #include <petsc/private/dmpleximpl.h>
13 
14 int main(int argc, char **argv) {
15   DM                dm, dmU, dmA, dmS, dmUA, dmUA2, *dmList;
16   Vec               X, U, A, S, UA, UA2;
17   IS                isU, isA, isS, isUA;
18   PetscSection      section;
19   const PetscInt    fieldU = 0;
20   const PetscInt    fieldA = 2;
21   const PetscInt    fieldS = 1;
22   const PetscInt    fieldUA[2] = {0, 2};
23   char              ifilename[PETSC_MAX_PATH_LEN], ofilename[PETSC_MAX_PATH_LEN];
24   int               exoid = -1;
25   IS                csIS;
26   const PetscInt   *csID;
27   PetscInt         *pStartDepth, *pEndDepth;
28   PetscInt          order = 1;
29   PetscInt          sdim, d, pStart, pEnd, p, numCS, set;
30   PetscMPIInt       rank, size;
31   PetscViewer       viewer;
32   PetscErrorCode    ierr;
33 
34   ierr = PetscInitialize(&argc, &argv,NULL, help);if (ierr) return ierr;
35   ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRMPI(ierr);
36   ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRMPI(ierr);
37   ierr = PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FEM Layout Options", "ex26");CHKERRQ(ierr);
38   ierr = PetscOptionsString("-i", "Filename to read", "ex26", ifilename, ifilename, sizeof(ifilename), NULL);CHKERRQ(ierr);
39   ierr = PetscOptionsString("-o", "Filename to write", "ex26", ofilename, ofilename, sizeof(ofilename), NULL);CHKERRQ(ierr);
40   ierr = PetscOptionsBoundedInt("-order", "FEM polynomial order", "ex26", order, &order, NULL,1);CHKERRQ(ierr);
41   ierr = PetscOptionsEnd();CHKERRQ(ierr);
42   if ((order > 2) || (order < 1)) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported polynomial order %D not in [1, 2]", order);
43 
44   /* Read the mesh from a file in any supported format */
45   ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename, PETSC_TRUE, &dm);CHKERRQ(ierr);
46   ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
47   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
48   ierr = DMGetDimension(dm, &sdim);CHKERRQ(ierr);
49 
50   /* Create the exodus result file */
51   {
52     PetscInt      numstep = 3, step;
53     char         *nodalVarName[4];
54     char         *zonalVarName[6];
55     int          *truthtable;
56     PetscInt      numNodalVar, numZonalVar, i;
57 
58     /* enable exodus debugging information */
59     ex_opts(EX_VERBOSE|EX_DEBUG);
60     /* Create the exodus file */
61     ierr = PetscViewerExodusIIOpen(PETSC_COMM_WORLD,ofilename,FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
62     /* The long way would be */
63     /*
64       ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
65       ierr = PetscViewerSetType(viewer,PETSCVIEWEREXODUSII);CHKERRQ(ierr);
66       ierr = PetscViewerFileSetMode(viewer,FILE_MODE_APPEND);CHKERRQ(ierr);
67       ierr = PetscViewerFileSetName(viewer,ofilename);CHKERRQ(ierr);
68     */
69     /* set the mesh order */
70     ierr = PetscViewerExodusIISetOrder(viewer,order);CHKERRQ(ierr);
71     ierr = PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
72     /*
73       Notice how the exodus file is actually NOT open at this point (exoid is -1)
74       Since we are overwritting the file (mode is FILE_MODE_WRITE), we are going to have to
75       write the geometry (the DM), which can only be done on a brand new file.
76     */
77 
78     /* Save the geometry to the file, erasing all previous content */
79     ierr = DMView(dm,viewer);CHKERRQ(ierr);
80     ierr = PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
81     /*
82       Note how the exodus file is now open
83     */
84 
85     /* "Format" the exodus result file, i.e. allocate space for nodal and zonal variables */
86     switch (sdim) {
87     case 2:
88       numNodalVar = 3;
89       nodalVarName[0] = (char *) "U_x";
90       nodalVarName[1] = (char *) "U_y";
91       nodalVarName[2] = (char *) "Alpha";
92       numZonalVar = 3;
93       zonalVarName[0] = (char *) "Sigma_11";
94       zonalVarName[1] = (char *) "Sigma_22";
95       zonalVarName[2] = (char *) "Sigma_12";
96       break;
97     case 3:
98       numNodalVar = 4;
99       nodalVarName[0] = (char *) "U_x";
100       nodalVarName[1] = (char *) "U_y";
101       nodalVarName[2] = (char *) "U_z";
102       nodalVarName[3] = (char *) "Alpha";
103       numZonalVar = 6;
104       zonalVarName[0] = (char *) "Sigma_11";
105       zonalVarName[1] = (char *) "Sigma_22";
106       zonalVarName[2] = (char *) "Sigma_33";
107       zonalVarName[3] = (char *) "Sigma_23";
108       zonalVarName[4] = (char *) "Sigma_13";
109       zonalVarName[5] = (char *) "Sigma_12";
110       break;
111     default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No layout for dimension %D", sdim);
112     }
113     ierr = PetscViewerExodusIIGetId(viewer,&exoid);CHKERRQ(ierr);
114     PetscStackCallStandard(ex_put_variable_param,(exoid, EX_ELEM_BLOCK, numZonalVar));
115     PetscStackCallStandard(ex_put_variable_names,(exoid, EX_ELEM_BLOCK, numZonalVar, zonalVarName));
116     PetscStackCallStandard(ex_put_variable_param,(exoid, EX_NODAL, numNodalVar));
117     PetscStackCallStandard(ex_put_variable_names,(exoid, EX_NODAL, numNodalVar, nodalVarName));
118     numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK);
119 
120     /*
121       An exodusII truth table specifies which fields are saved at which time step
122       It speeds up I/O but reserving space for fieldsin the file ahead of time.
123     */
124     ierr = PetscMalloc1(numZonalVar * numCS, &truthtable);CHKERRQ(ierr);
125     for (i = 0; i < numZonalVar * numCS; ++i) truthtable[i] = 1;
126     PetscStackCallStandard(ex_put_truth_table,(exoid, EX_ELEM_BLOCK, numCS, numZonalVar, truthtable));
127     ierr = PetscFree(truthtable);CHKERRQ(ierr);
128 
129     /* Writing time step information in the file. Note that this is currently broken in the exodus library for netcdf4 (HDF5-based) files */
130     for (step = 0; step < numstep; ++step) {
131       PetscReal time = step;
132       PetscStackCallStandard(ex_put_time,(exoid, step+1, &time));
133     }
134   }
135 
136   /* Create the main section containing all fields */
137   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &section);CHKERRQ(ierr);
138   ierr = PetscSectionSetNumFields(section, 3);CHKERRQ(ierr);
139   ierr = PetscSectionSetFieldName(section, fieldU, "U");CHKERRQ(ierr);
140   ierr = PetscSectionSetFieldName(section, fieldA, "Alpha");CHKERRQ(ierr);
141   ierr = PetscSectionSetFieldName(section, fieldS, "Sigma");CHKERRQ(ierr);
142   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
143   ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr);
144   ierr = PetscMalloc2(sdim+1, &pStartDepth, sdim+1, &pEndDepth);CHKERRQ(ierr);
145   for (d = 0; d <= sdim; ++d) {ierr = DMPlexGetDepthStratum(dm, d, &pStartDepth[d], &pEndDepth[d]);CHKERRQ(ierr);}
146   /* Vector field U, Scalar field Alpha, Tensor field Sigma */
147   ierr = PetscSectionSetFieldComponents(section, fieldU, sdim);CHKERRQ(ierr);
148   ierr = PetscSectionSetFieldComponents(section, fieldA, 1);CHKERRQ(ierr);
149   ierr = PetscSectionSetFieldComponents(section, fieldS, sdim*(sdim+1)/2);CHKERRQ(ierr);
150 
151   /* Going through cell sets then cells, and setting up storage for the sections */
152   ierr = DMGetLabelSize(dm, "Cell Sets", &numCS);CHKERRQ(ierr);
153   ierr = DMGetLabelIdIS(dm, "Cell Sets", &csIS);CHKERRQ(ierr);
154   if (csIS) {ierr = ISGetIndices(csIS, &csID);CHKERRQ(ierr);}
155   for (set = 0; set < numCS; set++) {
156     IS                cellIS;
157     const PetscInt   *cellID;
158     PetscInt          numCells, cell, closureSize, *closureA = NULL;
159 
160     ierr = DMGetStratumSize(dm, "Cell Sets", csID[set], &numCells);CHKERRQ(ierr);
161     ierr = DMGetStratumIS(dm, "Cell Sets", csID[set], &cellIS);CHKERRQ(ierr);
162     if (numCells > 0) {
163       /* dof layout ordered by increasing height in the DAG: cell, face, edge, vertex */
164       PetscInt          dofUP1Tri[]  = {2, 0, 0};
165       PetscInt          dofAP1Tri[]  = {1, 0, 0};
166       PetscInt          dofUP2Tri[]  = {2, 2, 0};
167       PetscInt          dofAP2Tri[]  = {1, 1, 0};
168       PetscInt          dofUP1Quad[] = {2, 0, 0};
169       PetscInt          dofAP1Quad[] = {1, 0, 0};
170       PetscInt          dofUP2Quad[] = {2, 2, 2};
171       PetscInt          dofAP2Quad[] = {1, 1, 1};
172       PetscInt          dofS2D[]     = {0, 0, 3};
173       PetscInt          dofUP1Tet[]  = {3, 0, 0, 0};
174       PetscInt          dofAP1Tet[]  = {1, 0, 0, 0};
175       PetscInt          dofUP2Tet[]  = {3, 3, 0, 0};
176       PetscInt          dofAP2Tet[]  = {1, 1, 0, 0};
177       PetscInt          dofUP1Hex[]  = {3, 0, 0, 0};
178       PetscInt          dofAP1Hex[]  = {1, 0, 0, 0};
179       PetscInt          dofUP2Hex[]  = {3, 3, 3, 3};
180       PetscInt          dofAP2Hex[]  = {1, 1, 1, 1};
181       PetscInt          dofS3D[]     = {0, 0, 0, 6};
182       PetscInt         *dofU, *dofA, *dofS;
183 
184       switch (sdim) {
185       case 2: dofS = dofS2D;break;
186       case 3: dofS = dofS3D;break;
187       default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No layout for dimension %D", sdim);
188       }
189 
190       /* Identify cell type based on closure size only. This works for Tri/Tet/Quad/Hex meshes
191          It will not be enough to identify more exotic elements like pyramid or prisms...  */
192       ierr = ISGetIndices(cellIS, &cellID);CHKERRQ(ierr);
193       ierr = DMPlexGetTransitiveClosure(dm, cellID[0], PETSC_TRUE, &closureSize, &closureA);CHKERRQ(ierr);
194       switch (closureSize) {
195         case 7: /* Tri */
196         if (order == 1) {
197           dofU = dofUP1Tri;
198           dofA = dofAP1Tri;
199         } else {
200           dofU = dofUP2Tri;
201           dofA = dofAP2Tri;
202         }
203         break;
204         case 9: /* Quad */
205         if (order == 1) {
206           dofU = dofUP1Quad;
207           dofA = dofAP1Quad;
208         } else {
209           dofU = dofUP2Quad;
210           dofA = dofAP2Quad;
211         }
212         break;
213         case 15: /* Tet */
214         if (order == 1) {
215           dofU = dofUP1Tet;
216           dofA = dofAP1Tet;
217         } else {
218           dofU = dofUP2Tet;
219           dofA = dofAP2Tet;
220         }
221         break;
222         case 27: /* Hex */
223         if (order == 1) {
224           dofU = dofUP1Hex;
225           dofA = dofAP1Hex;
226         } else {
227           dofU = dofUP2Hex;
228           dofA = dofAP2Hex;
229         }
230         break;
231         default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Unknown element with closure size %D\n", closureSize);
232       }
233       ierr = DMPlexRestoreTransitiveClosure(dm, cellID[0], PETSC_TRUE, &closureSize, &closureA);CHKERRQ(ierr);
234 
235       for (cell = 0; cell < numCells; cell++) {
236         PetscInt *closure = NULL;
237 
238         ierr = DMPlexGetTransitiveClosure(dm, cellID[cell], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
239         for (p = 0; p < closureSize; ++p) {
240           /* Find depth of p */
241           for (d = 0; d <= sdim; ++d) {
242             if ((closure[2*p] >= pStartDepth[d]) && (closure[2*p] < pEndDepth[d])) {
243               ierr = PetscSectionSetDof(section, closure[2*p], dofU[d]+dofA[d]+dofS[d]);CHKERRQ(ierr);
244               ierr = PetscSectionSetFieldDof(section, closure[2*p], fieldU, dofU[d]);CHKERRQ(ierr);
245               ierr = PetscSectionSetFieldDof(section, closure[2*p], fieldA, dofA[d]);CHKERRQ(ierr);
246               ierr = PetscSectionSetFieldDof(section, closure[2*p], fieldS, dofS[d]);CHKERRQ(ierr);
247             }
248           }
249         }
250         ierr = DMPlexRestoreTransitiveClosure(dm, cellID[cell], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
251       }
252       ierr = ISRestoreIndices(cellIS, &cellID);CHKERRQ(ierr);
253       ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
254     }
255   }
256   if (csIS) {ierr = ISRestoreIndices(csIS, &csID);CHKERRQ(ierr);}
257   ierr = ISDestroy(&csIS);CHKERRQ(ierr);
258   ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
259   ierr = DMSetLocalSection(dm, section);CHKERRQ(ierr);
260   ierr = PetscObjectViewFromOptions((PetscObject) section, NULL, "-dm_section_view");CHKERRQ(ierr);
261   ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
262 
263   {
264     DM               pdm;
265     PetscSF          migrationSF;
266     PetscInt         ovlp = 0;
267     PetscPartitioner part;
268 
269     ierr = DMSetUseNatural(dm,PETSC_TRUE);CHKERRQ(ierr);
270     ierr = DMPlexGetPartitioner(dm,&part);CHKERRQ(ierr);
271     ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
272     ierr = DMPlexDistribute(dm,ovlp,&migrationSF,&pdm);CHKERRQ(ierr);
273     if (pdm) {
274       ierr = DMPlexSetMigrationSF(pdm,migrationSF);CHKERRQ(ierr);
275       ierr = PetscSFDestroy(&migrationSF);CHKERRQ(ierr);
276       ierr = DMDestroy(&dm);CHKERRQ(ierr);
277       dm = pdm;
278       ierr = DMViewFromOptions(dm,NULL,"-dm_view");CHKERRQ(ierr);
279     }
280   }
281 
282   /* Get DM and IS for each field of dm */
283   ierr = DMCreateSubDM(dm, 1, &fieldU, &isU,  &dmU);CHKERRQ(ierr);
284   ierr = DMCreateSubDM(dm, 1, &fieldA, &isA,  &dmA);CHKERRQ(ierr);
285   ierr = DMCreateSubDM(dm, 1, &fieldS, &isS,  &dmS);CHKERRQ(ierr);
286   ierr = DMCreateSubDM(dm, 2, fieldUA, &isUA, &dmUA);CHKERRQ(ierr);
287 
288   ierr = PetscMalloc1(2,&dmList);CHKERRQ(ierr);
289   dmList[0] = dmU;
290   dmList[1] = dmA;
291   /* We temporarily disable dmU->useNatural to test that we can reconstruct the
292      NaturaltoGlobal SF from any of the dm in dms
293   */
294   dmU->useNatural = PETSC_FALSE;
295   ierr = DMCreateSuperDM(dmList,2,NULL,&dmUA2);CHKERRQ(ierr);
296   dmU->useNatural = PETSC_TRUE;
297   ierr = PetscFree(dmList);CHKERRQ(ierr);
298 
299   ierr = DMGetGlobalVector(dm,   &X);CHKERRQ(ierr);
300   ierr = DMGetGlobalVector(dmU,  &U);CHKERRQ(ierr);
301   ierr = DMGetGlobalVector(dmA,  &A);CHKERRQ(ierr);
302   ierr = DMGetGlobalVector(dmS,  &S);CHKERRQ(ierr);
303   ierr = DMGetGlobalVector(dmUA, &UA);CHKERRQ(ierr);
304   ierr = DMGetGlobalVector(dmUA2, &UA2);CHKERRQ(ierr);
305 
306   ierr = PetscObjectSetName((PetscObject) U,  "U");CHKERRQ(ierr);
307   ierr = PetscObjectSetName((PetscObject) A,  "Alpha");CHKERRQ(ierr);
308   ierr = PetscObjectSetName((PetscObject) S,  "Sigma");CHKERRQ(ierr);
309   ierr = PetscObjectSetName((PetscObject) UA, "UAlpha");CHKERRQ(ierr);
310   ierr = PetscObjectSetName((PetscObject) UA2, "UAlpha2");CHKERRQ(ierr);
311   ierr = VecSet(X, -111.);CHKERRQ(ierr);
312 
313   /* Setting u to [x,y,z]  and alpha to x^2+y^2+z^2 by writing in UAlpha then restricting to U and Alpha */
314   {
315     PetscSection sectionUA;
316     Vec          UALoc;
317     PetscSection coordSection;
318     Vec          coord;
319     PetscScalar *cval, *xyz;
320     PetscInt     clSize, i, j;
321 
322     ierr = DMGetLocalSection(dmUA, &sectionUA);CHKERRQ(ierr);
323     ierr = DMGetLocalVector(dmUA, &UALoc);CHKERRQ(ierr);
324     ierr = VecGetArray(UALoc, &cval);CHKERRQ(ierr);
325     ierr = DMGetCoordinateSection(dmUA, &coordSection);CHKERRQ(ierr);
326     ierr = DMGetCoordinatesLocal(dmUA, &coord);CHKERRQ(ierr);
327     ierr = DMPlexGetChart(dmUA, &pStart, &pEnd);CHKERRQ(ierr);
328 
329     for (p = pStart; p < pEnd; ++p) {
330       PetscInt dofUA, offUA;
331 
332       ierr = PetscSectionGetDof(sectionUA, p, &dofUA);CHKERRQ(ierr);
333       if (dofUA > 0) {
334         xyz=NULL;
335         ierr = PetscSectionGetOffset(sectionUA, p, &offUA);CHKERRQ(ierr);
336         ierr = DMPlexVecGetClosure(dmUA, coordSection, coord, p, &clSize, &xyz);CHKERRQ(ierr);
337         cval[offUA+sdim] = 0.;
338         for (i = 0; i < sdim; ++i) {
339           cval[offUA+i] = 0;
340           for (j = 0; j < clSize/sdim; ++j) {
341             cval[offUA+i] += xyz[j*sdim+i];
342           }
343           cval[offUA+i] = cval[offUA+i] * sdim / clSize;
344           cval[offUA+sdim] += PetscSqr(cval[offUA+i]);
345         }
346         ierr = DMPlexVecRestoreClosure(dmUA, coordSection, coord, p, &clSize, &xyz);CHKERRQ(ierr);
347       }
348     }
349     ierr = VecRestoreArray(UALoc, &cval);CHKERRQ(ierr);
350     ierr = DMLocalToGlobalBegin(dmUA, UALoc, INSERT_VALUES, UA);CHKERRQ(ierr);
351     ierr = DMLocalToGlobalEnd(dmUA, UALoc, INSERT_VALUES, UA);CHKERRQ(ierr);
352     ierr = DMRestoreLocalVector(dmUA, &UALoc);CHKERRQ(ierr);
353 
354     /* Update X */
355     ierr = VecISCopy(X, isUA, SCATTER_FORWARD, UA);CHKERRQ(ierr);
356     ierr = VecViewFromOptions(UA, NULL, "-ua_vec_view");CHKERRQ(ierr);
357 
358     /* Restrict to U and Alpha */
359     ierr = VecISCopy(X, isU, SCATTER_REVERSE, U);CHKERRQ(ierr);
360     ierr = VecISCopy(X, isA, SCATTER_REVERSE, A);CHKERRQ(ierr);
361 
362     /* restrict to UA2 */
363     ierr = VecISCopy(X, isUA, SCATTER_REVERSE, UA2);CHKERRQ(ierr);
364     ierr = VecViewFromOptions(UA2, NULL, "-ua2_vec_view");CHKERRQ(ierr);
365   }
366 
367   {
368     Vec          tmpVec;
369     PetscSection coordSection;
370     Vec          coord;
371     PetscReal    norm;
372     PetscReal    time = 1.234;
373 
374     /* Writing nodal variables to ExodusII file */
375     ierr = DMSetOutputSequenceNumber(dmU,0,time);CHKERRQ(ierr);
376     ierr = DMSetOutputSequenceNumber(dmA,0,time);CHKERRQ(ierr);
377 
378     ierr = VecView(U, viewer);CHKERRQ(ierr);
379     ierr = VecView(A, viewer);CHKERRQ(ierr);
380 
381     /* Saving U and Alpha in one shot.
382        For this, we need to cheat and change the Vec's name
383        Note that in the end we write variables one component at a time,
384        so that there is no real values in doing this
385     */
386 
387     ierr = DMSetOutputSequenceNumber(dmUA,1,time);CHKERRQ(ierr);
388     ierr = DMGetGlobalVector(dmUA, &tmpVec);CHKERRQ(ierr);
389     ierr = VecCopy(UA, tmpVec);CHKERRQ(ierr);
390     ierr = PetscObjectSetName((PetscObject) tmpVec, "U");CHKERRQ(ierr);
391     ierr = VecView(tmpVec, viewer);CHKERRQ(ierr);
392     /* Reading nodal variables in Exodus file */
393     ierr = VecSet(tmpVec, -1000.0);CHKERRQ(ierr);
394     ierr = VecLoad(tmpVec, viewer);CHKERRQ(ierr);
395     ierr = VecAXPY(UA, -1.0, tmpVec);CHKERRQ(ierr);
396     ierr = VecNorm(UA, NORM_INFINITY, &norm);CHKERRQ(ierr);
397     if (norm > PETSC_SQRT_MACHINE_EPSILON) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "UAlpha ||Vin - Vout|| = %g\n", (double) norm);
398     ierr = DMRestoreGlobalVector(dmUA, &tmpVec);CHKERRQ(ierr);
399 
400     /* same thing with the UA2 Vec obtained from the superDM */
401     ierr = DMGetGlobalVector(dmUA2, &tmpVec);CHKERRQ(ierr);
402     ierr = VecCopy(UA2, tmpVec);CHKERRQ(ierr);
403     ierr = PetscObjectSetName((PetscObject) tmpVec, "U");CHKERRQ(ierr);
404     ierr = DMSetOutputSequenceNumber(dmUA2,2,time);CHKERRQ(ierr);
405     ierr = VecView(tmpVec, viewer);CHKERRQ(ierr);
406     /* Reading nodal variables in Exodus file */
407     ierr = VecSet(tmpVec, -1000.0);CHKERRQ(ierr);
408     ierr = VecLoad(tmpVec,viewer);CHKERRQ(ierr);
409     ierr = VecAXPY(UA2, -1.0, tmpVec);CHKERRQ(ierr);
410     ierr = VecNorm(UA2, NORM_INFINITY, &norm);CHKERRQ(ierr);
411     if (norm > PETSC_SQRT_MACHINE_EPSILON) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "UAlpha2 ||Vin - Vout|| = %g\n", (double) norm);
412     ierr = DMRestoreGlobalVector(dmUA2, &tmpVec);CHKERRQ(ierr);
413 
414     /* Building and saving Sigma
415        We set sigma_0 = rank (to see partitioning)
416               sigma_1 = cell set ID
417               sigma_2 = x_coordinate of the cell center of mass
418     */
419     ierr = DMGetCoordinateSection(dmS, &coordSection);CHKERRQ(ierr);
420     ierr = DMGetCoordinatesLocal(dmS, &coord);CHKERRQ(ierr);
421     ierr = DMGetLabelIdIS(dmS, "Cell Sets", &csIS);CHKERRQ(ierr);
422     ierr = DMGetLabelSize(dmS, "Cell Sets", &numCS);CHKERRQ(ierr);
423     ierr = ISGetIndices(csIS, &csID);CHKERRQ(ierr);
424     for (set = 0; set < numCS; ++set) {
425       /* We know that all cells in a cell set have the same type, so we can dimension cval and xyz once for each cell set */
426       IS              cellIS;
427       const PetscInt *cellID;
428       PetscInt        numCells, cell;
429       PetscScalar    *cval = NULL, *xyz  = NULL;
430       PetscInt        clSize, cdimCoord, c;
431 
432       ierr = DMGetStratumIS(dmS, "Cell Sets", csID[set], &cellIS);CHKERRQ(ierr);
433       ierr = ISGetIndices(cellIS, &cellID);CHKERRQ(ierr);
434       ierr = ISGetSize(cellIS, &numCells);CHKERRQ(ierr);
435       for (cell = 0; cell < numCells; cell++) {
436         ierr = DMPlexVecGetClosure(dmS, NULL, S, cellID[cell], &clSize, &cval);CHKERRQ(ierr);
437         ierr = DMPlexVecGetClosure(dmS, coordSection, coord, cellID[cell], &cdimCoord, &xyz);CHKERRQ(ierr);
438         cval[0] = rank;
439         cval[1] = csID[set];
440         cval[2] = 0.;
441         for (c = 0; c < cdimCoord/sdim; c++) cval[2] += xyz[c*sdim];
442         cval[2] = cval[2] * sdim / cdimCoord;
443         ierr = DMPlexVecSetClosure(dmS, NULL, S, cellID[cell], cval, INSERT_ALL_VALUES);CHKERRQ(ierr);
444       }
445       ierr = DMPlexVecRestoreClosure(dmS, NULL, S, cellID[0], &clSize, &cval);CHKERRQ(ierr);
446       ierr = DMPlexVecRestoreClosure(dmS, coordSection, coord, cellID[0], NULL, &xyz);CHKERRQ(ierr);
447       ierr = ISRestoreIndices(cellIS, &cellID);CHKERRQ(ierr);
448       ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
449     }
450     ierr = ISRestoreIndices(csIS, &csID);CHKERRQ(ierr);
451     ierr = ISDestroy(&csIS);CHKERRQ(ierr);
452     ierr = VecViewFromOptions(S, NULL, "-s_vec_view");CHKERRQ(ierr);
453 
454     /* Writing zonal variables in Exodus file */
455     ierr = DMSetOutputSequenceNumber(dmS,0,time);CHKERRQ(ierr);
456     ierr = VecView(S,viewer);CHKERRQ(ierr);
457 
458     /* Reading zonal variables in Exodus file */
459     ierr = DMGetGlobalVector(dmS, &tmpVec);CHKERRQ(ierr);
460     ierr = VecSet(tmpVec, -1000.0);CHKERRQ(ierr);
461     ierr = PetscObjectSetName((PetscObject) tmpVec, "Sigma");CHKERRQ(ierr);
462     ierr = VecLoad(tmpVec,viewer);CHKERRQ(ierr);
463     ierr = VecAXPY(S, -1.0, tmpVec);CHKERRQ(ierr);
464     ierr = VecNorm(S, NORM_INFINITY, &norm);CHKERRQ(ierr);
465     if (norm > PETSC_SQRT_MACHINE_EPSILON) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "Sigma ||Vin - Vout|| = %g\n", (double) norm);
466     ierr = DMRestoreGlobalVector(dmS, &tmpVec);CHKERRQ(ierr);
467   }
468   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
469 
470   ierr = DMRestoreGlobalVector(dmUA2, &UA2);CHKERRQ(ierr);
471   ierr = DMRestoreGlobalVector(dmUA, &UA);CHKERRQ(ierr);
472   ierr = DMRestoreGlobalVector(dmS,  &S);CHKERRQ(ierr);
473   ierr = DMRestoreGlobalVector(dmA,  &A);CHKERRQ(ierr);
474   ierr = DMRestoreGlobalVector(dmU,  &U);CHKERRQ(ierr);
475   ierr = DMRestoreGlobalVector(dm,   &X);CHKERRQ(ierr);
476   ierr = DMDestroy(&dmU);CHKERRQ(ierr); ierr = ISDestroy(&isU);CHKERRQ(ierr);
477   ierr = DMDestroy(&dmA);CHKERRQ(ierr); ierr = ISDestroy(&isA);CHKERRQ(ierr);
478   ierr = DMDestroy(&dmS);CHKERRQ(ierr); ierr = ISDestroy(&isS);CHKERRQ(ierr);
479   ierr = DMDestroy(&dmUA);CHKERRQ(ierr);ierr = ISDestroy(&isUA);CHKERRQ(ierr);
480   ierr = DMDestroy(&dmUA2);CHKERRQ(ierr);
481   ierr = DMDestroy(&dm);CHKERRQ(ierr);
482   ierr = PetscFree2(pStartDepth, pEndDepth);CHKERRQ(ierr);
483   ierr = PetscFinalize();
484   return ierr;
485 }
486 
487 /*TEST
488 
489   build:
490     requires: exodusii pnetcdf !complex
491   # 2D seq
492   test:
493     suffix: 0
494     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 1
495     #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints
496   test:
497     suffix: 1
498     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 1
499   test:
500     suffix: 2
501     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 1
502     #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints
503   test:
504     suffix: 3
505     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 2
506   test:
507     suffix: 4
508     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 2
509   test:
510     suffix: 5
511     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 2
512 
513   # 2D par
514   test:
515     suffix: 6
516     nsize: 2
517     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 1
518     #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints
519   test:
520     suffix: 7
521     nsize: 2
522     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 1
523   test:
524     suffix: 8
525     nsize: 2
526     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 1
527     #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: invalid dimension ID or name
528   test:
529     suffix: 9
530     nsize: 2
531     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 2
532   test:
533     suffix: 10
534     nsize: 2
535     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 2
536   test:
537     # Something is now broken with parallel read/write for wahtever shape H is
538     TODO: broken
539     suffix: 11
540     nsize: 2
541     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 2
542 
543   #3d seq
544   test:
545     suffix: 12
546     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 1
547   test:
548     suffix: 13
549     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 1
550     #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints
551   test:
552     suffix: 14
553     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 2
554   test:
555     suffix: 15
556     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 2
557     #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints
558   #3d par
559   test:
560     suffix: 16
561     nsize: 2
562     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 1
563   test:
564     suffix: 17
565     nsize: 2
566     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 1
567     #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints
568   test:
569     suffix: 18
570     nsize: 2
571     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 2
572   test:
573     suffix: 19
574     nsize: 2
575     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -dm_section_view -petscpartitioner_type simple -order 2
576     #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: One or more variable sizes violate format constraints
577 
578 TEST*/
579