Actual source code: ex10.c
petsc-3.6.4 2016-04-12
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: {
22: DM da2D;
23: PetscInt i,j,ixs, ixm, iys, iym;;
24: PetscViewer H5viewer;
25: PetscScalar xm = -1.0, xp=1.0;
26: PetscScalar ym = -1.0, yp=1.0;
27: PetscScalar value = 1.0,dx,dy;
28: PetscInt Nx = 40, Ny=40;
29: Vec gauss,input;
30: PetscScalar **gauss_ptr;
31: PetscReal norm;
32: const char *vecname;
34: dx=(xp-xm)/(Nx-1);
35: dy=(yp-ym)/(Ny-1);
37: /* Initialize the Petsc context */
38: PetscInitialize(&argc,&argv,(char*)0,help);
40: /* Build of the DMDA */
41: DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,Nx,Ny,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da2D);
43: /* Set the coordinates */
44: DMDASetUniformCoordinates(da2D, 0.0, 1.0, 0.0, 1.0, 0.0, 0.0);
46: /* Declare gauss as a DMDA component */
47: DMCreateGlobalVector(da2D,&gauss);
48: PetscObjectSetName((PetscObject) gauss, "pressure");
50: /* Initialize vector gauss with a constant value (=1) */
51: VecSet(gauss,value);
53: /* Get the coordinates of the corners for each process */
54: DMDAGetCorners(da2D, &ixs, &iys, 0, &ixm, &iym, 0);
56: /* Build the gaussian profile (exp(-x^2-y^2)) */
57: DMDAVecGetArray(da2D,gauss,&gauss_ptr);
58: for (j=iys; j<iys+iym; j++) {
59: for (i=ixs; i<ixs+ixm; i++) {
60: gauss_ptr[j][i]=PetscExpScalar(-(xm+i*dx)*(xm+i*dx)-(ym+j*dy)*(ym+j*dy));
61: }
62: }
63: DMDAVecRestoreArray(da2D,gauss,&gauss_ptr);
65: /* Create the HDF5 viewer */
66: PetscViewerHDF5Open(PETSC_COMM_WORLD,"gauss.h5",FILE_MODE_WRITE,&H5viewer);
67: PetscViewerSetFromOptions(H5viewer);
69: /* Write the H5 file */
70: VecView(gauss,H5viewer);
72: /* Close the viewer */
73: PetscViewerDestroy(&H5viewer);
75: VecDuplicate(gauss,&input);
76: PetscObjectGetName((PetscObject)gauss,&vecname);
77: PetscObjectSetName((PetscObject)input,vecname);
79: /* Create the HDF5 viewer for reading */
80: PetscViewerHDF5Open(PETSC_COMM_WORLD,"gauss.h5",FILE_MODE_READ,&H5viewer);
81: PetscViewerSetFromOptions(H5viewer);
82: VecLoad(input,H5viewer);
83: PetscViewerDestroy(&H5viewer);
85: VecAXPY(input,-1.0,gauss);
86: VecNorm(input,NORM_2,&norm);
87: if (norm > 1.e-6) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Vec read in does not match vector written out");
89: VecDestroy(&input);
90: VecDestroy(&gauss);
91: DMDestroy(&da2D);
92: PetscFinalize();
93: return 0;
94: }