Actual source code: ex10.c

  1: /*
  2:    Demonstrates using the HDF5 viewer with a DMDA Vec
  3:  - create a global vector containing a gauss profile (exp(-x^2-y^2))
  4:  - write the global vector in a hdf5 file

  6:    The resulting file gauss.h5 can be viewed with Visit (an open source visualization package)
  7:    Or with some versions of MATLAB with data=hdfread('gauss.h5','pressure'); mesh(data);

  9:    The file storage of the vector is independent of the number of processes used.
 10:  */

 12: #include <petscdm.h>
 13: #include <petscdmda.h>
 14: #include <petscsys.h>
 15: #include <petscviewerhdf5.h>

 17: static char help[] = "Test to write HDF5 file from PETSc DMDA Vec.\n\n";

 19: int main(int argc,char **argv)
 20: {
 21:   DM             da2D;
 22:   PetscInt       i,j,ixs, ixm, iys, iym;
 23:   PetscViewer    H5viewer;
 24:   PetscScalar    xm    = -1.0, xp=1.0;
 25:   PetscScalar    ym    = -1.0, yp=1.0;
 26:   PetscScalar    value = 1.0,dx,dy;
 27:   PetscInt       Nx    = 40, Ny=40;
 28:   Vec            gauss,input;
 29:   PetscScalar    **gauss_ptr;
 30:   PetscReal      norm;
 31:   const char     *vecname;

 33:   dx=(xp-xm)/(Nx-1);
 34:   dy=(yp-ym)/(Ny-1);

 36:   /* Initialize the Petsc context */
 37:   PetscInitialize(&argc,&argv,(char*)0,help);
 38:   DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,Nx,Ny,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da2D);
 39:   DMSetFromOptions(da2D);
 40:   DMSetUp(da2D);

 42:   /* Set the coordinates */
 43:   DMDASetUniformCoordinates(da2D, 0.0, 1.0, 0.0, 1.0, 0.0, 0.0);

 45:   /* Declare gauss as a DMDA component */
 46:   DMCreateGlobalVector(da2D,&gauss);
 47:   PetscObjectSetName((PetscObject) gauss, "pressure");

 49:   /* Initialize vector gauss with a constant value (=1) */
 50:   VecSet(gauss,value);

 52:   /* Get the coordinates of the corners for each process */
 53:   DMDAGetCorners(da2D, &ixs, &iys, 0, &ixm, &iym, 0);

 55:   /* Build the gaussian profile (exp(-x^2-y^2)) */
 56:   DMDAVecGetArray(da2D,gauss,&gauss_ptr);
 57:   for (j=iys; j<iys+iym; j++) {
 58:     for (i=ixs; i<ixs+ixm; i++) {
 59:       gauss_ptr[j][i]=PetscExpScalar(-(xm+i*dx)*(xm+i*dx)-(ym+j*dy)*(ym+j*dy));
 60:     }
 61:   }
 62:   DMDAVecRestoreArray(da2D,gauss,&gauss_ptr);

 64:   /* Create the HDF5 viewer */
 65:   PetscViewerHDF5Open(PETSC_COMM_WORLD,"gauss.h5",FILE_MODE_WRITE,&H5viewer);
 66:   PetscViewerSetFromOptions(H5viewer);

 68:   /* Write the H5 file */
 69:   VecView(gauss,H5viewer);

 71:   /* Close the viewer */
 72:   PetscViewerDestroy(&H5viewer);

 74:   VecDuplicate(gauss,&input);
 75:   PetscObjectGetName((PetscObject)gauss,&vecname);
 76:   PetscObjectSetName((PetscObject)input,vecname);

 78:   /* Create the HDF5 viewer for reading */
 79:   PetscViewerHDF5Open(PETSC_COMM_WORLD,"gauss.h5",FILE_MODE_READ,&H5viewer);
 80:   PetscViewerSetFromOptions(H5viewer);
 81:   VecLoad(input,H5viewer);
 82:   PetscViewerDestroy(&H5viewer);

 84:   VecAXPY(input,-1.0,gauss);
 85:   VecNorm(input,NORM_2,&norm);

 88:   VecDestroy(&input);
 89:   VecDestroy(&gauss);
 90:   DMDestroy(&da2D);
 91:   PetscFinalize();
 92:   return 0;
 93: }

 95: /*TEST

 97:       build:
 98:          requires: hdf5 !defined(PETSC_USE_CXXCOMPLEX)

100:       test:
101:          nsize: 4

103:       test:
104:          nsize: 4
105:          suffix: 2
106:          args: -viewer_hdf5_base_dimension2
107:          output_file: output/ex10_1.out

109:       test:
110:          nsize: 4
111:          suffix: 3
112:          args: -viewer_hdf5_sp_output
113:          output_file: output/ex10_1.out

115: TEST*/