Reads in a atomic transitions table from file and puts it into a AtomicTransitions type
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer(kind=HID_T), | intent(in) | :: | fid | HDF5 file ID |
||
character(len=*), | intent(in) | :: | grp | HDF5 group to read from |
||
real(kind=Float64), | intent(in), | dimension(2) | :: | b_amu | Atomic masses of "beam" species (beam ion and thermal ion) |
|
real(kind=Float64), | intent(in) | :: | t_amu | Atomic mass of "target" species (thermal ion) |
||
type(AtomicTransitions), | intent(inout) | :: | rates | Atomic transitions |
subroutine read_atomic_transitions(fid, grp, b_amu, t_amu, rates)
!+ Reads in a atomic transitions table from file
!+ and puts it into a [[AtomicTransitions]] type
integer(HID_T), intent(in) :: fid
!+ HDF5 file ID
character(len=*), intent(in) :: grp
!+ HDF5 group to read from
real(Float64), dimension(2), intent(in) :: b_amu
!+ Atomic masses of "beam" species (beam ion and thermal ion)
real(Float64), intent(in) :: t_amu
!+ Atomic mass of "target" species (thermal ion)
type(AtomicTransitions), intent(inout) :: rates
!+ Atomic transitions
integer(HSIZE_T), dimension(2) :: dim2
integer(HSIZE_T), dimension(4) :: dim4
integer(HSIZE_T), dimension(5) :: dim5
logical :: path_valid
integer :: i, j, n, n_max, m_max, error
integer :: n_bt_amu, tt_ind, bt_ind, drank
real(Float64) :: emin,emax,tmin,tmax,rmin
real(Float64) :: bt_min, tt_min, tt_dum, bt_dum
real(Float64), dimension(2) :: bt_amu, tt_amu
real(Float64), dimension(:,:), allocatable :: dummy2
real(Float64), dimension(:,:,:,:), allocatable :: dummy4
real(Float64), dimension(:,:,:,:,:), allocatable :: dummy5
call h5ltpath_valid_f(fid, grp, .True., path_valid, error)
if(.not.path_valid) then
if(inputs%verbose.ge.0) then
write(*,'(a,a)') 'READ_ATOMIC_TRANSITIONS: Unknown atomic interaction: ', trim(grp)
endif
stop
endif
call h5ltread_dataset_int_scalar_f(fid, grp//"/n_bt_amu", n_bt_amu, error)
allocate(dummy2(2, n_bt_amu))
dim2 = [2, n_bt_amu]
call h5ltread_dataset_double_f(fid, grp//"/bt_amu", dummy2, dim2, error)
call h5ltread_dataset_int_scalar_f(fid, grp//"/n_max", n_max, error)
call h5ltread_dataset_int_scalar_f(fid, grp//"/m_max", m_max, error)
call h5ltread_dataset_int_scalar_f(fid, grp//"/nenergy", rates%nenergy, error)
call h5ltread_dataset_double_scalar_f(fid, grp//"/emin", emin, error)
call h5ltread_dataset_double_scalar_f(fid, grp//"/emax", emax, error)
call h5ltread_dataset_double_scalar_f(fid, grp//"/dlogE", rates%dlogE, error)
call h5ltread_dataset_int_scalar_f(fid, grp//"/ntemp", rates%ntemp, error)
call h5ltread_dataset_double_scalar_f(fid, grp//"/tmin", tmin, error)
call h5ltread_dataset_double_scalar_f(fid, grp//"/tmax", tmax, error)
call h5ltread_dataset_double_scalar_f(fid, grp//"/dlogT", rates%dlogT, error)
rates%logemin = log10(emin)
rates%logemax = log10(emax)
rates%logtmin = log10(tmin)
rates%logtmax = log10(tmax)
bt_ind = 1
tt_ind = 1
bt_amu = [b_amu(1), t_amu]
tt_amu = [b_amu(2), t_amu]
bt_min = norm2(bt_amu - dummy2(:,1))
tt_min = norm2(tt_amu - dummy2(:,1))
do i=2, n_bt_amu
bt_dum = norm2(bt_amu - dummy2(:,i))
tt_dum = norm2(tt_amu - dummy2(:,i))
if(bt_dum.lt.bt_min) then
bt_min = bt_dum
bt_ind = i
endif
if(tt_dum.lt.tt_min) then
tt_min = tt_dum
tt_ind = i
endif
enddo
rates%ab(1) = dummy2(1,bt_ind)
rates%ab(2) = dummy2(1,tt_ind)
deallocate(dummy2)
allocate(rates%log_pop(&
rates%m_max, &
rates%n_max, &
rates%nenergy, &
rates%ntemp, 2))
allocate(rates%log_depop(&
rates%n_max, &
rates%nenergy, &
rates%ntemp, 2))
rates%log_pop = 0.d0
rates%log_depop = 0.d0
!!Read CX
call h5ltpath_valid_f(fid, grp//"/cx", .True., path_valid, error)
if(path_valid) then
call h5ltget_dataset_ndims_f(fid, grp//"/cx", drank, error)
if(drank.eq.4) then
allocate(dummy4(n_max, &
rates%nenergy, &
rates%ntemp, n_bt_amu))
dim4 = [n_max, rates%nenergy, rates%ntemp,n_bt_amu]
call h5ltread_dataset_double_f(fid, grp//"/cx", dummy4, dim4, error)
do j=1,rates%ntemp
do i=1,rates%nenergy
do n=1,rates%n_max
rates%log_depop(n,i,j,1) = dummy4(n,i,j,bt_ind)
rates%log_depop(n,i,j,2) = dummy4(n,i,j,tt_ind)
enddo
enddo
enddo
deallocate(dummy4)
endif
if(drank.eq.5) then
allocate(dummy5(n_max, m_max, &
rates%nenergy, &
rates%ntemp, n_bt_amu))
dim5 = [n_max, m_max, rates%nenergy, rates%ntemp,n_bt_amu]
call h5ltread_dataset_double_f(fid, grp//"/cx", dummy5, dim5, error)
do j=1,rates%ntemp
do i=1,rates%nenergy
do n=1,rates%n_max
rates%log_depop(n,i,j,1) = sum(dummy5(n,:,i,j,bt_ind))
rates%log_depop(n,i,j,2) = sum(dummy5(n,:,i,j,tt_ind))
enddo
enddo
enddo
deallocate(dummy5)
endif
endif
!!Read ionization
call h5ltpath_valid_f(fid, grp//"/ionization", .True., path_valid, error)
if(path_valid) then
allocate(dummy4(n_max, &
rates%nenergy, &
rates%ntemp, n_bt_amu))
dim4 = [n_max, rates%nenergy, rates%ntemp,n_bt_amu]
call h5ltread_dataset_double_f(fid, grp//"/ionization", dummy4, dim4, error)
do j=1,rates%ntemp
do i=1,rates%nenergy
do n=1,rates%n_max
rates%log_depop(n,i,j,1) = rates%log_depop(n,i,j,1) + &
dummy4(n,i,j,bt_ind)
rates%log_depop(n,i,j,2) = rates%log_depop(n,i,j,2) + &
dummy4(n,i,j,tt_ind)
enddo
enddo
enddo
deallocate(dummy4)
endif
!!Read excitation
call h5ltpath_valid_f(fid, grp//"/excitation", .True., path_valid, error)
if(path_valid) then
allocate(dummy5(n_max, m_max,&
rates%nenergy, &
rates%ntemp, n_bt_amu))
dim5 = [n_max, m_max, rates%nenergy, rates%ntemp,n_bt_amu]
call h5ltread_dataset_double_f(fid, grp//"/excitation", dummy5, dim5, error)
do j=1,rates%ntemp
do i=1,rates%nenergy
rates%log_pop(:,:,i,j,1) = transpose(dummy5(1:nlevs,1:nlevs,i,j,bt_ind))
rates%log_pop(:,:,i,j,2) = transpose(dummy5(1:nlevs,1:nlevs,i,j,tt_ind))
do n=1,rates%n_max
rates%log_depop(n,i,j,1) = rates%log_depop(n,i,j,1) + &
sum(dummy5(n,:,i,j,bt_ind))
rates%log_depop(n,i,j,2) = rates%log_depop(n,i,j,2) + &
sum(dummy5(n,:,i,j,tt_ind))
enddo
enddo
enddo
deallocate(dummy5)
endif
rmin = minval(rates%log_depop, rates%log_depop.gt.0.d0)
where (rates%log_depop.le.0.d0)
rates%log_depop = 0.9*rmin
end where
rates%minlog_depop = log10(rmin)
rates%log_depop = log10(rates%log_depop)
rmin = minval(rates%log_pop, rates%log_pop.gt.0.d0)
where (rates%log_pop.le.0.d0)
rates%log_pop = 0.9*rmin
end where
rates%minlog_pop = log10(rmin)
rates%log_pop = log10(rates%log_pop)
end subroutine read_atomic_transitions