1*49c89c76SBlaise Bourdin static char help[] = "Test PetscViewer_ExodusII\n\n"; 2*49c89c76SBlaise Bourdin 3*49c89c76SBlaise Bourdin #include <petsc.h> 4*49c89c76SBlaise Bourdin #include <exodusII.h> 5*49c89c76SBlaise Bourdin 6*49c89c76SBlaise Bourdin int main(int argc, char **argv) 7*49c89c76SBlaise Bourdin { 8*49c89c76SBlaise Bourdin DM dm, pdm; 9*49c89c76SBlaise Bourdin PetscInt ovlp = 0; 10*49c89c76SBlaise Bourdin char ifilename[PETSC_MAX_PATH_LEN], ofilename[PETSC_MAX_PATH_LEN]; 11*49c89c76SBlaise Bourdin int numZVars, numNVars; 12*49c89c76SBlaise Bourdin int nNodalVar = 4; 13*49c89c76SBlaise Bourdin int nZonalVar = 3; 14*49c89c76SBlaise Bourdin int order = 1; 15*49c89c76SBlaise Bourdin PetscViewer viewer; 16*49c89c76SBlaise Bourdin int index = -1; 17*49c89c76SBlaise Bourdin const char *nodalVarName[4] = {"U_x", "U_y", "Alpha", "Beta"}; 18*49c89c76SBlaise Bourdin const char *zonalVarName[3] = {"Sigma_11", "Sigma_12", "Sigma_22"}; 19*49c89c76SBlaise Bourdin const char *testNames[3] = {"U", "Sigma", "Gamma"}; 20*49c89c76SBlaise Bourdin const char *varName = NULL; 21*49c89c76SBlaise Bourdin 22*49c89c76SBlaise Bourdin PetscFunctionBeginUser; 23*49c89c76SBlaise Bourdin PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 24*49c89c76SBlaise Bourdin PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "PetscViewer_ExodusII test", "ex95"); 25*49c89c76SBlaise Bourdin PetscCall(PetscOptionsString("-i", "Filename to read", "ex95", ifilename, ifilename, sizeof(ifilename), NULL)); 26*49c89c76SBlaise Bourdin PetscCall(PetscOptionsString("-o", "Filename to write", "ex95", ofilename, ofilename, sizeof(ofilename), NULL)); 27*49c89c76SBlaise Bourdin PetscOptionsEnd(); 28*49c89c76SBlaise Bourdin 29*49c89c76SBlaise Bourdin #ifdef PETSC_USE_DEBUG 30*49c89c76SBlaise Bourdin PetscCallExternal(ex_opts, EX_VERBOSE + EX_DEBUG); 31*49c89c76SBlaise Bourdin #endif 32*49c89c76SBlaise Bourdin 33*49c89c76SBlaise Bourdin PetscCall(DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename, NULL, PETSC_TRUE, &dm)); 34*49c89c76SBlaise Bourdin PetscCall(DMPlexDistributeSetDefault(dm, PETSC_FALSE)); 35*49c89c76SBlaise Bourdin PetscCall(DMSetFromOptions(dm)); 36*49c89c76SBlaise Bourdin PetscCall(PetscObjectSetName((PetscObject)dm, "ex95")); 37*49c89c76SBlaise Bourdin PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 38*49c89c76SBlaise Bourdin 39*49c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIOpen(PETSC_COMM_WORLD, ofilename, FILE_MODE_WRITE, &viewer)); 40*49c89c76SBlaise Bourdin 41*49c89c76SBlaise Bourdin /* Save the geometry to the file, erasing all previous content */ 42*49c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIISetOrder(viewer, order)); 43*49c89c76SBlaise Bourdin PetscCall(DMView(dm, viewer)); 44*49c89c76SBlaise Bourdin PetscCall(PetscViewerView(viewer, PETSC_VIEWER_STDOUT_WORLD)); 45*49c89c76SBlaise Bourdin PetscCall(PetscViewerFlush(viewer)); 46*49c89c76SBlaise Bourdin 47*49c89c76SBlaise Bourdin PetscCall(DMPlexDistribute(dm, ovlp, NULL, &pdm)); 48*49c89c76SBlaise Bourdin if (!pdm) pdm = dm; 49*49c89c76SBlaise Bourdin 50*49c89c76SBlaise Bourdin /* Testing Variable Number*/ 51*49c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIISetZonalVariable(viewer, nZonalVar)); 52*49c89c76SBlaise Bourdin nZonalVar = -1; 53*49c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetZonalVariable(viewer, &nZonalVar)); 54*49c89c76SBlaise Bourdin PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of zonal variables: %d\n", nZonalVar)); 55*49c89c76SBlaise Bourdin 56*49c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIISetNodalVariable(viewer, nNodalVar)); 57*49c89c76SBlaise Bourdin nNodalVar = -1; 58*49c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetNodalVariable(viewer, &nNodalVar)); 59*49c89c76SBlaise Bourdin PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of nodal variables: %d\n", nNodalVar)); 60*49c89c76SBlaise Bourdin PetscCall(PetscViewerView(viewer, PETSC_VIEWER_STDOUT_WORLD)); 61*49c89c76SBlaise Bourdin 62*49c89c76SBlaise Bourdin /* 63*49c89c76SBlaise Bourdin Test of PetscViewerExodusIISet[Nodal/Zonal]VariableName 64*49c89c76SBlaise Bourdin */ 65*49c89c76SBlaise Bourdin PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Testing PetscViewerExodusIISet[Nodal/Zonal]VariableName\n")); 66*49c89c76SBlaise Bourdin for (int i = 0; i < nNodalVar; i++) { PetscCall(PetscViewerExodusIISetNodalVariableName(viewer, i, nodalVarName[i])); } 67*49c89c76SBlaise Bourdin for (int i = 0; i < nZonalVar; i++) { PetscCall(PetscViewerExodusIISetZonalVariableName(viewer, i, zonalVarName[i])); } 68*49c89c76SBlaise Bourdin PetscCall(PetscViewerView(viewer, PETSC_VIEWER_STDOUT_WORLD)); 69*49c89c76SBlaise Bourdin PetscCall(PetscViewerDestroy(&viewer)); 70*49c89c76SBlaise Bourdin 71*49c89c76SBlaise Bourdin /* 72*49c89c76SBlaise Bourdin Test of PetscViewerExodusIIGet[Nodal/Zonal]VariableName 73*49c89c76SBlaise Bourdin */ 74*49c89c76SBlaise Bourdin PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n\nReopenning the output file in Read-only mode\n")); 75*49c89c76SBlaise Bourdin PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Testing PetscViewerExodusIIGet[Nodal/Zonal]VariableName\n")); 76*49c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIOpen(PETSC_COMM_WORLD, ofilename, FILE_MODE_APPEND, &viewer)); 77*49c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIISetOrder(viewer, order)); 78*49c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetZonalVariable(viewer, &numZVars)); 79*49c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetNodalVariable(viewer, &numNVars)); 80*49c89c76SBlaise Bourdin 81*49c89c76SBlaise Bourdin for (int i = 0; i < numNVars; i++) { 82*49c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetNodalVariableName(viewer, i, &varName)); 83*49c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetNodalVariableIndex(viewer, varName, &index)); 84*49c89c76SBlaise Bourdin PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Nodal variable %d: %s, index in file %d\n", i, varName, index)); 85*49c89c76SBlaise Bourdin } 86*49c89c76SBlaise Bourdin for (int i = 0; i < 3; i++) { 87*49c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetNodalVariableIndex(viewer, testNames[i], &index)); 88*49c89c76SBlaise Bourdin PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Nodal variable %d: %s, index in file %d\n", i, testNames[i], index)); 89*49c89c76SBlaise Bourdin } 90*49c89c76SBlaise Bourdin PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n")); 91*49c89c76SBlaise Bourdin 92*49c89c76SBlaise Bourdin for (int i = 0; i < numZVars; i++) { 93*49c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetZonalVariableName(viewer, i, &varName)); 94*49c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetZonalVariableIndex(viewer, varName, &index)); 95*49c89c76SBlaise Bourdin PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Zonal variable %d: %s, index in file %d\n", i, varName, index)); 96*49c89c76SBlaise Bourdin } 97*49c89c76SBlaise Bourdin for (int i = 0; i < 3; i++) { 98*49c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetZonalVariableIndex(viewer, testNames[i], &index)); 99*49c89c76SBlaise Bourdin PetscCall(PetscPrintf(PETSC_COMM_WORLD, " Zonal variable %d: %s, index in file %d\n", i, testNames[i], index)); 100*49c89c76SBlaise Bourdin } 101*49c89c76SBlaise Bourdin 102*49c89c76SBlaise Bourdin PetscCall(PetscViewerDestroy(&viewer)); 103*49c89c76SBlaise Bourdin PetscCall(DMDestroy(&dm)); 104*49c89c76SBlaise Bourdin PetscCall(PetscFinalize()); 105*49c89c76SBlaise Bourdin return 0; 106*49c89c76SBlaise Bourdin } 107*49c89c76SBlaise Bourdin 108*49c89c76SBlaise Bourdin /*TEST 109*49c89c76SBlaise Bourdin 110*49c89c76SBlaise Bourdin build: 111*49c89c76SBlaise Bourdin requires: exodusii pnetcdf !complex 112*49c89c76SBlaise Bourdin test: 113*49c89c76SBlaise Bourdin suffix: 0 114*49c89c76SBlaise Bourdin nsize: 1 115*49c89c76SBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -o FourSquareT-large_out.exo 116*49c89c76SBlaise Bourdin 117*49c89c76SBlaise Bourdin TEST*/ 118