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 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; 33c4762a1bSJed Brown 34327415f7SBarry Smith PetscFunctionBeginUser; 359566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 38d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FEM Layout Options", "ex26"); 399566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-i", "Filename to read", "ex26", ifilename, ifilename, sizeof(ifilename), NULL)); 409566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-o", "Filename to write", "ex26", ofilename, ofilename, sizeof(ofilename), NULL)); 419566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-order", "FEM polynomial order", "ex26", order, &order, NULL, 1)); 42d0609cedSBarry Smith PetscOptionsEnd(); 4363a3b9bcSJacob Faibussowitsch PetscCheck((order >= 1) && (order <= 2), PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported polynomial order %" PetscInt_FMT " not in [1, 2]", order); 44c4762a1bSJed Brown 456823f3c5SBlaise Bourdin /* Read the mesh from a file in any supported format */ 469566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename, NULL, PETSC_TRUE, &dm)); 479566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeSetDefault(dm, PETSC_FALSE)); 489566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 499566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 509566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &sdim)); 51c4762a1bSJed Brown 526823f3c5SBlaise Bourdin /* Create the exodus result file */ 536823f3c5SBlaise Bourdin { 546823f3c5SBlaise Bourdin PetscInt numstep = 3, step; 556823f3c5SBlaise Bourdin char *nodalVarName[4]; 566823f3c5SBlaise Bourdin char *zonalVarName[6]; 576823f3c5SBlaise Bourdin int *truthtable; 586823f3c5SBlaise Bourdin PetscInt numNodalVar, numZonalVar, i; 596823f3c5SBlaise Bourdin 60a5b23f4aSJose E. Roman /* enable exodus debugging information */ 616823f3c5SBlaise Bourdin ex_opts(EX_VERBOSE | EX_DEBUG); 626823f3c5SBlaise Bourdin /* Create the exodus file */ 639566063dSJacob Faibussowitsch PetscCall(PetscViewerExodusIIOpen(PETSC_COMM_WORLD, ofilename, FILE_MODE_WRITE, &viewer)); 646823f3c5SBlaise Bourdin /* The long way would be */ 656823f3c5SBlaise Bourdin /* 669566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PETSC_COMM_WORLD,&viewer)); 679566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(viewer,PETSCVIEWEREXODUSII)); 689566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(viewer,FILE_MODE_APPEND)); 699566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(viewer,ofilename)); 706823f3c5SBlaise Bourdin */ 716823f3c5SBlaise Bourdin /* set the mesh order */ 729566063dSJacob Faibussowitsch PetscCall(PetscViewerExodusIISetOrder(viewer, order)); 739566063dSJacob Faibussowitsch PetscCall(PetscViewerView(viewer, PETSC_VIEWER_STDOUT_WORLD)); 746823f3c5SBlaise Bourdin /* 756823f3c5SBlaise Bourdin Notice how the exodus file is actually NOT open at this point (exoid is -1) 76*d5b43468SJose E. Roman Since we are overwriting the file (mode is FILE_MODE_WRITE), we are going to have to 776823f3c5SBlaise Bourdin write the geometry (the DM), which can only be done on a brand new file. 786823f3c5SBlaise Bourdin */ 796823f3c5SBlaise Bourdin 806823f3c5SBlaise Bourdin /* Save the geometry to the file, erasing all previous content */ 819566063dSJacob Faibussowitsch PetscCall(DMView(dm, viewer)); 829566063dSJacob Faibussowitsch PetscCall(PetscViewerView(viewer, PETSC_VIEWER_STDOUT_WORLD)); 836823f3c5SBlaise Bourdin /* 846823f3c5SBlaise Bourdin Note how the exodus file is now open 856823f3c5SBlaise Bourdin */ 866823f3c5SBlaise Bourdin 876823f3c5SBlaise Bourdin /* "Format" the exodus result file, i.e. allocate space for nodal and zonal variables */ 886823f3c5SBlaise Bourdin switch (sdim) { 896823f3c5SBlaise Bourdin case 2: 906823f3c5SBlaise Bourdin numNodalVar = 3; 916823f3c5SBlaise Bourdin nodalVarName[0] = (char *)"U_x"; 926823f3c5SBlaise Bourdin nodalVarName[1] = (char *)"U_y"; 936823f3c5SBlaise Bourdin nodalVarName[2] = (char *)"Alpha"; 946823f3c5SBlaise Bourdin numZonalVar = 3; 956823f3c5SBlaise Bourdin zonalVarName[0] = (char *)"Sigma_11"; 966823f3c5SBlaise Bourdin zonalVarName[1] = (char *)"Sigma_22"; 976823f3c5SBlaise Bourdin zonalVarName[2] = (char *)"Sigma_12"; 986823f3c5SBlaise Bourdin break; 996823f3c5SBlaise Bourdin case 3: 1006823f3c5SBlaise Bourdin numNodalVar = 4; 1016823f3c5SBlaise Bourdin nodalVarName[0] = (char *)"U_x"; 1026823f3c5SBlaise Bourdin nodalVarName[1] = (char *)"U_y"; 1036823f3c5SBlaise Bourdin nodalVarName[2] = (char *)"U_z"; 1046823f3c5SBlaise Bourdin nodalVarName[3] = (char *)"Alpha"; 1056823f3c5SBlaise Bourdin numZonalVar = 6; 1066823f3c5SBlaise Bourdin zonalVarName[0] = (char *)"Sigma_11"; 1076823f3c5SBlaise Bourdin zonalVarName[1] = (char *)"Sigma_22"; 1086823f3c5SBlaise Bourdin zonalVarName[2] = (char *)"Sigma_33"; 1096823f3c5SBlaise Bourdin zonalVarName[3] = (char *)"Sigma_23"; 1106823f3c5SBlaise Bourdin zonalVarName[4] = (char *)"Sigma_13"; 1116823f3c5SBlaise Bourdin zonalVarName[5] = (char *)"Sigma_12"; 1126823f3c5SBlaise Bourdin break; 113d71ae5a4SJacob Faibussowitsch default: 114d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No layout for dimension %" PetscInt_FMT, sdim); 1156823f3c5SBlaise Bourdin } 1169566063dSJacob Faibussowitsch PetscCall(PetscViewerExodusIIGetId(viewer, &exoid)); 117792fecdfSBarry Smith PetscCallExternal(ex_put_variable_param, exoid, EX_ELEM_BLOCK, numZonalVar); 118792fecdfSBarry Smith PetscCallExternal(ex_put_variable_names, exoid, EX_ELEM_BLOCK, numZonalVar, zonalVarName); 119792fecdfSBarry Smith PetscCallExternal(ex_put_variable_param, exoid, EX_NODAL, numNodalVar); 120792fecdfSBarry Smith PetscCallExternal(ex_put_variable_names, exoid, EX_NODAL, numNodalVar, nodalVarName); 1216823f3c5SBlaise Bourdin numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK); 1226823f3c5SBlaise Bourdin 1236823f3c5SBlaise Bourdin /* 1246823f3c5SBlaise Bourdin An exodusII truth table specifies which fields are saved at which time step 1256823f3c5SBlaise Bourdin It speeds up I/O but reserving space for fieldsin the file ahead of time. 1266823f3c5SBlaise Bourdin */ 1279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numZonalVar * numCS, &truthtable)); 1286823f3c5SBlaise Bourdin for (i = 0; i < numZonalVar * numCS; ++i) truthtable[i] = 1; 129792fecdfSBarry Smith PetscCallExternal(ex_put_truth_table, exoid, EX_ELEM_BLOCK, numCS, numZonalVar, truthtable); 1309566063dSJacob Faibussowitsch PetscCall(PetscFree(truthtable)); 1316823f3c5SBlaise Bourdin 1326823f3c5SBlaise Bourdin /* Writing time step information in the file. Note that this is currently broken in the exodus library for netcdf4 (HDF5-based) files */ 1336823f3c5SBlaise Bourdin for (step = 0; step < numstep; ++step) { 1346823f3c5SBlaise Bourdin PetscReal time = step; 135792fecdfSBarry Smith PetscCallExternal(ex_put_time, exoid, step + 1, &time); 1366823f3c5SBlaise Bourdin } 1376823f3c5SBlaise Bourdin } 1386823f3c5SBlaise Bourdin 1396823f3c5SBlaise Bourdin /* Create the main section containing all fields */ 1409566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ion)); 1419566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(section, 3)); 1429566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(section, fieldU, "U")); 1439566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(section, fieldA, "Alpha")); 1449566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldName(section, fieldS, "Sigma")); 1459566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 1469566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(section, pStart, pEnd)); 1479566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(sdim + 1, &pStartDepth, sdim + 1, &pEndDepth)); 14848a46eb9SPierre Jolivet for (d = 0; d <= sdim; ++d) PetscCall(DMPlexGetDepthStratum(dm, d, &pStartDepth[d], &pEndDepth[d])); 149c4762a1bSJed Brown /* Vector field U, Scalar field Alpha, Tensor field Sigma */ 1509566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(section, fieldU, sdim)); 1519566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(section, fieldA, 1)); 1529566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(section, fieldS, sdim * (sdim + 1) / 2)); 153c4762a1bSJed Brown 154c4762a1bSJed Brown /* Going through cell sets then cells, and setting up storage for the sections */ 1559566063dSJacob Faibussowitsch PetscCall(DMGetLabelSize(dm, "Cell Sets", &numCS)); 1569566063dSJacob Faibussowitsch PetscCall(DMGetLabelIdIS(dm, "Cell Sets", &csIS)); 15748a46eb9SPierre Jolivet if (csIS) PetscCall(ISGetIndices(csIS, &csID)); 158c4762a1bSJed Brown for (set = 0; set < numCS; set++) { 159c4762a1bSJed Brown IS cellIS; 160c4762a1bSJed Brown const PetscInt *cellID; 161c4762a1bSJed Brown PetscInt numCells, cell, closureSize, *closureA = NULL; 162c4762a1bSJed Brown 1639566063dSJacob Faibussowitsch PetscCall(DMGetStratumSize(dm, "Cell Sets", csID[set], &numCells)); 1649566063dSJacob Faibussowitsch PetscCall(DMGetStratumIS(dm, "Cell Sets", csID[set], &cellIS)); 165c4762a1bSJed Brown if (numCells > 0) { 166c4762a1bSJed Brown /* dof layout ordered by increasing height in the DAG: cell, face, edge, vertex */ 167c4762a1bSJed Brown PetscInt dofUP1Tri[] = {2, 0, 0}; 168c4762a1bSJed Brown PetscInt dofAP1Tri[] = {1, 0, 0}; 169c4762a1bSJed Brown PetscInt dofUP2Tri[] = {2, 2, 0}; 170c4762a1bSJed Brown PetscInt dofAP2Tri[] = {1, 1, 0}; 171c4762a1bSJed Brown PetscInt dofUP1Quad[] = {2, 0, 0}; 172c4762a1bSJed Brown PetscInt dofAP1Quad[] = {1, 0, 0}; 173c4762a1bSJed Brown PetscInt dofUP2Quad[] = {2, 2, 2}; 174c4762a1bSJed Brown PetscInt dofAP2Quad[] = {1, 1, 1}; 175c4762a1bSJed Brown PetscInt dofS2D[] = {0, 0, 3}; 176c4762a1bSJed Brown PetscInt dofUP1Tet[] = {3, 0, 0, 0}; 177c4762a1bSJed Brown PetscInt dofAP1Tet[] = {1, 0, 0, 0}; 178c4762a1bSJed Brown PetscInt dofUP2Tet[] = {3, 3, 0, 0}; 179c4762a1bSJed Brown PetscInt dofAP2Tet[] = {1, 1, 0, 0}; 180c4762a1bSJed Brown PetscInt dofUP1Hex[] = {3, 0, 0, 0}; 181c4762a1bSJed Brown PetscInt dofAP1Hex[] = {1, 0, 0, 0}; 182c4762a1bSJed Brown PetscInt dofUP2Hex[] = {3, 3, 3, 3}; 183c4762a1bSJed Brown PetscInt dofAP2Hex[] = {1, 1, 1, 1}; 184c4762a1bSJed Brown PetscInt dofS3D[] = {0, 0, 0, 6}; 185c4762a1bSJed Brown PetscInt *dofU, *dofA, *dofS; 186c4762a1bSJed Brown 187c4762a1bSJed Brown switch (sdim) { 188d71ae5a4SJacob Faibussowitsch case 2: 189d71ae5a4SJacob Faibussowitsch dofS = dofS2D; 190d71ae5a4SJacob Faibussowitsch break; 191d71ae5a4SJacob Faibussowitsch case 3: 192d71ae5a4SJacob Faibussowitsch dofS = dofS3D; 193d71ae5a4SJacob Faibussowitsch break; 194d71ae5a4SJacob Faibussowitsch default: 195d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No layout for dimension %" PetscInt_FMT, sdim); 196c4762a1bSJed Brown } 197c4762a1bSJed Brown 198c4762a1bSJed Brown /* Identify cell type based on closure size only. This works for Tri/Tet/Quad/Hex meshes 199c4762a1bSJed Brown It will not be enough to identify more exotic elements like pyramid or prisms... */ 2009566063dSJacob Faibussowitsch PetscCall(ISGetIndices(cellIS, &cellID)); 2019566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, cellID[0], PETSC_TRUE, &closureSize, &closureA)); 202c4762a1bSJed Brown switch (closureSize) { 203c4762a1bSJed Brown case 7: /* Tri */ 204c4762a1bSJed Brown if (order == 1) { 205c4762a1bSJed Brown dofU = dofUP1Tri; 206c4762a1bSJed Brown dofA = dofAP1Tri; 207c4762a1bSJed Brown } else { 208c4762a1bSJed Brown dofU = dofUP2Tri; 209c4762a1bSJed Brown dofA = dofAP2Tri; 210c4762a1bSJed Brown } 211c4762a1bSJed Brown break; 212c4762a1bSJed Brown case 9: /* Quad */ 213c4762a1bSJed Brown if (order == 1) { 214c4762a1bSJed Brown dofU = dofUP1Quad; 215c4762a1bSJed Brown dofA = dofAP1Quad; 216c4762a1bSJed Brown } else { 217c4762a1bSJed Brown dofU = dofUP2Quad; 218c4762a1bSJed Brown dofA = dofAP2Quad; 219c4762a1bSJed Brown } 220c4762a1bSJed Brown break; 221c4762a1bSJed Brown case 15: /* Tet */ 222c4762a1bSJed Brown if (order == 1) { 223c4762a1bSJed Brown dofU = dofUP1Tet; 224c4762a1bSJed Brown dofA = dofAP1Tet; 225c4762a1bSJed Brown } else { 226c4762a1bSJed Brown dofU = dofUP2Tet; 227c4762a1bSJed Brown dofA = dofAP2Tet; 228c4762a1bSJed Brown } 229c4762a1bSJed Brown break; 230c4762a1bSJed Brown case 27: /* Hex */ 231c4762a1bSJed Brown if (order == 1) { 232c4762a1bSJed Brown dofU = dofUP1Hex; 233c4762a1bSJed Brown dofA = dofAP1Hex; 234c4762a1bSJed Brown } else { 235c4762a1bSJed Brown dofU = dofUP2Hex; 236c4762a1bSJed Brown dofA = dofAP2Hex; 237c4762a1bSJed Brown } 238c4762a1bSJed Brown break; 239d71ae5a4SJacob Faibussowitsch default: 240d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Unknown element with closure size %" PetscInt_FMT, closureSize); 241c4762a1bSJed Brown } 2429566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, cellID[0], PETSC_TRUE, &closureSize, &closureA)); 243c4762a1bSJed Brown 244c4762a1bSJed Brown for (cell = 0; cell < numCells; cell++) { 245c4762a1bSJed Brown PetscInt *closure = NULL; 246c4762a1bSJed Brown 2479566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, cellID[cell], PETSC_TRUE, &closureSize, &closure)); 248c4762a1bSJed Brown for (p = 0; p < closureSize; ++p) { 249c4762a1bSJed Brown /* Find depth of p */ 250c4762a1bSJed Brown for (d = 0; d <= sdim; ++d) { 251c4762a1bSJed Brown if ((closure[2 * p] >= pStartDepth[d]) && (closure[2 * p] < pEndDepth[d])) { 2529566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(section, closure[2 * p], dofU[d] + dofA[d] + dofS[d])); 2539566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(section, closure[2 * p], fieldU, dofU[d])); 2549566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(section, closure[2 * p], fieldA, dofA[d])); 2559566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(section, closure[2 * p], fieldS, dofS[d])); 256c4762a1bSJed Brown } 257c4762a1bSJed Brown } 258c4762a1bSJed Brown } 2599566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, cellID[cell], PETSC_TRUE, &closureSize, &closure)); 260c4762a1bSJed Brown } 2619566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(cellIS, &cellID)); 2629566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cellIS)); 263c4762a1bSJed Brown } 264c4762a1bSJed Brown } 26548a46eb9SPierre Jolivet if (csIS) PetscCall(ISRestoreIndices(csIS, &csID)); 2669566063dSJacob Faibussowitsch PetscCall(ISDestroy(&csIS)); 2679566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(section)); 2689566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(dm, section)); 2699566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)section, NULL, "-dm_section_view")); 2709566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(§ion)); 271c4762a1bSJed Brown 272c4762a1bSJed Brown { 273c4762a1bSJed Brown PetscSF migrationSF; 274c4762a1bSJed Brown PetscInt ovlp = 0; 275c4762a1bSJed Brown PetscPartitioner part; 276c4762a1bSJed Brown 2779566063dSJacob Faibussowitsch PetscCall(DMSetUseNatural(dm, PETSC_TRUE)); 2789566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(dm, &part)); 2799566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetFromOptions(part)); 2809566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(dm, ovlp, &migrationSF, &pdm)); 281e2739ba6SAlexis Marboeuf if (!pdm) pdm = dm; 282e2739ba6SAlexis Marboeuf /* Set the migrationSF is mandatory since useNatural = PETSC_TRUE */ 283e2739ba6SAlexis Marboeuf if (migrationSF) { 2849566063dSJacob Faibussowitsch PetscCall(DMPlexSetMigrationSF(pdm, migrationSF)); 2859566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&migrationSF)); 286c4762a1bSJed Brown } 287e2739ba6SAlexis Marboeuf PetscCall(DMViewFromOptions(pdm, NULL, "-dm_view")); 288c4762a1bSJed Brown } 289c4762a1bSJed Brown 290c4762a1bSJed Brown /* Get DM and IS for each field of dm */ 291e2739ba6SAlexis Marboeuf PetscCall(DMCreateSubDM(pdm, 1, &fieldU, &isU, &dmU)); 292e2739ba6SAlexis Marboeuf PetscCall(DMCreateSubDM(pdm, 1, &fieldA, &isA, &dmA)); 293e2739ba6SAlexis Marboeuf PetscCall(DMCreateSubDM(pdm, 1, &fieldS, &isS, &dmS)); 294e2739ba6SAlexis Marboeuf PetscCall(DMCreateSubDM(pdm, 2, fieldUA, &isUA, &dmUA)); 295c4762a1bSJed Brown 2969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2, &dmList)); 297c4762a1bSJed Brown dmList[0] = dmU; 298c4762a1bSJed Brown dmList[1] = dmA; 299e2739ba6SAlexis Marboeuf 3009566063dSJacob Faibussowitsch PetscCall(DMCreateSuperDM(dmList, 2, NULL, &dmUA2)); 3019566063dSJacob Faibussowitsch PetscCall(PetscFree(dmList)); 302c4762a1bSJed Brown 303e2739ba6SAlexis Marboeuf PetscCall(DMGetGlobalVector(pdm, &X)); 3049566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmU, &U)); 3059566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmA, &A)); 3069566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmS, &S)); 3079566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmUA, &UA)); 3089566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmUA2, &UA2)); 309c4762a1bSJed Brown 3109566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)U, "U")); 3119566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)A, "Alpha")); 3129566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)S, "Sigma")); 3139566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)UA, "UAlpha")); 3149566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)UA2, "UAlpha2")); 3159566063dSJacob Faibussowitsch PetscCall(VecSet(X, -111.)); 316c4762a1bSJed Brown 317c4762a1bSJed 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 */ 318c4762a1bSJed Brown { 319c4762a1bSJed Brown PetscSection sectionUA; 320c4762a1bSJed Brown Vec UALoc; 321c4762a1bSJed Brown PetscSection coordSection; 322c4762a1bSJed Brown Vec coord; 323c4762a1bSJed Brown PetscScalar *cval, *xyz; 3246823f3c5SBlaise Bourdin PetscInt clSize, i, j; 325c4762a1bSJed Brown 3269566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dmUA, §ionUA)); 3279566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dmUA, &UALoc)); 3289566063dSJacob Faibussowitsch PetscCall(VecGetArray(UALoc, &cval)); 3299566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dmUA, &coordSection)); 3309566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmUA, &coord)); 3319566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dmUA, &pStart, &pEnd)); 3326823f3c5SBlaise Bourdin 333c4762a1bSJed Brown for (p = pStart; p < pEnd; ++p) { 334c4762a1bSJed Brown PetscInt dofUA, offUA; 335c4762a1bSJed Brown 3369566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(sectionUA, p, &dofUA)); 337c4762a1bSJed Brown if (dofUA > 0) { 3386823f3c5SBlaise Bourdin xyz = NULL; 3399566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionUA, p, &offUA)); 3409566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmUA, coordSection, coord, p, &clSize, &xyz)); 341c4762a1bSJed Brown cval[offUA + sdim] = 0.; 342c4762a1bSJed Brown for (i = 0; i < sdim; ++i) { 343c4762a1bSJed Brown cval[offUA + i] = 0; 344ad540459SPierre Jolivet for (j = 0; j < clSize / sdim; ++j) cval[offUA + i] += xyz[j * sdim + i]; 345c4762a1bSJed Brown cval[offUA + i] = cval[offUA + i] * sdim / clSize; 346c4762a1bSJed Brown cval[offUA + sdim] += PetscSqr(cval[offUA + i]); 347c4762a1bSJed Brown } 3489566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dmUA, coordSection, coord, p, &clSize, &xyz)); 349c4762a1bSJed Brown } 350c4762a1bSJed Brown } 3519566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(UALoc, &cval)); 3529566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dmUA, UALoc, INSERT_VALUES, UA)); 3539566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dmUA, UALoc, INSERT_VALUES, UA)); 3549566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dmUA, &UALoc)); 3556823f3c5SBlaise Bourdin 356c4762a1bSJed Brown /* Update X */ 3579566063dSJacob Faibussowitsch PetscCall(VecISCopy(X, isUA, SCATTER_FORWARD, UA)); 3589566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(UA, NULL, "-ua_vec_view")); 3596823f3c5SBlaise Bourdin 360c4762a1bSJed Brown /* Restrict to U and Alpha */ 3619566063dSJacob Faibussowitsch PetscCall(VecISCopy(X, isU, SCATTER_REVERSE, U)); 3629566063dSJacob Faibussowitsch PetscCall(VecISCopy(X, isA, SCATTER_REVERSE, A)); 3636823f3c5SBlaise Bourdin 364c4762a1bSJed Brown /* restrict to UA2 */ 3659566063dSJacob Faibussowitsch PetscCall(VecISCopy(X, isUA, SCATTER_REVERSE, UA2)); 3669566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(UA2, NULL, "-ua2_vec_view")); 367c4762a1bSJed Brown } 368c4762a1bSJed Brown 369c4762a1bSJed Brown { 370c4762a1bSJed Brown Vec tmpVec; 371c4762a1bSJed Brown PetscSection coordSection; 372c4762a1bSJed Brown Vec coord; 373c4762a1bSJed Brown PetscReal norm; 3746823f3c5SBlaise Bourdin PetscReal time = 1.234; 375c4762a1bSJed Brown 376c4762a1bSJed Brown /* Writing nodal variables to ExodusII file */ 3779566063dSJacob Faibussowitsch PetscCall(DMSetOutputSequenceNumber(dmU, 0, time)); 3789566063dSJacob Faibussowitsch PetscCall(DMSetOutputSequenceNumber(dmA, 0, time)); 3796823f3c5SBlaise Bourdin 3809566063dSJacob Faibussowitsch PetscCall(VecView(U, viewer)); 3819566063dSJacob Faibussowitsch PetscCall(VecView(A, viewer)); 3826823f3c5SBlaise Bourdin 383c4762a1bSJed Brown /* Saving U and Alpha in one shot. 384c4762a1bSJed Brown For this, we need to cheat and change the Vec's name 3856823f3c5SBlaise Bourdin Note that in the end we write variables one component at a time, 3866823f3c5SBlaise Bourdin so that there is no real values in doing this 3876823f3c5SBlaise Bourdin */ 3886823f3c5SBlaise Bourdin 3899566063dSJacob Faibussowitsch PetscCall(DMSetOutputSequenceNumber(dmUA, 1, time)); 3909566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmUA, &tmpVec)); 3919566063dSJacob Faibussowitsch PetscCall(VecCopy(UA, tmpVec)); 3929566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)tmpVec, "U")); 3939566063dSJacob Faibussowitsch PetscCall(VecView(tmpVec, viewer)); 394c4762a1bSJed Brown /* Reading nodal variables in Exodus file */ 3959566063dSJacob Faibussowitsch PetscCall(VecSet(tmpVec, -1000.0)); 3969566063dSJacob Faibussowitsch PetscCall(VecLoad(tmpVec, viewer)); 3979566063dSJacob Faibussowitsch PetscCall(VecAXPY(UA, -1.0, tmpVec)); 3989566063dSJacob Faibussowitsch PetscCall(VecNorm(UA, NORM_INFINITY, &norm)); 39908401ef6SPierre Jolivet PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "UAlpha ||Vin - Vout|| = %g", (double)norm); 4009566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmUA, &tmpVec)); 401c4762a1bSJed Brown 402c4762a1bSJed Brown /* same thing with the UA2 Vec obtained from the superDM */ 4039566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmUA2, &tmpVec)); 4049566063dSJacob Faibussowitsch PetscCall(VecCopy(UA2, tmpVec)); 4059566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)tmpVec, "U")); 4069566063dSJacob Faibussowitsch PetscCall(DMSetOutputSequenceNumber(dmUA2, 2, time)); 4079566063dSJacob Faibussowitsch PetscCall(VecView(tmpVec, viewer)); 408c4762a1bSJed Brown /* Reading nodal variables in Exodus file */ 4099566063dSJacob Faibussowitsch PetscCall(VecSet(tmpVec, -1000.0)); 4109566063dSJacob Faibussowitsch PetscCall(VecLoad(tmpVec, viewer)); 4119566063dSJacob Faibussowitsch PetscCall(VecAXPY(UA2, -1.0, tmpVec)); 4129566063dSJacob Faibussowitsch PetscCall(VecNorm(UA2, NORM_INFINITY, &norm)); 41308401ef6SPierre Jolivet PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "UAlpha2 ||Vin - Vout|| = %g", (double)norm); 4149566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmUA2, &tmpVec)); 415c4762a1bSJed Brown 416c4762a1bSJed Brown /* Building and saving Sigma 417c4762a1bSJed Brown We set sigma_0 = rank (to see partitioning) 418c4762a1bSJed Brown sigma_1 = cell set ID 4196823f3c5SBlaise Bourdin sigma_2 = x_coordinate of the cell center of mass 4206823f3c5SBlaise Bourdin */ 4219566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(dmS, &coordSection)); 4229566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmS, &coord)); 4239566063dSJacob Faibussowitsch PetscCall(DMGetLabelIdIS(dmS, "Cell Sets", &csIS)); 4249566063dSJacob Faibussowitsch PetscCall(DMGetLabelSize(dmS, "Cell Sets", &numCS)); 4259566063dSJacob Faibussowitsch PetscCall(ISGetIndices(csIS, &csID)); 426c4762a1bSJed Brown for (set = 0; set < numCS; ++set) { 427c4762a1bSJed 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 */ 428c4762a1bSJed Brown IS cellIS; 429c4762a1bSJed Brown const PetscInt *cellID; 430c4762a1bSJed Brown PetscInt numCells, cell; 431c4762a1bSJed Brown PetscScalar *cval = NULL, *xyz = NULL; 432c4762a1bSJed Brown PetscInt clSize, cdimCoord, c; 433c4762a1bSJed Brown 4349566063dSJacob Faibussowitsch PetscCall(DMGetStratumIS(dmS, "Cell Sets", csID[set], &cellIS)); 4359566063dSJacob Faibussowitsch PetscCall(ISGetIndices(cellIS, &cellID)); 4369566063dSJacob Faibussowitsch PetscCall(ISGetSize(cellIS, &numCells)); 437c4762a1bSJed Brown for (cell = 0; cell < numCells; cell++) { 4389566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmS, NULL, S, cellID[cell], &clSize, &cval)); 4399566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(dmS, coordSection, coord, cellID[cell], &cdimCoord, &xyz)); 440c4762a1bSJed Brown cval[0] = rank; 441c4762a1bSJed Brown cval[1] = csID[set]; 442c4762a1bSJed Brown cval[2] = 0.; 443c4762a1bSJed Brown for (c = 0; c < cdimCoord / sdim; c++) cval[2] += xyz[c * sdim]; 444c4762a1bSJed Brown cval[2] = cval[2] * sdim / cdimCoord; 4459566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetClosure(dmS, NULL, S, cellID[cell], cval, INSERT_ALL_VALUES)); 446c4762a1bSJed Brown } 4479566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dmS, NULL, S, cellID[0], &clSize, &cval)); 4489566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(dmS, coordSection, coord, cellID[0], NULL, &xyz)); 4499566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(cellIS, &cellID)); 4509566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cellIS)); 451c4762a1bSJed Brown } 4529566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(csIS, &csID)); 4539566063dSJacob Faibussowitsch PetscCall(ISDestroy(&csIS)); 4549566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(S, NULL, "-s_vec_view")); 4556823f3c5SBlaise Bourdin 456c4762a1bSJed Brown /* Writing zonal variables in Exodus file */ 4579566063dSJacob Faibussowitsch PetscCall(DMSetOutputSequenceNumber(dmS, 0, time)); 4589566063dSJacob Faibussowitsch PetscCall(VecView(S, viewer)); 4596823f3c5SBlaise Bourdin 460c4762a1bSJed Brown /* Reading zonal variables in Exodus file */ 4619566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmS, &tmpVec)); 4629566063dSJacob Faibussowitsch PetscCall(VecSet(tmpVec, -1000.0)); 4639566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)tmpVec, "Sigma")); 4649566063dSJacob Faibussowitsch PetscCall(VecLoad(tmpVec, viewer)); 4659566063dSJacob Faibussowitsch PetscCall(VecAXPY(S, -1.0, tmpVec)); 4669566063dSJacob Faibussowitsch PetscCall(VecNorm(S, NORM_INFINITY, &norm)); 46708401ef6SPierre Jolivet PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Sigma ||Vin - Vout|| = %g", (double)norm); 4689566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmS, &tmpVec)); 469c4762a1bSJed Brown } 4709566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 471c4762a1bSJed Brown 4729566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmUA2, &UA2)); 4739566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmUA, &UA)); 4749566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmS, &S)); 4759566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmA, &A)); 4769566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmU, &U)); 477e2739ba6SAlexis Marboeuf PetscCall(DMRestoreGlobalVector(pdm, &X)); 4789371c9d4SSatish Balay PetscCall(DMDestroy(&dmU)); 4799371c9d4SSatish Balay PetscCall(ISDestroy(&isU)); 4809371c9d4SSatish Balay PetscCall(DMDestroy(&dmA)); 4819371c9d4SSatish Balay PetscCall(ISDestroy(&isA)); 4829371c9d4SSatish Balay PetscCall(DMDestroy(&dmS)); 4839371c9d4SSatish Balay PetscCall(ISDestroy(&isS)); 4849371c9d4SSatish Balay PetscCall(DMDestroy(&dmUA)); 4859371c9d4SSatish Balay PetscCall(ISDestroy(&isUA)); 4869566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmUA2)); 487e2739ba6SAlexis Marboeuf PetscCall(DMDestroy(&pdm)); 488e2739ba6SAlexis Marboeuf if (size > 1) PetscCall(DMDestroy(&dm)); 4899566063dSJacob Faibussowitsch PetscCall(PetscFree2(pStartDepth, pEndDepth)); 4909566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 491b122ec5aSJacob Faibussowitsch return 0; 492c4762a1bSJed Brown } 493c4762a1bSJed Brown 494c4762a1bSJed Brown /*TEST 495c4762a1bSJed Brown 496c4762a1bSJed Brown build: 497c4762a1bSJed Brown requires: exodusii pnetcdf !complex 498e2739ba6SAlexis Marboeuf # 2D seq 499e2739ba6SAlexis Marboeuf test: 500e2739ba6SAlexis Marboeuf suffix: 0 5012e8d78feSBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -petscpartitioner_type simple -order 1 502e2739ba6SAlexis 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 503e2739ba6SAlexis Marboeuf test: 504e2739ba6SAlexis Marboeuf suffix: 1 5052e8d78feSBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -petscpartitioner_type simple -order 1 506e2739ba6SAlexis Marboeuf test: 507e2739ba6SAlexis Marboeuf suffix: 2 5082e8d78feSBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -petscpartitioner_type simple -order 1 509e2739ba6SAlexis 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 510e2739ba6SAlexis Marboeuf test: 511e2739ba6SAlexis Marboeuf suffix: 3 5122e8d78feSBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -petscpartitioner_type simple -order 2 513e2739ba6SAlexis Marboeuf test: 514e2739ba6SAlexis Marboeuf suffix: 4 5152e8d78feSBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -petscpartitioner_type simple -order 2 516e2739ba6SAlexis Marboeuf test: 517e2739ba6SAlexis Marboeuf suffix: 5 5182e8d78feSBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -petscpartitioner_type simple -order 2 519c4762a1bSJed Brown 520e2739ba6SAlexis Marboeuf # 2D par 521e2739ba6SAlexis Marboeuf test: 522e2739ba6SAlexis Marboeuf suffix: 6 523e2739ba6SAlexis Marboeuf nsize: 2 5242e8d78feSBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -petscpartitioner_type simple -order 1 525e2739ba6SAlexis 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 526e2739ba6SAlexis Marboeuf test: 527e2739ba6SAlexis Marboeuf suffix: 7 528e2739ba6SAlexis Marboeuf nsize: 2 5292e8d78feSBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -petscpartitioner_type simple -order 1 530e2739ba6SAlexis Marboeuf test: 531e2739ba6SAlexis Marboeuf suffix: 8 532e2739ba6SAlexis Marboeuf nsize: 2 5332e8d78feSBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -petscpartitioner_type simple -order 1 534e2739ba6SAlexis Marboeuf #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: invalid dimension ID or name 535e2739ba6SAlexis Marboeuf test: 536e2739ba6SAlexis Marboeuf suffix: 9 537e2739ba6SAlexis Marboeuf nsize: 2 5382e8d78feSBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo -dm_view -petscpartitioner_type simple -order 2 539e2739ba6SAlexis Marboeuf test: 540e2739ba6SAlexis Marboeuf suffix: 10 541e2739ba6SAlexis Marboeuf nsize: 2 5422e8d78feSBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -o FourSquareQ-large_out.exo -dm_view -petscpartitioner_type simple -order 2 543e2739ba6SAlexis Marboeuf test: 544e2739ba6SAlexis Marboeuf # Something is now broken with parallel read/write for wahtever shape H is 545e2739ba6SAlexis Marboeuf TODO: broken 546e2739ba6SAlexis Marboeuf suffix: 11 547e2739ba6SAlexis Marboeuf nsize: 2 5482e8d78feSBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -o FourSquareH-large_out.exo -dm_view -petscpartitioner_type simple -order 2 549c4762a1bSJed Brown 550e2739ba6SAlexis Marboeuf #3d seq 551e2739ba6SAlexis Marboeuf test: 552e2739ba6SAlexis Marboeuf suffix: 12 5532e8d78feSBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -petscpartitioner_type simple -order 1 554e2739ba6SAlexis Marboeuf test: 555e2739ba6SAlexis Marboeuf suffix: 13 5562e8d78feSBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -petscpartitioner_type simple -order 1 557e2739ba6SAlexis 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 558e2739ba6SAlexis Marboeuf test: 559e2739ba6SAlexis Marboeuf suffix: 14 5602e8d78feSBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -petscpartitioner_type simple -order 2 561e2739ba6SAlexis Marboeuf test: 562e2739ba6SAlexis Marboeuf suffix: 15 5632e8d78feSBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -petscpartitioner_type simple -order 2 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 #3d par 566e2739ba6SAlexis Marboeuf test: 567e2739ba6SAlexis Marboeuf suffix: 16 568e2739ba6SAlexis Marboeuf nsize: 2 5692e8d78feSBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -petscpartitioner_type simple -order 1 570e2739ba6SAlexis Marboeuf test: 571e2739ba6SAlexis Marboeuf suffix: 17 572e2739ba6SAlexis Marboeuf nsize: 2 5732e8d78feSBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -petscpartitioner_type simple -order 1 574e2739ba6SAlexis 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 575e2739ba6SAlexis Marboeuf test: 576e2739ba6SAlexis Marboeuf suffix: 18 577e2739ba6SAlexis Marboeuf nsize: 2 5782e8d78feSBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -o FourBrickHex-large_out.exo -dm_view -petscpartitioner_type simple -order 2 579e2739ba6SAlexis Marboeuf test: 580e2739ba6SAlexis Marboeuf suffix: 19 581e2739ba6SAlexis Marboeuf nsize: 2 5822e8d78feSBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -o FourBrickTet-large_out.exo -dm_view -petscpartitioner_type simple -order 2 583e2739ba6SAlexis 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 584c4762a1bSJed Brown 585c4762a1bSJed Brown TEST*/ 586