xref: /petsc/src/dm/impls/plex/tests/ex26.c (revision 792fecdfe9134cce4d631112660ddd34f063bc17)
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 
339566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv,NULL, help));
349566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
359566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
36d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FEM Layout Options", "ex26");
379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-i", "Filename to read", "ex26", ifilename, ifilename, sizeof(ifilename), NULL));
389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-o", "Filename to write", "ex26", ofilename, ofilename, sizeof(ofilename), NULL));
399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-order", "FEM polynomial order", "ex26", order, &order, NULL,1));
40d0609cedSBarry Smith   PetscOptionsEnd();
4163a3b9bcSJacob Faibussowitsch   PetscCheck((order >= 1) && (order <= 2),PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported polynomial order %" PetscInt_FMT " not in [1, 2]", order);
42c4762a1bSJed Brown 
436823f3c5SBlaise Bourdin   /* Read the mesh from a file in any supported format */
449566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename, NULL, PETSC_TRUE, &dm));
459566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeSetDefault(dm, PETSC_FALSE));
469566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
479566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
489566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &sdim));
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 */
619566063dSJacob Faibussowitsch     PetscCall(PetscViewerExodusIIOpen(PETSC_COMM_WORLD,ofilename,FILE_MODE_WRITE,&viewer));
626823f3c5SBlaise Bourdin     /* The long way would be */
636823f3c5SBlaise Bourdin     /*
649566063dSJacob Faibussowitsch       PetscCall(PetscViewerCreate(PETSC_COMM_WORLD,&viewer));
659566063dSJacob Faibussowitsch       PetscCall(PetscViewerSetType(viewer,PETSCVIEWEREXODUSII));
669566063dSJacob Faibussowitsch       PetscCall(PetscViewerFileSetMode(viewer,FILE_MODE_APPEND));
679566063dSJacob Faibussowitsch       PetscCall(PetscViewerFileSetName(viewer,ofilename));
686823f3c5SBlaise Bourdin     */
696823f3c5SBlaise Bourdin     /* set the mesh order */
709566063dSJacob Faibussowitsch     PetscCall(PetscViewerExodusIISetOrder(viewer,order));
719566063dSJacob Faibussowitsch     PetscCall(PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD));
726823f3c5SBlaise Bourdin     /*
736823f3c5SBlaise Bourdin       Notice how the exodus file is actually NOT open at this point (exoid is -1)
746aad120cSJose E. Roman       Since we are overwriting 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 */
799566063dSJacob Faibussowitsch     PetscCall(DMView(dm,viewer));
809566063dSJacob Faibussowitsch     PetscCall(PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD));
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;
11163a3b9bcSJacob Faibussowitsch     default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No layout for dimension %" PetscInt_FMT, sdim);
1126823f3c5SBlaise Bourdin     }
1139566063dSJacob Faibussowitsch     PetscCall(PetscViewerExodusIIGetId(viewer,&exoid));
114*792fecdfSBarry Smith     PetscCallExternal(ex_put_variable_param,exoid, EX_ELEM_BLOCK, numZonalVar);
115*792fecdfSBarry Smith     PetscCallExternal(ex_put_variable_names,exoid, EX_ELEM_BLOCK, numZonalVar, zonalVarName);
116*792fecdfSBarry Smith     PetscCallExternal(ex_put_variable_param,exoid, EX_NODAL, numNodalVar);
117*792fecdfSBarry Smith     PetscCallExternal(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     */
1249566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numZonalVar * numCS, &truthtable));
1256823f3c5SBlaise Bourdin     for (i = 0; i < numZonalVar * numCS; ++i) truthtable[i] = 1;
126*792fecdfSBarry Smith     PetscCallExternal(ex_put_truth_table,exoid, EX_ELEM_BLOCK, numCS, numZonalVar, truthtable);
1279566063dSJacob Faibussowitsch     PetscCall(PetscFree(truthtable));
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;
132*792fecdfSBarry Smith       PetscCallExternal(ex_put_time,exoid, step+1, &time);
1336823f3c5SBlaise Bourdin     }
1346823f3c5SBlaise Bourdin   }
1356823f3c5SBlaise Bourdin 
1366823f3c5SBlaise Bourdin   /* Create the main section containing all fields */
1379566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) dm), &section));
1389566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetNumFields(section, 3));
1399566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFieldName(section, fieldU, "U"));
1409566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFieldName(section, fieldA, "Alpha"));
1419566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFieldName(section, fieldS, "Sigma"));
1429566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
1439566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(section, pStart, pEnd));
1449566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(sdim+1, &pStartDepth, sdim+1, &pEndDepth));
1459566063dSJacob Faibussowitsch   for (d = 0; d <= sdim; ++d) PetscCall(DMPlexGetDepthStratum(dm, d, &pStartDepth[d], &pEndDepth[d]));
146c4762a1bSJed Brown   /* Vector field U, Scalar field Alpha, Tensor field Sigma */
1479566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFieldComponents(section, fieldU, sdim));
1489566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFieldComponents(section, fieldA, 1));
1499566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFieldComponents(section, fieldS, sdim*(sdim+1)/2));
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   /* Going through cell sets then cells, and setting up storage for the sections */
1529566063dSJacob Faibussowitsch   PetscCall(DMGetLabelSize(dm, "Cell Sets", &numCS));
1539566063dSJacob Faibussowitsch   PetscCall(DMGetLabelIdIS(dm, "Cell Sets", &csIS));
1549566063dSJacob Faibussowitsch   if (csIS) PetscCall(ISGetIndices(csIS, &csID));
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 
1609566063dSJacob Faibussowitsch     PetscCall(DMGetStratumSize(dm, "Cell Sets", csID[set], &numCells));
1619566063dSJacob Faibussowitsch     PetscCall(DMGetStratumIS(dm, "Cell Sets", csID[set], &cellIS));
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;
18763a3b9bcSJacob Faibussowitsch       default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No layout for dimension %" PetscInt_FMT, 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...  */
1929566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(cellIS, &cellID));
1939566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(dm, cellID[0], PETSC_TRUE, &closureSize, &closureA));
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;
23163a3b9bcSJacob Faibussowitsch         default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Unknown element with closure size %" PetscInt_FMT, closureSize);
232c4762a1bSJed Brown       }
2339566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(dm, cellID[0], PETSC_TRUE, &closureSize, &closureA));
234c4762a1bSJed Brown 
235c4762a1bSJed Brown       for (cell = 0; cell < numCells; cell++) {
236c4762a1bSJed Brown         PetscInt *closure = NULL;
237c4762a1bSJed Brown 
2389566063dSJacob Faibussowitsch         PetscCall(DMPlexGetTransitiveClosure(dm, cellID[cell], PETSC_TRUE, &closureSize, &closure));
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])) {
2439566063dSJacob Faibussowitsch               PetscCall(PetscSectionSetDof(section, closure[2*p], dofU[d]+dofA[d]+dofS[d]));
2449566063dSJacob Faibussowitsch               PetscCall(PetscSectionSetFieldDof(section, closure[2*p], fieldU, dofU[d]));
2459566063dSJacob Faibussowitsch               PetscCall(PetscSectionSetFieldDof(section, closure[2*p], fieldA, dofA[d]));
2469566063dSJacob Faibussowitsch               PetscCall(PetscSectionSetFieldDof(section, closure[2*p], fieldS, dofS[d]));
247c4762a1bSJed Brown             }
248c4762a1bSJed Brown           }
249c4762a1bSJed Brown         }
2509566063dSJacob Faibussowitsch         PetscCall(DMPlexRestoreTransitiveClosure(dm, cellID[cell], PETSC_TRUE, &closureSize, &closure));
251c4762a1bSJed Brown       }
2529566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(cellIS, &cellID));
2539566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&cellIS));
254c4762a1bSJed Brown     }
255c4762a1bSJed Brown   }
2569566063dSJacob Faibussowitsch   if (csIS) PetscCall(ISRestoreIndices(csIS, &csID));
2579566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&csIS));
2589566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(section));
2599566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(dm, section));
2609566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject) section, NULL, "-dm_section_view"));
2619566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&section));
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   {
264c4762a1bSJed Brown     DM               pdm;
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));
273c4762a1bSJed Brown     if (pdm) {
2749566063dSJacob Faibussowitsch       PetscCall(DMPlexSetMigrationSF(pdm,migrationSF));
2759566063dSJacob Faibussowitsch       PetscCall(PetscSFDestroy(&migrationSF));
2769566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&dm));
277c4762a1bSJed Brown       dm = pdm;
2789566063dSJacob Faibussowitsch       PetscCall(DMViewFromOptions(dm,NULL,"-dm_view"));
279c4762a1bSJed Brown     }
280c4762a1bSJed Brown   }
281c4762a1bSJed Brown 
282c4762a1bSJed Brown   /* Get DM and IS for each field of dm */
2839566063dSJacob Faibussowitsch   PetscCall(DMCreateSubDM(dm, 1, &fieldU, &isU,  &dmU));
2849566063dSJacob Faibussowitsch   PetscCall(DMCreateSubDM(dm, 1, &fieldA, &isA,  &dmA));
2859566063dSJacob Faibussowitsch   PetscCall(DMCreateSubDM(dm, 1, &fieldS, &isS,  &dmS));
2869566063dSJacob Faibussowitsch   PetscCall(DMCreateSubDM(dm, 2, fieldUA, &isUA, &dmUA));
287c4762a1bSJed Brown 
2889566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2,&dmList));
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;
2959566063dSJacob Faibussowitsch   PetscCall(DMCreateSuperDM(dmList,2,NULL,&dmUA2));
296c4762a1bSJed Brown   dmU->useNatural = PETSC_TRUE;
2979566063dSJacob Faibussowitsch   PetscCall(PetscFree(dmList));
298c4762a1bSJed Brown 
2999566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm,   &X));
3009566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dmU,  &U));
3019566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dmA,  &A));
3029566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dmS,  &S));
3039566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dmUA, &UA));
3049566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dmUA2, &UA2));
305c4762a1bSJed Brown 
3069566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) U,  "U"));
3079566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) A,  "Alpha"));
3089566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) S,  "Sigma"));
3099566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) UA, "UAlpha"));
3109566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) UA2, "UAlpha2"));
3119566063dSJacob Faibussowitsch   PetscCall(VecSet(X, -111.));
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 
3229566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dmUA, &sectionUA));
3239566063dSJacob Faibussowitsch     PetscCall(DMGetLocalVector(dmUA, &UALoc));
3249566063dSJacob Faibussowitsch     PetscCall(VecGetArray(UALoc, &cval));
3259566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateSection(dmUA, &coordSection));
3269566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dmUA, &coord));
3279566063dSJacob Faibussowitsch     PetscCall(DMPlexGetChart(dmUA, &pStart, &pEnd));
3286823f3c5SBlaise Bourdin 
329c4762a1bSJed Brown     for (p = pStart; p < pEnd; ++p) {
330c4762a1bSJed Brown       PetscInt dofUA, offUA;
331c4762a1bSJed Brown 
3329566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(sectionUA, p, &dofUA));
333c4762a1bSJed Brown       if (dofUA > 0) {
3346823f3c5SBlaise Bourdin         xyz=NULL;
3359566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionUA, p, &offUA));
3369566063dSJacob Faibussowitsch         PetscCall(DMPlexVecGetClosure(dmUA, coordSection, coord, p, &clSize, &xyz));
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         }
3469566063dSJacob Faibussowitsch         PetscCall(DMPlexVecRestoreClosure(dmUA, coordSection, coord, p, &clSize, &xyz));
347c4762a1bSJed Brown       }
348c4762a1bSJed Brown     }
3499566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(UALoc, &cval));
3509566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalBegin(dmUA, UALoc, INSERT_VALUES, UA));
3519566063dSJacob Faibussowitsch     PetscCall(DMLocalToGlobalEnd(dmUA, UALoc, INSERT_VALUES, UA));
3529566063dSJacob Faibussowitsch     PetscCall(DMRestoreLocalVector(dmUA, &UALoc));
3536823f3c5SBlaise Bourdin 
354c4762a1bSJed Brown     /* Update X */
3559566063dSJacob Faibussowitsch     PetscCall(VecISCopy(X, isUA, SCATTER_FORWARD, UA));
3569566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(UA, NULL, "-ua_vec_view"));
3576823f3c5SBlaise Bourdin 
358c4762a1bSJed Brown     /* Restrict to U and Alpha */
3599566063dSJacob Faibussowitsch     PetscCall(VecISCopy(X, isU, SCATTER_REVERSE, U));
3609566063dSJacob Faibussowitsch     PetscCall(VecISCopy(X, isA, SCATTER_REVERSE, A));
3616823f3c5SBlaise Bourdin 
362c4762a1bSJed Brown     /* restrict to UA2 */
3639566063dSJacob Faibussowitsch     PetscCall(VecISCopy(X, isUA, SCATTER_REVERSE, UA2));
3649566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(UA2, NULL, "-ua2_vec_view"));
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 */
3759566063dSJacob Faibussowitsch     PetscCall(DMSetOutputSequenceNumber(dmU,0,time));
3769566063dSJacob Faibussowitsch     PetscCall(DMSetOutputSequenceNumber(dmA,0,time));
3776823f3c5SBlaise Bourdin 
3789566063dSJacob Faibussowitsch     PetscCall(VecView(U, viewer));
3799566063dSJacob Faibussowitsch     PetscCall(VecView(A, viewer));
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 
3879566063dSJacob Faibussowitsch     PetscCall(DMSetOutputSequenceNumber(dmUA,1,time));
3889566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(dmUA, &tmpVec));
3899566063dSJacob Faibussowitsch     PetscCall(VecCopy(UA, tmpVec));
3909566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) tmpVec, "U"));
3919566063dSJacob Faibussowitsch     PetscCall(VecView(tmpVec, viewer));
392c4762a1bSJed Brown     /* Reading nodal variables in Exodus file */
3939566063dSJacob Faibussowitsch     PetscCall(VecSet(tmpVec, -1000.0));
3949566063dSJacob Faibussowitsch     PetscCall(VecLoad(tmpVec, viewer));
3959566063dSJacob Faibussowitsch     PetscCall(VecAXPY(UA, -1.0, tmpVec));
3969566063dSJacob Faibussowitsch     PetscCall(VecNorm(UA, NORM_INFINITY, &norm));
39708401ef6SPierre Jolivet     PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON,PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "UAlpha ||Vin - Vout|| = %g", (double) norm);
3989566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(dmUA, &tmpVec));
399c4762a1bSJed Brown 
400c4762a1bSJed Brown     /* same thing with the UA2 Vec obtained from the superDM */
4019566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(dmUA2, &tmpVec));
4029566063dSJacob Faibussowitsch     PetscCall(VecCopy(UA2, tmpVec));
4039566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) tmpVec, "U"));
4049566063dSJacob Faibussowitsch     PetscCall(DMSetOutputSequenceNumber(dmUA2,2,time));
4059566063dSJacob Faibussowitsch     PetscCall(VecView(tmpVec, viewer));
406c4762a1bSJed Brown     /* Reading nodal variables in Exodus file */
4079566063dSJacob Faibussowitsch     PetscCall(VecSet(tmpVec, -1000.0));
4089566063dSJacob Faibussowitsch     PetscCall(VecLoad(tmpVec,viewer));
4099566063dSJacob Faibussowitsch     PetscCall(VecAXPY(UA2, -1.0, tmpVec));
4109566063dSJacob Faibussowitsch     PetscCall(VecNorm(UA2, NORM_INFINITY, &norm));
41108401ef6SPierre Jolivet     PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON,PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "UAlpha2 ||Vin - Vout|| = %g", (double) norm);
4129566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(dmUA2, &tmpVec));
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     */
4199566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateSection(dmS, &coordSection));
4209566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dmS, &coord));
4219566063dSJacob Faibussowitsch     PetscCall(DMGetLabelIdIS(dmS, "Cell Sets", &csIS));
4229566063dSJacob Faibussowitsch     PetscCall(DMGetLabelSize(dmS, "Cell Sets", &numCS));
4239566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(csIS, &csID));
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 
4329566063dSJacob Faibussowitsch       PetscCall(DMGetStratumIS(dmS, "Cell Sets", csID[set], &cellIS));
4339566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(cellIS, &cellID));
4349566063dSJacob Faibussowitsch       PetscCall(ISGetSize(cellIS, &numCells));
435c4762a1bSJed Brown       for (cell = 0; cell < numCells; cell++) {
4369566063dSJacob Faibussowitsch         PetscCall(DMPlexVecGetClosure(dmS, NULL, S, cellID[cell], &clSize, &cval));
4379566063dSJacob Faibussowitsch         PetscCall(DMPlexVecGetClosure(dmS, coordSection, coord, cellID[cell], &cdimCoord, &xyz));
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;
4439566063dSJacob Faibussowitsch         PetscCall(DMPlexVecSetClosure(dmS, NULL, S, cellID[cell], cval, INSERT_ALL_VALUES));
444c4762a1bSJed Brown       }
4459566063dSJacob Faibussowitsch       PetscCall(DMPlexVecRestoreClosure(dmS, NULL, S, cellID[0], &clSize, &cval));
4469566063dSJacob Faibussowitsch       PetscCall(DMPlexVecRestoreClosure(dmS, coordSection, coord, cellID[0], NULL, &xyz));
4479566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(cellIS, &cellID));
4489566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&cellIS));
449c4762a1bSJed Brown     }
4509566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(csIS, &csID));
4519566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&csIS));
4529566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(S, NULL, "-s_vec_view"));
4536823f3c5SBlaise Bourdin 
454c4762a1bSJed Brown     /* Writing zonal variables in Exodus file */
4559566063dSJacob Faibussowitsch     PetscCall(DMSetOutputSequenceNumber(dmS,0,time));
4569566063dSJacob Faibussowitsch     PetscCall(VecView(S,viewer));
4576823f3c5SBlaise Bourdin 
458c4762a1bSJed Brown     /* Reading zonal variables in Exodus file */
4599566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(dmS, &tmpVec));
4609566063dSJacob Faibussowitsch     PetscCall(VecSet(tmpVec, -1000.0));
4619566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) tmpVec, "Sigma"));
4629566063dSJacob Faibussowitsch     PetscCall(VecLoad(tmpVec,viewer));
4639566063dSJacob Faibussowitsch     PetscCall(VecAXPY(S, -1.0, tmpVec));
4649566063dSJacob Faibussowitsch     PetscCall(VecNorm(S, NORM_INFINITY, &norm));
46508401ef6SPierre Jolivet     PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON,PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "Sigma ||Vin - Vout|| = %g", (double) norm);
4669566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(dmS, &tmpVec));
467c4762a1bSJed Brown   }
4689566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
469c4762a1bSJed Brown 
4709566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dmUA2, &UA2));
4719566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dmUA, &UA));
4729566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dmS,  &S));
4739566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dmA,  &A));
4749566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dmU,  &U));
4759566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm,   &X));
4769566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmU)); PetscCall(ISDestroy(&isU));
4779566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmA)); PetscCall(ISDestroy(&isA));
4789566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmS)); PetscCall(ISDestroy(&isS));
4799566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmUA));PetscCall(ISDestroy(&isUA));
4809566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmUA2));
4819566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
4829566063dSJacob Faibussowitsch   PetscCall(PetscFree2(pStartDepth, pEndDepth));
4839566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
484b122ec5aSJacob Faibussowitsch   return 0;
485c4762a1bSJed Brown }
486c4762a1bSJed Brown 
487c4762a1bSJed Brown /*TEST
488c4762a1bSJed Brown 
489c4762a1bSJed Brown   build:
490c4762a1bSJed Brown     requires: exodusii pnetcdf !complex
491445019dbSMatthew G. Knepley   ### # 2D seq
492445019dbSMatthew G. Knepley   ### test:
493445019dbSMatthew G. Knepley   ###   suffix: 0
494445019dbSMatthew G. Knepley   ###   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
495445019dbSMatthew G. Knepley   ###   #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
496445019dbSMatthew G. Knepley   ### test:
497445019dbSMatthew G. Knepley   ###   suffix: 1
498445019dbSMatthew G. Knepley   ###   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
499445019dbSMatthew G. Knepley   ### test:
500445019dbSMatthew G. Knepley   ###   suffix: 2
501445019dbSMatthew G. Knepley   ###   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
502445019dbSMatthew G. Knepley   ###   #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
503445019dbSMatthew G. Knepley   ### test:
504445019dbSMatthew G. Knepley   ###   suffix: 3
505445019dbSMatthew G. Knepley   ###   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
506445019dbSMatthew G. Knepley   ### test:
507445019dbSMatthew G. Knepley   ###   suffix: 4
508445019dbSMatthew G. Knepley   ###   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
509445019dbSMatthew G. Knepley   ### test:
510445019dbSMatthew G. Knepley   ###   suffix: 5
511445019dbSMatthew G. Knepley   ###   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 
513445019dbSMatthew G. Knepley   ### # 2D par
514445019dbSMatthew G. Knepley   ### test:
515445019dbSMatthew G. Knepley   ###   suffix: 6
516445019dbSMatthew G. Knepley   ###   nsize: 2
517445019dbSMatthew G. Knepley   ###   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
518445019dbSMatthew G. Knepley   ###   #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
519445019dbSMatthew G. Knepley   ### test:
520445019dbSMatthew G. Knepley   ###   suffix: 7
521445019dbSMatthew G. Knepley   ###   nsize: 2
522445019dbSMatthew G. Knepley   ###   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
523445019dbSMatthew G. Knepley   ### test:
524445019dbSMatthew G. Knepley   ###   suffix: 8
525445019dbSMatthew G. Knepley   ###   nsize: 2
526445019dbSMatthew G. Knepley   ###   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
527445019dbSMatthew G. Knepley   ###   #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: invalid dimension ID or name
528445019dbSMatthew G. Knepley   ### test:
529445019dbSMatthew G. Knepley   ###   suffix: 9
530445019dbSMatthew G. Knepley   ###   nsize: 2
531445019dbSMatthew G. Knepley   ###   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
532445019dbSMatthew G. Knepley   ### test:
533445019dbSMatthew G. Knepley   ###   suffix: 10
534445019dbSMatthew G. Knepley   ###   nsize: 2
535445019dbSMatthew G. Knepley   ###   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
536445019dbSMatthew G. Knepley   ### test:
537445019dbSMatthew G. Knepley   ###   # Something is now broken with parallel read/write for wahtever shape H is
538445019dbSMatthew G. Knepley   ###   TODO: broken
539445019dbSMatthew G. Knepley   ###   suffix: 11
540445019dbSMatthew G. Knepley   ###   nsize: 2
541445019dbSMatthew G. Knepley   ###   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 
543445019dbSMatthew G. Knepley   ### #3d seq
544445019dbSMatthew G. Knepley   ### test:
545445019dbSMatthew G. Knepley   ###   suffix: 12
546445019dbSMatthew G. Knepley   ###   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
547445019dbSMatthew G. Knepley   ### test:
548445019dbSMatthew G. Knepley   ###   suffix: 13
549445019dbSMatthew G. Knepley   ###   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
550445019dbSMatthew G. Knepley   ###   #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
551445019dbSMatthew G. Knepley   ### test:
552445019dbSMatthew G. Knepley   ###   suffix: 14
553445019dbSMatthew G. Knepley   ###   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
554445019dbSMatthew G. Knepley   ### test:
555445019dbSMatthew G. Knepley   ###   suffix: 15
556445019dbSMatthew G. Knepley   ###   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
557445019dbSMatthew G. Knepley   ###   #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
558445019dbSMatthew G. Knepley   ### #3d par
559445019dbSMatthew G. Knepley   ### test:
560445019dbSMatthew G. Knepley   ###   suffix: 16
561445019dbSMatthew G. Knepley   ###   nsize: 2
562445019dbSMatthew G. Knepley   ###   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
563445019dbSMatthew G. Knepley   ### test:
564445019dbSMatthew G. Knepley   ###   suffix: 17
565445019dbSMatthew G. Knepley   ###   nsize: 2
566445019dbSMatthew G. Knepley   ###   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
567445019dbSMatthew G. Knepley   ###   #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
568445019dbSMatthew G. Knepley   ### test:
569445019dbSMatthew G. Knepley   ###   suffix: 18
570445019dbSMatthew G. Knepley   ###   nsize: 2
571445019dbSMatthew G. Knepley   ###   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
572445019dbSMatthew G. Knepley   ### test:
573445019dbSMatthew G. Knepley   ###   suffix: 19
574445019dbSMatthew G. Knepley   ###   nsize: 2
575445019dbSMatthew G. Knepley   ###   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
576445019dbSMatthew G. Knepley   ###   #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