read_atomic_cross Subroutine

public subroutine read_atomic_cross(fid, grp, cross)

Reads in a cross section table from file and puts it into a AtomicCrossSection type

Arguments

Type IntentOptional AttributesName
integer(kind=HID_T), intent(in) :: fid

HDF5 file ID

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

HDF5 group to read from

type(AtomicCrossSection), intent(inout) :: cross

Atomic cross section


Calls

proc~~read_atomic_cross~~CallsGraph proc~read_atomic_cross read_atomic_cross proc~h5ltread_dataset_int_scalar_f h5ltread_dataset_int_scalar_f proc~read_atomic_cross->proc~h5ltread_dataset_int_scalar_f h5ltpath_valid_f h5ltpath_valid_f proc~read_atomic_cross->h5ltpath_valid_f h5ltread_dataset_double_f h5ltread_dataset_double_f proc~read_atomic_cross->h5ltread_dataset_double_f proc~h5ltread_dataset_double_scalar_f h5ltread_dataset_double_scalar_f proc~read_atomic_cross->proc~h5ltread_dataset_double_scalar_f h5ltread_dataset_int_f h5ltread_dataset_int_f proc~h5ltread_dataset_int_scalar_f->h5ltread_dataset_int_f proc~h5ltread_dataset_double_scalar_f->h5ltread_dataset_double_f

Called by

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

Contents

Source Code


Source Code

subroutine read_atomic_cross(fid, grp, cross)
    !+ Reads in a cross section table from file
    !+ and puts it into a [[AtomicCrossSection]] type
    integer(HID_T), intent(in)              :: fid
        !+ HDF5 file ID
    character(len=*), intent(in)            :: grp
        !+ HDF5 group to read from
    type(AtomicCrossSection), intent(inout) :: cross
        !+ Atomic cross section

    integer(HSIZE_T), dimension(3) :: dim3
    real(Float64) :: emin, emax, rmin
    integer :: i, n_max, m_max, error
    real(Float64), dimension(:,:,:), allocatable :: dummy3
    logical :: path_valid

    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_CROSS: Unknown atomic interaction: ', trim(grp)
        endif
        stop
    endif

    call h5ltread_dataset_int_scalar_f(fid, grp//"/nenergy", cross%nenergy, 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_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", cross%dlogE, error)
    cross%logemin = log10(emin)
    cross%logemax = log10(emax)

    allocate(dummy3(n_max, m_max, cross%nenergy))

    allocate(cross%log_cross(cross%m_max,cross%n_max, cross%nenergy))

    dim3 = [n_max, m_max, cross%nenergy]
    call h5ltread_dataset_double_f(fid,grp//"/cx", dummy3, dim3, error)
    rmin = minval(dummy3,dummy3.gt.0.d0)
    where (dummy3.le.0.0)
        dummy3 = 0.9*rmin
    end where
    cross%minlog_cross = log10(rmin)
    do i=1, cross%nenergy
        cross%log_cross(:,:,i) = log10(transpose(dummy3(1:nlevs,1:nlevs,i)))
    enddo
    deallocate(dummy3)

end subroutine read_atomic_cross