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