xref: /petsc/src/dm/impls/plex/tests/ex26.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
1c4762a1bSJed Brown static char help[] = "Test FEM layout with DM and ExodusII storage\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown   In order to see the vectors which are being tested, use
5c4762a1bSJed Brown 
6c4762a1bSJed Brown      -ua_vec_view -s_vec_view
7c4762a1bSJed Brown */
8c4762a1bSJed Brown 
9c4762a1bSJed Brown #include <petsc.h>
10c4762a1bSJed Brown #include <exodusII.h>
11c4762a1bSJed Brown 
12c4762a1bSJed Brown #include <petsc/private/dmpleximpl.h>
13c4762a1bSJed Brown 
14c4762a1bSJed Brown int main(int argc, char **argv) {
15c4762a1bSJed Brown   DM                dm, dmU, dmA, dmS, dmUA, dmUA2, *dmList;
16c4762a1bSJed Brown   Vec               X, U, A, S, UA, UA2;
17c4762a1bSJed Brown   IS                isU, isA, isS, isUA;
18c4762a1bSJed Brown   PetscSection      section;
19c4762a1bSJed Brown   const PetscInt    fieldU = 0;
20c4762a1bSJed Brown   const PetscInt    fieldA = 2;
21c4762a1bSJed Brown   const PetscInt    fieldS = 1;
22c4762a1bSJed Brown   const PetscInt    fieldUA[2] = {0, 2};
23c4762a1bSJed Brown   char              ifilename[PETSC_MAX_PATH_LEN], ofilename[PETSC_MAX_PATH_LEN];
24c4762a1bSJed Brown   int               exoid = -1;
25c4762a1bSJed Brown   IS                csIS;
26c4762a1bSJed Brown   const PetscInt   *csID;
27c4762a1bSJed Brown   PetscInt         *pStartDepth, *pEndDepth;
28c4762a1bSJed Brown   PetscInt          order = 1;
29c4762a1bSJed Brown   PetscInt          sdim, d, pStart, pEnd, p, numCS, set;
30c4762a1bSJed Brown   PetscMPIInt       rank, size;
316823f3c5SBlaise Bourdin   PetscViewer       viewer;
32c4762a1bSJed Brown   PetscErrorCode    ierr;
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv,NULL, help);if (ierr) return ierr;
35ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRMPI(ierr);
36ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRMPI(ierr);
37c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FEM Layout Options", "ex26");CHKERRQ(ierr);
38c4762a1bSJed Brown   ierr = PetscOptionsString("-i", "Filename to read", "ex26", ifilename, ifilename, sizeof(ifilename), NULL);CHKERRQ(ierr);
39c4762a1bSJed Brown   ierr = PetscOptionsString("-o", "Filename to write", "ex26", ofilename, ofilename, sizeof(ofilename), NULL);CHKERRQ(ierr);
40c4762a1bSJed Brown   ierr = PetscOptionsBoundedInt("-order", "FEM polynomial order", "ex26", order, &order, NULL,1);CHKERRQ(ierr);
41c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
42*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse((order > 2) || (order < 1),PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported polynomial order %D not in [1, 2]", order);
43c4762a1bSJed Brown 
446823f3c5SBlaise Bourdin   /* Read the mesh from a file in any supported format */
45ed5e4e85SVaclav Hapla   ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename, NULL, PETSC_TRUE, &dm);CHKERRQ(ierr);
46c4762a1bSJed Brown   ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
47c4762a1bSJed Brown   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
486823f3c5SBlaise Bourdin   ierr = DMGetDimension(dm, &sdim);CHKERRQ(ierr);
49c4762a1bSJed Brown 
506823f3c5SBlaise Bourdin   /* Create the exodus result file */
516823f3c5SBlaise Bourdin   {
526823f3c5SBlaise Bourdin     PetscInt      numstep = 3, step;
536823f3c5SBlaise Bourdin     char         *nodalVarName[4];
546823f3c5SBlaise Bourdin     char         *zonalVarName[6];
556823f3c5SBlaise Bourdin     int          *truthtable;
566823f3c5SBlaise Bourdin     PetscInt      numNodalVar, numZonalVar, i;
576823f3c5SBlaise Bourdin 
58a5b23f4aSJose E. Roman     /* enable exodus debugging information */
596823f3c5SBlaise Bourdin     ex_opts(EX_VERBOSE|EX_DEBUG);
606823f3c5SBlaise Bourdin     /* Create the exodus file */
616823f3c5SBlaise Bourdin     ierr = PetscViewerExodusIIOpen(PETSC_COMM_WORLD,ofilename,FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
626823f3c5SBlaise Bourdin     /* The long way would be */
636823f3c5SBlaise Bourdin     /*
646823f3c5SBlaise Bourdin       ierr = PetscViewerCreate(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr);
656823f3c5SBlaise Bourdin       ierr = PetscViewerSetType(viewer,PETSCVIEWEREXODUSII);CHKERRQ(ierr);
666823f3c5SBlaise Bourdin       ierr = PetscViewerFileSetMode(viewer,FILE_MODE_APPEND);CHKERRQ(ierr);
676823f3c5SBlaise Bourdin       ierr = PetscViewerFileSetName(viewer,ofilename);CHKERRQ(ierr);
686823f3c5SBlaise Bourdin     */
696823f3c5SBlaise Bourdin     /* set the mesh order */
706823f3c5SBlaise Bourdin     ierr = PetscViewerExodusIISetOrder(viewer,order);CHKERRQ(ierr);
711e1ea65dSPierre Jolivet     ierr = PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
726823f3c5SBlaise Bourdin     /*
736823f3c5SBlaise Bourdin       Notice how the exodus file is actually NOT open at this point (exoid is -1)
746823f3c5SBlaise Bourdin       Since we are overwritting the file (mode is FILE_MODE_WRITE), we are going to have to
756823f3c5SBlaise Bourdin       write the geometry (the DM), which can only be done on a brand new file.
766823f3c5SBlaise Bourdin     */
776823f3c5SBlaise Bourdin 
786823f3c5SBlaise Bourdin     /* Save the geometry to the file, erasing all previous content */
796823f3c5SBlaise Bourdin     ierr = DMView(dm,viewer);CHKERRQ(ierr);
801e1ea65dSPierre Jolivet     ierr = PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
816823f3c5SBlaise Bourdin     /*
826823f3c5SBlaise Bourdin       Note how the exodus file is now open
836823f3c5SBlaise Bourdin     */
846823f3c5SBlaise Bourdin 
856823f3c5SBlaise Bourdin     /* "Format" the exodus result file, i.e. allocate space for nodal and zonal variables */
866823f3c5SBlaise Bourdin     switch (sdim) {
876823f3c5SBlaise Bourdin     case 2:
886823f3c5SBlaise Bourdin       numNodalVar = 3;
896823f3c5SBlaise Bourdin       nodalVarName[0] = (char *) "U_x";
906823f3c5SBlaise Bourdin       nodalVarName[1] = (char *) "U_y";
916823f3c5SBlaise Bourdin       nodalVarName[2] = (char *) "Alpha";
926823f3c5SBlaise Bourdin       numZonalVar = 3;
936823f3c5SBlaise Bourdin       zonalVarName[0] = (char *) "Sigma_11";
946823f3c5SBlaise Bourdin       zonalVarName[1] = (char *) "Sigma_22";
956823f3c5SBlaise Bourdin       zonalVarName[2] = (char *) "Sigma_12";
966823f3c5SBlaise Bourdin       break;
976823f3c5SBlaise Bourdin     case 3:
986823f3c5SBlaise Bourdin       numNodalVar = 4;
996823f3c5SBlaise Bourdin       nodalVarName[0] = (char *) "U_x";
1006823f3c5SBlaise Bourdin       nodalVarName[1] = (char *) "U_y";
1016823f3c5SBlaise Bourdin       nodalVarName[2] = (char *) "U_z";
1026823f3c5SBlaise Bourdin       nodalVarName[3] = (char *) "Alpha";
1036823f3c5SBlaise Bourdin       numZonalVar = 6;
1046823f3c5SBlaise Bourdin       zonalVarName[0] = (char *) "Sigma_11";
1056823f3c5SBlaise Bourdin       zonalVarName[1] = (char *) "Sigma_22";
1066823f3c5SBlaise Bourdin       zonalVarName[2] = (char *) "Sigma_33";
1076823f3c5SBlaise Bourdin       zonalVarName[3] = (char *) "Sigma_23";
1086823f3c5SBlaise Bourdin       zonalVarName[4] = (char *) "Sigma_13";
1096823f3c5SBlaise Bourdin       zonalVarName[5] = (char *) "Sigma_12";
1106823f3c5SBlaise Bourdin       break;
11198921bdaSJacob Faibussowitsch     default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No layout for dimension %D", sdim);
1126823f3c5SBlaise Bourdin     }
1136823f3c5SBlaise Bourdin     ierr = PetscViewerExodusIIGetId(viewer,&exoid);CHKERRQ(ierr);
114a74df02fSJacob Faibussowitsch     PetscStackCallStandard(ex_put_variable_param,exoid, EX_ELEM_BLOCK, numZonalVar);
115a74df02fSJacob Faibussowitsch     PetscStackCallStandard(ex_put_variable_names,exoid, EX_ELEM_BLOCK, numZonalVar, zonalVarName);
116a74df02fSJacob Faibussowitsch     PetscStackCallStandard(ex_put_variable_param,exoid, EX_NODAL, numNodalVar);
117a74df02fSJacob Faibussowitsch     PetscStackCallStandard(ex_put_variable_names,exoid, EX_NODAL, numNodalVar, nodalVarName);
1186823f3c5SBlaise Bourdin     numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK);
1196823f3c5SBlaise Bourdin 
1206823f3c5SBlaise Bourdin     /*
1216823f3c5SBlaise Bourdin       An exodusII truth table specifies which fields are saved at which time step
1226823f3c5SBlaise Bourdin       It speeds up I/O but reserving space for fieldsin the file ahead of time.
1236823f3c5SBlaise Bourdin     */
1246823f3c5SBlaise Bourdin     ierr = PetscMalloc1(numZonalVar * numCS, &truthtable);CHKERRQ(ierr);
1256823f3c5SBlaise Bourdin     for (i = 0; i < numZonalVar * numCS; ++i) truthtable[i] = 1;
126a74df02fSJacob Faibussowitsch     PetscStackCallStandard(ex_put_truth_table,exoid, EX_ELEM_BLOCK, numCS, numZonalVar, truthtable);
1276823f3c5SBlaise Bourdin     ierr = PetscFree(truthtable);CHKERRQ(ierr);
1286823f3c5SBlaise Bourdin 
1296823f3c5SBlaise Bourdin     /* Writing time step information in the file. Note that this is currently broken in the exodus library for netcdf4 (HDF5-based) files */
1306823f3c5SBlaise Bourdin     for (step = 0; step < numstep; ++step) {
1316823f3c5SBlaise Bourdin       PetscReal time = step;
132a74df02fSJacob Faibussowitsch       PetscStackCallStandard(ex_put_time,exoid, step+1, &time);
1336823f3c5SBlaise Bourdin     }
1346823f3c5SBlaise Bourdin   }
1356823f3c5SBlaise Bourdin 
1366823f3c5SBlaise Bourdin   /* Create the main section containing all fields */
137c4762a1bSJed Brown   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &section);CHKERRQ(ierr);
138c4762a1bSJed Brown   ierr = PetscSectionSetNumFields(section, 3);CHKERRQ(ierr);
139c4762a1bSJed Brown   ierr = PetscSectionSetFieldName(section, fieldU, "U");CHKERRQ(ierr);
140c4762a1bSJed Brown   ierr = PetscSectionSetFieldName(section, fieldA, "Alpha");CHKERRQ(ierr);
141c4762a1bSJed Brown   ierr = PetscSectionSetFieldName(section, fieldS, "Sigma");CHKERRQ(ierr);
142c4762a1bSJed Brown   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
143c4762a1bSJed Brown   ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr);
144c4762a1bSJed Brown   ierr = PetscMalloc2(sdim+1, &pStartDepth, sdim+1, &pEndDepth);CHKERRQ(ierr);
145c4762a1bSJed Brown   for (d = 0; d <= sdim; ++d) {ierr = DMPlexGetDepthStratum(dm, d, &pStartDepth[d], &pEndDepth[d]);CHKERRQ(ierr);}
146c4762a1bSJed Brown   /* Vector field U, Scalar field Alpha, Tensor field Sigma */
147c4762a1bSJed Brown   ierr = PetscSectionSetFieldComponents(section, fieldU, sdim);CHKERRQ(ierr);
148c4762a1bSJed Brown   ierr = PetscSectionSetFieldComponents(section, fieldA, 1);CHKERRQ(ierr);
149c4762a1bSJed Brown   ierr = PetscSectionSetFieldComponents(section, fieldS, sdim*(sdim+1)/2);CHKERRQ(ierr);
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   /* Going through cell sets then cells, and setting up storage for the sections */
1526823f3c5SBlaise Bourdin   ierr = DMGetLabelSize(dm, "Cell Sets", &numCS);CHKERRQ(ierr);
153c4762a1bSJed Brown   ierr = DMGetLabelIdIS(dm, "Cell Sets", &csIS);CHKERRQ(ierr);
154c4762a1bSJed Brown   if (csIS) {ierr = ISGetIndices(csIS, &csID);CHKERRQ(ierr);}
155c4762a1bSJed Brown   for (set = 0; set < numCS; set++) {
156c4762a1bSJed Brown     IS                cellIS;
157c4762a1bSJed Brown     const PetscInt   *cellID;
158c4762a1bSJed Brown     PetscInt          numCells, cell, closureSize, *closureA = NULL;
159c4762a1bSJed Brown 
160c4762a1bSJed Brown     ierr = DMGetStratumSize(dm, "Cell Sets", csID[set], &numCells);CHKERRQ(ierr);
161c4762a1bSJed Brown     ierr = DMGetStratumIS(dm, "Cell Sets", csID[set], &cellIS);CHKERRQ(ierr);
162c4762a1bSJed Brown     if (numCells > 0) {
163c4762a1bSJed Brown       /* dof layout ordered by increasing height in the DAG: cell, face, edge, vertex */
164c4762a1bSJed Brown       PetscInt          dofUP1Tri[]  = {2, 0, 0};
165c4762a1bSJed Brown       PetscInt          dofAP1Tri[]  = {1, 0, 0};
166c4762a1bSJed Brown       PetscInt          dofUP2Tri[]  = {2, 2, 0};
167c4762a1bSJed Brown       PetscInt          dofAP2Tri[]  = {1, 1, 0};
168c4762a1bSJed Brown       PetscInt          dofUP1Quad[] = {2, 0, 0};
169c4762a1bSJed Brown       PetscInt          dofAP1Quad[] = {1, 0, 0};
170c4762a1bSJed Brown       PetscInt          dofUP2Quad[] = {2, 2, 2};
171c4762a1bSJed Brown       PetscInt          dofAP2Quad[] = {1, 1, 1};
172c4762a1bSJed Brown       PetscInt          dofS2D[]     = {0, 0, 3};
173c4762a1bSJed Brown       PetscInt          dofUP1Tet[]  = {3, 0, 0, 0};
174c4762a1bSJed Brown       PetscInt          dofAP1Tet[]  = {1, 0, 0, 0};
175c4762a1bSJed Brown       PetscInt          dofUP2Tet[]  = {3, 3, 0, 0};
176c4762a1bSJed Brown       PetscInt          dofAP2Tet[]  = {1, 1, 0, 0};
177c4762a1bSJed Brown       PetscInt          dofUP1Hex[]  = {3, 0, 0, 0};
178c4762a1bSJed Brown       PetscInt          dofAP1Hex[]  = {1, 0, 0, 0};
179c4762a1bSJed Brown       PetscInt          dofUP2Hex[]  = {3, 3, 3, 3};
180c4762a1bSJed Brown       PetscInt          dofAP2Hex[]  = {1, 1, 1, 1};
181c4762a1bSJed Brown       PetscInt          dofS3D[]     = {0, 0, 0, 6};
182c4762a1bSJed Brown       PetscInt         *dofU, *dofA, *dofS;
183c4762a1bSJed Brown 
184c4762a1bSJed Brown       switch (sdim) {
185c4762a1bSJed Brown       case 2: dofS = dofS2D;break;
186c4762a1bSJed Brown       case 3: dofS = dofS3D;break;
18798921bdaSJacob Faibussowitsch       default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No layout for dimension %D", sdim);
188c4762a1bSJed Brown       }
189c4762a1bSJed Brown 
190c4762a1bSJed Brown       /* Identify cell type based on closure size only. This works for Tri/Tet/Quad/Hex meshes
191c4762a1bSJed Brown          It will not be enough to identify more exotic elements like pyramid or prisms...  */
192c4762a1bSJed Brown       ierr = ISGetIndices(cellIS, &cellID);CHKERRQ(ierr);
193c4762a1bSJed Brown       ierr = DMPlexGetTransitiveClosure(dm, cellID[0], PETSC_TRUE, &closureSize, &closureA);CHKERRQ(ierr);
194c4762a1bSJed Brown       switch (closureSize) {
195c4762a1bSJed Brown         case 7: /* Tri */
196c4762a1bSJed Brown         if (order == 1) {
197c4762a1bSJed Brown           dofU = dofUP1Tri;
198c4762a1bSJed Brown           dofA = dofAP1Tri;
199c4762a1bSJed Brown         } else {
200c4762a1bSJed Brown           dofU = dofUP2Tri;
201c4762a1bSJed Brown           dofA = dofAP2Tri;
202c4762a1bSJed Brown         }
203c4762a1bSJed Brown         break;
204c4762a1bSJed Brown         case 9: /* Quad */
205c4762a1bSJed Brown         if (order == 1) {
206c4762a1bSJed Brown           dofU = dofUP1Quad;
207c4762a1bSJed Brown           dofA = dofAP1Quad;
208c4762a1bSJed Brown         } else {
209c4762a1bSJed Brown           dofU = dofUP2Quad;
210c4762a1bSJed Brown           dofA = dofAP2Quad;
211c4762a1bSJed Brown         }
212c4762a1bSJed Brown         break;
213c4762a1bSJed Brown         case 15: /* Tet */
214c4762a1bSJed Brown         if (order == 1) {
215c4762a1bSJed Brown           dofU = dofUP1Tet;
216c4762a1bSJed Brown           dofA = dofAP1Tet;
217c4762a1bSJed Brown         } else {
218c4762a1bSJed Brown           dofU = dofUP2Tet;
219c4762a1bSJed Brown           dofA = dofAP2Tet;
220c4762a1bSJed Brown         }
221c4762a1bSJed Brown         break;
222c4762a1bSJed Brown         case 27: /* Hex */
223c4762a1bSJed Brown         if (order == 1) {
224c4762a1bSJed Brown           dofU = dofUP1Hex;
225c4762a1bSJed Brown           dofA = dofAP1Hex;
226c4762a1bSJed Brown         } else {
227c4762a1bSJed Brown           dofU = dofUP2Hex;
228c4762a1bSJed Brown           dofA = dofAP2Hex;
229c4762a1bSJed Brown         }
230c4762a1bSJed Brown         break;
23198921bdaSJacob Faibussowitsch         default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Unknown element with closure size %D", closureSize);
232c4762a1bSJed Brown       }
233c4762a1bSJed Brown       ierr = DMPlexRestoreTransitiveClosure(dm, cellID[0], PETSC_TRUE, &closureSize, &closureA);CHKERRQ(ierr);
234c4762a1bSJed Brown 
235c4762a1bSJed Brown       for (cell = 0; cell < numCells; cell++) {
236c4762a1bSJed Brown         PetscInt *closure = NULL;
237c4762a1bSJed Brown 
238c4762a1bSJed Brown         ierr = DMPlexGetTransitiveClosure(dm, cellID[cell], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
239c4762a1bSJed Brown         for (p = 0; p < closureSize; ++p) {
240c4762a1bSJed Brown           /* Find depth of p */
241c4762a1bSJed Brown           for (d = 0; d <= sdim; ++d) {
242c4762a1bSJed Brown             if ((closure[2*p] >= pStartDepth[d]) && (closure[2*p] < pEndDepth[d])) {
243c4762a1bSJed Brown               ierr = PetscSectionSetDof(section, closure[2*p], dofU[d]+dofA[d]+dofS[d]);CHKERRQ(ierr);
244c4762a1bSJed Brown               ierr = PetscSectionSetFieldDof(section, closure[2*p], fieldU, dofU[d]);CHKERRQ(ierr);
245c4762a1bSJed Brown               ierr = PetscSectionSetFieldDof(section, closure[2*p], fieldA, dofA[d]);CHKERRQ(ierr);
246c4762a1bSJed Brown               ierr = PetscSectionSetFieldDof(section, closure[2*p], fieldS, dofS[d]);CHKERRQ(ierr);
247c4762a1bSJed Brown             }
248c4762a1bSJed Brown           }
249c4762a1bSJed Brown         }
250c4762a1bSJed Brown         ierr = DMPlexRestoreTransitiveClosure(dm, cellID[cell], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
251c4762a1bSJed Brown       }
252c4762a1bSJed Brown       ierr = ISRestoreIndices(cellIS, &cellID);CHKERRQ(ierr);
253c4762a1bSJed Brown       ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
254c4762a1bSJed Brown     }
255c4762a1bSJed Brown   }
256c4762a1bSJed Brown   if (csIS) {ierr = ISRestoreIndices(csIS, &csID);CHKERRQ(ierr);}
257c4762a1bSJed Brown   ierr = ISDestroy(&csIS);CHKERRQ(ierr);
258c4762a1bSJed Brown   ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
259c4762a1bSJed Brown   ierr = DMSetLocalSection(dm, section);CHKERRQ(ierr);
260c4762a1bSJed Brown   ierr = PetscObjectViewFromOptions((PetscObject) section, NULL, "-dm_section_view");CHKERRQ(ierr);
261c4762a1bSJed Brown   ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   {
264c4762a1bSJed Brown     DM               pdm;
265c4762a1bSJed Brown     PetscSF          migrationSF;
266c4762a1bSJed Brown     PetscInt         ovlp = 0;
267c4762a1bSJed Brown     PetscPartitioner part;
268c4762a1bSJed Brown 
269c4762a1bSJed Brown     ierr = DMSetUseNatural(dm,PETSC_TRUE);CHKERRQ(ierr);
270c4762a1bSJed Brown     ierr = DMPlexGetPartitioner(dm,&part);CHKERRQ(ierr);
271c4762a1bSJed Brown     ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
272c4762a1bSJed Brown     ierr = DMPlexDistribute(dm,ovlp,&migrationSF,&pdm);CHKERRQ(ierr);
273c4762a1bSJed Brown     if (pdm) {
274c4762a1bSJed Brown       ierr = DMPlexSetMigrationSF(pdm,migrationSF);CHKERRQ(ierr);
275c4762a1bSJed Brown       ierr = PetscSFDestroy(&migrationSF);CHKERRQ(ierr);
276c4762a1bSJed Brown       ierr = DMDestroy(&dm);CHKERRQ(ierr);
277c4762a1bSJed Brown       dm = pdm;
278c4762a1bSJed Brown       ierr = DMViewFromOptions(dm,NULL,"-dm_view");CHKERRQ(ierr);
279c4762a1bSJed Brown     }
280c4762a1bSJed Brown   }
281c4762a1bSJed Brown 
282c4762a1bSJed Brown   /* Get DM and IS for each field of dm */
283c4762a1bSJed Brown   ierr = DMCreateSubDM(dm, 1, &fieldU, &isU,  &dmU);CHKERRQ(ierr);
284c4762a1bSJed Brown   ierr = DMCreateSubDM(dm, 1, &fieldA, &isA,  &dmA);CHKERRQ(ierr);
285c4762a1bSJed Brown   ierr = DMCreateSubDM(dm, 1, &fieldS, &isS,  &dmS);CHKERRQ(ierr);
286c4762a1bSJed Brown   ierr = DMCreateSubDM(dm, 2, fieldUA, &isUA, &dmUA);CHKERRQ(ierr);
287c4762a1bSJed Brown 
288c4762a1bSJed Brown   ierr = PetscMalloc1(2,&dmList);CHKERRQ(ierr);
289c4762a1bSJed Brown   dmList[0] = dmU;
290c4762a1bSJed Brown   dmList[1] = dmA;
291c4762a1bSJed Brown   /* We temporarily disable dmU->useNatural to test that we can reconstruct the
292c4762a1bSJed Brown      NaturaltoGlobal SF from any of the dm in dms
293c4762a1bSJed Brown   */
294c4762a1bSJed Brown   dmU->useNatural = PETSC_FALSE;
295c4762a1bSJed Brown   ierr = DMCreateSuperDM(dmList,2,NULL,&dmUA2);CHKERRQ(ierr);
296c4762a1bSJed Brown   dmU->useNatural = PETSC_TRUE;
297c4762a1bSJed Brown   ierr = PetscFree(dmList);CHKERRQ(ierr);
298c4762a1bSJed Brown 
299c4762a1bSJed Brown   ierr = DMGetGlobalVector(dm,   &X);CHKERRQ(ierr);
300c4762a1bSJed Brown   ierr = DMGetGlobalVector(dmU,  &U);CHKERRQ(ierr);
301c4762a1bSJed Brown   ierr = DMGetGlobalVector(dmA,  &A);CHKERRQ(ierr);
302c4762a1bSJed Brown   ierr = DMGetGlobalVector(dmS,  &S);CHKERRQ(ierr);
303c4762a1bSJed Brown   ierr = DMGetGlobalVector(dmUA, &UA);CHKERRQ(ierr);
304c4762a1bSJed Brown   ierr = DMGetGlobalVector(dmUA2, &UA2);CHKERRQ(ierr);
305c4762a1bSJed Brown 
306c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) U,  "U");CHKERRQ(ierr);
307c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) A,  "Alpha");CHKERRQ(ierr);
308c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) S,  "Sigma");CHKERRQ(ierr);
309c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) UA, "UAlpha");CHKERRQ(ierr);
310c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) UA2, "UAlpha2");CHKERRQ(ierr);
311c4762a1bSJed Brown   ierr = VecSet(X, -111.);CHKERRQ(ierr);
312c4762a1bSJed Brown 
313c4762a1bSJed Brown   /* 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 */
314c4762a1bSJed Brown   {
315c4762a1bSJed Brown     PetscSection sectionUA;
316c4762a1bSJed Brown     Vec          UALoc;
317c4762a1bSJed Brown     PetscSection coordSection;
318c4762a1bSJed Brown     Vec          coord;
319c4762a1bSJed Brown     PetscScalar *cval, *xyz;
3206823f3c5SBlaise Bourdin     PetscInt     clSize, i, j;
321c4762a1bSJed Brown 
322c4762a1bSJed Brown     ierr = DMGetLocalSection(dmUA, &sectionUA);CHKERRQ(ierr);
323c4762a1bSJed Brown     ierr = DMGetLocalVector(dmUA, &UALoc);CHKERRQ(ierr);
324c4762a1bSJed Brown     ierr = VecGetArray(UALoc, &cval);CHKERRQ(ierr);
325c4762a1bSJed Brown     ierr = DMGetCoordinateSection(dmUA, &coordSection);CHKERRQ(ierr);
326c4762a1bSJed Brown     ierr = DMGetCoordinatesLocal(dmUA, &coord);CHKERRQ(ierr);
327c4762a1bSJed Brown     ierr = DMPlexGetChart(dmUA, &pStart, &pEnd);CHKERRQ(ierr);
3286823f3c5SBlaise Bourdin 
329c4762a1bSJed Brown     for (p = pStart; p < pEnd; ++p) {
330c4762a1bSJed Brown       PetscInt dofUA, offUA;
331c4762a1bSJed Brown 
332c4762a1bSJed Brown       ierr = PetscSectionGetDof(sectionUA, p, &dofUA);CHKERRQ(ierr);
333c4762a1bSJed Brown       if (dofUA > 0) {
3346823f3c5SBlaise Bourdin         xyz=NULL;
335c4762a1bSJed Brown         ierr = PetscSectionGetOffset(sectionUA, p, &offUA);CHKERRQ(ierr);
336c4762a1bSJed Brown         ierr = DMPlexVecGetClosure(dmUA, coordSection, coord, p, &clSize, &xyz);CHKERRQ(ierr);
337c4762a1bSJed Brown         cval[offUA+sdim] = 0.;
338c4762a1bSJed Brown         for (i = 0; i < sdim; ++i) {
339c4762a1bSJed Brown           cval[offUA+i] = 0;
340c4762a1bSJed Brown           for (j = 0; j < clSize/sdim; ++j) {
341c4762a1bSJed Brown             cval[offUA+i] += xyz[j*sdim+i];
342c4762a1bSJed Brown           }
343c4762a1bSJed Brown           cval[offUA+i] = cval[offUA+i] * sdim / clSize;
344c4762a1bSJed Brown           cval[offUA+sdim] += PetscSqr(cval[offUA+i]);
345c4762a1bSJed Brown         }
3466823f3c5SBlaise Bourdin         ierr = DMPlexVecRestoreClosure(dmUA, coordSection, coord, p, &clSize, &xyz);CHKERRQ(ierr);
347c4762a1bSJed Brown       }
348c4762a1bSJed Brown     }
349c4762a1bSJed Brown     ierr = VecRestoreArray(UALoc, &cval);CHKERRQ(ierr);
350c4762a1bSJed Brown     ierr = DMLocalToGlobalBegin(dmUA, UALoc, INSERT_VALUES, UA);CHKERRQ(ierr);
351c4762a1bSJed Brown     ierr = DMLocalToGlobalEnd(dmUA, UALoc, INSERT_VALUES, UA);CHKERRQ(ierr);
352c4762a1bSJed Brown     ierr = DMRestoreLocalVector(dmUA, &UALoc);CHKERRQ(ierr);
3536823f3c5SBlaise Bourdin 
354c4762a1bSJed Brown     /* Update X */
355c4762a1bSJed Brown     ierr = VecISCopy(X, isUA, SCATTER_FORWARD, UA);CHKERRQ(ierr);
356c4762a1bSJed Brown     ierr = VecViewFromOptions(UA, NULL, "-ua_vec_view");CHKERRQ(ierr);
3576823f3c5SBlaise Bourdin 
358c4762a1bSJed Brown     /* Restrict to U and Alpha */
359c4762a1bSJed Brown     ierr = VecISCopy(X, isU, SCATTER_REVERSE, U);CHKERRQ(ierr);
360c4762a1bSJed Brown     ierr = VecISCopy(X, isA, SCATTER_REVERSE, A);CHKERRQ(ierr);
3616823f3c5SBlaise Bourdin 
362c4762a1bSJed Brown     /* restrict to UA2 */
363c4762a1bSJed Brown     ierr = VecISCopy(X, isUA, SCATTER_REVERSE, UA2);CHKERRQ(ierr);
364c4762a1bSJed Brown     ierr = VecViewFromOptions(UA2, NULL, "-ua2_vec_view");CHKERRQ(ierr);
365c4762a1bSJed Brown   }
366c4762a1bSJed Brown 
367c4762a1bSJed Brown   {
368c4762a1bSJed Brown     Vec          tmpVec;
369c4762a1bSJed Brown     PetscSection coordSection;
370c4762a1bSJed Brown     Vec          coord;
371c4762a1bSJed Brown     PetscReal    norm;
3726823f3c5SBlaise Bourdin     PetscReal    time = 1.234;
373c4762a1bSJed Brown 
374c4762a1bSJed Brown     /* Writing nodal variables to ExodusII file */
3756823f3c5SBlaise Bourdin     ierr = DMSetOutputSequenceNumber(dmU,0,time);CHKERRQ(ierr);
3766823f3c5SBlaise Bourdin     ierr = DMSetOutputSequenceNumber(dmA,0,time);CHKERRQ(ierr);
3776823f3c5SBlaise Bourdin 
3786823f3c5SBlaise Bourdin     ierr = VecView(U, viewer);CHKERRQ(ierr);
3796823f3c5SBlaise Bourdin     ierr = VecView(A, viewer);CHKERRQ(ierr);
3806823f3c5SBlaise Bourdin 
381c4762a1bSJed Brown     /* Saving U and Alpha in one shot.
382c4762a1bSJed Brown        For this, we need to cheat and change the Vec's name
3836823f3c5SBlaise Bourdin        Note that in the end we write variables one component at a time,
3846823f3c5SBlaise Bourdin        so that there is no real values in doing this
3856823f3c5SBlaise Bourdin     */
3866823f3c5SBlaise Bourdin 
3876823f3c5SBlaise Bourdin     ierr = DMSetOutputSequenceNumber(dmUA,1,time);CHKERRQ(ierr);
388c4762a1bSJed Brown     ierr = DMGetGlobalVector(dmUA, &tmpVec);CHKERRQ(ierr);
389c4762a1bSJed Brown     ierr = VecCopy(UA, tmpVec);CHKERRQ(ierr);
390c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject) tmpVec, "U");CHKERRQ(ierr);
3916823f3c5SBlaise Bourdin     ierr = VecView(tmpVec, viewer);CHKERRQ(ierr);
392c4762a1bSJed Brown     /* Reading nodal variables in Exodus file */
393c4762a1bSJed Brown     ierr = VecSet(tmpVec, -1000.0);CHKERRQ(ierr);
3946823f3c5SBlaise Bourdin     ierr = VecLoad(tmpVec, viewer);CHKERRQ(ierr);
395c4762a1bSJed Brown     ierr = VecAXPY(UA, -1.0, tmpVec);CHKERRQ(ierr);
396c4762a1bSJed Brown     ierr = VecNorm(UA, NORM_INFINITY, &norm);CHKERRQ(ierr);
397*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(norm > PETSC_SQRT_MACHINE_EPSILON,PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "UAlpha ||Vin - Vout|| = %g", (double) norm);
398c4762a1bSJed Brown     ierr = DMRestoreGlobalVector(dmUA, &tmpVec);CHKERRQ(ierr);
399c4762a1bSJed Brown 
400c4762a1bSJed Brown     /* same thing with the UA2 Vec obtained from the superDM */
401c4762a1bSJed Brown     ierr = DMGetGlobalVector(dmUA2, &tmpVec);CHKERRQ(ierr);
402c4762a1bSJed Brown     ierr = VecCopy(UA2, tmpVec);CHKERRQ(ierr);
403c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject) tmpVec, "U");CHKERRQ(ierr);
4046823f3c5SBlaise Bourdin     ierr = DMSetOutputSequenceNumber(dmUA2,2,time);CHKERRQ(ierr);
4056823f3c5SBlaise Bourdin     ierr = VecView(tmpVec, viewer);CHKERRQ(ierr);
406c4762a1bSJed Brown     /* Reading nodal variables in Exodus file */
407c4762a1bSJed Brown     ierr = VecSet(tmpVec, -1000.0);CHKERRQ(ierr);
4086823f3c5SBlaise Bourdin     ierr = VecLoad(tmpVec,viewer);CHKERRQ(ierr);
409c4762a1bSJed Brown     ierr = VecAXPY(UA2, -1.0, tmpVec);CHKERRQ(ierr);
410c4762a1bSJed Brown     ierr = VecNorm(UA2, NORM_INFINITY, &norm);CHKERRQ(ierr);
411*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(norm > PETSC_SQRT_MACHINE_EPSILON,PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "UAlpha2 ||Vin - Vout|| = %g", (double) norm);
412c4762a1bSJed Brown     ierr = DMRestoreGlobalVector(dmUA2, &tmpVec);CHKERRQ(ierr);
413c4762a1bSJed Brown 
414c4762a1bSJed Brown     /* Building and saving Sigma
415c4762a1bSJed Brown        We set sigma_0 = rank (to see partitioning)
416c4762a1bSJed Brown               sigma_1 = cell set ID
4176823f3c5SBlaise Bourdin               sigma_2 = x_coordinate of the cell center of mass
4186823f3c5SBlaise Bourdin     */
419c4762a1bSJed Brown     ierr = DMGetCoordinateSection(dmS, &coordSection);CHKERRQ(ierr);
420c4762a1bSJed Brown     ierr = DMGetCoordinatesLocal(dmS, &coord);CHKERRQ(ierr);
421c4762a1bSJed Brown     ierr = DMGetLabelIdIS(dmS, "Cell Sets", &csIS);CHKERRQ(ierr);
4226823f3c5SBlaise Bourdin     ierr = DMGetLabelSize(dmS, "Cell Sets", &numCS);CHKERRQ(ierr);
423c4762a1bSJed Brown     ierr = ISGetIndices(csIS, &csID);CHKERRQ(ierr);
424c4762a1bSJed Brown     for (set = 0; set < numCS; ++set) {
425c4762a1bSJed Brown       /* 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 */
426c4762a1bSJed Brown       IS              cellIS;
427c4762a1bSJed Brown       const PetscInt *cellID;
428c4762a1bSJed Brown       PetscInt        numCells, cell;
429c4762a1bSJed Brown       PetscScalar    *cval = NULL, *xyz  = NULL;
430c4762a1bSJed Brown       PetscInt        clSize, cdimCoord, c;
431c4762a1bSJed Brown 
432c4762a1bSJed Brown       ierr = DMGetStratumIS(dmS, "Cell Sets", csID[set], &cellIS);CHKERRQ(ierr);
433c4762a1bSJed Brown       ierr = ISGetIndices(cellIS, &cellID);CHKERRQ(ierr);
434c4762a1bSJed Brown       ierr = ISGetSize(cellIS, &numCells);CHKERRQ(ierr);
435c4762a1bSJed Brown       for (cell = 0; cell < numCells; cell++) {
436c4762a1bSJed Brown         ierr = DMPlexVecGetClosure(dmS, NULL, S, cellID[cell], &clSize, &cval);CHKERRQ(ierr);
437c4762a1bSJed Brown         ierr = DMPlexVecGetClosure(dmS, coordSection, coord, cellID[cell], &cdimCoord, &xyz);CHKERRQ(ierr);
438c4762a1bSJed Brown         cval[0] = rank;
439c4762a1bSJed Brown         cval[1] = csID[set];
440c4762a1bSJed Brown         cval[2] = 0.;
441c4762a1bSJed Brown         for (c = 0; c < cdimCoord/sdim; c++) cval[2] += xyz[c*sdim];
442c4762a1bSJed Brown         cval[2] = cval[2] * sdim / cdimCoord;
443c4762a1bSJed Brown         ierr = DMPlexVecSetClosure(dmS, NULL, S, cellID[cell], cval, INSERT_ALL_VALUES);CHKERRQ(ierr);
444c4762a1bSJed Brown       }
445c4762a1bSJed Brown       ierr = DMPlexVecRestoreClosure(dmS, NULL, S, cellID[0], &clSize, &cval);CHKERRQ(ierr);
446c4762a1bSJed Brown       ierr = DMPlexVecRestoreClosure(dmS, coordSection, coord, cellID[0], NULL, &xyz);CHKERRQ(ierr);
447c4762a1bSJed Brown       ierr = ISRestoreIndices(cellIS, &cellID);CHKERRQ(ierr);
448c4762a1bSJed Brown       ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
449c4762a1bSJed Brown     }
450c4762a1bSJed Brown     ierr = ISRestoreIndices(csIS, &csID);CHKERRQ(ierr);
451c4762a1bSJed Brown     ierr = ISDestroy(&csIS);CHKERRQ(ierr);
452c4762a1bSJed Brown     ierr = VecViewFromOptions(S, NULL, "-s_vec_view");CHKERRQ(ierr);
4536823f3c5SBlaise Bourdin 
454c4762a1bSJed Brown     /* Writing zonal variables in Exodus file */
4556823f3c5SBlaise Bourdin     ierr = DMSetOutputSequenceNumber(dmS,0,time);CHKERRQ(ierr);
4566823f3c5SBlaise Bourdin     ierr = VecView(S,viewer);CHKERRQ(ierr);
4576823f3c5SBlaise Bourdin 
458c4762a1bSJed Brown     /* Reading zonal variables in Exodus file */
459c4762a1bSJed Brown     ierr = DMGetGlobalVector(dmS, &tmpVec);CHKERRQ(ierr);
460c4762a1bSJed Brown     ierr = VecSet(tmpVec, -1000.0);CHKERRQ(ierr);
461c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject) tmpVec, "Sigma");CHKERRQ(ierr);
4626823f3c5SBlaise Bourdin     ierr = VecLoad(tmpVec,viewer);CHKERRQ(ierr);
463c4762a1bSJed Brown     ierr = VecAXPY(S, -1.0, tmpVec);CHKERRQ(ierr);
4641e1ea65dSPierre Jolivet     ierr = VecNorm(S, NORM_INFINITY, &norm);CHKERRQ(ierr);
465*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(norm > PETSC_SQRT_MACHINE_EPSILON,PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "Sigma ||Vin - Vout|| = %g", (double) norm);
466c4762a1bSJed Brown     ierr = DMRestoreGlobalVector(dmS, &tmpVec);CHKERRQ(ierr);
467c4762a1bSJed Brown   }
4686823f3c5SBlaise Bourdin   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
469c4762a1bSJed Brown 
470c4762a1bSJed Brown   ierr = DMRestoreGlobalVector(dmUA2, &UA2);CHKERRQ(ierr);
471c4762a1bSJed Brown   ierr = DMRestoreGlobalVector(dmUA, &UA);CHKERRQ(ierr);
472c4762a1bSJed Brown   ierr = DMRestoreGlobalVector(dmS,  &S);CHKERRQ(ierr);
473c4762a1bSJed Brown   ierr = DMRestoreGlobalVector(dmA,  &A);CHKERRQ(ierr);
474c4762a1bSJed Brown   ierr = DMRestoreGlobalVector(dmU,  &U);CHKERRQ(ierr);
475c4762a1bSJed Brown   ierr = DMRestoreGlobalVector(dm,   &X);CHKERRQ(ierr);
476c4762a1bSJed Brown   ierr = DMDestroy(&dmU);CHKERRQ(ierr); ierr = ISDestroy(&isU);CHKERRQ(ierr);
477c4762a1bSJed Brown   ierr = DMDestroy(&dmA);CHKERRQ(ierr); ierr = ISDestroy(&isA);CHKERRQ(ierr);
478c4762a1bSJed Brown   ierr = DMDestroy(&dmS);CHKERRQ(ierr); ierr = ISDestroy(&isS);CHKERRQ(ierr);
479c4762a1bSJed Brown   ierr = DMDestroy(&dmUA);CHKERRQ(ierr);ierr = ISDestroy(&isUA);CHKERRQ(ierr);
480c4762a1bSJed Brown   ierr = DMDestroy(&dmUA2);CHKERRQ(ierr);
481c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
482c4762a1bSJed Brown   ierr = PetscFree2(pStartDepth, pEndDepth);CHKERRQ(ierr);
483c4762a1bSJed Brown   ierr = PetscFinalize();
484c4762a1bSJed Brown   return ierr;
485c4762a1bSJed Brown }
486c4762a1bSJed Brown 
487c4762a1bSJed Brown /*TEST
488c4762a1bSJed Brown 
489c4762a1bSJed Brown   build:
490c4762a1bSJed Brown     requires: exodusii pnetcdf !complex
491c4762a1bSJed Brown   # 2D seq
492c4762a1bSJed Brown   test:
493c4762a1bSJed Brown     suffix: 0
494c4762a1bSJed Brown     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
495c4762a1bSJed Brown     #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
496c4762a1bSJed Brown   test:
497c4762a1bSJed Brown     suffix: 1
498c4762a1bSJed Brown     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
499c4762a1bSJed Brown   test:
500c4762a1bSJed Brown     suffix: 2
501c4762a1bSJed Brown     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
502c4762a1bSJed Brown     #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
503c4762a1bSJed Brown   test:
504c4762a1bSJed Brown     suffix: 3
505c4762a1bSJed Brown     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
506c4762a1bSJed Brown   test:
507c4762a1bSJed Brown     suffix: 4
508c4762a1bSJed Brown     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
509c4762a1bSJed Brown   test:
510c4762a1bSJed Brown     suffix: 5
511c4762a1bSJed Brown     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
512c4762a1bSJed Brown 
513c4762a1bSJed Brown   # 2D par
514c4762a1bSJed Brown   test:
515c4762a1bSJed Brown     suffix: 6
516c4762a1bSJed Brown     nsize: 2
517c4762a1bSJed Brown     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
518c4762a1bSJed Brown     #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
519c4762a1bSJed Brown   test:
520c4762a1bSJed Brown     suffix: 7
521c4762a1bSJed Brown     nsize: 2
522c4762a1bSJed Brown     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
523c4762a1bSJed Brown   test:
524c4762a1bSJed Brown     suffix: 8
525c4762a1bSJed Brown     nsize: 2
526c4762a1bSJed Brown     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
527c4762a1bSJed Brown     #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: invalid dimension ID or name
528c4762a1bSJed Brown   test:
529c4762a1bSJed Brown     suffix: 9
530c4762a1bSJed Brown     nsize: 2
531c4762a1bSJed Brown     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
532c4762a1bSJed Brown   test:
533c4762a1bSJed Brown     suffix: 10
534c4762a1bSJed Brown     nsize: 2
535c4762a1bSJed Brown     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
536c4762a1bSJed Brown   test:
53754fcfd0cSMatthew G. Knepley     # Something is now broken with parallel read/write for wahtever shape H is
53854fcfd0cSMatthew G. Knepley     TODO: broken
539c4762a1bSJed Brown     suffix: 11
540c4762a1bSJed Brown     nsize: 2
541c4762a1bSJed Brown     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
542c4762a1bSJed Brown 
543c4762a1bSJed Brown   #3d seq
544c4762a1bSJed Brown   test:
545c4762a1bSJed Brown     suffix: 12
546c4762a1bSJed Brown     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
547c4762a1bSJed Brown   test:
548c4762a1bSJed Brown     suffix: 13
549c4762a1bSJed Brown     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
550c4762a1bSJed Brown     #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
551c4762a1bSJed Brown   test:
552c4762a1bSJed Brown     suffix: 14
553c4762a1bSJed Brown     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
554c4762a1bSJed Brown   test:
555c4762a1bSJed Brown     suffix: 15
556c4762a1bSJed Brown     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
557c4762a1bSJed Brown     #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
558c4762a1bSJed Brown   #3d par
559c4762a1bSJed Brown   test:
560c4762a1bSJed Brown     suffix: 16
561c4762a1bSJed Brown     nsize: 2
562c4762a1bSJed Brown     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
563c4762a1bSJed Brown   test:
564c4762a1bSJed Brown     suffix: 17
565c4762a1bSJed Brown     nsize: 2
566c4762a1bSJed Brown     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
567c4762a1bSJed Brown     #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
568c4762a1bSJed Brown   test:
569c4762a1bSJed Brown     suffix: 18
570c4762a1bSJed Brown     nsize: 2
571c4762a1bSJed Brown     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
572c4762a1bSJed Brown   test:
573c4762a1bSJed Brown     suffix: 19
574c4762a1bSJed Brown     nsize: 2
575c4762a1bSJed Brown     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
576c4762a1bSJed Brown     #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
577c4762a1bSJed Brown 
578c4762a1bSJed Brown TEST*/
579