Something went wrong on our end
Select Git revision
modreadmeta.f90
-
Roy Fabrice authored
Some parameters were moved from type_parameter_pfof_snap to type_parameter_halofinder They were also added to type_parameter_conecreator_part
Roy Fabrice authoredSome parameters were moved from type_parameter_pfof_snap to type_parameter_halofinder They were also added to type_parameter_conecreator_part
modreadmeta.f90 20.46 KiB
module modreadmeta
use iso_fortran_env, only : OUTPUT_UNIT
use modconstant, only : type_parameter_halofinder,&
type_parameter_pfof,&
type_parameter_pfof_snap,&
type_parameter_pfof_cone,&
type_parameter_conecreator,&
type_parameter_conecreator_part,&
type_parameter_conecreator_grav,&
type_parameter_psod_snap,&
type_info_ramses,&
type_info_cone,&
type_info_cone_part,&
type_info_cone_grav,&
type_common_metadata,&
NAME_PFOF_SNAP,&
NAME_PFOF_CONE,&
NAME_PSOD_SNAP,&
NAME_CONECREATOR_PART,&
NAME_CONECREATOR_GRAV,&
GIT_VERSION
use modhdf5, only : hid_t,&
H5_STR_LEN,&
Hdf5_create_group,&
Hdf5_open_group,&
Hdf5_close_group,&
Hdf5_read_attr,&
Hdf5_read_data
use modmpicommons, only : procid
implicit none
private
public :: Read_meta_common, &
Read_meta_conecreator_parameter, &
Read_meta_halofinder_parameter, &
Read_meta_info_cone, &
Read_meta_info_ramses
contains
!=======================================================================
!> Read metadata common to every hdf5 ouput files of the pfof toolbox
subroutine Read_meta_common(file_id, common)
integer(kind=hid_t), intent(in) :: file_id
type(type_common_metadata), intent(out) :: common
character(len=H5_STR_LEN) :: groupname
character(len=H5_STR_LEN) :: name
integer(kind=hid_t) :: meta_id
integer(kind=4) :: length
#ifdef DEBUG
write(OUTPUT_UNIT,*) 'Enter Read_meta_common on process ', procid
#endif
groupname = 'metadata'
call Hdf5_open_group(file_id, groupname, meta_id)
length = len(common%created_by)
name = 'created_by'
call Hdf5_read_attr(meta_id, name, length, common%created_by)
name = 'version'
call Hdf5_read_attr(meta_id, name, common%version)
length = len(common%simulation_code)
name = 'simulation_code'
call Hdf5_read_attr(meta_id, name, length, common%simulation_code)
length = len(common%particle_type)
name = 'particle_type'
call Hdf5_read_attr(meta_id, name, length, common%particle_type)
name = 'constant_mass'
call Hdf5_read_attr(meta_id, name, common%constant_mass)
length = len(common%units)
name = 'units'
call Hdf5_read_attr(meta_id, name, length, common%units)
name = 'npart_simulation'
call Hdf5_read_data(meta_id, name, common%npart_simulation)
call Hdf5_close_group(meta_id)
#ifdef DEBUG
write(OUTPUT_UNIT,*) 'Exit Read_meta_common on process ', procid
#endif
end subroutine Read_meta_common
!=======================================================================
!> Read conecreator input parameters as metadata
subroutine Read_meta_conecreator_parameter(file_id, param_cone)
integer(kind=hid_t), intent(in) :: file_id
class(type_parameter_conecreator), intent(out) :: param_cone
character(len=H5_STR_LEN) :: groupname
character(len=H5_STR_LEN) :: aname
integer(kind=hid_t) :: gr_id
integer(kind=hid_t) :: gr_param_id
integer(kind=hid_t) :: gr_input_id
integer(kind=hid_t) :: gr_output_id
integer(kind=4) :: tmpint4
character(len=16) :: c_b_name
#ifdef DEBUG
write(OUTPUT_UNIT,*) 'Enter Read_meta_conecreator_parameter on process ', procid
#endif
select type(param_cone)
type is(type_parameter_conecreator_part)
c_b_name = NAME_CONECREATOR_PART
type is(type_parameter_conecreator_grav)
c_b_name = NAME_CONECREATOR_GRAV
end select
groupname = 'metadata'
call Hdf5_open_group(file_id, groupname, gr_id)
groupname = c_b_name//'_parameters'
call Hdf5_open_group(gr_id, groupname, gr_param_id)
groupname = 'input_parameters'
call Hdf5_open_group(gr_param_id, groupname, gr_input_id)
aname = 'input_path'
call Hdf5_read_attr(gr_input_id, aname, len(param_cone%input_path), &
param_cone%input_path)
aname = 'cone_input_file'
call Hdf5_read_attr(gr_input_id, aname, len(param_cone%cone_input_file), &
param_cone%cone_input_file)
aname = 'info_cone_input_file'
call Hdf5_read_attr(gr_input_id, aname, len(param_cone%info_cone_input_file), &
param_cone%info_cone_input_file)
aname = 'info_ramses_input_file'
call Hdf5_read_attr(gr_input_id, aname, len(param_cone%info_ramses_input_file), &
param_cone%info_ramses_input_file)
aname = 'nfile'
call Hdf5_read_attr(gr_input_id, aname, param_cone%nfile)
aname = 'first_file'
call Hdf5_read_attr(gr_input_id, aname, param_cone%first_file)
aname = 'last_file'
call Hdf5_read_attr(gr_input_id, aname, param_cone%last_file)
aname = 'cone_max_radius'
call Hdf5_read_attr(gr_input_id, aname, param_cone%cone_max_radius)
select type(param_cone)
type is (type_parameter_conecreator_part)
aname = 'do_read_ramses_part_id'
call Hdf5_read_attr(gr_input_id, aname, tmpint4)
param_cone%do_read_ramses_part_id = .false.
if(tmpint4==1) param_cone%do_read_ramses_part_id = .true.
aname = 'star'
call Hdf5_read_attr(gr_input_id, aname, tmpint4)
param_cone%star = .false.
if(tmpint4==1) param_cone%star = .true.
aname = 'metal'
call Hdf5_read_attr(gr_input_id, aname, tmpint4)
param_cone%metal = .false.
if(tmpint4==1) param_cone%metal = .true.
aname = 'do_skip_star'
call Hdf5_read_attr(gr_input_id, aname, tmpint4)
param_cone%do_skip_star = .false.
if(tmpint4==1) param_cone%do_skip_star = .true.
aname = 'do_skip_metal'
call Hdf5_read_attr(gr_input_id, aname, tmpint4)
param_cone%do_skip_metal = .false.
if(tmpint4==1) param_cone%do_skip_metal = .true.
aname = 'do_skip_mass'
call Hdf5_read_attr(gr_input_id, aname, tmpint4)
param_cone%do_skip_mass = .false.
if(tmpint4==1) param_cone%do_skip_mass = .true.
aname = 'do_read_potential'
call Hdf5_read_attr(gr_input_id, aname, tmpint4)
param_cone%do_read_potential = .false.
if(tmpint4==1) param_cone%do_read_potential = .true.
aname = 'do_read_gravitational_field'
call Hdf5_read_attr(gr_input_id, aname, tmpint4)
param_cone%do_read_gravitational_field = .false.
if(tmpint4==1) param_cone%do_read_gravitational_field = .true.
type is (type_parameter_conecreator_grav)
aname = 'nlevel'
call Hdf5_read_attr(gr_input_id, aname, param_cone%nlevel)
aname = 'levelmin'
call Hdf5_read_attr(gr_input_id, aname, param_cone%levelmin)
end select
call Hdf5_close_group(gr_input_id)
groupname = 'output_parameters'
call Hdf5_open_group(gr_param_id, groupname, gr_output_id)
aname = 'simulation_name'
call Hdf5_read_attr(gr_output_id, aname, len(param_cone%simulation_name), &
param_cone%simulation_name)
aname = 'cube_size'
call Hdf5_read_attr(gr_output_id, aname, param_cone%cube_size)
call Hdf5_close_group(gr_output_id)
call Hdf5_close_group(gr_param_id)
call Hdf5_close_group(gr_id)
#ifdef DEBUG
write(OUTPUT_UNIT,*) 'Exit Read_meta_conecreator_parameter on process ', procid
#endif
end subroutine Read_meta_conecreator_parameter
!=======================================================================
!> Read pfof input parameters as metadata
subroutine Read_meta_halofinder_parameter(file_id, param_pfof)
integer(kind=hid_t), intent(in) :: file_id
class(type_parameter_pfof), intent(out) :: param_pfof
character(len=H5_STR_LEN) :: groupname
character(len=H5_STR_LEN) :: aname
integer(kind=hid_t) :: gr_id
integer(kind=hid_t) :: gr_pfof_id
integer(kind=hid_t) :: gr_input_id
integer(kind=hid_t) :: gr_fof_id
integer(kind=hid_t) :: gr_output_id
integer(kind=4) :: tmpint4
character(len=16) :: c_b_name
#ifdef DEBUG
write(OUTPUT_UNIT,*) 'Enter Read_meta_halofinder_parameter on process ', procid
#endif
select type (param_pfof)
type is (type_parameter_pfof_snap)
c_b_name = NAME_PFOF_SNAP
type is (type_parameter_pfof_cone)
c_b_name = NAME_PFOF_CONE
end select
groupname = 'metadata'
call Hdf5_open_group(file_id,groupname, gr_id)
! pfof parameters:
groupname = trim(c_b_name)//'_parameters'
call Hdf5_open_group(gr_id, groupname, gr_pfof_id)
groupname = 'input_parameters'
call Hdf5_open_group(gr_pfof_id, groupname, gr_input_id)
aname = 'input_path'
call Hdf5_read_attr(gr_input_id, aname, len(param_pfof%input_path), &
param_pfof%input_path)
aname = 'part_input_file'
call Hdf5_read_attr(gr_input_id, aname, len(param_pfof%part_input_file),&
param_pfof%part_input_file)
aname = 'star'
call Hdf5_read_attr(gr_input_id, aname, tmpint4)
param_pfof%star = .false.
if(tmpint4==1) param_pfof%star = .true.
aname = 'metal'
call Hdf5_read_attr(gr_input_id, aname, tmpint4)
param_pfof%metal = .false.
if(tmpint4==1) param_pfof%metal = .true.
aname = 'do_skip_star'
call Hdf5_read_attr(gr_input_id, aname, tmpint4)
param_pfof%do_skip_star = .false.
if(tmpint4==1) param_pfof%do_skip_star = .true.
aname = 'do_skip_metal'
call Hdf5_read_attr(gr_input_id, aname, tmpint4)
param_pfof%do_skip_metal = .false.
if(tmpint4==1) param_pfof%do_skip_metal = .true.
aname = 'do_skip_mass'
call Hdf5_read_attr(gr_input_id, aname, tmpint4)
param_pfof%do_skip_mass = .false.
if(tmpint4==1) param_pfof%do_skip_mass = .true.
aname = 'do_read_potential'
call Hdf5_read_attr(gr_input_id, aname, tmpint4)
param_pfof%do_read_potential = .false.
if(tmpint4==1) param_pfof%do_read_potential = .true.
aname = 'do_read_gravitational_field'
call Hdf5_read_attr(gr_input_id, aname, tmpint4)
param_pfof%do_read_gravitational_field = .false.
if(tmpint4==1) param_pfof%do_read_gravitational_field = .true.
select type (param_pfof)
type is (type_parameter_pfof_snap)
aname = 'info_input_file'
call Hdf5_read_attr(gr_input_id, aname, len(param_pfof%info_input_file), &
param_pfof%info_input_file)
aname = 'grpsize'
call Hdf5_read_attr(gr_input_id, aname, param_pfof%grpsize)
aname = 'code_index'
call Hdf5_read_attr(gr_input_id, aname, len(param_pfof%code_index), param_pfof%code_index)
aname='do_read_from_cube'
call Hdf5_read_attr(gr_input_id, aname, tmpint4)
param_pfof%do_read_from_cube = .false.
if(tmpint4==1) param_pfof%do_read_from_cube = .true.
aname = 'gatherread_factor'
call Hdf5_read_attr(gr_input_id, aname, param_pfof%gatherread_factor)
type is (type_parameter_pfof_cone)
aname = 'shell_first_id'
call Hdf5_read_attr(gr_input_id, aname, param_pfof%shell_first_id)
aname = 'shell_last_id'
call Hdf5_read_attr(gr_input_id, aname, param_pfof%shell_last_id)
end select
call Hdf5_close_group(gr_input_id)
groupname = 'fof_parameters'
call Hdf5_open_group(gr_pfof_id, groupname, gr_fof_id)
! Write percolation parameter b as attribute
aname = 'percolation_length'
call Hdf5_read_attr(gr_fof_id, aname, param_pfof%percolation_length)
! Write minimum halo mass Mmin as attribute
aname = 'npart_halo_min'
call Hdf5_read_attr(gr_fof_id, aname, param_pfof%mmin)
aname = 'npart_halo_max'
call Hdf5_read_attr(gr_fof_id, aname, param_pfof%mmax)
! Write doUnbinding as attribute (not implemented yet)
aname = 'do_unbinding'
call Hdf5_read_attr(gr_fof_id, aname, tmpint4)
param_pfof%do_unbinding = .false.
if(tmpint4==1) param_pfof%do_unbinding = .true.
! Write doSubHalo as attribute (not implemented yet)
aname = 'do_subhalo'
call Hdf5_read_attr(gr_fof_id, aname, tmpint4)
param_pfof%do_subhalo = .false.
if(tmpint4==1) param_pfof%do_subhalo = .true.
select type (param_pfof)
type is (type_parameter_pfof_snap)
aname = 'do_fof'
call Hdf5_read_attr(gr_fof_id, aname, tmpint4)
param_pfof%do_fof = .false.
if(tmpint4==1) param_pfof%do_fof = .true.
end select
call Hdf5_close_group(gr_fof_id)
groupname = 'output_parameters'
call Hdf5_open_group(gr_pfof_id, groupname, gr_output_id)
! Write the simulation name as attribute
aname = 'simulation_name'
call Hdf5_read_attr(gr_output_id, aname, len(param_pfof%simulation_name), &
param_pfof%simulation_name)
! Write do_timings as attribute
aname = 'do_timings'
call Hdf5_read_attr(gr_fof_id, aname, tmpint4)
param_pfof%do_timings = .false.
if(tmpint4==1) param_pfof%do_timings = .true.
select type (param_pfof)
type is (type_parameter_pfof_snap)
! write the snapshot number: if halo computed from a lightcone then outputNB = -1
aname = 'snapshot'
call Hdf5_read_attr(gr_output_id, aname, param_pfof%snapshot)
aname = 'do_write_cube'
call Hdf5_read_attr(gr_output_id, aname, tmpint4)
param_pfof%do_write_cube = .false.
if(tmpint4==1) param_pfof%do_write_cube = .true.
aname = 'do_sort_cube'
call Hdf5_read_attr(gr_output_id, aname, tmpint4)
param_pfof%do_sort_cube = .false.
if(tmpint4==1) param_pfof%do_sort_cube = .true.
aname = 'gatherwrite_factor'
call Hdf5_read_attr(gr_output_id, aname, param_pfof%gatherwrite_factor)
end select
call Hdf5_close_group(gr_output_id)
call Hdf5_close_group(gr_pfof_id)
call Hdf5_close_group(gr_id)
#ifdef DEBUG
write(OUTPUT_UNIT,*) 'Exit Read_meta_halofinder_parameter on process ', procid
#endif
end subroutine Read_meta_halofinder_parameter
!=======================================================================
!> Read light cone info as metadata
subroutine Read_meta_info_cone(file_id, infocone, islast)
integer(kind=hid_t), intent(in) :: file_id
class(type_info_cone), intent(out) :: infocone
logical(kind=4), intent(in) :: islast
character(len=H5_STR_LEN) :: groupname
character(len=H5_STR_LEN) :: aname
integer(kind=hid_t) :: gr_id
integer(kind=hid_t) :: gr_cone_id
#ifdef DEBUG
write(OUTPUT_UNIT,*) 'Enter Read_meta_info_cone on process ', procid
#endif
groupname = 'metadata'
call Hdf5_open_group(file_id, groupname, gr_id)
groupname = 'cone_info'
if(islast) then
groupname = 'cone_info_last'
end if
call Hdf5_open_group(gr_id, groupname, gr_cone_id)
aname= "ncpu"
call Hdf5_read_attr(gr_cone_id, aname, infocone%ncpu)
aname = 'nstride'
call Hdf5_read_attr(gr_cone_id, aname, infocone%ncpu)
aname = 'nstep_coarse'
call Hdf5_read_attr(gr_cone_id, aname, infocone%nstep_coarse)
aname = 'cone_id'
call Hdf5_read_attr(gr_cone_id, aname, infocone%cone_id)
aname = 'nglobalfile'
call Hdf5_read_attr(gr_cone_id, aname, infocone%nglobalfile)
aname = 'isfullsky'
call Hdf5_read_attr(gr_cone_id, aname, infocone%isfullsky)
aname = 'aexp'
call Hdf5_read_attr(gr_cone_id, aname, infocone%aexp)
aname = 'observer_x'
call Hdf5_read_attr(gr_cone_id, aname, infocone%observer_x)
aname = 'observer_y'
call Hdf5_read_attr(gr_cone_id, aname, infocone%observer_y)
aname = 'observer_z'
call Hdf5_read_attr(gr_cone_id, aname, infocone%observer_z)
aname = 'observer_rds'
call Hdf5_read_attr(gr_cone_id, aname, infocone%observer_rds)
aname = 'cone_zlim'
call Hdf5_read_attr(gr_cone_id, aname, infocone%cone_zlim)
aname = 'amax'
call Hdf5_read_attr(gr_cone_id, aname, infocone%amax)
aname = 'amin'
call Hdf5_read_attr(gr_cone_id, aname, infocone%amin)
aname = 'zmax'
call Hdf5_read_attr(gr_cone_id, aname, infocone%zmax)
aname = 'zmin'
call Hdf5_read_attr(gr_cone_id, aname, infocone%zmin)
aname = 'dmax'
call Hdf5_read_attr(gr_cone_id, aname, infocone%dmax)
aname = 'dmin'
call Hdf5_read_attr(gr_cone_id, aname, infocone%dmin)
aname = 'dtol'
call Hdf5_read_attr(gr_cone_id, aname, infocone%dtol)
aname = 'thetay'
call Hdf5_read_attr(gr_cone_id, aname, infocone%thetay)
aname = 'thetaz'
call Hdf5_read_attr(gr_cone_id, aname, infocone%thetaz)
aname = 'theta'
call Hdf5_read_attr(gr_cone_id, aname, infocone%theta)
aname = 'phi'
call Hdf5_read_attr(gr_cone_id, aname, infocone%phi)
aname = 'aendconem2'
call Hdf5_read_attr(gr_cone_id, aname, infocone%aendconem2)
aname = 'aendconem1'
call Hdf5_read_attr(gr_cone_id, aname, infocone%aendconem1)
aname = 'aendcone'
call Hdf5_read_attr(gr_cone_id, aname, infocone%aendcone)
aname = 'zendconem2'
call Hdf5_read_attr(gr_cone_id, aname, infocone%zendconem2)
aname = 'zendconem1'
call Hdf5_read_attr(gr_cone_id, aname, infocone%zendconem1)
aname = 'zendcone'
call Hdf5_read_attr(gr_cone_id, aname, infocone%zendcone)
aname = 'dendconem2'
call Hdf5_read_attr(gr_cone_id, aname, infocone%dendconem2)
aname = 'dendconem1'
call Hdf5_read_attr(gr_cone_id, aname, infocone%dendconem1)
aname = 'dendcone'
call Hdf5_read_attr(gr_cone_id, aname, infocone%dendcone)
aname = 'future'
call Hdf5_read_attr(gr_cone_id, aname, infocone%future)
select type(infocone)
type is(type_info_cone_part)
aname = 'npart'
call Hdf5_read_data(gr_cone_id, aname, infocone%npart)
aname = 'aexpold'
call Hdf5_read_attr(gr_cone_id, aname, infocone%aexpold)
aname = 'zexpold'
call Hdf5_read_attr(gr_cone_id, aname, infocone%zexpold)
aname = 'zexp'
call Hdf5_read_attr(gr_cone_id, aname, infocone%zexp)
aname = 'dexpold'
call Hdf5_read_attr(gr_cone_id, aname, infocone%dexpold)
aname = 'dexp'
call Hdf5_read_attr(gr_cone_id, aname, infocone%dexp)
type is(type_info_cone_grav)
aname = 'nglobalcell'
call Hdf5_read_data(gr_cone_id, aname, infocone%nglobalcell)
aname = 'nlevel'
call Hdf5_read_attr(gr_cone_id, aname, infocone%nlevel)
aname = 'levelmin'
call Hdf5_read_attr(gr_cone_id, aname, infocone%levelmin)
aname = 'levelmax'
call Hdf5_read_attr(gr_cone_id, aname, infocone%levelmax)
end select
call Hdf5_close_group(gr_cone_id)
call Hdf5_close_group(gr_id)
#ifdef DEBUG
write(OUTPUT_UNIT,*) 'Exit Read_meta_info_cone on process ', procid
#endif
end subroutine Read_meta_info_cone
!=======================================================================
!> Read ramses info as metadata
subroutine Read_meta_info_ramses(file_id, inforamses, islast)
integer(kind=hid_t), intent(in) :: file_id
type(type_info_ramses), intent(out) :: inforamses
logical(kind=4), intent(in) :: islast
character(len=H5_STR_LEN) :: groupname
character(len=H5_STR_LEN) :: aname
integer(kind=hid_t) :: gr_id
integer(kind=hid_t) :: gr_ramses_id
#ifdef DEBUG
write(OUTPUT_UNIT,*) 'Enter Read_meta_info_ramses on process ', procid
#endif
groupname = 'metadata'
call Hdf5_open_group(file_id, groupname, gr_id)
! Ramses Info Metadata
groupname = 'ramses_info'
if(islast) then
groupname = 'ramses_info_last'
end if
call Hdf5_open_group(gr_id,groupname,gr_ramses_id)
aname = 'ncpu'
call Hdf5_read_attr(gr_ramses_id,aname,inforamses%ncpu)
aname = 'ndim'
call Hdf5_read_attr(gr_ramses_id,aname,inforamses%ndim)
aname = 'levelmin'
call Hdf5_read_attr(gr_ramses_id,aname,inforamses%levelmin)
aname = 'levelmax'
call Hdf5_read_attr(gr_ramses_id,aname,inforamses%levelmax)
aname = 'ngridmax'
call Hdf5_read_attr(gr_ramses_id,aname,inforamses%ngridmax)
aname = 'nstep_coarse'
call Hdf5_read_attr(gr_ramses_id,aname,inforamses%nstep_coarse)
aname = 'boxlen'
call Hdf5_read_attr(gr_ramses_id,aname,inforamses%boxlen)
aname = 'time'
call Hdf5_read_attr(gr_ramses_id,aname,inforamses%time)
aname = 'aexp'
call Hdf5_read_attr(gr_ramses_id,aname,inforamses%aexp)
aname = 'h0'
call Hdf5_read_attr(gr_ramses_id,aname,inforamses%h0)
aname = 'omega_m'
call Hdf5_read_attr(gr_ramses_id,aname,inforamses%omega_m)
aname = 'omega_l'
call Hdf5_read_attr(gr_ramses_id,aname,inforamses%omega_l)
aname = 'omega_k'
call Hdf5_read_attr(gr_ramses_id,aname,inforamses%omega_k)
aname = 'omega_b'
call Hdf5_read_attr(gr_ramses_id,aname,inforamses%omega_b)
aname = 'unit_l'
call Hdf5_read_attr(gr_ramses_id,aname,inforamses%unit_l)
aname = 'unit_d'
call Hdf5_read_attr(gr_ramses_id,aname,inforamses%unit_d)
aname = 'unit_t'
call Hdf5_read_attr(gr_ramses_id,aname,inforamses%unit_t)
call Hdf5_close_group(gr_ramses_id)
call Hdf5_close_group(gr_id)
#ifdef DEBUG
write(OUTPUT_UNIT,*) 'Exit Read_meta_info_ramses on process ', procid
#endif
end subroutine Read_meta_info_ramses
end module modreadmeta