xref: /petsc/src/sys/tests/ex80f.F90 (revision a968899dcf54030bc4915260052f186034d09103)
1*a968899dSTapashree Pradhan!
2*a968899dSTapashree Pradhan! PETSc Program to test HDF5 viewer and HDF5 attribute I/O
3*a968899dSTapashree Pradhan!
4*a968899dSTapashree Pradhanprogram main
5*a968899dSTapashree Pradhan#include <petsc/finclude/petscsys.h>
6*a968899dSTapashree Pradhan#include <petsc/finclude/petscvec.h>
7*a968899dSTapashree Pradhan  use petscsys
8*a968899dSTapashree Pradhan  use petscvec
9*a968899dSTapashree Pradhan  implicit none
10*a968899dSTapashree Pradhan
11*a968899dSTapashree Pradhan  PetscViewer :: viewer
12*a968899dSTapashree Pradhan  PetscErrorCode :: ierr
13*a968899dSTapashree Pradhan  Vec :: x
14*a968899dSTapashree Pradhan  PetscReal, parameter :: one = 1.0
15*a968899dSTapashree Pradhan  PetscInt :: ival = 42
16*a968899dSTapashree Pradhan  PetscReal :: rval = 3.14
17*a968899dSTapashree Pradhan  ! initialize PETSc
18*a968899dSTapashree Pradhan  PetscCallA(PetscInitialize(ierr))
19*a968899dSTapashree Pradhan  ! create and write a vector
20*a968899dSTapashree Pradhan  PetscCallA(VecCreate(PETSC_COMM_WORLD,x,ierr))
21*a968899dSTapashree Pradhan  PetscCallA(PetscObjectSetName(x,"vec",ierr))
22*a968899dSTapashree Pradhan  PetscCallA(VecSetSizes(x,3,PETSC_DETERMINE,ierr))
23*a968899dSTapashree Pradhan  PetscCallA(VecSetType(x,VECSTANDARD,ierr))
24*a968899dSTapashree Pradhan  PetscCallA(VecSet(x,one,ierr))
25*a968899dSTapashree Pradhan  PetscCallA(PetscViewerCreate(PETSC_COMM_WORLD,viewer,ierr))
26*a968899dSTapashree Pradhan  PetscCallA(PetscViewerSetType(viewer, PETSCVIEWERHDF5, ierr))
27*a968899dSTapashree Pradhan  PetscCallA(PetscViewerFileSetMode(viewer, FILE_MODE_WRITE, ierr))
28*a968899dSTapashree Pradhan  PetscCallA(PetscViewerFileSetName(viewer, "ex80f.hdf5", ierr))
29*a968899dSTapashree Pradhan  PetscCallA(VecView(x,viewer,ierr))
30*a968899dSTapashree Pradhan  PetscCallA(PetscViewerHDF5WriteAttribute(viewer, "vec", "int_attribute", ival, ierr))
31*a968899dSTapashree Pradhan  PetscCallA(PetscViewerHDF5WriteAttribute(viewer, "vec", "float_attribute", rval, ierr))
32*a968899dSTapashree Pradhan  PetscCallA(PetscViewerDestroy(viewer, ierr))
33*a968899dSTapashree Pradhan  PetscCallA(VecDestroy(x, ierr))
34*a968899dSTapashree Pradhan  PetscCallA(PetscFinalize(ierr))
35*a968899dSTapashree Pradhanend program main
36*a968899dSTapashree Pradhan!/*TEST
37*a968899dSTapashree Pradhan!  build:
38*a968899dSTapashree Pradhan!    requires: hdf5
39*a968899dSTapashree Pradhan!TEST*/