Actual source code: ex9.c

petsc-3.3-p7 2013-05-11
  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 <petscdmda.h>


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

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

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

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

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

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

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