Actual source code: ex9.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  1: static char help[] = "Demonstrates HD5 vector input/ouput\n\n";

  3: /*T
  4:    Concepts: viewers
  5:    Concepts: HDF5
  6:    Processors: n
  7: T*/
  8: #include <petscsys.h>
  9: #include <petscdm.h>
 10: #include <petscdmda.h>
 11: #include <petscviewerhdf5.h>

 15: int main(int argc,char **argv)
 16: {
 18:   PetscViewer    viewer;
 19:   DM             da;
 20:   Vec            global,local,global2;
 21:   PetscMPIInt    rank;
 22:   PetscBool      flg;

 24:   /*
 25:     Every PETSc routine should begin with the PetscInitialize() routine.
 26:     argc, argv - These command line arguments are taken to extract the options
 27:                  supplied to PETSc and options supplied to MPI.
 28:     help       - When PETSc executable is invoked with the option -help,
 29:                  it prints the various options that can be applied at
 30:                  runtime.  The user can use the "help" variable place
 31:                  additional help messages in this printout.
 32:   */
 33:   PetscInitialize(&argc,&argv,(char*)0,help);

 35:   /* Create a DMDA and an associated vector */
 36:   DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,100,90,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,NULL,&da);
 37:   DMCreateGlobalVector(da,&global);
 38:   DMCreateLocalVector(da,&local);
 39:   VecSet(global,-1.0);
 40:   DMGlobalToLocalBegin(da,global,INSERT_VALUES,local);
 41:   DMGlobalToLocalEnd(da,global,INSERT_VALUES,local);
 42:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 43:   VecScale(local,rank+1);
 44:   DMLocalToGlobalBegin(da,local,ADD_VALUES,global);
 45:   DMLocalToGlobalEnd(da,local,ADD_VALUES,global);

 47:   /*
 48:      Write output file with PetscViewerHDF5 viewer.

 50:   */
 51:   PetscViewerHDF5Open(PETSC_COMM_WORLD,"hdf5output",FILE_MODE_WRITE,&viewer);
 52:   VecView(global,viewer);
 53:   PetscViewerDestroy(&viewer);
 54:   MPI_Barrier(PETSC_COMM_WORLD);

 56:   VecDuplicate(global,&global2);
 57:   VecCopy(global,global2);
 58:   PetscViewerHDF5Open(PETSC_COMM_WORLD,"hdf5output",FILE_MODE_READ,&viewer);
 59:   VecLoad(global,viewer);
 60:   PetscViewerDestroy(&viewer);
 61:   VecEqual(global,global2,&flg);
 62:   if (flg) {
 63:     PetscPrintf(PETSC_COMM_WORLD,"Vectors are equal\n");
 64:   } else {
 65:     PetscPrintf(PETSC_COMM_WORLD,"Vectors are not equal\n");
 66:   }

 68:   /* clean up and exit */
 69:   DMDestroy(&da);
 70:   VecDestroy(&local);
 71:   VecDestroy(&global);
 72:   VecDestroy(&global2);
 73:   PetscFinalize();
 74:   return 0;
 75: }