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