xref: /petsc/src/dm/impls/plex/tests/ex26.c (revision 2e8d78fe20bce725752f1fc408b7d282cf1a9b39)
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) {
15e2739ba6SAlexis Marboeuf   DM                dm, pdm, 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 
33327415f7SBarry Smith   PetscFunctionBeginUser;
349566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv,NULL, help));
359566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
369566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
37d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FEM Layout Options", "ex26");
389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-i", "Filename to read", "ex26", ifilename, ifilename, sizeof(ifilename), NULL));
399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-o", "Filename to write", "ex26", ofilename, ofilename, sizeof(ofilename), NULL));
409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-order", "FEM polynomial order", "ex26", order, &order, NULL,1));
41d0609cedSBarry Smith   PetscOptionsEnd();
4263a3b9bcSJacob Faibussowitsch   PetscCheck((order >= 1) && (order <= 2),PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported polynomial order %" PetscInt_FMT " not in [1, 2]", order);
43c4762a1bSJed Brown 
446823f3c5SBlaise Bourdin   /* Read the mesh from a file in any supported format */
459566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename, NULL, PETSC_TRUE, &dm));
469566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeSetDefault(dm, PETSC_FALSE));
479566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
489566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
499566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &sdim));
50c4762a1bSJed Brown 
516823f3c5SBlaise Bourdin   /* Create the exodus result file */
526823f3c5SBlaise Bourdin   {
536823f3c5SBlaise Bourdin     PetscInt      numstep = 3, step;
546823f3c5SBlaise Bourdin     char         *nodalVarName[4];
556823f3c5SBlaise Bourdin     char         *zonalVarName[6];
566823f3c5SBlaise Bourdin     int          *truthtable;
576823f3c5SBlaise Bourdin     PetscInt      numNodalVar, numZonalVar, i;
586823f3c5SBlaise Bourdin 
59a5b23f4aSJose E. Roman     /* enable exodus debugging information */
606823f3c5SBlaise Bourdin     ex_opts(EX_VERBOSE|EX_DEBUG);
616823f3c5SBlaise Bourdin     /* Create the exodus file */
629566063dSJacob Faibussowitsch     PetscCall(PetscViewerExodusIIOpen(PETSC_COMM_WORLD,ofilename,FILE_MODE_WRITE,&viewer));
636823f3c5SBlaise Bourdin     /* The long way would be */
646823f3c5SBlaise Bourdin     /*
659566063dSJacob Faibussowitsch       PetscCall(PetscViewerCreate(PETSC_COMM_WORLD,&viewer));
669566063dSJacob Faibussowitsch       PetscCall(PetscViewerSetType(viewer,PETSCVIEWEREXODUSII));
679566063dSJacob Faibussowitsch       PetscCall(PetscViewerFileSetMode(viewer,FILE_MODE_APPEND));
689566063dSJacob Faibussowitsch       PetscCall(PetscViewerFileSetName(viewer,ofilename));
696823f3c5SBlaise Bourdin     */
706823f3c5SBlaise Bourdin     /* set the mesh order */
719566063dSJacob Faibussowitsch     PetscCall(PetscViewerExodusIISetOrder(viewer,order));
729566063dSJacob Faibussowitsch     PetscCall(PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD));
736823f3c5SBlaise Bourdin     /*
746823f3c5SBlaise Bourdin       Notice how the exodus file is actually NOT open at this point (exoid is -1)
75e2739ba6SAlexis Marboeuf       Since we are overwritting the file (mode is FILE_MODE_WRITE), we are going to have to
766823f3c5SBlaise Bourdin       write the geometry (the DM), which can only be done on a brand new file.
776823f3c5SBlaise Bourdin     */
786823f3c5SBlaise Bourdin 
796823f3c5SBlaise Bourdin     /* Save the geometry to the file, erasing all previous content */
809566063dSJacob Faibussowitsch     PetscCall(DMView(dm,viewer));
819566063dSJacob Faibussowitsch     PetscCall(PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD));
826823f3c5SBlaise Bourdin     /*
836823f3c5SBlaise Bourdin       Note how the exodus file is now open
846823f3c5SBlaise Bourdin     */
856823f3c5SBlaise Bourdin 
866823f3c5SBlaise Bourdin     /* "Format" the exodus result file, i.e. allocate space for nodal and zonal variables */
876823f3c5SBlaise Bourdin     switch (sdim) {
886823f3c5SBlaise Bourdin     case 2:
896823f3c5SBlaise Bourdin       numNodalVar = 3;
906823f3c5SBlaise Bourdin       nodalVarName[0] = (char *) "U_x";
916823f3c5SBlaise Bourdin       nodalVarName[1] = (char *) "U_y";
926823f3c5SBlaise Bourdin       nodalVarName[2] = (char *) "Alpha";
936823f3c5SBlaise Bourdin       numZonalVar = 3;
946823f3c5SBlaise Bourdin       zonalVarName[0] = (char *) "Sigma_11";
956823f3c5SBlaise Bourdin       zonalVarName[1] = (char *) "Sigma_22";
966823f3c5SBlaise Bourdin       zonalVarName[2] = (char *) "Sigma_12";
976823f3c5SBlaise Bourdin       break;
986823f3c5SBlaise Bourdin     case 3:
996823f3c5SBlaise Bourdin       numNodalVar = 4;
1006823f3c5SBlaise Bourdin       nodalVarName[0] = (char *) "U_x";
1016823f3c5SBlaise Bourdin       nodalVarName[1] = (char *) "U_y";
1026823f3c5SBlaise Bourdin       nodalVarName[2] = (char *) "U_z";
1036823f3c5SBlaise Bourdin       nodalVarName[3] = (char *) "Alpha";
1046823f3c5SBlaise Bourdin       numZonalVar = 6;
1056823f3c5SBlaise Bourdin       zonalVarName[0] = (char *) "Sigma_11";
1066823f3c5SBlaise Bourdin       zonalVarName[1] = (char *) "Sigma_22";
1076823f3c5SBlaise Bourdin       zonalVarName[2] = (char *) "Sigma_33";
1086823f3c5SBlaise Bourdin       zonalVarName[3] = (char *) "Sigma_23";
1096823f3c5SBlaise Bourdin       zonalVarName[4] = (char *) "Sigma_13";
1106823f3c5SBlaise Bourdin       zonalVarName[5] = (char *) "Sigma_12";
1116823f3c5SBlaise Bourdin       break;
11263a3b9bcSJacob Faibussowitsch     default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No layout for dimension %" PetscInt_FMT, sdim);
1136823f3c5SBlaise Bourdin     }
1149566063dSJacob Faibussowitsch     PetscCall(PetscViewerExodusIIGetId(viewer,&exoid));
115792fecdfSBarry Smith     PetscCallExternal(ex_put_variable_param,exoid, EX_ELEM_BLOCK, numZonalVar);
116792fecdfSBarry Smith     PetscCallExternal(ex_put_variable_names,exoid, EX_ELEM_BLOCK, numZonalVar, zonalVarName);
117792fecdfSBarry Smith     PetscCallExternal(ex_put_variable_param,exoid, EX_NODAL, numNodalVar);
118792fecdfSBarry Smith     PetscCallExternal(ex_put_variable_names,exoid, EX_NODAL, numNodalVar, nodalVarName);
1196823f3c5SBlaise Bourdin     numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK);
1206823f3c5SBlaise Bourdin 
1216823f3c5SBlaise Bourdin     /*
1226823f3c5SBlaise Bourdin       An exodusII truth table specifies which fields are saved at which time step
1236823f3c5SBlaise Bourdin       It speeds up I/O but reserving space for fieldsin the file ahead of time.
1246823f3c5SBlaise Bourdin     */
1259566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numZonalVar * numCS, &truthtable));
1266823f3c5SBlaise Bourdin     for (i = 0; i < numZonalVar * numCS; ++i) truthtable[i] = 1;
127792fecdfSBarry Smith     PetscCallExternal(ex_put_truth_table,exoid, EX_ELEM_BLOCK, numCS, numZonalVar, truthtable);
1289566063dSJacob Faibussowitsch     PetscCall(PetscFree(truthtable));
1296823f3c5SBlaise Bourdin 
1306823f3c5SBlaise Bourdin     /* Writing time step information in the file. Note that this is currently broken in the exodus library for netcdf4 (HDF5-based) files */
1316823f3c5SBlaise Bourdin     for (step = 0; step < numstep; ++step) {
1326823f3c5SBlaise Bourdin       PetscReal time = step;
133792fecdfSBarry Smith       PetscCallExternal(ex_put_time,exoid, step+1, &time);
1346823f3c5SBlaise Bourdin     }
1356823f3c5SBlaise Bourdin   }
1366823f3c5SBlaise Bourdin 
1376823f3c5SBlaise Bourdin   /* Create the main section containing all fields */
1389566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) dm), &section));
1399566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetNumFields(section, 3));
1409566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFieldName(section, fieldU, "U"));
1419566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFieldName(section, fieldA, "Alpha"));
1429566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFieldName(section, fieldS, "Sigma"));
1439566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
1449566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(section, pStart, pEnd));
1459566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(sdim+1, &pStartDepth, sdim+1, &pEndDepth));
146e2739ba6SAlexis Marboeuf   for (d = 0; d <= sdim; ++d) {PetscCall(DMPlexGetDepthStratum(dm, d, &pStartDepth[d], &pEndDepth[d]));}
147c4762a1bSJed Brown   /* Vector field U, Scalar field Alpha, Tensor field Sigma */
1489566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFieldComponents(section, fieldU, sdim));
1499566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFieldComponents(section, fieldA, 1));
1509566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFieldComponents(section, fieldS, sdim*(sdim+1)/2));
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   /* Going through cell sets then cells, and setting up storage for the sections */
1539566063dSJacob Faibussowitsch   PetscCall(DMGetLabelSize(dm, "Cell Sets", &numCS));
1549566063dSJacob Faibussowitsch   PetscCall(DMGetLabelIdIS(dm, "Cell Sets", &csIS));
155e2739ba6SAlexis Marboeuf   if (csIS) {PetscCall(ISGetIndices(csIS, &csID));}
156c4762a1bSJed Brown   for (set = 0; set < numCS; set++) {
157c4762a1bSJed Brown     IS                cellIS;
158c4762a1bSJed Brown     const PetscInt   *cellID;
159c4762a1bSJed Brown     PetscInt          numCells, cell, closureSize, *closureA = NULL;
160c4762a1bSJed Brown 
1619566063dSJacob Faibussowitsch     PetscCall(DMGetStratumSize(dm, "Cell Sets", csID[set], &numCells));
1629566063dSJacob Faibussowitsch     PetscCall(DMGetStratumIS(dm, "Cell Sets", csID[set], &cellIS));
163c4762a1bSJed Brown     if (numCells > 0) {
164c4762a1bSJed Brown       /* dof layout ordered by increasing height in the DAG: cell, face, edge, vertex */
165c4762a1bSJed Brown       PetscInt          dofUP1Tri[]  = {2, 0, 0};
166c4762a1bSJed Brown       PetscInt          dofAP1Tri[]  = {1, 0, 0};
167c4762a1bSJed Brown       PetscInt          dofUP2Tri[]  = {2, 2, 0};
168c4762a1bSJed Brown       PetscInt          dofAP2Tri[]  = {1, 1, 0};
169c4762a1bSJed Brown       PetscInt          dofUP1Quad[] = {2, 0, 0};
170c4762a1bSJed Brown       PetscInt          dofAP1Quad[] = {1, 0, 0};
171c4762a1bSJed Brown       PetscInt          dofUP2Quad[] = {2, 2, 2};
172c4762a1bSJed Brown       PetscInt          dofAP2Quad[] = {1, 1, 1};
173c4762a1bSJed Brown       PetscInt          dofS2D[]     = {0, 0, 3};
174c4762a1bSJed Brown       PetscInt          dofUP1Tet[]  = {3, 0, 0, 0};
175c4762a1bSJed Brown       PetscInt          dofAP1Tet[]  = {1, 0, 0, 0};
176c4762a1bSJed Brown       PetscInt          dofUP2Tet[]  = {3, 3, 0, 0};
177c4762a1bSJed Brown       PetscInt          dofAP2Tet[]  = {1, 1, 0, 0};
178c4762a1bSJed Brown       PetscInt          dofUP1Hex[]  = {3, 0, 0, 0};
179c4762a1bSJed Brown       PetscInt          dofAP1Hex[]  = {1, 0, 0, 0};
180c4762a1bSJed Brown       PetscInt          dofUP2Hex[]  = {3, 3, 3, 3};
181c4762a1bSJed Brown       PetscInt          dofAP2Hex[]  = {1, 1, 1, 1};
182c4762a1bSJed Brown       PetscInt          dofS3D[]     = {0, 0, 0, 6};
183c4762a1bSJed Brown       PetscInt         *dofU, *dofA, *dofS;
184c4762a1bSJed Brown 
185c4762a1bSJed Brown       switch (sdim) {
186c4762a1bSJed Brown       case 2: dofS = dofS2D;break;
187c4762a1bSJed Brown       case 3: dofS = dofS3D;break;
18863a3b9bcSJacob Faibussowitsch       default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No layout for dimension %" PetscInt_FMT, sdim);
189c4762a1bSJed Brown       }
190c4762a1bSJed Brown 
191c4762a1bSJed Brown       /* Identify cell type based on closure size only. This works for Tri/Tet/Quad/Hex meshes
192c4762a1bSJed Brown          It will not be enough to identify more exotic elements like pyramid or prisms...  */
1939566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(cellIS, &cellID));
1949566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(dm, cellID[0], PETSC_TRUE, &closureSize, &closureA));
195c4762a1bSJed Brown       switch (closureSize) {
196c4762a1bSJed Brown         case 7: /* Tri */
197c4762a1bSJed Brown         if (order == 1) {
198c4762a1bSJed Brown           dofU = dofUP1Tri;
199c4762a1bSJed Brown           dofA = dofAP1Tri;
200c4762a1bSJed Brown         } else {
201c4762a1bSJed Brown           dofU = dofUP2Tri;
202c4762a1bSJed Brown           dofA = dofAP2Tri;
203c4762a1bSJed Brown         }
204c4762a1bSJed Brown         break;
205c4762a1bSJed Brown         case 9: /* Quad */
206c4762a1bSJed Brown         if (order == 1) {
207c4762a1bSJed Brown           dofU = dofUP1Quad;
208c4762a1bSJed Brown           dofA = dofAP1Quad;
209c4762a1bSJed Brown         } else {
210c4762a1bSJed Brown           dofU = dofUP2Quad;
211c4762a1bSJed Brown           dofA = dofAP2Quad;
212c4762a1bSJed Brown         }
213c4762a1bSJed Brown         break;
214c4762a1bSJed Brown         case 15: /* Tet */
215c4762a1bSJed Brown         if (order == 1) {
216c4762a1bSJed Brown           dofU = dofUP1Tet;
217c4762a1bSJed Brown           dofA = dofAP1Tet;
218c4762a1bSJed Brown         } else {
219c4762a1bSJed Brown           dofU = dofUP2Tet;
220c4762a1bSJed Brown           dofA = dofAP2Tet;
221c4762a1bSJed Brown         }
222c4762a1bSJed Brown         break;
223c4762a1bSJed Brown         case 27: /* Hex */
224c4762a1bSJed Brown         if (order == 1) {
225c4762a1bSJed Brown           dofU = dofUP1Hex;
226c4762a1bSJed Brown           dofA = dofAP1Hex;
227c4762a1bSJed Brown         } else {
228c4762a1bSJed Brown           dofU = dofUP2Hex;
229c4762a1bSJed Brown           dofA = dofAP2Hex;
230c4762a1bSJed Brown         }
231c4762a1bSJed Brown         break;
23263a3b9bcSJacob Faibussowitsch         default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Unknown element with closure size %" PetscInt_FMT, closureSize);
233c4762a1bSJed Brown       }
2349566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(dm, cellID[0], PETSC_TRUE, &closureSize, &closureA));
235c4762a1bSJed Brown 
236c4762a1bSJed Brown       for (cell = 0; cell < numCells; cell++) {
237c4762a1bSJed Brown         PetscInt *closure = NULL;
238c4762a1bSJed Brown 
2399566063dSJacob Faibussowitsch         PetscCall(DMPlexGetTransitiveClosure(dm, cellID[cell], PETSC_TRUE, &closureSize, &closure));
240c4762a1bSJed Brown         for (p = 0; p < closureSize; ++p) {
241c4762a1bSJed Brown           /* Find depth of p */
242c4762a1bSJed Brown           for (d = 0; d <= sdim; ++d) {
243c4762a1bSJed Brown             if ((closure[2*p] >= pStartDepth[d]) && (closure[2*p] < pEndDepth[d])) {
2449566063dSJacob Faibussowitsch               PetscCall(PetscSectionSetDof(section, closure[2*p], dofU[d]+dofA[d]+dofS[d]));
2459566063dSJacob Faibussowitsch               PetscCall(PetscSectionSetFieldDof(section, closure[2*p], fieldU, dofU[d]));
2469566063dSJacob Faibussowitsch               PetscCall(PetscSectionSetFieldDof(section, closure[2*p], fieldA, dofA[d]));
2479566063dSJacob Faibussowitsch               PetscCall(PetscSectionSetFieldDof(section, closure[2*p], fieldS, dofS[d]));
248c4762a1bSJed Brown             }
249c4762a1bSJed Brown           }
250c4762a1bSJed Brown         }
2519566063dSJacob Faibussowitsch         PetscCall(DMPlexRestoreTransitiveClosure(dm, cellID[cell], PETSC_TRUE, &closureSize, &closure));
252c4762a1bSJed Brown       }
2539566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(cellIS, &cellID));
2549566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&cellIS));
255c4762a1bSJed Brown     }
256c4762a1bSJed Brown   }
257e2739ba6SAlexis Marboeuf   if (csIS) {PetscCall(ISRestoreIndices(csIS, &csID));}
2589566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&csIS));
2599566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(section));
2609566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(dm, section));
2619566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject) section, NULL, "-dm_section_view"));
2629566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&section));
263c4762a1bSJed Brown 
264c4762a1bSJed Brown   {
265c4762a1bSJed Brown     PetscSF          migrationSF;
266c4762a1bSJed Brown     PetscInt         ovlp = 0;
267c4762a1bSJed Brown     PetscPartitioner part;
268c4762a1bSJed Brown 
2699566063dSJacob Faibussowitsch     PetscCall(DMSetUseNatural(dm,PETSC_TRUE));
2709566063dSJacob Faibussowitsch     PetscCall(DMPlexGetPartitioner(dm,&part));
2719566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerSetFromOptions(part));
2729566063dSJacob Faibussowitsch     PetscCall(DMPlexDistribute(dm,ovlp,&migrationSF,&pdm));
273e2739ba6SAlexis Marboeuf     if (!pdm) pdm = dm;
274e2739ba6SAlexis Marboeuf     /* Set the migrationSF is mandatory since useNatural = PETSC_TRUE */
275e2739ba6SAlexis Marboeuf     if (migrationSF) {
2769566063dSJacob Faibussowitsch       PetscCall(DMPlexSetMigrationSF(pdm, migrationSF));
2779566063dSJacob Faibussowitsch       PetscCall(PetscSFDestroy(&migrationSF));
278c4762a1bSJed Brown     }
279e2739ba6SAlexis Marboeuf     PetscCall(DMViewFromOptions(pdm, NULL, "-dm_view"));
280c4762a1bSJed Brown   }
281c4762a1bSJed Brown 
282c4762a1bSJed Brown   /* Get DM and IS for each field of dm */
283e2739ba6SAlexis Marboeuf   PetscCall(DMCreateSubDM(pdm, 1, &fieldU, &isU,  &dmU));
284e2739ba6SAlexis Marboeuf   PetscCall(DMCreateSubDM(pdm, 1, &fieldA, &isA,  &dmA));
285e2739ba6SAlexis Marboeuf   PetscCall(DMCreateSubDM(pdm, 1, &fieldS, &isS,  &dmS));
286e2739ba6SAlexis Marboeuf   PetscCall(DMCreateSubDM(pdm, 2, fieldUA, &isUA, &dmUA));
287c4762a1bSJed Brown 
2889566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2,&dmList));
289c4762a1bSJed Brown   dmList[0] = dmU;
290c4762a1bSJed Brown   dmList[1] = dmA;
291e2739ba6SAlexis Marboeuf 
2929566063dSJacob Faibussowitsch   PetscCall(DMCreateSuperDM(dmList,2,NULL,&dmUA2));
2939566063dSJacob Faibussowitsch   PetscCall(PetscFree(dmList));
294c4762a1bSJed Brown 
295e2739ba6SAlexis Marboeuf   PetscCall(DMGetGlobalVector(pdm,  &X));
2969566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dmU,  &U));
2979566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dmA,  &A));
2989566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dmS,  &S));
2999566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dmUA, &UA));
3009566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dmUA2, &UA2));
301c4762a1bSJed Brown 
3029566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) U,  "U"));
3039566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) A,  "Alpha"));
3049566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) S,  "Sigma"));
3059566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) UA, "UAlpha"));
3069566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) UA2, "UAlpha2"));
3079566063dSJacob Faibussowitsch   PetscCall(VecSet(X, -111.));
308c4762a1bSJed Brown 
309c4762a1bSJed 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 */
310c4762a1bSJed Brown   {
311c4762a1bSJed Brown     PetscSection sectionUA;
312c4762a1bSJed Brown     Vec          UALoc;
313c4762a1bSJed Brown     PetscSection coordSection;
314c4762a1bSJed Brown     Vec          coord;
315c4762a1bSJed Brown     PetscScalar *cval, *xyz;
3166823f3c5SBlaise Bourdin     PetscInt     clSize, i, j;
317c4762a1bSJed Brown 
3189566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dmUA, &sectionUA));
3199566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dmUA, &UALoc));
3209566063dSJacob Faibussowitsch     PetscCall(VecGetArray(UALoc, &cval));
3219566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateSection(dmUA, &coordSection));
3229566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dmUA, &coord));
3239566063dSJacob Faibussowitsch     PetscCall(DMPlexGetChart(dmUA, &pStart, &pEnd));
3246823f3c5SBlaise Bourdin 
325c4762a1bSJed Brown     for (p = pStart; p < pEnd; ++p) {
326c4762a1bSJed Brown       PetscInt dofUA, offUA;
327c4762a1bSJed Brown 
3289566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(sectionUA, p, &dofUA));
329c4762a1bSJed Brown       if (dofUA > 0) {
3306823f3c5SBlaise Bourdin         xyz=NULL;
3319566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionUA, p, &offUA));
3329566063dSJacob Faibussowitsch         PetscCall(DMPlexVecGetClosure(dmUA, coordSection, coord, p, &clSize, &xyz));
333c4762a1bSJed Brown         cval[offUA+sdim] = 0.;
334c4762a1bSJed Brown         for (i = 0; i < sdim; ++i) {
335c4762a1bSJed Brown           cval[offUA+i] = 0;
336c4762a1bSJed Brown           for (j = 0; j < clSize/sdim; ++j) {
337c4762a1bSJed Brown             cval[offUA+i] += xyz[j*sdim+i];
338c4762a1bSJed Brown           }
339c4762a1bSJed Brown           cval[offUA+i] = cval[offUA+i] * sdim / clSize;
340c4762a1bSJed Brown           cval[offUA+sdim] += PetscSqr(cval[offUA+i]);
341c4762a1bSJed Brown         }
3429566063dSJacob Faibussowitsch         PetscCall(DMPlexVecRestoreClosure(dmUA, coordSection, coord, p, &clSize, &xyz));
343c4762a1bSJed Brown       }
344c4762a1bSJed Brown     }
3459566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(UALoc, &cval));
3469566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(dmUA, UALoc, INSERT_VALUES, UA));
3479566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(dmUA, UALoc, INSERT_VALUES, UA));
3489566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dmUA, &UALoc));
3496823f3c5SBlaise Bourdin 
350c4762a1bSJed Brown     /* Update X */
3519566063dSJacob Faibussowitsch     PetscCall(VecISCopy(X, isUA, SCATTER_FORWARD, UA));
3529566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(UA, NULL, "-ua_vec_view"));
3536823f3c5SBlaise Bourdin 
354c4762a1bSJed Brown     /* Restrict to U and Alpha */
3559566063dSJacob Faibussowitsch     PetscCall(VecISCopy(X, isU, SCATTER_REVERSE, U));
3569566063dSJacob Faibussowitsch     PetscCall(VecISCopy(X, isA, SCATTER_REVERSE, A));
3576823f3c5SBlaise Bourdin 
358c4762a1bSJed Brown     /* restrict to UA2 */
3599566063dSJacob Faibussowitsch     PetscCall(VecISCopy(X, isUA, SCATTER_REVERSE, UA2));
3609566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(UA2, NULL, "-ua2_vec_view"));
361c4762a1bSJed Brown   }
362c4762a1bSJed Brown 
363c4762a1bSJed Brown   {
364c4762a1bSJed Brown     Vec          tmpVec;
365c4762a1bSJed Brown     PetscSection coordSection;
366c4762a1bSJed Brown     Vec          coord;
367c4762a1bSJed Brown     PetscReal    norm;
3686823f3c5SBlaise Bourdin     PetscReal    time = 1.234;
369c4762a1bSJed Brown 
370c4762a1bSJed Brown     /* Writing nodal variables to ExodusII file */
3719566063dSJacob Faibussowitsch     PetscCall(DMSetOutputSequenceNumber(dmU,0,time));
3729566063dSJacob Faibussowitsch     PetscCall(DMSetOutputSequenceNumber(dmA,0,time));
3736823f3c5SBlaise Bourdin 
3749566063dSJacob Faibussowitsch     PetscCall(VecView(U, viewer));
3759566063dSJacob Faibussowitsch     PetscCall(VecView(A, viewer));
3766823f3c5SBlaise Bourdin 
377c4762a1bSJed Brown     /* Saving U and Alpha in one shot.
378c4762a1bSJed Brown        For this, we need to cheat and change the Vec's name
3796823f3c5SBlaise Bourdin        Note that in the end we write variables one component at a time,
3806823f3c5SBlaise Bourdin        so that there is no real values in doing this
3816823f3c5SBlaise Bourdin     */
3826823f3c5SBlaise Bourdin 
3839566063dSJacob Faibussowitsch     PetscCall(DMSetOutputSequenceNumber(dmUA,1,time));
3849566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(dmUA, &tmpVec));
3859566063dSJacob Faibussowitsch     PetscCall(VecCopy(UA, tmpVec));
3869566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) tmpVec, "U"));
3879566063dSJacob Faibussowitsch     PetscCall(VecView(tmpVec, viewer));
388c4762a1bSJed Brown     /* Reading nodal variables in Exodus file */
3899566063dSJacob Faibussowitsch     PetscCall(VecSet(tmpVec, -1000.0));
3909566063dSJacob Faibussowitsch     PetscCall(VecLoad(tmpVec, viewer));
3919566063dSJacob Faibussowitsch     PetscCall(VecAXPY(UA, -1.0, tmpVec));
3929566063dSJacob Faibussowitsch     PetscCall(VecNorm(UA, NORM_INFINITY, &norm));
39308401ef6SPierre Jolivet     PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON,PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "UAlpha ||Vin - Vout|| = %g", (double) norm);
3949566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(dmUA, &tmpVec));
395c4762a1bSJed Brown 
396c4762a1bSJed Brown     /* same thing with the UA2 Vec obtained from the superDM */
3979566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(dmUA2, &tmpVec));
3989566063dSJacob Faibussowitsch     PetscCall(VecCopy(UA2, tmpVec));
3999566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) tmpVec, "U"));
4009566063dSJacob Faibussowitsch     PetscCall(DMSetOutputSequenceNumber(dmUA2,2,time));
4019566063dSJacob Faibussowitsch     PetscCall(VecView(tmpVec, viewer));
402c4762a1bSJed Brown     /* Reading nodal variables in Exodus file */
4039566063dSJacob Faibussowitsch     PetscCall(VecSet(tmpVec, -1000.0));
4049566063dSJacob Faibussowitsch     PetscCall(VecLoad(tmpVec,viewer));
4059566063dSJacob Faibussowitsch     PetscCall(VecAXPY(UA2, -1.0, tmpVec));
4069566063dSJacob Faibussowitsch     PetscCall(VecNorm(UA2, NORM_INFINITY, &norm));
40708401ef6SPierre Jolivet     PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON,PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "UAlpha2 ||Vin - Vout|| = %g", (double) norm);
4089566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(dmUA2, &tmpVec));
409c4762a1bSJed Brown 
410c4762a1bSJed Brown     /* Building and saving Sigma
411c4762a1bSJed Brown        We set sigma_0 = rank (to see partitioning)
412c4762a1bSJed Brown               sigma_1 = cell set ID
4136823f3c5SBlaise Bourdin               sigma_2 = x_coordinate of the cell center of mass
4146823f3c5SBlaise Bourdin     */
4159566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateSection(dmS, &coordSection));
4169566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dmS, &coord));
4179566063dSJacob Faibussowitsch     PetscCall(DMGetLabelIdIS(dmS, "Cell Sets", &csIS));
4189566063dSJacob Faibussowitsch     PetscCall(DMGetLabelSize(dmS, "Cell Sets", &numCS));
4199566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(csIS, &csID));
420c4762a1bSJed Brown     for (set = 0; set < numCS; ++set) {
421c4762a1bSJed 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 */
422c4762a1bSJed Brown       IS              cellIS;
423c4762a1bSJed Brown       const PetscInt *cellID;
424c4762a1bSJed Brown       PetscInt        numCells, cell;
425c4762a1bSJed Brown       PetscScalar    *cval = NULL, *xyz  = NULL;
426c4762a1bSJed Brown       PetscInt        clSize, cdimCoord, c;
427c4762a1bSJed Brown 
4289566063dSJacob Faibussowitsch       PetscCall(DMGetStratumIS(dmS, "Cell Sets", csID[set], &cellIS));
4299566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(cellIS, &cellID));
4309566063dSJacob Faibussowitsch       PetscCall(ISGetSize(cellIS, &numCells));
431c4762a1bSJed Brown       for (cell = 0; cell < numCells; cell++) {
4329566063dSJacob Faibussowitsch         PetscCall(DMPlexVecGetClosure(dmS, NULL, S, cellID[cell], &clSize, &cval));
4339566063dSJacob Faibussowitsch         PetscCall(DMPlexVecGetClosure(dmS, coordSection, coord, cellID[cell], &cdimCoord, &xyz));
434c4762a1bSJed Brown         cval[0] = rank;
435c4762a1bSJed Brown         cval[1] = csID[set];
436c4762a1bSJed Brown         cval[2] = 0.;
437c4762a1bSJed Brown         for (c = 0; c < cdimCoord/sdim; c++) cval[2] += xyz[c*sdim];
438c4762a1bSJed Brown         cval[2] = cval[2] * sdim / cdimCoord;
4399566063dSJacob Faibussowitsch         PetscCall(DMPlexVecSetClosure(dmS, NULL, S, cellID[cell], cval, INSERT_ALL_VALUES));
440c4762a1bSJed Brown       }
4419566063dSJacob Faibussowitsch       PetscCall(DMPlexVecRestoreClosure(dmS, NULL, S, cellID[0], &clSize, &cval));
4429566063dSJacob Faibussowitsch       PetscCall(DMPlexVecRestoreClosure(dmS, coordSection, coord, cellID[0], NULL, &xyz));
4439566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(cellIS, &cellID));
4449566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&cellIS));
445c4762a1bSJed Brown     }
4469566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(csIS, &csID));
4479566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&csIS));
4489566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(S, NULL, "-s_vec_view"));
4496823f3c5SBlaise Bourdin 
450c4762a1bSJed Brown     /* Writing zonal variables in Exodus file */
4519566063dSJacob Faibussowitsch     PetscCall(DMSetOutputSequenceNumber(dmS,0,time));
4529566063dSJacob Faibussowitsch     PetscCall(VecView(S,viewer));
4536823f3c5SBlaise Bourdin 
454c4762a1bSJed Brown     /* Reading zonal variables in Exodus file */
4559566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(dmS, &tmpVec));
4569566063dSJacob Faibussowitsch     PetscCall(VecSet(tmpVec, -1000.0));
4579566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) tmpVec, "Sigma"));
4589566063dSJacob Faibussowitsch     PetscCall(VecLoad(tmpVec,viewer));
4599566063dSJacob Faibussowitsch     PetscCall(VecAXPY(S, -1.0, tmpVec));
4609566063dSJacob Faibussowitsch     PetscCall(VecNorm(S, NORM_INFINITY, &norm));
46108401ef6SPierre Jolivet     PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON,PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "Sigma ||Vin - Vout|| = %g", (double) norm);
4629566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(dmS, &tmpVec));
463c4762a1bSJed Brown   }
4649566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
465c4762a1bSJed Brown 
4669566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dmUA2, &UA2));
4679566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dmUA, &UA));
4689566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dmS,  &S));
4699566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dmA,  &A));
4709566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dmU,  &U));
471e2739ba6SAlexis Marboeuf   PetscCall(DMRestoreGlobalVector(pdm,   &X));
4729566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmU));PetscCall(ISDestroy(&isU));
4739566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmA));PetscCall(ISDestroy(&isA));
4749566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmS));PetscCall(ISDestroy(&isS));
4759566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmUA));PetscCall(ISDestroy(&isUA));
4769566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmUA2));
477e2739ba6SAlexis Marboeuf   PetscCall(DMDestroy(&pdm));
478e2739ba6SAlexis Marboeuf   if (size > 1) PetscCall(DMDestroy(&dm));
4799566063dSJacob Faibussowitsch   PetscCall(PetscFree2(pStartDepth, pEndDepth));
4809566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
481b122ec5aSJacob Faibussowitsch   return 0;
482c4762a1bSJed Brown }
483c4762a1bSJed Brown 
484c4762a1bSJed Brown /*TEST
485c4762a1bSJed Brown 
486c4762a1bSJed Brown   build:
487c4762a1bSJed Brown     requires: exodusii pnetcdf !complex
488e2739ba6SAlexis Marboeuf   # 2D seq
489e2739ba6SAlexis Marboeuf   test:
490e2739ba6SAlexis Marboeuf     suffix: 0
491*2e8d78feSBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -petscpartitioner_type simple -order 1
492e2739ba6SAlexis Marboeuf     #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
493e2739ba6SAlexis Marboeuf   test:
494e2739ba6SAlexis Marboeuf     suffix: 1
495*2e8d78feSBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -petscpartitioner_type simple -order 1
496e2739ba6SAlexis Marboeuf   test:
497e2739ba6SAlexis Marboeuf     suffix: 2
498*2e8d78feSBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -petscpartitioner_type simple -order 1
499e2739ba6SAlexis Marboeuf     #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
500e2739ba6SAlexis Marboeuf   test:
501e2739ba6SAlexis Marboeuf     suffix: 3
502*2e8d78feSBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -petscpartitioner_type simple -order 2
503e2739ba6SAlexis Marboeuf   test:
504e2739ba6SAlexis Marboeuf     suffix: 4
505*2e8d78feSBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -petscpartitioner_type simple -order 2
506e2739ba6SAlexis Marboeuf   test:
507e2739ba6SAlexis Marboeuf     suffix: 5
508*2e8d78feSBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -petscpartitioner_type simple -order 2
509c4762a1bSJed Brown 
510e2739ba6SAlexis Marboeuf   # 2D par
511e2739ba6SAlexis Marboeuf   test:
512e2739ba6SAlexis Marboeuf     suffix: 6
513e2739ba6SAlexis Marboeuf     nsize: 2
514*2e8d78feSBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -petscpartitioner_type simple -order 1
515e2739ba6SAlexis Marboeuf     #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
516e2739ba6SAlexis Marboeuf   test:
517e2739ba6SAlexis Marboeuf     suffix: 7
518e2739ba6SAlexis Marboeuf     nsize: 2
519*2e8d78feSBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -petscpartitioner_type simple -order 1
520e2739ba6SAlexis Marboeuf   test:
521e2739ba6SAlexis Marboeuf     suffix: 8
522e2739ba6SAlexis Marboeuf     nsize: 2
523*2e8d78feSBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -petscpartitioner_type simple -order 1
524e2739ba6SAlexis Marboeuf     #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: invalid dimension ID or name
525e2739ba6SAlexis Marboeuf   test:
526e2739ba6SAlexis Marboeuf     suffix: 9
527e2739ba6SAlexis Marboeuf     nsize: 2
528*2e8d78feSBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -petscpartitioner_type simple -order 2
529e2739ba6SAlexis Marboeuf   test:
530e2739ba6SAlexis Marboeuf     suffix: 10
531e2739ba6SAlexis Marboeuf     nsize: 2
532*2e8d78feSBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -petscpartitioner_type simple -order 2
533e2739ba6SAlexis Marboeuf   test:
534e2739ba6SAlexis Marboeuf     # Something is now broken with parallel read/write for wahtever shape H is
535e2739ba6SAlexis Marboeuf     TODO: broken
536e2739ba6SAlexis Marboeuf     suffix: 11
537e2739ba6SAlexis Marboeuf     nsize: 2
538*2e8d78feSBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -petscpartitioner_type simple -order 2
539c4762a1bSJed Brown 
540e2739ba6SAlexis Marboeuf   #3d seq
541e2739ba6SAlexis Marboeuf   test:
542e2739ba6SAlexis Marboeuf     suffix: 12
543*2e8d78feSBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -petscpartitioner_type simple -order 1
544e2739ba6SAlexis Marboeuf   test:
545e2739ba6SAlexis Marboeuf     suffix: 13
546*2e8d78feSBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -petscpartitioner_type simple -order 1
547e2739ba6SAlexis Marboeuf     #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
548e2739ba6SAlexis Marboeuf   test:
549e2739ba6SAlexis Marboeuf     suffix: 14
550*2e8d78feSBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -petscpartitioner_type simple -order 2
551e2739ba6SAlexis Marboeuf   test:
552e2739ba6SAlexis Marboeuf     suffix: 15
553*2e8d78feSBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -petscpartitioner_type simple -order 2
554e2739ba6SAlexis Marboeuf     #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
555e2739ba6SAlexis Marboeuf   #3d par
556e2739ba6SAlexis Marboeuf   test:
557e2739ba6SAlexis Marboeuf     suffix: 16
558e2739ba6SAlexis Marboeuf     nsize: 2
559*2e8d78feSBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -petscpartitioner_type simple -order 1
560e2739ba6SAlexis Marboeuf   test:
561e2739ba6SAlexis Marboeuf     suffix: 17
562e2739ba6SAlexis Marboeuf     nsize: 2
563*2e8d78feSBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -petscpartitioner_type simple -order 1
564e2739ba6SAlexis Marboeuf     #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
565e2739ba6SAlexis Marboeuf   test:
566e2739ba6SAlexis Marboeuf     suffix: 18
567e2739ba6SAlexis Marboeuf     nsize: 2
568*2e8d78feSBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -petscpartitioner_type simple -order 2
569e2739ba6SAlexis Marboeuf   test:
570e2739ba6SAlexis Marboeuf     suffix: 19
571e2739ba6SAlexis Marboeuf     nsize: 2
572*2e8d78feSBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -petscpartitioner_type simple -order 2
573e2739ba6SAlexis Marboeuf     #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
574c4762a1bSJed Brown 
575c4762a1bSJed Brown TEST*/
576