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
main(int argc,char ** argv)14d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
15d71ae5a4SJacob Faibussowitsch {
16e2739ba6SAlexis Marboeuf DM dm, pdm, dmU, dmA, dmS, dmUA, dmUA2, *dmList;
17c4762a1bSJed Brown Vec X, U, A, S, UA, UA2;
18c4762a1bSJed Brown IS isU, isA, isS, isUA;
19c4762a1bSJed Brown PetscSection section;
20c4762a1bSJed Brown const PetscInt fieldU = 0;
21c4762a1bSJed Brown const PetscInt fieldA = 2;
22c4762a1bSJed Brown const PetscInt fieldS = 1;
23c4762a1bSJed Brown const PetscInt fieldUA[2] = {0, 2};
24c4762a1bSJed Brown char ifilename[PETSC_MAX_PATH_LEN], ofilename[PETSC_MAX_PATH_LEN];
25c4762a1bSJed Brown int exoid = -1;
26c4762a1bSJed Brown IS csIS;
27c4762a1bSJed Brown const PetscInt *csID;
28c4762a1bSJed Brown PetscInt *pStartDepth, *pEndDepth;
29c4762a1bSJed Brown PetscInt order = 1;
30c4762a1bSJed Brown PetscInt sdim, d, pStart, pEnd, p, numCS, set;
31c4762a1bSJed Brown PetscMPIInt rank, size;
326823f3c5SBlaise Bourdin PetscViewer viewer;
3349c89c76SBlaise Bourdin const char *nodalVarName2D[4] = {"U_x", "U_y", "Alpha", "Beta"};
3449c89c76SBlaise Bourdin const char *zonalVarName2D[3] = {"Sigma_11", "Sigma_22", "Sigma_12"};
3549c89c76SBlaise Bourdin const char *nodalVarName3D[5] = {"U_x", "U_y", "U_z", "Alpha", "Beta"};
3649c89c76SBlaise Bourdin const char *zonalVarName3D[6] = {"Sigma_11", "Sigma_22", "Sigma_33", "Sigma_23", "Sigma_13", "Sigma_12"};
37c4762a1bSJed Brown
38327415f7SBarry Smith PetscFunctionBeginUser;
399566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
42d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FEM Layout Options", "ex26");
439566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-i", "Filename to read", "ex26", ifilename, ifilename, sizeof(ifilename), NULL));
449566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-o", "Filename to write", "ex26", ofilename, ofilename, sizeof(ofilename), NULL));
459566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-order", "FEM polynomial order", "ex26", order, &order, NULL, 1));
46d0609cedSBarry Smith PetscOptionsEnd();
4763a3b9bcSJacob Faibussowitsch PetscCheck((order >= 1) && (order <= 2), PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported polynomial order %" PetscInt_FMT " not in [1, 2]", order);
48c4762a1bSJed Brown
496823f3c5SBlaise Bourdin /* Read the mesh from a file in any supported format */
509566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename, NULL, PETSC_TRUE, &dm));
519566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(dm, PETSC_FALSE));
529566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm));
539566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
549566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &sdim));
55c4762a1bSJed Brown
566823f3c5SBlaise Bourdin /* Create the exodus result file */
576823f3c5SBlaise Bourdin {
586823f3c5SBlaise Bourdin PetscInt numstep = 3, step;
596823f3c5SBlaise Bourdin int *truthtable;
60*0a5cf5d8SBlaise Bourdin int numNodalVar, numZonalVar, i;
616823f3c5SBlaise Bourdin
62a5b23f4aSJose E. Roman /* enable exodus debugging information */
636823f3c5SBlaise Bourdin ex_opts(EX_VERBOSE | EX_DEBUG);
646823f3c5SBlaise Bourdin /* Create the exodus file */
659566063dSJacob Faibussowitsch PetscCall(PetscViewerExodusIIOpen(PETSC_COMM_WORLD, ofilename, FILE_MODE_WRITE, &viewer));
666823f3c5SBlaise Bourdin /* The long way would be */
676823f3c5SBlaise Bourdin /*
689566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PETSC_COMM_WORLD,&viewer));
699566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(viewer,PETSCVIEWEREXODUSII));
709566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(viewer,FILE_MODE_APPEND));
719566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(viewer,ofilename));
726823f3c5SBlaise Bourdin */
736823f3c5SBlaise Bourdin /* set the mesh order */
749566063dSJacob Faibussowitsch PetscCall(PetscViewerExodusIISetOrder(viewer, order));
759566063dSJacob Faibussowitsch PetscCall(PetscViewerView(viewer, PETSC_VIEWER_STDOUT_WORLD));
766823f3c5SBlaise Bourdin /*
776823f3c5SBlaise Bourdin Notice how the exodus file is actually NOT open at this point (exoid is -1)
78d5b43468SJose E. Roman Since we are overwriting the file (mode is FILE_MODE_WRITE), we are going to have to
796823f3c5SBlaise Bourdin write the geometry (the DM), which can only be done on a brand new file.
806823f3c5SBlaise Bourdin */
816823f3c5SBlaise Bourdin
826823f3c5SBlaise Bourdin /* Save the geometry to the file, erasing all previous content */
839566063dSJacob Faibussowitsch PetscCall(DMView(dm, viewer));
849566063dSJacob Faibussowitsch PetscCall(PetscViewerView(viewer, PETSC_VIEWER_STDOUT_WORLD));
856823f3c5SBlaise Bourdin /*
866823f3c5SBlaise Bourdin Note how the exodus file is now open
876823f3c5SBlaise Bourdin */
8849c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetId(viewer, &exoid));
896823f3c5SBlaise Bourdin
906823f3c5SBlaise Bourdin /* "Format" the exodus result file, i.e. allocate space for nodal and zonal variables */
916823f3c5SBlaise Bourdin switch (sdim) {
926823f3c5SBlaise Bourdin case 2:
9349c89c76SBlaise Bourdin numNodalVar = 4;
946823f3c5SBlaise Bourdin numZonalVar = 3;
9549c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIISetZonalVariable(viewer, numZonalVar));
9649c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIISetZonalVariableNames(viewer, zonalVarName2D));
9749c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIISetNodalVariable(viewer, numNodalVar));
9849c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIISetNodalVariableNames(viewer, nodalVarName2D));
996823f3c5SBlaise Bourdin break;
1006823f3c5SBlaise Bourdin case 3:
10149c89c76SBlaise Bourdin numNodalVar = 5;
1026823f3c5SBlaise Bourdin numZonalVar = 6;
10349c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIISetZonalVariable(viewer, numZonalVar));
10449c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIISetZonalVariableNames(viewer, zonalVarName3D));
10549c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIISetNodalVariable(viewer, numNodalVar));
10649c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIISetNodalVariableNames(viewer, nodalVarName3D));
1076823f3c5SBlaise Bourdin break;
108d71ae5a4SJacob Faibussowitsch default:
109d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No layout for dimension %" PetscInt_FMT, sdim);
1106823f3c5SBlaise Bourdin }
11149c89c76SBlaise Bourdin PetscCall(PetscViewerView(viewer, PETSC_VIEWER_STDOUT_WORLD));
1126823f3c5SBlaise Bourdin
1136823f3c5SBlaise Bourdin /*
1146823f3c5SBlaise Bourdin An exodusII truth table specifies which fields are saved at which time step
1156823f3c5SBlaise Bourdin It speeds up I/O but reserving space for fields in the file ahead of time.
1166823f3c5SBlaise Bourdin */
11749c89c76SBlaise Bourdin numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK);
1189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numZonalVar * numCS, &truthtable));
1196823f3c5SBlaise Bourdin for (i = 0; i < numZonalVar * numCS; ++i) truthtable[i] = 1;
120792fecdfSBarry Smith PetscCallExternal(ex_put_truth_table, exoid, EX_ELEM_BLOCK, numCS, numZonalVar, truthtable);
1219566063dSJacob Faibussowitsch PetscCall(PetscFree(truthtable));
1226823f3c5SBlaise Bourdin
1236823f3c5SBlaise Bourdin /* Writing time step information in the file. Note that this is currently broken in the exodus library for netcdf4 (HDF5-based) files */
1246823f3c5SBlaise Bourdin for (step = 0; step < numstep; ++step) {
1256823f3c5SBlaise Bourdin PetscReal time = step;
126792fecdfSBarry Smith PetscCallExternal(ex_put_time, exoid, step + 1, &time);
1276823f3c5SBlaise Bourdin }
1286823f3c5SBlaise Bourdin }
1296823f3c5SBlaise Bourdin
1306823f3c5SBlaise Bourdin /* Create the main section containing all fields */
1319566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ion));
1329566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(section, 3));
1339566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(section, fieldU, "U"));
1349566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(section, fieldA, "Alpha"));
1359566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(section, fieldS, "Sigma"));
1369566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
1379566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(section, pStart, pEnd));
1389566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(sdim + 1, &pStartDepth, sdim + 1, &pEndDepth));
13948a46eb9SPierre Jolivet for (d = 0; d <= sdim; ++d) PetscCall(DMPlexGetDepthStratum(dm, d, &pStartDepth[d], &pEndDepth[d]));
140c4762a1bSJed Brown /* Vector field U, Scalar field Alpha, Tensor field Sigma */
1419566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(section, fieldU, sdim));
1429566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(section, fieldA, 1));
1439566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(section, fieldS, sdim * (sdim + 1) / 2));
144c4762a1bSJed Brown
145c4762a1bSJed Brown /* Going through cell sets then cells, and setting up storage for the sections */
1469566063dSJacob Faibussowitsch PetscCall(DMGetLabelSize(dm, "Cell Sets", &numCS));
1479566063dSJacob Faibussowitsch PetscCall(DMGetLabelIdIS(dm, "Cell Sets", &csIS));
14848a46eb9SPierre Jolivet if (csIS) PetscCall(ISGetIndices(csIS, &csID));
149c4762a1bSJed Brown for (set = 0; set < numCS; set++) {
150c4762a1bSJed Brown IS cellIS;
151c4762a1bSJed Brown const PetscInt *cellID;
152c4762a1bSJed Brown PetscInt numCells, cell, closureSize, *closureA = NULL;
153c4762a1bSJed Brown
1549566063dSJacob Faibussowitsch PetscCall(DMGetStratumSize(dm, "Cell Sets", csID[set], &numCells));
1559566063dSJacob Faibussowitsch PetscCall(DMGetStratumIS(dm, "Cell Sets", csID[set], &cellIS));
156c4762a1bSJed Brown if (numCells > 0) {
157c4762a1bSJed Brown /* dof layout ordered by increasing height in the DAG: cell, face, edge, vertex */
158c4762a1bSJed Brown PetscInt dofUP1Tri[] = {2, 0, 0};
159c4762a1bSJed Brown PetscInt dofAP1Tri[] = {1, 0, 0};
160c4762a1bSJed Brown PetscInt dofUP2Tri[] = {2, 2, 0};
161c4762a1bSJed Brown PetscInt dofAP2Tri[] = {1, 1, 0};
162c4762a1bSJed Brown PetscInt dofUP1Quad[] = {2, 0, 0};
163c4762a1bSJed Brown PetscInt dofAP1Quad[] = {1, 0, 0};
164c4762a1bSJed Brown PetscInt dofUP2Quad[] = {2, 2, 2};
165c4762a1bSJed Brown PetscInt dofAP2Quad[] = {1, 1, 1};
166c4762a1bSJed Brown PetscInt dofS2D[] = {0, 0, 3};
167c4762a1bSJed Brown PetscInt dofUP1Tet[] = {3, 0, 0, 0};
168c4762a1bSJed Brown PetscInt dofAP1Tet[] = {1, 0, 0, 0};
169c4762a1bSJed Brown PetscInt dofUP2Tet[] = {3, 3, 0, 0};
170c4762a1bSJed Brown PetscInt dofAP2Tet[] = {1, 1, 0, 0};
171c4762a1bSJed Brown PetscInt dofUP1Hex[] = {3, 0, 0, 0};
172c4762a1bSJed Brown PetscInt dofAP1Hex[] = {1, 0, 0, 0};
173c4762a1bSJed Brown PetscInt dofUP2Hex[] = {3, 3, 3, 3};
174c4762a1bSJed Brown PetscInt dofAP2Hex[] = {1, 1, 1, 1};
175c4762a1bSJed Brown PetscInt dofS3D[] = {0, 0, 0, 6};
176c4762a1bSJed Brown PetscInt *dofU, *dofA, *dofS;
177c4762a1bSJed Brown
178c4762a1bSJed Brown switch (sdim) {
179d71ae5a4SJacob Faibussowitsch case 2:
180d71ae5a4SJacob Faibussowitsch dofS = dofS2D;
181d71ae5a4SJacob Faibussowitsch break;
182d71ae5a4SJacob Faibussowitsch case 3:
183d71ae5a4SJacob Faibussowitsch dofS = dofS3D;
184d71ae5a4SJacob Faibussowitsch break;
185d71ae5a4SJacob Faibussowitsch default:
186d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No layout for dimension %" PetscInt_FMT, sdim);
187c4762a1bSJed Brown }
188c4762a1bSJed Brown
189c4762a1bSJed Brown /* Identify cell type based on closure size only. This works for Tri/Tet/Quad/Hex meshes
190c4762a1bSJed Brown It will not be enough to identify more exotic elements like pyramid or prisms... */
1919566063dSJacob Faibussowitsch PetscCall(ISGetIndices(cellIS, &cellID));
1929566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, cellID[0], PETSC_TRUE, &closureSize, &closureA));
193c4762a1bSJed Brown switch (closureSize) {
194c4762a1bSJed Brown case 7: /* Tri */
195c4762a1bSJed Brown if (order == 1) {
196c4762a1bSJed Brown dofU = dofUP1Tri;
197c4762a1bSJed Brown dofA = dofAP1Tri;
198c4762a1bSJed Brown } else {
199c4762a1bSJed Brown dofU = dofUP2Tri;
200c4762a1bSJed Brown dofA = dofAP2Tri;
201c4762a1bSJed Brown }
202c4762a1bSJed Brown break;
203c4762a1bSJed Brown case 9: /* Quad */
204c4762a1bSJed Brown if (order == 1) {
205c4762a1bSJed Brown dofU = dofUP1Quad;
206c4762a1bSJed Brown dofA = dofAP1Quad;
207c4762a1bSJed Brown } else {
208c4762a1bSJed Brown dofU = dofUP2Quad;
209c4762a1bSJed Brown dofA = dofAP2Quad;
210c4762a1bSJed Brown }
211c4762a1bSJed Brown break;
212c4762a1bSJed Brown case 15: /* Tet */
213c4762a1bSJed Brown if (order == 1) {
214c4762a1bSJed Brown dofU = dofUP1Tet;
215c4762a1bSJed Brown dofA = dofAP1Tet;
216c4762a1bSJed Brown } else {
217c4762a1bSJed Brown dofU = dofUP2Tet;
218c4762a1bSJed Brown dofA = dofAP2Tet;
219c4762a1bSJed Brown }
220c4762a1bSJed Brown break;
221c4762a1bSJed Brown case 27: /* Hex */
222c4762a1bSJed Brown if (order == 1) {
223c4762a1bSJed Brown dofU = dofUP1Hex;
224c4762a1bSJed Brown dofA = dofAP1Hex;
225c4762a1bSJed Brown } else {
226c4762a1bSJed Brown dofU = dofUP2Hex;
227c4762a1bSJed Brown dofA = dofAP2Hex;
228c4762a1bSJed Brown }
229c4762a1bSJed Brown break;
230d71ae5a4SJacob Faibussowitsch default:
231d71ae5a4SJacob Faibussowitsch 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 }
25648a46eb9SPierre Jolivet 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(§ion));
262c4762a1bSJed Brown
263c4762a1bSJed Brown {
264c4762a1bSJed Brown PetscSF migrationSF;
265c4762a1bSJed Brown PetscInt ovlp = 0;
266c4762a1bSJed Brown PetscPartitioner part;
267c4762a1bSJed Brown
2689566063dSJacob Faibussowitsch PetscCall(DMSetUseNatural(dm, PETSC_TRUE));
2699566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(dm, &part));
2709566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetFromOptions(part));
2719566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(dm, ovlp, &migrationSF, &pdm));
272e2739ba6SAlexis Marboeuf if (!pdm) pdm = dm;
273e2739ba6SAlexis Marboeuf /* Set the migrationSF is mandatory since useNatural = PETSC_TRUE */
274e2739ba6SAlexis Marboeuf if (migrationSF) {
2759566063dSJacob Faibussowitsch PetscCall(DMPlexSetMigrationSF(pdm, migrationSF));
2769566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&migrationSF));
277c4762a1bSJed Brown }
278e2739ba6SAlexis Marboeuf PetscCall(DMViewFromOptions(pdm, NULL, "-dm_view"));
279c4762a1bSJed Brown }
280c4762a1bSJed Brown
281c4762a1bSJed Brown /* Get DM and IS for each field of dm */
282e2739ba6SAlexis Marboeuf PetscCall(DMCreateSubDM(pdm, 1, &fieldU, &isU, &dmU));
283e2739ba6SAlexis Marboeuf PetscCall(DMCreateSubDM(pdm, 1, &fieldA, &isA, &dmA));
284e2739ba6SAlexis Marboeuf PetscCall(DMCreateSubDM(pdm, 1, &fieldS, &isS, &dmS));
285e2739ba6SAlexis Marboeuf PetscCall(DMCreateSubDM(pdm, 2, fieldUA, &isUA, &dmUA));
286c4762a1bSJed Brown
2879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2, &dmList));
288c4762a1bSJed Brown dmList[0] = dmU;
289c4762a1bSJed Brown dmList[1] = dmA;
290e2739ba6SAlexis Marboeuf
2919566063dSJacob Faibussowitsch PetscCall(DMCreateSuperDM(dmList, 2, NULL, &dmUA2));
2929566063dSJacob Faibussowitsch PetscCall(PetscFree(dmList));
293c4762a1bSJed Brown
294e2739ba6SAlexis Marboeuf PetscCall(DMGetGlobalVector(pdm, &X));
2959566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmU, &U));
2969566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmA, &A));
2979566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmS, &S));
2989566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmUA, &UA));
2999566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmUA2, &UA2));
300c4762a1bSJed Brown
3019566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)U, "U"));
3029566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)A, "Alpha"));
3039566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)S, "Sigma"));
3049566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)UA, "UAlpha"));
3059566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)UA2, "UAlpha2"));
3069566063dSJacob Faibussowitsch PetscCall(VecSet(X, -111.));
307c4762a1bSJed Brown
308c4762a1bSJed 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 */
309c4762a1bSJed Brown {
310c4762a1bSJed Brown PetscSection sectionUA;
311c4762a1bSJed Brown Vec UALoc;
312c4762a1bSJed Brown PetscSection coordSection;
313c4762a1bSJed Brown Vec coord;
314c4762a1bSJed Brown PetscScalar *cval, *xyz;
3156823f3c5SBlaise Bourdin PetscInt clSize, i, j;
316c4762a1bSJed Brown
3179566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dmUA, §ionUA));
3189566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dmUA, &UALoc));
3199566063dSJacob Faibussowitsch PetscCall(VecGetArray(UALoc, &cval));
3209566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dmUA, &coordSection));
3219566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmUA, &coord));
3229566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dmUA, &pStart, &pEnd));
3236823f3c5SBlaise Bourdin
324c4762a1bSJed Brown for (p = pStart; p < pEnd; ++p) {
325c4762a1bSJed Brown PetscInt dofUA, offUA;
326c4762a1bSJed Brown
3279566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(sectionUA, p, &dofUA));
328c4762a1bSJed Brown if (dofUA > 0) {
3296823f3c5SBlaise Bourdin xyz = NULL;
3309566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionUA, p, &offUA));
3319566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmUA, coordSection, coord, p, &clSize, &xyz));
332c4762a1bSJed Brown cval[offUA + sdim] = 0.;
333c4762a1bSJed Brown for (i = 0; i < sdim; ++i) {
334c4762a1bSJed Brown cval[offUA + i] = 0;
335ad540459SPierre Jolivet for (j = 0; j < clSize / sdim; ++j) cval[offUA + i] += xyz[j * sdim + i];
336c4762a1bSJed Brown cval[offUA + i] = cval[offUA + i] * sdim / clSize;
337c4762a1bSJed Brown cval[offUA + sdim] += PetscSqr(cval[offUA + i]);
338c4762a1bSJed Brown }
3399566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dmUA, coordSection, coord, p, &clSize, &xyz));
340c4762a1bSJed Brown }
341c4762a1bSJed Brown }
3429566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(UALoc, &cval));
3439566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dmUA, UALoc, INSERT_VALUES, UA));
3449566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dmUA, UALoc, INSERT_VALUES, UA));
3459566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dmUA, &UALoc));
3466823f3c5SBlaise Bourdin
347c4762a1bSJed Brown /* Update X */
3489566063dSJacob Faibussowitsch PetscCall(VecISCopy(X, isUA, SCATTER_FORWARD, UA));
3499566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(UA, NULL, "-ua_vec_view"));
3506823f3c5SBlaise Bourdin
351c4762a1bSJed Brown /* Restrict to U and Alpha */
3529566063dSJacob Faibussowitsch PetscCall(VecISCopy(X, isU, SCATTER_REVERSE, U));
3539566063dSJacob Faibussowitsch PetscCall(VecISCopy(X, isA, SCATTER_REVERSE, A));
3546823f3c5SBlaise Bourdin
355c4762a1bSJed Brown /* restrict to UA2 */
3569566063dSJacob Faibussowitsch PetscCall(VecISCopy(X, isUA, SCATTER_REVERSE, UA2));
3579566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(UA2, NULL, "-ua2_vec_view"));
358c4762a1bSJed Brown }
359c4762a1bSJed Brown
360c4762a1bSJed Brown {
361c4762a1bSJed Brown Vec tmpVec;
362c4762a1bSJed Brown PetscSection coordSection;
363c4762a1bSJed Brown Vec coord;
364c4762a1bSJed Brown PetscReal norm;
3656823f3c5SBlaise Bourdin PetscReal time = 1.234;
366c4762a1bSJed Brown
367c4762a1bSJed Brown /* Writing nodal variables to ExodusII file */
3689566063dSJacob Faibussowitsch PetscCall(DMSetOutputSequenceNumber(dmU, 0, time));
3699566063dSJacob Faibussowitsch PetscCall(DMSetOutputSequenceNumber(dmA, 0, time));
3709566063dSJacob Faibussowitsch PetscCall(VecView(U, viewer));
3719566063dSJacob Faibussowitsch PetscCall(VecView(A, viewer));
3726823f3c5SBlaise Bourdin
373c4762a1bSJed Brown /* Saving U and Alpha in one shot.
374c4762a1bSJed Brown For this, we need to cheat and change the Vec's name
3756823f3c5SBlaise Bourdin Note that in the end we write variables one component at a time,
3766823f3c5SBlaise Bourdin so that there is no real values in doing this
3776823f3c5SBlaise Bourdin */
3786823f3c5SBlaise Bourdin
3799566063dSJacob Faibussowitsch PetscCall(DMSetOutputSequenceNumber(dmUA, 1, time));
3809566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmUA, &tmpVec));
3819566063dSJacob Faibussowitsch PetscCall(VecCopy(UA, tmpVec));
3829566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)tmpVec, "U"));
3839566063dSJacob Faibussowitsch PetscCall(VecView(tmpVec, viewer));
384c4762a1bSJed Brown /* Reading nodal variables in Exodus file */
3859566063dSJacob Faibussowitsch PetscCall(VecSet(tmpVec, -1000.0));
3869566063dSJacob Faibussowitsch PetscCall(VecLoad(tmpVec, viewer));
3879566063dSJacob Faibussowitsch PetscCall(VecAXPY(UA, -1.0, tmpVec));
3889566063dSJacob Faibussowitsch PetscCall(VecNorm(UA, NORM_INFINITY, &norm));
38908401ef6SPierre Jolivet PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "UAlpha ||Vin - Vout|| = %g", (double)norm);
3909566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmUA, &tmpVec));
391c4762a1bSJed Brown
392c4762a1bSJed Brown /* same thing with the UA2 Vec obtained from the superDM */
3939566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmUA2, &tmpVec));
3949566063dSJacob Faibussowitsch PetscCall(VecCopy(UA2, tmpVec));
3959566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)tmpVec, "U"));
3969566063dSJacob Faibussowitsch PetscCall(DMSetOutputSequenceNumber(dmUA2, 2, time));
3979566063dSJacob Faibussowitsch PetscCall(VecView(tmpVec, viewer));
398c4762a1bSJed Brown /* Reading nodal variables in Exodus file */
3999566063dSJacob Faibussowitsch PetscCall(VecSet(tmpVec, -1000.0));
4009566063dSJacob Faibussowitsch PetscCall(VecLoad(tmpVec, viewer));
4019566063dSJacob Faibussowitsch PetscCall(VecAXPY(UA2, -1.0, tmpVec));
4029566063dSJacob Faibussowitsch PetscCall(VecNorm(UA2, NORM_INFINITY, &norm));
40308401ef6SPierre Jolivet PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "UAlpha2 ||Vin - Vout|| = %g", (double)norm);
4049566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmUA2, &tmpVec));
405c4762a1bSJed Brown
406c4762a1bSJed Brown /* Building and saving Sigma
407c4762a1bSJed Brown We set sigma_0 = rank (to see partitioning)
408c4762a1bSJed Brown sigma_1 = cell set ID
4096823f3c5SBlaise Bourdin sigma_2 = x_coordinate of the cell center of mass
4106823f3c5SBlaise Bourdin */
4119566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dmS, &coordSection));
4129566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmS, &coord));
4139566063dSJacob Faibussowitsch PetscCall(DMGetLabelIdIS(dmS, "Cell Sets", &csIS));
4149566063dSJacob Faibussowitsch PetscCall(DMGetLabelSize(dmS, "Cell Sets", &numCS));
4159566063dSJacob Faibussowitsch PetscCall(ISGetIndices(csIS, &csID));
416c4762a1bSJed Brown for (set = 0; set < numCS; ++set) {
417c4762a1bSJed 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 */
418c4762a1bSJed Brown IS cellIS;
419c4762a1bSJed Brown const PetscInt *cellID;
420c4762a1bSJed Brown PetscInt numCells, cell;
421c4762a1bSJed Brown PetscScalar *cval = NULL, *xyz = NULL;
422c4762a1bSJed Brown PetscInt clSize, cdimCoord, c;
423c4762a1bSJed Brown
4249566063dSJacob Faibussowitsch PetscCall(DMGetStratumIS(dmS, "Cell Sets", csID[set], &cellIS));
4259566063dSJacob Faibussowitsch PetscCall(ISGetIndices(cellIS, &cellID));
4269566063dSJacob Faibussowitsch PetscCall(ISGetSize(cellIS, &numCells));
427c4762a1bSJed Brown for (cell = 0; cell < numCells; cell++) {
4289566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmS, NULL, S, cellID[cell], &clSize, &cval));
4299566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmS, coordSection, coord, cellID[cell], &cdimCoord, &xyz));
430c4762a1bSJed Brown cval[0] = rank;
431c4762a1bSJed Brown cval[1] = csID[set];
432c4762a1bSJed Brown cval[2] = 0.;
433c4762a1bSJed Brown for (c = 0; c < cdimCoord / sdim; c++) cval[2] += xyz[c * sdim];
434c4762a1bSJed Brown cval[2] = cval[2] * sdim / cdimCoord;
4359566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetClosure(dmS, NULL, S, cellID[cell], cval, INSERT_ALL_VALUES));
436c4762a1bSJed Brown }
4379566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dmS, NULL, S, cellID[0], &clSize, &cval));
4389566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dmS, coordSection, coord, cellID[0], NULL, &xyz));
4399566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(cellIS, &cellID));
4409566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cellIS));
441c4762a1bSJed Brown }
4429566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(csIS, &csID));
4439566063dSJacob Faibussowitsch PetscCall(ISDestroy(&csIS));
4449566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(S, NULL, "-s_vec_view"));
4456823f3c5SBlaise Bourdin
446c4762a1bSJed Brown /* Writing zonal variables in Exodus file */
4479566063dSJacob Faibussowitsch PetscCall(DMSetOutputSequenceNumber(dmS, 0, time));
4489566063dSJacob Faibussowitsch PetscCall(VecView(S, viewer));
4496823f3c5SBlaise Bourdin
450c4762a1bSJed Brown /* Reading zonal variables in Exodus file */
4519566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmS, &tmpVec));
4529566063dSJacob Faibussowitsch PetscCall(VecSet(tmpVec, -1000.0));
4539566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)tmpVec, "Sigma"));
4549566063dSJacob Faibussowitsch PetscCall(VecLoad(tmpVec, viewer));
4559566063dSJacob Faibussowitsch PetscCall(VecAXPY(S, -1.0, tmpVec));
4569566063dSJacob Faibussowitsch PetscCall(VecNorm(S, NORM_INFINITY, &norm));
45708401ef6SPierre Jolivet PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Sigma ||Vin - Vout|| = %g", (double)norm);
4589566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmS, &tmpVec));
459c4762a1bSJed Brown }
4609566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer));
461c4762a1bSJed Brown
4629566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmUA2, &UA2));
4639566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmUA, &UA));
4649566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmS, &S));
4659566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmA, &A));
4669566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmU, &U));
467e2739ba6SAlexis Marboeuf PetscCall(DMRestoreGlobalVector(pdm, &X));
4689371c9d4SSatish Balay PetscCall(DMDestroy(&dmU));
4699371c9d4SSatish Balay PetscCall(ISDestroy(&isU));
4709371c9d4SSatish Balay PetscCall(DMDestroy(&dmA));
4719371c9d4SSatish Balay PetscCall(ISDestroy(&isA));
4729371c9d4SSatish Balay PetscCall(DMDestroy(&dmS));
4739371c9d4SSatish Balay PetscCall(ISDestroy(&isS));
4749371c9d4SSatish Balay PetscCall(DMDestroy(&dmUA));
4759371c9d4SSatish Balay 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
4912e8d78feSBlaise 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
4952e8d78feSBlaise 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
4982e8d78feSBlaise 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
5022e8d78feSBlaise 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
5052e8d78feSBlaise 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
5082e8d78feSBlaise 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
5142e8d78feSBlaise 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
5192e8d78feSBlaise 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
5232e8d78feSBlaise 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
5282e8d78feSBlaise 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
5322e8d78feSBlaise 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:
5349c89aa79SPierre Jolivet # Something is now broken with parallel read/write for whatever shape H is
535e2739ba6SAlexis Marboeuf TODO: broken
536e2739ba6SAlexis Marboeuf suffix: 11
537e2739ba6SAlexis Marboeuf nsize: 2
5382e8d78feSBlaise 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
5432e8d78feSBlaise 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
5462e8d78feSBlaise 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
5502e8d78feSBlaise 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
5532e8d78feSBlaise 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
5592e8d78feSBlaise 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
5632e8d78feSBlaise 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
5682e8d78feSBlaise 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
5722e8d78feSBlaise 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