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; 31c4762a1bSJed Brown PetscErrorCode ierr; 32c4762a1bSJed Brown 33c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv,NULL, help);if (ierr) return ierr; 34c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRQ(ierr); 35c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRQ(ierr); 36c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FEM Layout Options", "ex26");CHKERRQ(ierr); 37c4762a1bSJed Brown ierr = PetscOptionsString("-i", "Filename to read", "ex26", ifilename, ifilename, sizeof(ifilename), NULL);CHKERRQ(ierr); 38c4762a1bSJed Brown ierr = PetscOptionsString("-o", "Filename to write", "ex26", ofilename, ofilename, sizeof(ofilename), NULL);CHKERRQ(ierr); 39c4762a1bSJed Brown ierr = PetscOptionsBoundedInt("-order", "FEM polynomial order", "ex26", order, &order, NULL,1);CHKERRQ(ierr); 40c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 41c4762a1bSJed Brown if ((order > 2) || (order < 1)) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported polynomial order %D not in [1, 2]", order); 42c4762a1bSJed Brown 43c4762a1bSJed Brown ex_opts(EX_VERBOSE+EX_DEBUG); 44c4762a1bSJed Brown ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename, PETSC_TRUE, &dm);CHKERRQ(ierr); 45c4762a1bSJed Brown ierr = DMSetFromOptions(dm);CHKERRQ(ierr); 46c4762a1bSJed Brown ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr); 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* Create the main section containning all fields */ 49c4762a1bSJed Brown ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), §ion);CHKERRQ(ierr); 50c4762a1bSJed Brown ierr = PetscSectionSetNumFields(section, 3);CHKERRQ(ierr); 51c4762a1bSJed Brown ierr = PetscSectionSetFieldName(section, fieldU, "U");CHKERRQ(ierr); 52c4762a1bSJed Brown ierr = PetscSectionSetFieldName(section, fieldA, "Alpha");CHKERRQ(ierr); 53c4762a1bSJed Brown ierr = PetscSectionSetFieldName(section, fieldS, "Sigma");CHKERRQ(ierr); 54c4762a1bSJed Brown ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 55c4762a1bSJed Brown ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr); 56c4762a1bSJed Brown ierr = DMGetDimension(dm, &sdim);CHKERRQ(ierr); 57c4762a1bSJed Brown ierr = PetscMalloc2(sdim+1, &pStartDepth, sdim+1, &pEndDepth);CHKERRQ(ierr); 58c4762a1bSJed Brown for (d = 0; d <= sdim; ++d) {ierr = DMPlexGetDepthStratum(dm, d, &pStartDepth[d], &pEndDepth[d]);CHKERRQ(ierr);} 59c4762a1bSJed Brown /* Vector field U, Scalar field Alpha, Tensor field Sigma */ 60c4762a1bSJed Brown ierr = PetscSectionSetFieldComponents(section, fieldU, sdim);CHKERRQ(ierr); 61c4762a1bSJed Brown ierr = PetscSectionSetFieldComponents(section, fieldA, 1);CHKERRQ(ierr); 62c4762a1bSJed Brown ierr = PetscSectionSetFieldComponents(section, fieldS, sdim*(sdim+1)/2);CHKERRQ(ierr); 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* Going through cell sets then cells, and setting up storage for the sections */ 65c4762a1bSJed Brown ierr = DMGetLabelSize(dm, "Cell Sets", &numCS); 66c4762a1bSJed Brown ierr = DMGetLabelIdIS(dm, "Cell Sets", &csIS);CHKERRQ(ierr); 67c4762a1bSJed Brown if (csIS) {ierr = ISGetIndices(csIS, &csID);CHKERRQ(ierr);} 68c4762a1bSJed Brown for (set = 0; set < numCS; set++) { 69c4762a1bSJed Brown IS cellIS; 70c4762a1bSJed Brown const PetscInt *cellID; 71c4762a1bSJed Brown PetscInt numCells, cell, closureSize, *closureA = NULL; 72c4762a1bSJed Brown 73c4762a1bSJed Brown ierr = DMGetStratumSize(dm, "Cell Sets", csID[set], &numCells);CHKERRQ(ierr); 74c4762a1bSJed Brown ierr = DMGetStratumIS(dm, "Cell Sets", csID[set], &cellIS);CHKERRQ(ierr); 75c4762a1bSJed Brown if (numCells > 0) { 76c4762a1bSJed Brown /* dof layout ordered by increasing height in the DAG: cell, face, edge, vertex */ 77c4762a1bSJed Brown PetscInt dofUP1Tri[] = {2, 0, 0}; 78c4762a1bSJed Brown PetscInt dofAP1Tri[] = {1, 0, 0}; 79c4762a1bSJed Brown PetscInt dofUP2Tri[] = {2, 2, 0}; 80c4762a1bSJed Brown PetscInt dofAP2Tri[] = {1, 1, 0}; 81c4762a1bSJed Brown PetscInt dofUP1Quad[] = {2, 0, 0}; 82c4762a1bSJed Brown PetscInt dofAP1Quad[] = {1, 0, 0}; 83c4762a1bSJed Brown PetscInt dofUP2Quad[] = {2, 2, 2}; 84c4762a1bSJed Brown PetscInt dofAP2Quad[] = {1, 1, 1}; 85c4762a1bSJed Brown PetscInt dofS2D[] = {0, 0, 3}; 86c4762a1bSJed Brown PetscInt dofUP1Tet[] = {3, 0, 0, 0}; 87c4762a1bSJed Brown PetscInt dofAP1Tet[] = {1, 0, 0, 0}; 88c4762a1bSJed Brown PetscInt dofUP2Tet[] = {3, 3, 0, 0}; 89c4762a1bSJed Brown PetscInt dofAP2Tet[] = {1, 1, 0, 0}; 90c4762a1bSJed Brown PetscInt dofUP1Hex[] = {3, 0, 0, 0}; 91c4762a1bSJed Brown PetscInt dofAP1Hex[] = {1, 0, 0, 0}; 92c4762a1bSJed Brown PetscInt dofUP2Hex[] = {3, 3, 3, 3}; 93c4762a1bSJed Brown PetscInt dofAP2Hex[] = {1, 1, 1, 1}; 94c4762a1bSJed Brown PetscInt dofS3D[] = {0, 0, 0, 6}; 95c4762a1bSJed Brown PetscInt *dofU, *dofA, *dofS; 96c4762a1bSJed Brown 97c4762a1bSJed Brown switch (sdim) { 98c4762a1bSJed Brown case 2: dofS = dofS2D;break; 99c4762a1bSJed Brown case 3: dofS = dofS3D;break; 100c4762a1bSJed Brown default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No layout for dimension %D", sdim); 101c4762a1bSJed Brown } 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* Identify cell type based on closure size only. This works for Tri/Tet/Quad/Hex meshes 104c4762a1bSJed Brown It will not be enough to identify more exotic elements like pyramid or prisms... */ 105c4762a1bSJed Brown ierr = ISGetIndices(cellIS, &cellID);CHKERRQ(ierr); 106c4762a1bSJed Brown ierr = DMPlexGetTransitiveClosure(dm, cellID[0], PETSC_TRUE, &closureSize, &closureA);CHKERRQ(ierr); 107c4762a1bSJed Brown switch (closureSize) { 108c4762a1bSJed Brown case 7: /* Tri */ 109c4762a1bSJed Brown if (order == 1) { 110c4762a1bSJed Brown dofU = dofUP1Tri; 111c4762a1bSJed Brown dofA = dofAP1Tri; 112c4762a1bSJed Brown } else { 113c4762a1bSJed Brown dofU = dofUP2Tri; 114c4762a1bSJed Brown dofA = dofAP2Tri; 115c4762a1bSJed Brown } 116c4762a1bSJed Brown break; 117c4762a1bSJed Brown case 9: /* Quad */ 118c4762a1bSJed Brown if (order == 1) { 119c4762a1bSJed Brown dofU = dofUP1Quad; 120c4762a1bSJed Brown dofA = dofAP1Quad; 121c4762a1bSJed Brown } else { 122c4762a1bSJed Brown dofU = dofUP2Quad; 123c4762a1bSJed Brown dofA = dofAP2Quad; 124c4762a1bSJed Brown } 125c4762a1bSJed Brown break; 126c4762a1bSJed Brown case 15: /* Tet */ 127c4762a1bSJed Brown if (order == 1) { 128c4762a1bSJed Brown dofU = dofUP1Tet; 129c4762a1bSJed Brown dofA = dofAP1Tet; 130c4762a1bSJed Brown } else { 131c4762a1bSJed Brown dofU = dofUP2Tet; 132c4762a1bSJed Brown dofA = dofAP2Tet; 133c4762a1bSJed Brown } 134c4762a1bSJed Brown break; 135c4762a1bSJed Brown case 27: /* Hex */ 136c4762a1bSJed Brown if (order == 1) { 137c4762a1bSJed Brown dofU = dofUP1Hex; 138c4762a1bSJed Brown dofA = dofAP1Hex; 139c4762a1bSJed Brown } else { 140c4762a1bSJed Brown dofU = dofUP2Hex; 141c4762a1bSJed Brown dofA = dofAP2Hex; 142c4762a1bSJed Brown } 143c4762a1bSJed Brown break; 144c4762a1bSJed Brown default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Unknown element with closure size %D\n", closureSize); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown ierr = DMPlexRestoreTransitiveClosure(dm, cellID[0], PETSC_TRUE, &closureSize, &closureA);CHKERRQ(ierr); 147c4762a1bSJed Brown 148c4762a1bSJed Brown for (cell = 0; cell < numCells; cell++) { 149c4762a1bSJed Brown PetscInt *closure = NULL; 150c4762a1bSJed Brown 151c4762a1bSJed Brown ierr = DMPlexGetTransitiveClosure(dm, cellID[cell], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 152c4762a1bSJed Brown for (p = 0; p < closureSize; ++p) { 153c4762a1bSJed Brown /* Find depth of p */ 154c4762a1bSJed Brown for (d = 0; d <= sdim; ++d) { 155c4762a1bSJed Brown if ((closure[2*p] >= pStartDepth[d]) && (closure[2*p] < pEndDepth[d])) { 156c4762a1bSJed Brown ierr = PetscSectionSetDof(section, closure[2*p], dofU[d]+dofA[d]+dofS[d]);CHKERRQ(ierr); 157c4762a1bSJed Brown ierr = PetscSectionSetFieldDof(section, closure[2*p], fieldU, dofU[d]);CHKERRQ(ierr); 158c4762a1bSJed Brown ierr = PetscSectionSetFieldDof(section, closure[2*p], fieldA, dofA[d]);CHKERRQ(ierr); 159c4762a1bSJed Brown ierr = PetscSectionSetFieldDof(section, closure[2*p], fieldS, dofS[d]);CHKERRQ(ierr); 160c4762a1bSJed Brown } 161c4762a1bSJed Brown } 162c4762a1bSJed Brown } 163c4762a1bSJed Brown ierr = DMPlexRestoreTransitiveClosure(dm, cellID[cell], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 164c4762a1bSJed Brown } 165c4762a1bSJed Brown ierr = ISRestoreIndices(cellIS, &cellID);CHKERRQ(ierr); 166c4762a1bSJed Brown ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 167c4762a1bSJed Brown } 168c4762a1bSJed Brown } 169c4762a1bSJed Brown if (csIS) {ierr = ISRestoreIndices(csIS, &csID);CHKERRQ(ierr);} 170c4762a1bSJed Brown ierr = ISDestroy(&csIS);CHKERRQ(ierr); 171c4762a1bSJed Brown ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 172c4762a1bSJed Brown ierr = DMSetLocalSection(dm, section);CHKERRQ(ierr); 173c4762a1bSJed Brown ierr = PetscObjectViewFromOptions((PetscObject) section, NULL, "-dm_section_view");CHKERRQ(ierr); 174c4762a1bSJed Brown ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 175c4762a1bSJed Brown 176c4762a1bSJed Brown { 177c4762a1bSJed Brown /* TODO: Replace with ExodusII viewer */ 178c4762a1bSJed Brown /* Create the exodus result file */ 179c4762a1bSJed Brown PetscInt numstep = 3, step; 180c4762a1bSJed Brown char *nodalVarName[4]; 181c4762a1bSJed Brown char *zonalVarName[6]; 182c4762a1bSJed Brown int *truthtable; 183c4762a1bSJed Brown PetscInt numNodalVar, numZonalVar, i; 184c4762a1bSJed Brown int CPU_word_size, IO_word_size, EXO_mode; 185c4762a1bSJed Brown 186c4762a1bSJed Brown ex_opts(EX_VERBOSE+EX_DEBUG); 187c4762a1bSJed Brown if (!rank) { 188c4762a1bSJed Brown CPU_word_size = sizeof(PetscReal); 189c4762a1bSJed Brown IO_word_size = sizeof(PetscReal); 190c4762a1bSJed Brown EXO_mode = EX_CLOBBER; 191c4762a1bSJed Brown #if defined(PETSC_USE_64BIT_INDICES) 192c4762a1bSJed Brown EXO_mode += EX_ALL_INT64_API; 193c4762a1bSJed Brown #endif 194c4762a1bSJed Brown exoid = ex_create(ofilename, EXO_mode, &CPU_word_size, &IO_word_size); 195c4762a1bSJed Brown if (exoid < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Unable to open exodus file %\n", ofilename); 196c4762a1bSJed Brown } 197c4762a1bSJed Brown ierr = DMPlexView_ExodusII_Internal(dm, exoid, order);CHKERRQ(ierr); 198c4762a1bSJed Brown 199c4762a1bSJed Brown if (!rank) { 200c4762a1bSJed Brown /* "Format" the exodus result file, i.e. allocate space for nodal and zonal variables */ 201c4762a1bSJed Brown switch (sdim) { 202c4762a1bSJed Brown case 2: 203c4762a1bSJed Brown numNodalVar = 3; 204c4762a1bSJed Brown nodalVarName[0] = (char *) "U_x"; 205c4762a1bSJed Brown nodalVarName[1] = (char *) "U_y"; 206c4762a1bSJed Brown nodalVarName[2] = (char *) "Alpha"; 207c4762a1bSJed Brown numZonalVar = 3; 208c4762a1bSJed Brown zonalVarName[0] = (char *) "Sigma_11"; 209c4762a1bSJed Brown zonalVarName[1] = (char *) "Sigma_22"; 210c4762a1bSJed Brown zonalVarName[2] = (char *) "Sigma_12"; 211c4762a1bSJed Brown break; 212c4762a1bSJed Brown case 3: 213c4762a1bSJed Brown numNodalVar = 4; 214c4762a1bSJed Brown nodalVarName[0] = (char *) "U_x"; 215c4762a1bSJed Brown nodalVarName[1] = (char *) "U_y"; 216c4762a1bSJed Brown nodalVarName[2] = (char *) "U_z"; 217c4762a1bSJed Brown nodalVarName[3] = (char *) "Alpha"; 218c4762a1bSJed Brown numZonalVar = 6; 219c4762a1bSJed Brown zonalVarName[0] = (char *) "Sigma_11"; 220c4762a1bSJed Brown zonalVarName[1] = (char *) "Sigma_22"; 221c4762a1bSJed Brown zonalVarName[2] = (char *) "Sigma_33"; 222c4762a1bSJed Brown zonalVarName[3] = (char *) "Sigma_23"; 223c4762a1bSJed Brown zonalVarName[4] = (char *) "Sigma_13"; 224c4762a1bSJed Brown zonalVarName[5] = (char *) "Sigma_12"; 225c4762a1bSJed Brown break; 226c4762a1bSJed Brown default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No layout for dimension %D", sdim); 227c4762a1bSJed Brown } 228c4762a1bSJed Brown ierr = ex_put_variable_param(exoid, EX_ELEM_BLOCK, numZonalVar);CHKERRQ(ierr); 229c4762a1bSJed Brown ierr = ex_put_variable_names(exoid, EX_ELEM_BLOCK, numZonalVar, zonalVarName);CHKERRQ(ierr); 230c4762a1bSJed Brown ierr = ex_put_variable_param(exoid, EX_NODAL, numNodalVar);CHKERRQ(ierr); 231c4762a1bSJed Brown ierr = ex_put_variable_names(exoid, EX_NODAL, numNodalVar, nodalVarName);CHKERRQ(ierr); 232c4762a1bSJed Brown numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK); 233c4762a1bSJed Brown ierr = PetscMalloc1(numZonalVar * numCS, &truthtable);CHKERRQ(ierr); 234c4762a1bSJed Brown for (i = 0; i < numZonalVar * numCS; ++i) truthtable[i] = 1; 235c4762a1bSJed Brown ierr = ex_put_truth_table(exoid, EX_ELEM_BLOCK, numCS, numZonalVar, truthtable);CHKERRQ(ierr); 236c4762a1bSJed Brown ierr = PetscFree(truthtable);CHKERRQ(ierr); 237c4762a1bSJed Brown /* Writing time step information in the file. Note that this is currently broken in the exodus library for netcdf4 (HDF5-based) files */ 238c4762a1bSJed Brown for (step = 0; step < numstep; ++step) { 239c4762a1bSJed Brown PetscReal time = step; 240c4762a1bSJed Brown ierr = ex_put_time(exoid, step+1, &time);CHKERRQ(ierr); 241c4762a1bSJed Brown } 242c4762a1bSJed Brown ierr = ex_close(exoid);CHKERRQ(ierr); 243c4762a1bSJed Brown } 244c4762a1bSJed Brown } 245c4762a1bSJed Brown 246c4762a1bSJed Brown { 247c4762a1bSJed Brown DM pdm; 248c4762a1bSJed Brown PetscSF migrationSF; 249c4762a1bSJed Brown PetscInt ovlp = 0; 250c4762a1bSJed Brown PetscPartitioner part; 251c4762a1bSJed Brown 252c4762a1bSJed Brown ierr = DMSetUseNatural(dm,PETSC_TRUE);CHKERRQ(ierr); 253c4762a1bSJed Brown ierr = DMPlexGetPartitioner(dm,&part);CHKERRQ(ierr); 254c4762a1bSJed Brown ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 255c4762a1bSJed Brown ierr = DMPlexDistribute(dm,ovlp,&migrationSF,&pdm);CHKERRQ(ierr); 256c4762a1bSJed Brown if (pdm) { 257c4762a1bSJed Brown ierr = DMPlexSetMigrationSF(pdm,migrationSF);CHKERRQ(ierr); 258c4762a1bSJed Brown ierr = PetscSFDestroy(&migrationSF);CHKERRQ(ierr); 259c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 260c4762a1bSJed Brown dm = pdm; 261c4762a1bSJed Brown ierr = DMViewFromOptions(dm,NULL,"-dm_view");CHKERRQ(ierr); 262c4762a1bSJed Brown } 263c4762a1bSJed Brown } 264c4762a1bSJed Brown 265c4762a1bSJed Brown { 266c4762a1bSJed Brown /* TODO Replace with ExodusII viewer */ 267c4762a1bSJed Brown /* Reopen the exodus result file on all processors */ 268c4762a1bSJed Brown MPI_Info mpi_info = MPI_INFO_NULL; 269c4762a1bSJed Brown int CPU_word_size, IO_word_size, EXO_mode; 270c4762a1bSJed Brown float EXO_version; 271c4762a1bSJed Brown 272c4762a1bSJed Brown EXO_mode = EX_WRITE; 273c4762a1bSJed Brown CPU_word_size = sizeof(PetscReal); 274c4762a1bSJed Brown IO_word_size = sizeof(PetscReal); 275c4762a1bSJed Brown #if defined(PETSC_USE_64BIT_INDICES) 276c4762a1bSJed Brown EXO_mode += EX_ALL_INT64_API; 277c4762a1bSJed Brown #endif 278c4762a1bSJed Brown exoid = ex_open_par(ofilename, EXO_mode, &CPU_word_size, &IO_word_size, &EXO_version, PetscObjectComm((PetscObject) dm), mpi_info); 279c4762a1bSJed Brown } 280c4762a1bSJed Brown 281c4762a1bSJed Brown /* Get DM and IS for each field of dm */ 282c4762a1bSJed Brown ierr = DMCreateSubDM(dm, 1, &fieldU, &isU, &dmU);CHKERRQ(ierr); 283c4762a1bSJed Brown ierr = DMCreateSubDM(dm, 1, &fieldA, &isA, &dmA);CHKERRQ(ierr); 284c4762a1bSJed Brown ierr = DMCreateSubDM(dm, 1, &fieldS, &isS, &dmS);CHKERRQ(ierr); 285c4762a1bSJed Brown ierr = DMCreateSubDM(dm, 2, fieldUA, &isUA, &dmUA);CHKERRQ(ierr); 286c4762a1bSJed Brown 287c4762a1bSJed Brown ierr = PetscMalloc1(2,&dmList);CHKERRQ(ierr); 288c4762a1bSJed Brown dmList[0] = dmU; 289c4762a1bSJed Brown dmList[1] = dmA; 290c4762a1bSJed Brown /* We temporarily disable dmU->useNatural to test that we can reconstruct the 291c4762a1bSJed Brown NaturaltoGlobal SF from any of the dm in dms 292c4762a1bSJed Brown */ 293c4762a1bSJed Brown dmU->useNatural = PETSC_FALSE; 294c4762a1bSJed Brown ierr = DMCreateSuperDM(dmList,2,NULL,&dmUA2);CHKERRQ(ierr); 295c4762a1bSJed Brown dmU->useNatural = PETSC_TRUE; 296c4762a1bSJed Brown ierr = PetscFree(dmList);CHKERRQ(ierr); 297c4762a1bSJed Brown 298c4762a1bSJed Brown ierr = DMGetGlobalVector(dm, &X);CHKERRQ(ierr); 299c4762a1bSJed Brown ierr = DMGetGlobalVector(dmU, &U);CHKERRQ(ierr); 300c4762a1bSJed Brown ierr = DMGetGlobalVector(dmA, &A);CHKERRQ(ierr); 301c4762a1bSJed Brown ierr = DMGetGlobalVector(dmS, &S);CHKERRQ(ierr); 302c4762a1bSJed Brown ierr = DMGetGlobalVector(dmUA, &UA);CHKERRQ(ierr); 303c4762a1bSJed Brown ierr = DMGetGlobalVector(dmUA2, &UA2);CHKERRQ(ierr); 304c4762a1bSJed Brown 305c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) U, "U");CHKERRQ(ierr); 306c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) A, "Alpha");CHKERRQ(ierr); 307c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) S, "Sigma");CHKERRQ(ierr); 308c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) UA, "UAlpha");CHKERRQ(ierr); 309c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) UA2, "UAlpha2");CHKERRQ(ierr); 310c4762a1bSJed Brown ierr = VecSet(X, -111.);CHKERRQ(ierr); 311c4762a1bSJed Brown 312c4762a1bSJed 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 */ 313c4762a1bSJed Brown { 314c4762a1bSJed Brown PetscSection sectionUA; 315c4762a1bSJed Brown Vec UALoc; 316c4762a1bSJed Brown PetscSection coordSection; 317c4762a1bSJed Brown Vec coord; 318c4762a1bSJed Brown PetscScalar *cval, *xyz; 319c4762a1bSJed Brown PetscInt cdimCoord = 24; 320c4762a1bSJed Brown 321c4762a1bSJed Brown ierr = DMGetLocalSection(dmUA, §ionUA);CHKERRQ(ierr); 322c4762a1bSJed Brown ierr = DMGetLocalVector(dmUA, &UALoc);CHKERRQ(ierr); 323c4762a1bSJed Brown ierr = VecGetArray(UALoc, &cval);CHKERRQ(ierr); 324c4762a1bSJed Brown ierr = DMGetCoordinateSection(dmUA, &coordSection);CHKERRQ(ierr); 325c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(dmUA, &coord);CHKERRQ(ierr); 326c4762a1bSJed Brown ierr = DMPlexGetChart(dmUA, &pStart, &pEnd);CHKERRQ(ierr); 327c4762a1bSJed Brown ierr = PetscMalloc1(cdimCoord, &xyz);CHKERRQ(ierr); 328c4762a1bSJed Brown /* I know that UA has 0 or sdim+1 dof at each point, since both U and Alpha use the same discretization 329c4762a1bSJed Brown The maximum possible size for the closure of the coordinate section is 8*3 at the cell center */ 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) { 335c4762a1bSJed Brown PetscInt clSize = cdimCoord, i, j; 336c4762a1bSJed Brown 337c4762a1bSJed Brown ierr = PetscSectionGetOffset(sectionUA, p, &offUA);CHKERRQ(ierr); 338c4762a1bSJed Brown ierr = DMPlexVecGetClosure(dmUA, coordSection, coord, p, &clSize, &xyz);CHKERRQ(ierr); 339c4762a1bSJed Brown cval[offUA+sdim] = 0.; 340c4762a1bSJed Brown for (i = 0; i < sdim; ++i) { 341c4762a1bSJed Brown cval[offUA+i] = 0; 342c4762a1bSJed Brown for (j = 0; j < clSize/sdim; ++j) { 343c4762a1bSJed Brown cval[offUA+i] += xyz[j*sdim+i]; 344c4762a1bSJed Brown } 345c4762a1bSJed Brown cval[offUA+i] = cval[offUA+i] * sdim / clSize; 346c4762a1bSJed Brown cval[offUA+sdim] += PetscSqr(cval[offUA+i]); 347c4762a1bSJed Brown } 348c4762a1bSJed Brown } 349c4762a1bSJed Brown } 350c4762a1bSJed Brown ierr = PetscFree(xyz);CHKERRQ(ierr); 351c4762a1bSJed Brown ierr = VecRestoreArray(UALoc, &cval);CHKERRQ(ierr); 352c4762a1bSJed Brown ierr = DMLocalToGlobalBegin(dmUA, UALoc, INSERT_VALUES, UA);CHKERRQ(ierr); 353c4762a1bSJed Brown ierr = DMLocalToGlobalEnd(dmUA, UALoc, INSERT_VALUES, UA);CHKERRQ(ierr); 354c4762a1bSJed Brown ierr = DMRestoreLocalVector(dmUA, &UALoc);CHKERRQ(ierr); 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); 358c4762a1bSJed Brown /* Restrict to U and Alpha */ 359c4762a1bSJed Brown ierr = VecISCopy(X, isU, SCATTER_REVERSE, U);CHKERRQ(ierr); 360c4762a1bSJed Brown ierr = VecISCopy(X, isA, SCATTER_REVERSE, A);CHKERRQ(ierr); 361c4762a1bSJed Brown /* restrict to UA2 */ 362c4762a1bSJed Brown ierr = VecISCopy(X, isUA, SCATTER_REVERSE, UA2);CHKERRQ(ierr); 363c4762a1bSJed Brown ierr = VecViewFromOptions(UA2, NULL, "-ua2_vec_view");CHKERRQ(ierr); 364c4762a1bSJed Brown } 365c4762a1bSJed Brown 366c4762a1bSJed Brown { 367c4762a1bSJed Brown Vec tmpVec; 368c4762a1bSJed Brown PetscSection coordSection; 369c4762a1bSJed Brown Vec coord; 370c4762a1bSJed Brown PetscReal norm; 371c4762a1bSJed Brown 372c4762a1bSJed Brown /* Writing nodal variables to ExodusII file */ 373c4762a1bSJed Brown ierr = VecViewPlex_ExodusII_Nodal_Internal(U, exoid, 1);CHKERRQ(ierr); 374c4762a1bSJed Brown ierr = VecViewPlex_ExodusII_Nodal_Internal(A, exoid, 1);CHKERRQ(ierr); 375c4762a1bSJed Brown /* Saving U and Alpha in one shot. 376c4762a1bSJed Brown For this, we need to cheat and change the Vec's name 377c4762a1bSJed Brown Note that in the end we write variables one component at a time, so that there is no real values in doing this */ 378c4762a1bSJed Brown ierr = DMGetGlobalVector(dmUA, &tmpVec);CHKERRQ(ierr); 379c4762a1bSJed Brown ierr = VecCopy(UA, tmpVec);CHKERRQ(ierr); 380c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) tmpVec, "U");CHKERRQ(ierr); 381c4762a1bSJed Brown ierr = VecViewPlex_ExodusII_Nodal_Internal(tmpVec, exoid, 2);CHKERRQ(ierr); 382c4762a1bSJed Brown /* Reading nodal variables in Exodus file */ 383c4762a1bSJed Brown ierr = VecSet(tmpVec, -1000.0);CHKERRQ(ierr); 384c4762a1bSJed Brown ierr = VecLoadPlex_ExodusII_Nodal_Internal(tmpVec, exoid, 2);CHKERRQ(ierr); 385c4762a1bSJed Brown ierr = VecAXPY(UA, -1.0, tmpVec);CHKERRQ(ierr); 386c4762a1bSJed Brown ierr = VecNorm(UA, NORM_INFINITY, &norm);CHKERRQ(ierr); 387c4762a1bSJed Brown if (norm > PETSC_SQRT_MACHINE_EPSILON) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "UAlpha ||Vin - Vout|| = %g\n", (double) norm); 388c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dmUA, &tmpVec);CHKERRQ(ierr); 389c4762a1bSJed Brown 390c4762a1bSJed Brown /* same thing with the UA2 Vec obtained from the superDM */ 391c4762a1bSJed Brown ierr = DMGetGlobalVector(dmUA2, &tmpVec);CHKERRQ(ierr); 392c4762a1bSJed Brown ierr = VecCopy(UA2, tmpVec);CHKERRQ(ierr); 393c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) tmpVec, "U");CHKERRQ(ierr); 394c4762a1bSJed Brown ierr = VecViewPlex_ExodusII_Nodal_Internal(tmpVec, exoid, 3);CHKERRQ(ierr); 395c4762a1bSJed Brown /* Reading nodal variables in Exodus file */ 396c4762a1bSJed Brown ierr = VecSet(tmpVec, -1000.0);CHKERRQ(ierr); 397c4762a1bSJed Brown ierr = VecLoadPlex_ExodusII_Nodal_Internal(tmpVec, exoid, 3);CHKERRQ(ierr); 398c4762a1bSJed Brown ierr = VecAXPY(UA2, -1.0, tmpVec);CHKERRQ(ierr); 399c4762a1bSJed Brown ierr = VecNorm(UA2, NORM_INFINITY, &norm);CHKERRQ(ierr); 400c4762a1bSJed Brown if (norm > PETSC_SQRT_MACHINE_EPSILON) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "UAlpha2 ||Vin - Vout|| = %g\n", (double) norm); 401c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dmUA2, &tmpVec);CHKERRQ(ierr); 402c4762a1bSJed Brown 403c4762a1bSJed Brown /* Building and saving Sigma 404c4762a1bSJed Brown We set sigma_0 = rank (to see partitioning) 405c4762a1bSJed Brown sigma_1 = cell set ID 406c4762a1bSJed Brown sigma_2 = x_coordinate of the cell center of mass */ 407c4762a1bSJed Brown ierr = DMGetCoordinateSection(dmS, &coordSection);CHKERRQ(ierr); 408c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(dmS, &coord);CHKERRQ(ierr); 409c4762a1bSJed Brown ierr = DMGetLabelIdIS(dmS, "Cell Sets", &csIS);CHKERRQ(ierr); 410c4762a1bSJed Brown ierr = DMGetLabelSize(dmS, "Cell Sets", &numCS); 411c4762a1bSJed Brown ierr = ISGetIndices(csIS, &csID);CHKERRQ(ierr); 412c4762a1bSJed Brown for (set = 0; set < numCS; ++set) { 413c4762a1bSJed 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 */ 414c4762a1bSJed Brown IS cellIS; 415c4762a1bSJed Brown const PetscInt *cellID; 416c4762a1bSJed Brown PetscInt numCells, cell; 417c4762a1bSJed Brown PetscScalar *cval = NULL, *xyz = NULL; 418c4762a1bSJed Brown PetscInt clSize, cdimCoord, c; 419c4762a1bSJed Brown 420c4762a1bSJed Brown ierr = DMGetStratumIS(dmS, "Cell Sets", csID[set], &cellIS);CHKERRQ(ierr); 421c4762a1bSJed Brown ierr = ISGetIndices(cellIS, &cellID);CHKERRQ(ierr); 422c4762a1bSJed Brown ierr = ISGetSize(cellIS, &numCells);CHKERRQ(ierr); 423c4762a1bSJed Brown for (cell = 0; cell < numCells; cell++) { 424c4762a1bSJed Brown ierr = DMPlexVecGetClosure(dmS, NULL, S, cellID[cell], &clSize, &cval);CHKERRQ(ierr); 425c4762a1bSJed Brown ierr = DMPlexVecGetClosure(dmS, coordSection, coord, cellID[cell], &cdimCoord, &xyz);CHKERRQ(ierr); 426c4762a1bSJed Brown cval[0] = rank; 427c4762a1bSJed Brown cval[1] = csID[set]; 428c4762a1bSJed Brown cval[2] = 0.; 429c4762a1bSJed Brown for (c = 0; c < cdimCoord/sdim; c++) cval[2] += xyz[c*sdim]; 430c4762a1bSJed Brown cval[2] = cval[2] * sdim / cdimCoord; 431c4762a1bSJed Brown ierr = DMPlexVecSetClosure(dmS, NULL, S, cellID[cell], cval, INSERT_ALL_VALUES);CHKERRQ(ierr); 432c4762a1bSJed Brown } 433c4762a1bSJed Brown ierr = DMPlexVecRestoreClosure(dmS, NULL, S, cellID[0], &clSize, &cval);CHKERRQ(ierr); 434c4762a1bSJed Brown ierr = DMPlexVecRestoreClosure(dmS, coordSection, coord, cellID[0], NULL, &xyz);CHKERRQ(ierr); 435c4762a1bSJed Brown ierr = ISRestoreIndices(cellIS, &cellID);CHKERRQ(ierr); 436c4762a1bSJed Brown ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 437c4762a1bSJed Brown } 438c4762a1bSJed Brown ierr = ISRestoreIndices(csIS, &csID);CHKERRQ(ierr); 439c4762a1bSJed Brown ierr = ISDestroy(&csIS);CHKERRQ(ierr); 440c4762a1bSJed Brown ierr = VecViewFromOptions(S, NULL, "-s_vec_view");CHKERRQ(ierr); 441c4762a1bSJed Brown /* Writing zonal variables in Exodus file */ 442c4762a1bSJed Brown ierr = VecViewPlex_ExodusII_Zonal_Internal(S, exoid, 1);CHKERRQ(ierr); 443c4762a1bSJed Brown /* Reading zonal variables in Exodus file */ 444c4762a1bSJed Brown ierr = DMGetGlobalVector(dmS, &tmpVec);CHKERRQ(ierr); 445c4762a1bSJed Brown ierr = VecSet(tmpVec, -1000.0);CHKERRQ(ierr); 446c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) tmpVec, "Sigma");CHKERRQ(ierr); 447c4762a1bSJed Brown ierr = VecLoadPlex_ExodusII_Zonal_Internal(tmpVec, exoid, 1);CHKERRQ(ierr); 448c4762a1bSJed Brown ierr = VecAXPY(S, -1.0, tmpVec);CHKERRQ(ierr); 449c4762a1bSJed Brown ierr = VecNorm(S, NORM_INFINITY, &norm); 450c4762a1bSJed Brown if (norm > PETSC_SQRT_MACHINE_EPSILON) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_PLIB, "Sigma ||Vin - Vout|| = %g\n", (double) norm); 451c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dmS, &tmpVec);CHKERRQ(ierr); 452c4762a1bSJed Brown } 453c4762a1bSJed Brown ierr = ex_close(exoid);CHKERRQ(ierr); 454c4762a1bSJed Brown 455c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dmUA2, &UA2);CHKERRQ(ierr); 456c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dmUA, &UA);CHKERRQ(ierr); 457c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dmS, &S);CHKERRQ(ierr); 458c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dmA, &A);CHKERRQ(ierr); 459c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dmU, &U);CHKERRQ(ierr); 460c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dm, &X);CHKERRQ(ierr); 461c4762a1bSJed Brown ierr = DMDestroy(&dmU);CHKERRQ(ierr); ierr = ISDestroy(&isU);CHKERRQ(ierr); 462c4762a1bSJed Brown ierr = DMDestroy(&dmA);CHKERRQ(ierr); ierr = ISDestroy(&isA);CHKERRQ(ierr); 463c4762a1bSJed Brown ierr = DMDestroy(&dmS);CHKERRQ(ierr); ierr = ISDestroy(&isS);CHKERRQ(ierr); 464c4762a1bSJed Brown ierr = DMDestroy(&dmUA);CHKERRQ(ierr);ierr = ISDestroy(&isUA);CHKERRQ(ierr); 465c4762a1bSJed Brown ierr = DMDestroy(&dmUA2);CHKERRQ(ierr); 466c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 467c4762a1bSJed Brown ierr = PetscFree2(pStartDepth, pEndDepth);CHKERRQ(ierr); 468c4762a1bSJed Brown ierr = PetscFinalize(); 469c4762a1bSJed Brown return ierr; 470c4762a1bSJed Brown } 471c4762a1bSJed Brown 472c4762a1bSJed Brown /*TEST 473c4762a1bSJed Brown 474c4762a1bSJed Brown build: 475c4762a1bSJed Brown requires: exodusii pnetcdf !complex 476c4762a1bSJed Brown # 2D seq 477c4762a1bSJed Brown test: 478c4762a1bSJed Brown suffix: 0 479c4762a1bSJed 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 480c4762a1bSJed 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 481c4762a1bSJed Brown test: 482c4762a1bSJed Brown suffix: 1 483c4762a1bSJed 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 484c4762a1bSJed Brown test: 485c4762a1bSJed Brown suffix: 2 486c4762a1bSJed 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 487c4762a1bSJed 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 488c4762a1bSJed Brown test: 489c4762a1bSJed Brown suffix: 3 490c4762a1bSJed 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 491c4762a1bSJed Brown test: 492c4762a1bSJed Brown suffix: 4 493c4762a1bSJed 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 494c4762a1bSJed Brown test: 495c4762a1bSJed Brown suffix: 5 496c4762a1bSJed 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 497c4762a1bSJed Brown 498c4762a1bSJed Brown # 2D par 499c4762a1bSJed Brown test: 500c4762a1bSJed Brown suffix: 6 501c4762a1bSJed Brown nsize: 2 502c4762a1bSJed 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 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: 7 506c4762a1bSJed Brown nsize: 2 507c4762a1bSJed 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 508c4762a1bSJed Brown test: 509c4762a1bSJed Brown suffix: 8 510c4762a1bSJed Brown nsize: 2 511c4762a1bSJed 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 512c4762a1bSJed Brown #TODO: bug in call to NetCDF failed to complete invalid type definition in file id 65536 NetCDF: invalid dimension ID or name 513c4762a1bSJed Brown test: 514c4762a1bSJed Brown suffix: 9 515c4762a1bSJed Brown nsize: 2 516c4762a1bSJed 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 517c4762a1bSJed Brown test: 518c4762a1bSJed Brown suffix: 10 519c4762a1bSJed Brown nsize: 2 520c4762a1bSJed 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 521c4762a1bSJed Brown test: 522*54fcfd0cSMatthew G. Knepley # Something is now broken with parallel read/write for wahtever shape H is 523*54fcfd0cSMatthew G. Knepley TODO: broken 524c4762a1bSJed Brown suffix: 11 525c4762a1bSJed Brown nsize: 2 526c4762a1bSJed 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 527c4762a1bSJed Brown 528c4762a1bSJed Brown #3d seq 529c4762a1bSJed Brown test: 530c4762a1bSJed Brown suffix: 12 531c4762a1bSJed 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 532c4762a1bSJed Brown test: 533c4762a1bSJed Brown suffix: 13 534c4762a1bSJed 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 535c4762a1bSJed 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 536c4762a1bSJed Brown test: 537c4762a1bSJed Brown suffix: 14 538c4762a1bSJed 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 539c4762a1bSJed Brown test: 540c4762a1bSJed Brown suffix: 15 541c4762a1bSJed 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 542c4762a1bSJed 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 543c4762a1bSJed Brown #3d par 544c4762a1bSJed Brown test: 545c4762a1bSJed Brown suffix: 16 546c4762a1bSJed Brown nsize: 2 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: 17 550c4762a1bSJed Brown nsize: 2 551c4762a1bSJed 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 552c4762a1bSJed 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 553c4762a1bSJed Brown test: 554c4762a1bSJed Brown suffix: 18 555c4762a1bSJed Brown nsize: 2 556c4762a1bSJed 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 557c4762a1bSJed Brown test: 558c4762a1bSJed Brown suffix: 19 559c4762a1bSJed Brown nsize: 2 560c4762a1bSJed 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 561c4762a1bSJed 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 562c4762a1bSJed Brown 563c4762a1bSJed Brown TEST*/ 564