-
Notifications
You must be signed in to change notification settings - Fork 16
HDF5
Bertrand Putigny edited this page Jan 21, 2015
·
14 revisions
If you installed ESTER with HDF5 (Hierarchical Data Format) support, writing/reading to/from a file ending with .hdf5 (-iand -o options) will write/read the file in HDF5 format.
You can visualize these files with a variety of tools such as
hdfview or h5dump, and you can easily read them from
several languages such as python, Fortran or IDL (see examples).
An ESTER's model written in HDF5 will have the following structure:
- a group named
starattached to the root of the file. - the parameters of the simulation are
attributesattached to this group- for instance
ndomainsornpts
- for instance
- the different fields are
datasetsattached to the/stargroup
-
/star/G: the stream function <math>\phi</math>. -
/star/N2: the squared Brunt-Väisälä frequency (in <math>rd^2 /s^2</math>). -
/star/R: radius of domain's boundary. -
/star/T: temperature. -
/star/X: hydrogen abundance. -
/star/Y: helium abundance. -
/star/Z: metal abundance. -
/star/nuc.eps: nuclear reaction. -
/star/p: pressure. -
/star/phi: gravitational potential. -
/star/phiex: gravitational potential in the external domain. -
/star/r: radius of collocation points. -
/star/rho: density. -
/star/th: colatitude. -
/star/w: angular velocity (<math>\omega</math>). -
/star/z: radial variable (<math>\zeta</math>).
You need to have the h5py package installed.
import h5py
f = h5py.File('star.hdf5', 'r') # open the file 'star.hdf5' ('r' = read only)
T = f['/star/T'][:] # read the temperature field
n = f['/star'].attrs['ndomains'] # read the number of domains
print "T (center): " + str(T[0][0])
print "T (equator): " + str(T[0][-1])You need to compile with the h5fc compiler wrapper.
program read_hdf5
use hdf5
implicit none
character*100 :: file_name = "star.hdf5"
integer status, error
integer(hsize_t) :: dims(2)
integer :: nr, nth
integer(HID_T) :: fid, gid, aid, did
double precision, allocatable :: T(:,:)
! init interface
call h5open_f(status)
! open the HDF5 file
call h5fopen_f(file_name, H5F_ACC_RDWR_F, fid, status)
! open the `star' group
call h5gopen_f(fid, "star", gid, error)
! open the `nr' attribute
call h5aopen_f(gid, "nr", aid, error)
! read the attribute
dims(1) = 1
dims(2) = 0
call h5aread_f(aid, H5T_NATIVE_INTEGER, nr, dims, error)
! close the attribute
call h5aclose_f(aid, error)
! open the `nth' attribute
call h5aopen_f(gid, "nth", aid, error)
! read the attribute
dims(1) = 1
dims(2) = 0
call h5aread_f(aid, H5T_NATIVE_INTEGER, nth, dims, error)
! close the attribute
call h5aclose_f(aid, error)
print *, "nr: ", nr
print *, "nth:", nth
! allocate memory for the temperature field
allocate(T(nr, nth))
! open the `T' dataset
call h5dopen_f(gid, "T", did, error)
! read the field
dims(1) = nr
dims(2) = nth
call h5dread_f(did, H5T_NATIVE_DOUBLE, T, dims, error)
print *, "T at the center: ", T(1, 1)
print *, "T at the equator:", T(nr, 1)
deallocate(T)
! close dataset, group and file
call h5dclose_f(did, error)
call h5gclose_f(gid, error)
call h5fclose_f(fid, error)
end programfile = 'star.hdf5'
fid = h5f_open(file)
did = h5d_open(fid, '/star/T')
T = h5d_read(did)
print, 'Temperature at the center: ', T[0, 0]
h5d_close, did
h5f_close, fid