1c4762a1bSJed Brown /*
2c4762a1bSJed Brown Demonstrates using the HDF5 viewer with a DMDA Vec
3c4762a1bSJed Brown - create a global vector containing a gauss profile (exp(-x^2-y^2))
4c4762a1bSJed Brown - write the global vector in a hdf5 file
5c4762a1bSJed Brown
6337bb527SBarry Smith The resulting file gauss.h5 can be viewed with Visit (an open-source visualization package)
7c4762a1bSJed Brown Or with some versions of MATLAB with data=hdfread('gauss.h5','pressure'); mesh(data);
8c4762a1bSJed Brown
9c4762a1bSJed Brown The file storage of the vector is independent of the number of processes used.
10c4762a1bSJed Brown */
11c4762a1bSJed Brown
12c4762a1bSJed Brown #include <petscdm.h>
13c4762a1bSJed Brown #include <petscdmda.h>
14c4762a1bSJed Brown #include <petscsys.h>
15c4762a1bSJed Brown #include <petscviewerhdf5.h>
16c4762a1bSJed Brown
17c4762a1bSJed Brown static char help[] = "Test to write HDF5 file from PETSc DMDA Vec.\n\n";
18c4762a1bSJed Brown
main(int argc,char ** argv)19d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
20d71ae5a4SJacob Faibussowitsch {
21c4762a1bSJed Brown DM da2D;
22c4762a1bSJed Brown PetscInt i, j, ixs, ixm, iys, iym;
23c4762a1bSJed Brown PetscViewer H5viewer;
24c4762a1bSJed Brown PetscScalar xm = -1.0, xp = 1.0;
25c4762a1bSJed Brown PetscScalar ym = -1.0, yp = 1.0;
26c4762a1bSJed Brown PetscScalar value = 1.0, dx, dy;
27c4762a1bSJed Brown PetscInt Nx = 40, Ny = 40;
28c4762a1bSJed Brown Vec gauss, input;
29c4762a1bSJed Brown PetscScalar **gauss_ptr;
30c4762a1bSJed Brown PetscReal norm;
31c4762a1bSJed Brown const char *vecname;
32c4762a1bSJed Brown
33c4762a1bSJed Brown dx = (xp - xm) / (Nx - 1);
34c4762a1bSJed Brown dy = (yp - ym) / (Ny - 1);
35c4762a1bSJed Brown
36f0b74427SPierre Jolivet /* Initialize the PETSc context */
37327415f7SBarry Smith PetscFunctionBeginUser;
38c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
399566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, Nx, Ny, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &da2D));
409566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da2D));
419566063dSJacob Faibussowitsch PetscCall(DMSetUp(da2D));
42c4762a1bSJed Brown
43c4762a1bSJed Brown /* Set the coordinates */
449566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da2D, 0.0, 1.0, 0.0, 1.0, 0.0, 0.0));
45c4762a1bSJed Brown
46c4762a1bSJed Brown /* Declare gauss as a DMDA component */
479566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da2D, &gauss));
489566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)gauss, "pressure"));
49c4762a1bSJed Brown
50c4762a1bSJed Brown /* Initialize vector gauss with a constant value (=1) */
519566063dSJacob Faibussowitsch PetscCall(VecSet(gauss, value));
52c4762a1bSJed Brown
53c4762a1bSJed Brown /* Get the coordinates of the corners for each process */
549566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da2D, &ixs, &iys, 0, &ixm, &iym, 0));
55c4762a1bSJed Brown
56c4762a1bSJed Brown /* Build the gaussian profile (exp(-x^2-y^2)) */
579566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da2D, gauss, &gauss_ptr));
58c4762a1bSJed Brown for (j = iys; j < iys + iym; j++) {
59ad540459SPierre Jolivet for (i = ixs; i < ixs + ixm; i++) gauss_ptr[j][i] = PetscExpScalar(-(xm + i * dx) * (xm + i * dx) - (ym + j * dy) * (ym + j * dy));
60c4762a1bSJed Brown }
619566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da2D, gauss, &gauss_ptr));
62c4762a1bSJed Brown
63c4762a1bSJed Brown /* Create the HDF5 viewer */
649566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD, "gauss.h5", FILE_MODE_WRITE, &H5viewer));
659566063dSJacob Faibussowitsch PetscCall(PetscViewerSetFromOptions(H5viewer));
66c4762a1bSJed Brown
67c4762a1bSJed Brown /* Write the H5 file */
689566063dSJacob Faibussowitsch PetscCall(VecView(gauss, H5viewer));
69c4762a1bSJed Brown
70c4762a1bSJed Brown /* Close the viewer */
719566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&H5viewer));
72c4762a1bSJed Brown
739566063dSJacob Faibussowitsch PetscCall(VecDuplicate(gauss, &input));
749566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)gauss, &vecname));
759566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)input, vecname));
76c4762a1bSJed Brown
77c4762a1bSJed Brown /* Create the HDF5 viewer for reading */
789566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Open(PETSC_COMM_WORLD, "gauss.h5", FILE_MODE_READ, &H5viewer));
799566063dSJacob Faibussowitsch PetscCall(PetscViewerSetFromOptions(H5viewer));
809566063dSJacob Faibussowitsch PetscCall(VecLoad(input, H5viewer));
819566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&H5viewer));
82c4762a1bSJed Brown
839566063dSJacob Faibussowitsch PetscCall(VecAXPY(input, -1.0, gauss));
849566063dSJacob Faibussowitsch PetscCall(VecNorm(input, NORM_2, &norm));
8508401ef6SPierre Jolivet PetscCheck(norm <= 1.e-6, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Vec read in does not match vector written out");
86c4762a1bSJed Brown
879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&input));
889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&gauss));
899566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da2D));
909566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
91b122ec5aSJacob Faibussowitsch return 0;
92c4762a1bSJed Brown }
93c4762a1bSJed Brown
94c4762a1bSJed Brown /*TEST
95c4762a1bSJed Brown
96c4762a1bSJed Brown build:
97dfd57a17SPierre Jolivet requires: hdf5 !defined(PETSC_USE_CXXCOMPLEX)
98c4762a1bSJed Brown
99c4762a1bSJed Brown test:
100c4762a1bSJed Brown nsize: 4
101*3886731fSPierre Jolivet output_file: output/empty.out
102c4762a1bSJed Brown
103c4762a1bSJed Brown test:
104c4762a1bSJed Brown nsize: 4
105c4762a1bSJed Brown suffix: 2
106c4762a1bSJed Brown args: -viewer_hdf5_base_dimension2
107*3886731fSPierre Jolivet output_file: output/empty.out
108c4762a1bSJed Brown
109c4762a1bSJed Brown test:
110c4762a1bSJed Brown nsize: 4
111c4762a1bSJed Brown suffix: 3
112c4762a1bSJed Brown args: -viewer_hdf5_sp_output
113*3886731fSPierre Jolivet output_file: output/empty.out
114c4762a1bSJed Brown
115c4762a1bSJed Brown TEST*/
116