1c4762a1bSJed Brown static char help[] = "Test FEM layout with DM and ExodusII storage\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown In order to see the vectors which are being tested, use 5c4762a1bSJed Brown 6c4762a1bSJed Brown -ua_vec_view -s_vec_view 7c4762a1bSJed Brown */ 8c4762a1bSJed Brown 9c4762a1bSJed Brown #include <petsc.h> 10c4762a1bSJed Brown #include <exodusII.h> 11c4762a1bSJed Brown 12c4762a1bSJed Brown #include <petsc/private/dmpleximpl.h> 13c4762a1bSJed Brown 14c4762a1bSJed Brown int main(int argc, char **argv) { 15c4762a1bSJed Brown DM dm, dmU, dmA, dmS, dmUA, dmUA2, *dmList; 16c4762a1bSJed Brown Vec X, U, A, S, UA, UA2; 17c4762a1bSJed Brown IS isU, isA, isS, isUA; 18c4762a1bSJed Brown PetscSection section; 19c4762a1bSJed Brown const PetscInt fieldU = 0; 20c4762a1bSJed Brown const PetscInt fieldA = 2; 21c4762a1bSJed Brown const PetscInt fieldS = 1; 22c4762a1bSJed Brown const PetscInt fieldUA[2] = {0, 2}; 23c4762a1bSJed Brown char ifilename[PETSC_MAX_PATH_LEN], ofilename[PETSC_MAX_PATH_LEN]; 24c4762a1bSJed Brown int exoid = -1; 25c4762a1bSJed Brown IS csIS; 26c4762a1bSJed Brown const PetscInt *csID; 27c4762a1bSJed Brown PetscInt *pStartDepth, *pEndDepth; 28c4762a1bSJed Brown PetscInt order = 1; 29c4762a1bSJed Brown PetscInt sdim, d, pStart, pEnd, p, numCS, set; 30c4762a1bSJed Brown PetscMPIInt rank, size; 316823f3c5SBlaise Bourdin PetscViewer viewer; 32c4762a1bSJed Brown PetscErrorCode ierr; 33c4762a1bSJed Brown 34*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv,NULL, help)); 355f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 365f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 37c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FEM Layout Options", "ex26");CHKERRQ(ierr); 385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsString("-i", "Filename to read", "ex26", ifilename, ifilename, sizeof(ifilename), NULL)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsString("-o", "Filename to write", "ex26", ofilename, ofilename, sizeof(ofilename), NULL)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBoundedInt("-order", "FEM polynomial order", "ex26", order, &order, NULL,1)); 41c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 422c71b3e2SJacob Faibussowitsch PetscCheckFalse((order > 2) || (order < 1),PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported polynomial order %D not in [1, 2]", order); 43c4762a1bSJed Brown 446823f3c5SBlaise Bourdin /* Read the mesh from a file in any supported format */ 455f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename, NULL, PETSC_TRUE, &dm)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistributeSetDefault(dm, PETSC_FALSE)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(dm)); 485f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(dm, NULL, "-dm_view")); 495f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(dm, &sdim)); 50c4762a1bSJed Brown 516823f3c5SBlaise Bourdin /* Create the exodus result file */ 526823f3c5SBlaise Bourdin { 536823f3c5SBlaise Bourdin PetscInt numstep = 3, step; 546823f3c5SBlaise Bourdin char *nodalVarName[4]; 556823f3c5SBlaise Bourdin char *zonalVarName[6]; 566823f3c5SBlaise Bourdin int *truthtable; 576823f3c5SBlaise Bourdin PetscInt numNodalVar, numZonalVar, i; 586823f3c5SBlaise Bourdin 59a5b23f4aSJose E. Roman /* enable exodus debugging information */ 606823f3c5SBlaise Bourdin ex_opts(EX_VERBOSE|EX_DEBUG); 616823f3c5SBlaise Bourdin /* Create the exodus file */ 625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerExodusIIOpen(PETSC_COMM_WORLD,ofilename,FILE_MODE_WRITE,&viewer)); 636823f3c5SBlaise Bourdin /* The long way would be */ 646823f3c5SBlaise Bourdin /* 655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerCreate(PETSC_COMM_WORLD,&viewer)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerSetType(viewer,PETSCVIEWEREXODUSII)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFileSetMode(viewer,FILE_MODE_APPEND)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFileSetName(viewer,ofilename)); 696823f3c5SBlaise Bourdin */ 706823f3c5SBlaise Bourdin /* set the mesh order */ 715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerExodusIISetOrder(viewer,order)); 725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD)); 736823f3c5SBlaise Bourdin /* 746823f3c5SBlaise Bourdin Notice how the exodus file is actually NOT open at this point (exoid is -1) 756823f3c5SBlaise Bourdin Since we are overwritting the file (mode is FILE_MODE_WRITE), we are going to have to 766823f3c5SBlaise Bourdin write the geometry (the DM), which can only be done on a brand new file. 776823f3c5SBlaise Bourdin */ 786823f3c5SBlaise Bourdin 796823f3c5SBlaise Bourdin /* Save the geometry to the file, erasing all previous content */ 805f80ce2aSJacob Faibussowitsch CHKERRQ(DMView(dm,viewer)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerView(viewer,PETSC_VIEWER_STDOUT_WORLD)); 826823f3c5SBlaise Bourdin /* 836823f3c5SBlaise Bourdin Note how the exodus file is now open 846823f3c5SBlaise Bourdin */ 856823f3c5SBlaise Bourdin 866823f3c5SBlaise Bourdin /* "Format" the exodus result file, i.e. allocate space for nodal and zonal variables */ 876823f3c5SBlaise Bourdin switch (sdim) { 886823f3c5SBlaise Bourdin case 2: 896823f3c5SBlaise Bourdin numNodalVar = 3; 906823f3c5SBlaise Bourdin nodalVarName[0] = (char *) "U_x"; 916823f3c5SBlaise Bourdin nodalVarName[1] = (char *) "U_y"; 926823f3c5SBlaise Bourdin nodalVarName[2] = (char *) "Alpha"; 936823f3c5SBlaise Bourdin numZonalVar = 3; 946823f3c5SBlaise Bourdin zonalVarName[0] = (char *) "Sigma_11"; 956823f3c5SBlaise Bourdin zonalVarName[1] = (char *) "Sigma_22"; 966823f3c5SBlaise Bourdin zonalVarName[2] = (char *) "Sigma_12"; 976823f3c5SBlaise Bourdin break; 986823f3c5SBlaise Bourdin case 3: 996823f3c5SBlaise Bourdin numNodalVar = 4; 1006823f3c5SBlaise Bourdin nodalVarName[0] = (char *) "U_x"; 1016823f3c5SBlaise Bourdin nodalVarName[1] = (char *) "U_y"; 1026823f3c5SBlaise Bourdin nodalVarName[2] = (char *) "U_z"; 1036823f3c5SBlaise Bourdin nodalVarName[3] = (char *) "Alpha"; 1046823f3c5SBlaise Bourdin numZonalVar = 6; 1056823f3c5SBlaise Bourdin zonalVarName[0] = (char *) "Sigma_11"; 1066823f3c5SBlaise Bourdin zonalVarName[1] = (char *) "Sigma_22"; 1076823f3c5SBlaise Bourdin zonalVarName[2] = (char *) "Sigma_33"; 1086823f3c5SBlaise Bourdin zonalVarName[3] = (char *) "Sigma_23"; 1096823f3c5SBlaise Bourdin zonalVarName[4] = (char *) "Sigma_13"; 1106823f3c5SBlaise Bourdin zonalVarName[5] = (char *) "Sigma_12"; 1116823f3c5SBlaise Bourdin break; 11298921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No layout for dimension %D", sdim); 1136823f3c5SBlaise Bourdin } 1145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerExodusIIGetId(viewer,&exoid)); 115a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_put_variable_param,exoid, EX_ELEM_BLOCK, numZonalVar); 116a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_put_variable_names,exoid, EX_ELEM_BLOCK, numZonalVar, zonalVarName); 117a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_put_variable_param,exoid, EX_NODAL, numNodalVar); 118a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_put_variable_names,exoid, EX_NODAL, numNodalVar, nodalVarName); 1196823f3c5SBlaise Bourdin numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK); 1206823f3c5SBlaise Bourdin 1216823f3c5SBlaise Bourdin /* 1226823f3c5SBlaise Bourdin An exodusII truth table specifies which fields are saved at which time step 1236823f3c5SBlaise Bourdin It speeds up I/O but reserving space for fieldsin the file ahead of time. 1246823f3c5SBlaise Bourdin */ 1255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(numZonalVar * numCS, &truthtable)); 1266823f3c5SBlaise Bourdin for (i = 0; i < numZonalVar * numCS; ++i) truthtable[i] = 1; 127a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_put_truth_table,exoid, EX_ELEM_BLOCK, numCS, numZonalVar, truthtable); 1285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(truthtable)); 1296823f3c5SBlaise Bourdin 1306823f3c5SBlaise Bourdin /* Writing time step information in the file. Note that this is currently broken in the exodus library for netcdf4 (HDF5-based) files */ 1316823f3c5SBlaise Bourdin for (step = 0; step < numstep; ++step) { 1326823f3c5SBlaise Bourdin PetscReal time = step; 133a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_put_time,exoid, step+1, &time); 1346823f3c5SBlaise Bourdin } 1356823f3c5SBlaise Bourdin } 1366823f3c5SBlaise Bourdin 1376823f3c5SBlaise Bourdin /* Create the main section containing all fields */ 1385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ion)); 1395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetNumFields(section, 3)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetFieldName(section, fieldU, "U")); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetFieldName(section, fieldA, "Alpha")); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetFieldName(section, fieldS, "Sigma")); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetChart(dm, &pStart, &pEnd)); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetChart(section, pStart, pEnd)); 1455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(sdim+1, &pStartDepth, sdim+1, &pEndDepth)); 1465f80ce2aSJacob Faibussowitsch for (d = 0; d <= sdim; ++d) CHKERRQ(DMPlexGetDepthStratum(dm, d, &pStartDepth[d], &pEndDepth[d])); 147c4762a1bSJed Brown /* Vector field U, Scalar field Alpha, Tensor field Sigma */ 1485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetFieldComponents(section, fieldU, sdim)); 1495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetFieldComponents(section, fieldA, 1)); 1505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetFieldComponents(section, fieldS, sdim*(sdim+1)/2)); 151c4762a1bSJed Brown 152c4762a1bSJed Brown /* Going through cell sets then cells, and setting up storage for the sections */ 1535f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabelSize(dm, "Cell Sets", &numCS)); 1545f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabelIdIS(dm, "Cell Sets", &csIS)); 1555f80ce2aSJacob Faibussowitsch if (csIS) CHKERRQ(ISGetIndices(csIS, &csID)); 156c4762a1bSJed Brown for (set = 0; set < numCS; set++) { 157c4762a1bSJed Brown IS cellIS; 158c4762a1bSJed Brown const PetscInt *cellID; 159c4762a1bSJed Brown PetscInt numCells, cell, closureSize, *closureA = NULL; 160c4762a1bSJed Brown 1615f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetStratumSize(dm, "Cell Sets", csID[set], &numCells)); 1625f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetStratumIS(dm, "Cell Sets", csID[set], &cellIS)); 163c4762a1bSJed Brown if (numCells > 0) { 164c4762a1bSJed Brown /* dof layout ordered by increasing height in the DAG: cell, face, edge, vertex */ 165c4762a1bSJed Brown PetscInt dofUP1Tri[] = {2, 0, 0}; 166c4762a1bSJed Brown PetscInt dofAP1Tri[] = {1, 0, 0}; 167c4762a1bSJed Brown PetscInt dofUP2Tri[] = {2, 2, 0}; 168c4762a1bSJed Brown PetscInt dofAP2Tri[] = {1, 1, 0}; 169c4762a1bSJed Brown PetscInt dofUP1Quad[] = {2, 0, 0}; 170c4762a1bSJed Brown PetscInt dofAP1Quad[] = {1, 0, 0}; 171c4762a1bSJed Brown PetscInt dofUP2Quad[] = {2, 2, 2}; 172c4762a1bSJed Brown PetscInt dofAP2Quad[] = {1, 1, 1}; 173c4762a1bSJed Brown PetscInt dofS2D[] = {0, 0, 3}; 174c4762a1bSJed Brown PetscInt dofUP1Tet[] = {3, 0, 0, 0}; 175c4762a1bSJed Brown PetscInt dofAP1Tet[] = {1, 0, 0, 0}; 176c4762a1bSJed Brown PetscInt dofUP2Tet[] = {3, 3, 0, 0}; 177c4762a1bSJed Brown PetscInt dofAP2Tet[] = {1, 1, 0, 0}; 178c4762a1bSJed Brown PetscInt dofUP1Hex[] = {3, 0, 0, 0}; 179c4762a1bSJed Brown PetscInt dofAP1Hex[] = {1, 0, 0, 0}; 180c4762a1bSJed Brown PetscInt dofUP2Hex[] = {3, 3, 3, 3}; 181c4762a1bSJed Brown PetscInt dofAP2Hex[] = {1, 1, 1, 1}; 182c4762a1bSJed Brown PetscInt dofS3D[] = {0, 0, 0, 6}; 183c4762a1bSJed Brown PetscInt *dofU, *dofA, *dofS; 184c4762a1bSJed Brown 185c4762a1bSJed Brown switch (sdim) { 186c4762a1bSJed Brown case 2: dofS = dofS2D;break; 187c4762a1bSJed Brown case 3: dofS = dofS3D;break; 18898921bdaSJacob Faibussowitsch default: SETERRQ(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... */ 1935f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(cellIS, &cellID)); 1945f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetTransitiveClosure(dm, cellID[0], PETSC_TRUE, &closureSize, &closureA)); 195c4762a1bSJed Brown switch (closureSize) { 196c4762a1bSJed Brown case 7: /* Tri */ 197c4762a1bSJed Brown if (order == 1) { 198c4762a1bSJed Brown dofU = dofUP1Tri; 199c4762a1bSJed Brown dofA = dofAP1Tri; 200c4762a1bSJed Brown } else { 201c4762a1bSJed Brown dofU = dofUP2Tri; 202c4762a1bSJed Brown dofA = dofAP2Tri; 203c4762a1bSJed Brown } 204c4762a1bSJed Brown break; 205c4762a1bSJed Brown case 9: /* Quad */ 206c4762a1bSJed Brown if (order == 1) { 207c4762a1bSJed Brown dofU = dofUP1Quad; 208c4762a1bSJed Brown dofA = dofAP1Quad; 209c4762a1bSJed Brown } else { 210c4762a1bSJed Brown dofU = dofUP2Quad; 211c4762a1bSJed Brown dofA = dofAP2Quad; 212c4762a1bSJed Brown } 213c4762a1bSJed Brown break; 214c4762a1bSJed Brown case 15: /* Tet */ 215c4762a1bSJed Brown if (order == 1) { 216c4762a1bSJed Brown dofU = dofUP1Tet; 217c4762a1bSJed Brown dofA = dofAP1Tet; 218c4762a1bSJed Brown } else { 219c4762a1bSJed Brown dofU = dofUP2Tet; 220c4762a1bSJed Brown dofA = dofAP2Tet; 221c4762a1bSJed Brown } 222c4762a1bSJed Brown break; 223c4762a1bSJed Brown case 27: /* Hex */ 224c4762a1bSJed Brown if (order == 1) { 225c4762a1bSJed Brown dofU = dofUP1Hex; 226c4762a1bSJed Brown dofA = dofAP1Hex; 227c4762a1bSJed Brown } else { 228c4762a1bSJed Brown dofU = dofUP2Hex; 229c4762a1bSJed Brown dofA = dofAP2Hex; 230c4762a1bSJed Brown } 231c4762a1bSJed Brown break; 23298921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Unknown element with closure size %D", closureSize); 233c4762a1bSJed Brown } 2345f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreTransitiveClosure(dm, cellID[0], PETSC_TRUE, &closureSize, &closureA)); 235c4762a1bSJed Brown 236c4762a1bSJed Brown for (cell = 0; cell < numCells; cell++) { 237c4762a1bSJed Brown PetscInt *closure = NULL; 238c4762a1bSJed Brown 2395f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetTransitiveClosure(dm, cellID[cell], PETSC_TRUE, &closureSize, &closure)); 240c4762a1bSJed Brown for (p = 0; p < closureSize; ++p) { 241c4762a1bSJed Brown /* Find depth of p */ 242c4762a1bSJed Brown for (d = 0; d <= sdim; ++d) { 243c4762a1bSJed Brown if ((closure[2*p] >= pStartDepth[d]) && (closure[2*p] < pEndDepth[d])) { 2445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetDof(section, closure[2*p], dofU[d]+dofA[d]+dofS[d])); 2455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetFieldDof(section, closure[2*p], fieldU, dofU[d])); 2465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetFieldDof(section, closure[2*p], fieldA, dofA[d])); 2475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetFieldDof(section, closure[2*p], fieldS, dofS[d])); 248c4762a1bSJed Brown } 249c4762a1bSJed Brown } 250c4762a1bSJed Brown } 2515f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreTransitiveClosure(dm, cellID[cell], PETSC_TRUE, &closureSize, &closure)); 252c4762a1bSJed Brown } 2535f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(cellIS, &cellID)); 2545f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&cellIS)); 255c4762a1bSJed Brown } 256c4762a1bSJed Brown } 2575f80ce2aSJacob Faibussowitsch if (csIS) CHKERRQ(ISRestoreIndices(csIS, &csID)); 2585f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&csIS)); 2595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetUp(section)); 2605f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetLocalSection(dm, section)); 2615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectViewFromOptions((PetscObject) section, NULL, "-dm_section_view")); 2625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionDestroy(§ion)); 263c4762a1bSJed Brown 264c4762a1bSJed Brown { 265c4762a1bSJed Brown DM pdm; 266c4762a1bSJed Brown PetscSF migrationSF; 267c4762a1bSJed Brown PetscInt ovlp = 0; 268c4762a1bSJed Brown PetscPartitioner part; 269c4762a1bSJed Brown 2705f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUseNatural(dm,PETSC_TRUE)); 2715f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetPartitioner(dm,&part)); 2725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPartitionerSetFromOptions(part)); 2735f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexDistribute(dm,ovlp,&migrationSF,&pdm)); 274c4762a1bSJed Brown if (pdm) { 2755f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetMigrationSF(pdm,migrationSF)); 2765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFDestroy(&migrationSF)); 2775f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 278c4762a1bSJed Brown dm = pdm; 2795f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(dm,NULL,"-dm_view")); 280c4762a1bSJed Brown } 281c4762a1bSJed Brown } 282c4762a1bSJed Brown 283c4762a1bSJed Brown /* Get DM and IS for each field of dm */ 2845f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateSubDM(dm, 1, &fieldU, &isU, &dmU)); 2855f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateSubDM(dm, 1, &fieldA, &isA, &dmA)); 2865f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateSubDM(dm, 1, &fieldS, &isS, &dmS)); 2875f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateSubDM(dm, 2, fieldUA, &isUA, &dmUA)); 288c4762a1bSJed Brown 2895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(2,&dmList)); 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; 2965f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateSuperDM(dmList,2,NULL,&dmUA2)); 297c4762a1bSJed Brown dmU->useNatural = PETSC_TRUE; 2985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(dmList)); 299c4762a1bSJed Brown 3005f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(dm, &X)); 3015f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(dmU, &U)); 3025f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(dmA, &A)); 3035f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(dmS, &S)); 3045f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(dmUA, &UA)); 3055f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(dmUA2, &UA2)); 306c4762a1bSJed Brown 3075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) U, "U")); 3085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) A, "Alpha")); 3095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) S, "Sigma")); 3105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) UA, "UAlpha")); 3115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) UA2, "UAlpha2")); 3125f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(X, -111.)); 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; 3216823f3c5SBlaise Bourdin PetscInt clSize, i, j; 322c4762a1bSJed Brown 3235f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalSection(dmUA, §ionUA)); 3245f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(dmUA, &UALoc)); 3255f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(UALoc, &cval)); 3265f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateSection(dmUA, &coordSection)); 3275f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocal(dmUA, &coord)); 3285f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetChart(dmUA, &pStart, &pEnd)); 3296823f3c5SBlaise Bourdin 330c4762a1bSJed Brown for (p = pStart; p < pEnd; ++p) { 331c4762a1bSJed Brown PetscInt dofUA, offUA; 332c4762a1bSJed Brown 3335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetDof(sectionUA, p, &dofUA)); 334c4762a1bSJed Brown if (dofUA > 0) { 3356823f3c5SBlaise Bourdin xyz=NULL; 3365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetOffset(sectionUA, p, &offUA)); 3375f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecGetClosure(dmUA, coordSection, coord, p, &clSize, &xyz)); 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 } 3475f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecRestoreClosure(dmUA, coordSection, coord, p, &clSize, &xyz)); 348c4762a1bSJed Brown } 349c4762a1bSJed Brown } 3505f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(UALoc, &cval)); 3515f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalBegin(dmUA, UALoc, INSERT_VALUES, UA)); 3525f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalEnd(dmUA, UALoc, INSERT_VALUES, UA)); 3535f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(dmUA, &UALoc)); 3546823f3c5SBlaise Bourdin 355c4762a1bSJed Brown /* Update X */ 3565f80ce2aSJacob Faibussowitsch CHKERRQ(VecISCopy(X, isUA, SCATTER_FORWARD, UA)); 3575f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(UA, NULL, "-ua_vec_view")); 3586823f3c5SBlaise Bourdin 359c4762a1bSJed Brown /* Restrict to U and Alpha */ 3605f80ce2aSJacob Faibussowitsch CHKERRQ(VecISCopy(X, isU, SCATTER_REVERSE, U)); 3615f80ce2aSJacob Faibussowitsch CHKERRQ(VecISCopy(X, isA, SCATTER_REVERSE, A)); 3626823f3c5SBlaise Bourdin 363c4762a1bSJed Brown /* restrict to UA2 */ 3645f80ce2aSJacob Faibussowitsch CHKERRQ(VecISCopy(X, isUA, SCATTER_REVERSE, UA2)); 3655f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(UA2, NULL, "-ua2_vec_view")); 366c4762a1bSJed Brown } 367c4762a1bSJed Brown 368c4762a1bSJed Brown { 369c4762a1bSJed Brown Vec tmpVec; 370c4762a1bSJed Brown PetscSection coordSection; 371c4762a1bSJed Brown Vec coord; 372c4762a1bSJed Brown PetscReal norm; 3736823f3c5SBlaise Bourdin PetscReal time = 1.234; 374c4762a1bSJed Brown 375c4762a1bSJed Brown /* Writing nodal variables to ExodusII file */ 3765f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetOutputSequenceNumber(dmU,0,time)); 3775f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetOutputSequenceNumber(dmA,0,time)); 3786823f3c5SBlaise Bourdin 3795f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(U, viewer)); 3805f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(A, viewer)); 3816823f3c5SBlaise Bourdin 382c4762a1bSJed Brown /* Saving U and Alpha in one shot. 383c4762a1bSJed Brown For this, we need to cheat and change the Vec's name 3846823f3c5SBlaise Bourdin Note that in the end we write variables one component at a time, 3856823f3c5SBlaise Bourdin so that there is no real values in doing this 3866823f3c5SBlaise Bourdin */ 3876823f3c5SBlaise Bourdin 3885f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetOutputSequenceNumber(dmUA,1,time)); 3895f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(dmUA, &tmpVec)); 3905f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(UA, tmpVec)); 3915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) tmpVec, "U")); 3925f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(tmpVec, viewer)); 393c4762a1bSJed Brown /* Reading nodal variables in Exodus file */ 3945f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(tmpVec, -1000.0)); 3955f80ce2aSJacob Faibussowitsch CHKERRQ(VecLoad(tmpVec, viewer)); 3965f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(UA, -1.0, tmpVec)); 3975f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(UA, NORM_INFINITY, &norm)); 3982c71b3e2SJacob Faibussowitsch PetscCheckFalse(norm > PETSC_SQRT_MACHINE_EPSILON,PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "UAlpha ||Vin - Vout|| = %g", (double) norm); 3995f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(dmUA, &tmpVec)); 400c4762a1bSJed Brown 401c4762a1bSJed Brown /* same thing with the UA2 Vec obtained from the superDM */ 4025f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(dmUA2, &tmpVec)); 4035f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(UA2, tmpVec)); 4045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) tmpVec, "U")); 4055f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetOutputSequenceNumber(dmUA2,2,time)); 4065f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(tmpVec, viewer)); 407c4762a1bSJed Brown /* Reading nodal variables in Exodus file */ 4085f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(tmpVec, -1000.0)); 4095f80ce2aSJacob Faibussowitsch CHKERRQ(VecLoad(tmpVec,viewer)); 4105f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(UA2, -1.0, tmpVec)); 4115f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(UA2, NORM_INFINITY, &norm)); 4122c71b3e2SJacob Faibussowitsch PetscCheckFalse(norm > PETSC_SQRT_MACHINE_EPSILON,PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "UAlpha2 ||Vin - Vout|| = %g", (double) norm); 4135f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(dmUA2, &tmpVec)); 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 4186823f3c5SBlaise Bourdin sigma_2 = x_coordinate of the cell center of mass 4196823f3c5SBlaise Bourdin */ 4205f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateSection(dmS, &coordSection)); 4215f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocal(dmS, &coord)); 4225f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabelIdIS(dmS, "Cell Sets", &csIS)); 4235f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabelSize(dmS, "Cell Sets", &numCS)); 4245f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(csIS, &csID)); 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 4335f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetStratumIS(dmS, "Cell Sets", csID[set], &cellIS)); 4345f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(cellIS, &cellID)); 4355f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetSize(cellIS, &numCells)); 436c4762a1bSJed Brown for (cell = 0; cell < numCells; cell++) { 4375f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecGetClosure(dmS, NULL, S, cellID[cell], &clSize, &cval)); 4385f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecGetClosure(dmS, coordSection, coord, cellID[cell], &cdimCoord, &xyz)); 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; 4445f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecSetClosure(dmS, NULL, S, cellID[cell], cval, INSERT_ALL_VALUES)); 445c4762a1bSJed Brown } 4465f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecRestoreClosure(dmS, NULL, S, cellID[0], &clSize, &cval)); 4475f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecRestoreClosure(dmS, coordSection, coord, cellID[0], NULL, &xyz)); 4485f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(cellIS, &cellID)); 4495f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&cellIS)); 450c4762a1bSJed Brown } 4515f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(csIS, &csID)); 4525f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&csIS)); 4535f80ce2aSJacob Faibussowitsch CHKERRQ(VecViewFromOptions(S, NULL, "-s_vec_view")); 4546823f3c5SBlaise Bourdin 455c4762a1bSJed Brown /* Writing zonal variables in Exodus file */ 4565f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetOutputSequenceNumber(dmS,0,time)); 4575f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(S,viewer)); 4586823f3c5SBlaise Bourdin 459c4762a1bSJed Brown /* Reading zonal variables in Exodus file */ 4605f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(dmS, &tmpVec)); 4615f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(tmpVec, -1000.0)); 4625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) tmpVec, "Sigma")); 4635f80ce2aSJacob Faibussowitsch CHKERRQ(VecLoad(tmpVec,viewer)); 4645f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(S, -1.0, tmpVec)); 4655f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(S, NORM_INFINITY, &norm)); 4662c71b3e2SJacob Faibussowitsch PetscCheckFalse(norm > PETSC_SQRT_MACHINE_EPSILON,PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "Sigma ||Vin - Vout|| = %g", (double) norm); 4675f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(dmS, &tmpVec)); 468c4762a1bSJed Brown } 4695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 470c4762a1bSJed Brown 4715f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(dmUA2, &UA2)); 4725f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(dmUA, &UA)); 4735f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(dmS, &S)); 4745f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(dmA, &A)); 4755f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(dmU, &U)); 4765f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(dm, &X)); 4775f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dmU)); CHKERRQ(ISDestroy(&isU)); 4785f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dmA)); CHKERRQ(ISDestroy(&isA)); 4795f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dmS)); CHKERRQ(ISDestroy(&isS)); 4805f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dmUA));CHKERRQ(ISDestroy(&isUA)); 4815f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dmUA2)); 4825f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 4835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(pStartDepth, pEndDepth)); 484*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 485*b122ec5aSJacob Faibussowitsch return 0; 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