Actual source code: ex47.c


  2: static char help[] = "Tests PetscViewerHDF5 VecView()/VecLoad() function.\n\n";

  4: #include <petscviewer.h>
  5: #include <petscviewerhdf5.h>
  6: #include <petscvec.h>

  8: int main(int argc,char **args)
  9: {
 10:   Vec            x,y;
 11:   PetscReal      norm,dnorm;
 12:   PetscViewer    H5viewer;
 13:   char           filename[PETSC_MAX_PATH_LEN];
 14:   PetscBool      flg;

 16:   PetscInitialize(&argc,&args,(char*)0,help);
 17:   PetscOptionsGetString(NULL,NULL,"-filename",filename,sizeof(filename),&flg);
 18:   if (!flg) PetscStrcpy(filename,"x.h5");
 19:   VecCreate(PETSC_COMM_WORLD,&x);
 20:   VecSetFromOptions(x);
 21:   VecSetSizes(x,11,PETSC_DETERMINE);
 22:   VecSet(x,22.3);

 24:   PetscViewerHDF5Open(PETSC_COMM_WORLD,filename,FILE_MODE_WRITE,&H5viewer);
 25:   PetscViewerSetFromOptions(H5viewer);

 27:   /* Write the Vec without one extra dimension for BS */
 28:   PetscViewerHDF5SetBaseDimension2(H5viewer, PETSC_FALSE);
 29:   PetscObjectSetName((PetscObject) x, "noBsDim");
 30:   VecView(x,H5viewer);

 32:   /* Write the Vec with one extra, 1-sized, dimension for BS */
 33:   PetscViewerHDF5SetBaseDimension2(H5viewer, PETSC_TRUE);
 34:   PetscObjectSetName((PetscObject) x, "bsDim");
 35:   VecView(x,H5viewer);

 37:   PetscViewerDestroy(&H5viewer);
 38:   VecDuplicate(x,&y);

 40:   /* Create the HDF5 viewer for reading */
 41:   PetscViewerHDF5Open(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&H5viewer);
 42:   PetscViewerSetFromOptions(H5viewer);

 44:   /* Load the Vec without the BS dim and compare */
 45:   PetscObjectSetName((PetscObject) y, "noBsDim");
 46:   VecLoad(y,H5viewer);
 47:   VecAXPY(y,-1.0,x);
 48:   VecNorm(y,NORM_2,&norm);
 49:   VecNorm(x,NORM_2,&dnorm);

 52:   /* Load the Vec with one extra, 1-sized, BS dim and compare */
 53:   PetscObjectSetName((PetscObject) y, "bsDim");
 54:   VecLoad(y,H5viewer);
 55:   VecAXPY(y,-1.0,x);
 56:   VecNorm(y,NORM_2,&norm);
 57:   VecNorm(x,NORM_2,&dnorm);

 60:   PetscViewerDestroy(&H5viewer);
 61:   VecDestroy(&y);
 62:   VecDestroy(&x);
 63:   PetscFinalize();
 64:   return 0;
 65: }

 67: /*TEST

 69:    build:
 70:      requires: hdf5

 72:    test:
 73:      requires: hdf5

 75:    test:
 76:      suffix: 2
 77:      nsize: 4

 79:    test:
 80:      suffix: 3
 81:      nsize: 4
 82:      args: -viewer_hdf5_base_dimension2

 84:    test:
 85:      suffix: 4
 86:      nsize: 4
 87:      args: -viewer_hdf5_sp_output

 89: TEST*/