read_atomic_transitions Subroutine

public subroutine read_atomic_transitions(fid, grp, rates)

Reads in a atomic transitions table from file and puts it into a AtomicTransitions type

Arguments

TypeIntentOptionalAttributesName
integer(kind=HID_T), intent(in) :: fid

HDF5 file ID

character(len=*), intent(in) :: grp

HDF5 group to read from

type(AtomicTransitions), intent(inout) :: rates

Atomic transitions


Calls

proc~~read_atomic_transitions~~CallsGraph proc~read_atomic_transitions read_atomic_transitions h5ltpath_valid_f h5ltpath_valid_f proc~read_atomic_transitions->h5ltpath_valid_f proc~h5ltread_dataset_double_scalar_f h5ltread_dataset_double_scalar_f proc~read_atomic_transitions->proc~h5ltread_dataset_double_scalar_f proc~h5ltread_dataset_int_scalar_f h5ltread_dataset_int_scalar_f proc~read_atomic_transitions->proc~h5ltread_dataset_int_scalar_f h5ltget_dataset_ndims_f h5ltget_dataset_ndims_f proc~read_atomic_transitions->h5ltget_dataset_ndims_f h5ltread_dataset_double_f h5ltread_dataset_double_f proc~read_atomic_transitions->h5ltread_dataset_double_f proc~h5ltread_dataset_double_scalar_f->h5ltread_dataset_double_f h5ltread_dataset_int_f h5ltread_dataset_int_f proc~h5ltread_dataset_int_scalar_f->h5ltread_dataset_int_f

Called by

proc~~read_atomic_transitions~~CalledByGraph proc~read_atomic_transitions read_atomic_transitions proc~read_tables read_tables proc~read_tables->proc~read_atomic_transitions program~fidasim fidasim program~fidasim->proc~read_tables

Contents


Source Code

subroutine read_atomic_transitions(fid, grp,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
    type(AtomicTransitions), intent(inout)  :: rates
        !+ Atomic transitions

    integer(HSIZE_T), dimension(2) :: dim2
    integer(HSIZE_T), dimension(3) :: dim3
    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 :: drank
    real(Float64) :: emin,emax,tmin,tmax,rmin
    real(Float64), dimension(:,:), allocatable :: dummy2
    real(Float64), dimension(:,:,:), allocatable :: dummy3
    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_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)

    allocate(rates%log_pop(&
                    rates%m_max, &
                    rates%n_max, &
                    rates%nenergy, &
                    rates%ntemp))
    allocate(rates%log_depop(&
                    rates%n_max, &
                    rates%nenergy, &
                    rates%ntemp))
    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.3) then
            allocate(dummy3(n_max, &
                           rates%nenergy, &
                           rates%ntemp))
            dim3 = [n_max, rates%nenergy, rates%ntemp]
            call h5ltread_dataset_double_f(fid, grp//"/cx", dummy3, dim3, error)
            do j=1,rates%ntemp
                do i=1,rates%nenergy
                    do n=1,rates%n_max
                        rates%log_depop(n,i,j) = dummy3(n,i,j)
                    enddo
                enddo
            enddo
            deallocate(dummy3)
        endif
        if(drank.eq.4) then
            allocate(dummy4(n_max, m_max, &
                           rates%nenergy, &
                           rates%ntemp))
            dim4 = [n_max, m_max, rates%nenergy, rates%ntemp]
            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) = sum(dummy4(n,:,i,j))
                    enddo
                enddo
            enddo
            deallocate(dummy4)
        endif
    endif

    !!Read ionization
    call h5ltpath_valid_f(fid, grp//"/ionization", .True., path_valid, error)
    if(path_valid) then
        allocate(dummy3(n_max, &
                       rates%nenergy, &
                       rates%ntemp))
        dim3 = [n_max, rates%nenergy, rates%ntemp]
        call h5ltread_dataset_double_f(fid, grp//"/ionization", dummy3, dim3, error)
        do j=1,rates%ntemp
            do i=1,rates%nenergy
                do n=1,rates%n_max
                    rates%log_depop(n,i,j) = rates%log_depop(n,i,j) + &
                                               dummy3(n,i,j)
                enddo
            enddo
        enddo
        deallocate(dummy3)
    endif

    !!Read excitation
    call h5ltpath_valid_f(fid, grp//"/excitation", .True., path_valid, error)
    if(path_valid) then
        allocate(dummy4(n_max, m_max,&
                       rates%nenergy, &
                       rates%ntemp))
        dim4 = [n_max, m_max, rates%nenergy, rates%ntemp]
        call h5ltread_dataset_double_f(fid, grp//"/excitation", dummy4, dim4, error)
        do j=1,rates%ntemp
            do i=1,rates%nenergy
                rates%log_pop(:,:,i,j) = transpose(dummy4(1:nlevs,1:nlevs,i,j))
                do n=1,rates%n_max
                    rates%log_depop(n,i,j) = rates%log_depop(n,i,j) + &
                                               sum(dummy4(n,:,i,j))
                enddo
            enddo
        enddo
        deallocate(dummy4)
    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