Actual source code: ex61view.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  1: static const char help[] = "Loads data generated by ex61 and VTK file suitable for Paraview or Visit.\n\n";

  3: /*

  5:       ./ex61gen  -random_seed <integer>                            :   creates ex61.random.<integer>
  6:       ./ex61 [-random_seed <integer>] -graphicsfile other options  :   creates ex61.data[.<integer>]
  7:       ./ex61view [-random_seed <integer>]                          :   creates ex61.data[.<integer>]_cv_%d.vtk  ex61.data[.<integer>]_eta_%d.vtk

  9:     where [ ] indicates the optional seed value;  <integer> the seed value provided and %d the timestep

 11:     The resulting .vtk files may be loaded into ParaView or ViSIT

 13: */

 15: #include <petscdmda.h>


 20: int main(int argc, char **argv)
 21: {
 23:   Vec            U,cv,eta;
 24:   DM             da,da2;
 25:   PetscViewer    viewer,view_vtk_cv,view_vtk_eta;
 26:   char           filename[PETSC_MAX_PATH_LEN],cv_filename[PETSC_MAX_PATH_LEN],eta_filename[PETSC_MAX_PATH_LEN];
 27:   PetscBool      flg,sflg = PETSC_FALSE;
 28:   PetscInt       i,n=10000;
 29:   PetscInt       seed;

 31:   PetscInitialize(&argc,&argv, (char*)0, help);
 32:   PetscOptionsSetValue("-viewer_binary_skip_info","true");
 33:   PetscOptionsGetString(NULL,"-f",filename,PETSC_MAX_PATH_LEN,&flg);
 34:   if (!flg) {
 35:     PetscOptionsGetInt(NULL,"-random_seed",&seed,&sflg);
 36:     if (!sflg) {
 37:       PetscStrcpy(filename,"ex61.data");
 38:     } else {
 39:       sprintf(filename,"ex61.data.%d",seed);
 40:     }
 41:   }

 43:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);

 45:   /* Get physics and time parameters */
 46:   DMCreate(PETSC_COMM_WORLD,&da);
 47:   DMLoad(da,viewer);
 48:   DMCreateGlobalVector(da,&U);
 49:   DMDAGetReducedDMDA(da,1,&da2);
 50:   DMCreateGlobalVector(da2,&cv);
 51:   DMCreateGlobalVector(da2,&eta);

 53:   for (i=0; i<n; i++) {
 54:     /* when this fails it simply means the file is finished */
 55:     VecLoad(U,viewer);
 56:     VecStrideGather(U,1,cv,INSERT_VALUES);
 57:     VecStrideGather(U,4,eta,INSERT_VALUES);
 58:     sprintf(cv_filename,"%s_cv_%d.vtk",filename,i);
 59:     sprintf(eta_filename,"%s_eta_%d.vtk",filename,i);
 60:     PetscViewerASCIIOpen(PETSC_COMM_WORLD,cv_filename,&view_vtk_cv);
 61:     PetscViewerASCIIOpen(PETSC_COMM_WORLD,eta_filename,&view_vtk_eta);
 62:     PetscViewerSetFormat(view_vtk_cv, PETSC_VIEWER_ASCII_VTK);
 63:     PetscViewerSetFormat(view_vtk_eta, PETSC_VIEWER_ASCII_VTK);
 64:     DMView(da2,view_vtk_cv);
 65:     DMView(da2,view_vtk_eta);
 66:     VecView(cv,view_vtk_cv);
 67:     VecView(eta,view_vtk_eta);
 68:     PetscViewerDestroy(&view_vtk_cv);
 69:     PetscViewerDestroy(&view_vtk_eta);
 70:   }
 71:   VecDestroy(&U);
 72:   VecDestroy(&cv);
 73:   VecDestroy(&eta);
 74:   DMDestroy(&da);
 75:   DMDestroy(&da2);
 76:   PetscFinalize();
 77:   return 0;
 78: }