write_bb_H_e Subroutine

public subroutine write_bb_H_e(id, namelist_file, n_max, m_max)

Write Hydrogen-Electron interaction cross sections to a HDF5 file

Arguments

Type IntentOptional AttributesName
integer(kind=HID_T), intent(inout) :: id

HDF5 file or group object id

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

Namelist file that contains settings

integer, intent(in) :: n_max

Number of initial atomic energy states to calculate

integer, intent(in) :: m_max

Number of final atomic energy states to calculate


Calls

proc~~write_bb_h_e~~CallsGraph proc~write_bb_h_e write_bb_H_e interface~h5ltmake_compressed_dataset_double_f h5ltmake_compressed_dataset_double_f proc~write_bb_h_e->interface~h5ltmake_compressed_dataset_double_f h5ltmake_dataset_int_f h5ltmake_dataset_int_f proc~write_bb_h_e->h5ltmake_dataset_int_f h5ltmake_dataset_double_f h5ltmake_dataset_double_f proc~write_bb_h_e->h5ltmake_dataset_double_f h5gclose_f h5gclose_f proc~write_bb_h_e->h5gclose_f h5ltset_attribute_string_f h5ltset_attribute_string_f proc~write_bb_h_e->h5ltset_attribute_string_f proc~e_excit e_excit proc~write_bb_h_e->proc~e_excit h5gcreate_f h5gcreate_f proc~write_bb_h_e->h5gcreate_f proc~e_ioniz e_ioniz proc~write_bb_h_e->proc~e_ioniz proc~h5ltmake_compressed_dataset_double_f_7 h5ltmake_compressed_dataset_double_f_7 interface~h5ltmake_compressed_dataset_double_f->proc~h5ltmake_compressed_dataset_double_f_7 proc~h5ltmake_compressed_dataset_double_f_6 h5ltmake_compressed_dataset_double_f_6 interface~h5ltmake_compressed_dataset_double_f->proc~h5ltmake_compressed_dataset_double_f_6 proc~h5ltmake_compressed_dataset_double_f_5 h5ltmake_compressed_dataset_double_f_5 interface~h5ltmake_compressed_dataset_double_f->proc~h5ltmake_compressed_dataset_double_f_5 proc~h5ltmake_compressed_dataset_double_f_4 h5ltmake_compressed_dataset_double_f_4 interface~h5ltmake_compressed_dataset_double_f->proc~h5ltmake_compressed_dataset_double_f_4 proc~h5ltmake_compressed_dataset_double_f_3 h5ltmake_compressed_dataset_double_f_3 interface~h5ltmake_compressed_dataset_double_f->proc~h5ltmake_compressed_dataset_double_f_3 proc~h5ltmake_compressed_dataset_double_f_2 h5ltmake_compressed_dataset_double_f_2 interface~h5ltmake_compressed_dataset_double_f->proc~h5ltmake_compressed_dataset_double_f_2 proc~h5ltmake_compressed_dataset_double_f_1 h5ltmake_compressed_dataset_double_f_1 interface~h5ltmake_compressed_dataset_double_f->proc~h5ltmake_compressed_dataset_double_f_1 proc~e_excit_n e_excit_n proc~e_excit->proc~e_excit_n proc~e_ioniz_n e_ioniz_n proc~e_ioniz->proc~e_ioniz_n proc~e_ioniz_2_janev e_ioniz_2_janev proc~e_ioniz_n->proc~e_ioniz_2_janev proc~e_ioniz_1_janev e_ioniz_1_janev proc~e_ioniz_n->proc~e_ioniz_1_janev proc~e_ioniz_3_janev e_ioniz_3_janev proc~e_ioniz_n->proc~e_ioniz_3_janev proc~e_excit_f e_excit_f proc~e_excit_n->proc~e_excit_f proc~e_excit_2_3_janev e_excit_2_3_janev proc~e_excit_n->proc~e_excit_2_3_janev proc~e_excit_1_janev e_excit_1_janev proc~e_excit_n->proc~e_excit_1_janev proc~h5ltmake_compressed_dataset_double_f_7->h5ltmake_dataset_double_f h5dclose_f h5dclose_f proc~h5ltmake_compressed_dataset_double_f_7->h5dclose_f h5sclose_f h5sclose_f proc~h5ltmake_compressed_dataset_double_f_7->h5sclose_f h5pset_chunk_f h5pset_chunk_f proc~h5ltmake_compressed_dataset_double_f_7->h5pset_chunk_f h5pcreate_f h5pcreate_f proc~h5ltmake_compressed_dataset_double_f_7->h5pcreate_f proc~chunk_size chunk_size proc~h5ltmake_compressed_dataset_double_f_7->proc~chunk_size h5pset_deflate_f h5pset_deflate_f proc~h5ltmake_compressed_dataset_double_f_7->h5pset_deflate_f h5pset_shuffle_f h5pset_shuffle_f proc~h5ltmake_compressed_dataset_double_f_7->h5pset_shuffle_f h5screate_simple_f h5screate_simple_f proc~h5ltmake_compressed_dataset_double_f_7->h5screate_simple_f h5dwrite_f h5dwrite_f proc~h5ltmake_compressed_dataset_double_f_7->h5dwrite_f h5pclose_f h5pclose_f proc~h5ltmake_compressed_dataset_double_f_7->h5pclose_f h5dcreate_f h5dcreate_f proc~h5ltmake_compressed_dataset_double_f_7->h5dcreate_f proc~h5ltmake_compressed_dataset_double_f_6->h5ltmake_dataset_double_f proc~h5ltmake_compressed_dataset_double_f_6->h5dclose_f proc~h5ltmake_compressed_dataset_double_f_6->h5sclose_f proc~h5ltmake_compressed_dataset_double_f_6->h5pset_chunk_f proc~h5ltmake_compressed_dataset_double_f_6->h5pcreate_f proc~h5ltmake_compressed_dataset_double_f_6->proc~chunk_size proc~h5ltmake_compressed_dataset_double_f_6->h5pset_deflate_f proc~h5ltmake_compressed_dataset_double_f_6->h5pset_shuffle_f proc~h5ltmake_compressed_dataset_double_f_6->h5screate_simple_f proc~h5ltmake_compressed_dataset_double_f_6->h5dwrite_f proc~h5ltmake_compressed_dataset_double_f_6->h5pclose_f proc~h5ltmake_compressed_dataset_double_f_6->h5dcreate_f proc~h5ltmake_compressed_dataset_double_f_5->h5ltmake_dataset_double_f proc~h5ltmake_compressed_dataset_double_f_5->h5dclose_f proc~h5ltmake_compressed_dataset_double_f_5->h5sclose_f proc~h5ltmake_compressed_dataset_double_f_5->h5pset_chunk_f proc~h5ltmake_compressed_dataset_double_f_5->h5pcreate_f proc~h5ltmake_compressed_dataset_double_f_5->proc~chunk_size proc~h5ltmake_compressed_dataset_double_f_5->h5pset_deflate_f proc~h5ltmake_compressed_dataset_double_f_5->h5pset_shuffle_f proc~h5ltmake_compressed_dataset_double_f_5->h5screate_simple_f proc~h5ltmake_compressed_dataset_double_f_5->h5dwrite_f proc~h5ltmake_compressed_dataset_double_f_5->h5pclose_f proc~h5ltmake_compressed_dataset_double_f_5->h5dcreate_f proc~h5ltmake_compressed_dataset_double_f_4->h5ltmake_dataset_double_f proc~h5ltmake_compressed_dataset_double_f_4->h5dclose_f proc~h5ltmake_compressed_dataset_double_f_4->h5sclose_f proc~h5ltmake_compressed_dataset_double_f_4->h5pset_chunk_f proc~h5ltmake_compressed_dataset_double_f_4->h5pcreate_f proc~h5ltmake_compressed_dataset_double_f_4->proc~chunk_size proc~h5ltmake_compressed_dataset_double_f_4->h5pset_deflate_f proc~h5ltmake_compressed_dataset_double_f_4->h5pset_shuffle_f proc~h5ltmake_compressed_dataset_double_f_4->h5screate_simple_f proc~h5ltmake_compressed_dataset_double_f_4->h5dwrite_f proc~h5ltmake_compressed_dataset_double_f_4->h5pclose_f proc~h5ltmake_compressed_dataset_double_f_4->h5dcreate_f proc~h5ltmake_compressed_dataset_double_f_3->h5ltmake_dataset_double_f proc~h5ltmake_compressed_dataset_double_f_3->h5dclose_f proc~h5ltmake_compressed_dataset_double_f_3->h5sclose_f proc~h5ltmake_compressed_dataset_double_f_3->h5pset_chunk_f proc~h5ltmake_compressed_dataset_double_f_3->h5pcreate_f proc~h5ltmake_compressed_dataset_double_f_3->proc~chunk_size proc~h5ltmake_compressed_dataset_double_f_3->h5pset_deflate_f proc~h5ltmake_compressed_dataset_double_f_3->h5pset_shuffle_f proc~h5ltmake_compressed_dataset_double_f_3->h5screate_simple_f proc~h5ltmake_compressed_dataset_double_f_3->h5dwrite_f proc~h5ltmake_compressed_dataset_double_f_3->h5pclose_f proc~h5ltmake_compressed_dataset_double_f_3->h5dcreate_f proc~h5ltmake_compressed_dataset_double_f_2->h5ltmake_dataset_double_f proc~h5ltmake_compressed_dataset_double_f_2->h5dclose_f proc~h5ltmake_compressed_dataset_double_f_2->h5sclose_f proc~h5ltmake_compressed_dataset_double_f_2->h5pset_chunk_f proc~h5ltmake_compressed_dataset_double_f_2->h5pcreate_f proc~h5ltmake_compressed_dataset_double_f_2->proc~chunk_size proc~h5ltmake_compressed_dataset_double_f_2->h5pset_deflate_f proc~h5ltmake_compressed_dataset_double_f_2->h5pset_shuffle_f proc~h5ltmake_compressed_dataset_double_f_2->h5screate_simple_f proc~h5ltmake_compressed_dataset_double_f_2->h5dwrite_f proc~h5ltmake_compressed_dataset_double_f_2->h5pclose_f proc~h5ltmake_compressed_dataset_double_f_2->h5dcreate_f proc~h5ltmake_compressed_dataset_double_f_1->h5ltmake_dataset_double_f proc~h5ltmake_compressed_dataset_double_f_1->h5dclose_f proc~h5ltmake_compressed_dataset_double_f_1->h5sclose_f proc~h5ltmake_compressed_dataset_double_f_1->h5pset_chunk_f proc~h5ltmake_compressed_dataset_double_f_1->h5pcreate_f proc~h5ltmake_compressed_dataset_double_f_1->proc~chunk_size proc~h5ltmake_compressed_dataset_double_f_1->h5pset_deflate_f proc~h5ltmake_compressed_dataset_double_f_1->h5pset_shuffle_f proc~h5ltmake_compressed_dataset_double_f_1->h5screate_simple_f proc~h5ltmake_compressed_dataset_double_f_1->h5dwrite_f proc~h5ltmake_compressed_dataset_double_f_1->h5pclose_f proc~h5ltmake_compressed_dataset_double_f_1->h5dcreate_f proc~e_excit_1_janev->proc~e_excit_f proc~e_excit_1_3_janev e_excit_1_3_janev proc~e_excit_1_janev->proc~e_excit_1_3_janev proc~e_excit_1_2_janev e_excit_1_2_janev proc~e_excit_1_janev->proc~e_excit_1_2_janev proc~e_excit_1_5_janev e_excit_1_5_janev proc~e_excit_1_janev->proc~e_excit_1_5_janev proc~e_excit_1_4_janev e_excit_1_4_janev proc~e_excit_1_janev->proc~e_excit_1_4_janev

Called by

proc~~write_bb_h_e~~CalledByGraph proc~write_bb_h_e write_bb_H_e program~generate_tables generate_tables program~generate_tables->proc~write_bb_h_e

Contents

Source Code


Source Code

subroutine write_bb_H_e(id, namelist_file, n_max, m_max)
    !+ Write Hydrogen-Electron interaction cross sections to a HDF5 file
    integer(HID_T), intent(inout) :: id
        !+ HDF5 file or group object id
    character(len=*), intent(in)  :: namelist_file
        !+ Namelist file that contains settings
    integer, intent(in)           :: n_max
        !+ Number of initial atomic energy states to calculate
    integer, intent(in)           :: m_max
        !+ Number of final atomic energy states to calculate

    real(Float64) :: emin
    real(Float64) :: emax
    integer :: nenergy

    real(Float64) :: eb
    real(Float64) :: dlogE
    real(Float64), dimension(:), allocatable :: ebarr

    real(Float64), dimension(:,:), allocatable :: ioniz
    real(Float64), dimension(:,:,:), allocatable :: excit

    integer(HID_T) :: gid
    integer(HSIZE_T), dimension(1) :: dim1
    integer(HSIZE_T), dimension(2) :: dim2
    integer(HSIZE_T), dimension(3) :: dim3

    integer :: i, cnt, error
    logical :: exis

    NAMELIST /H_e_cross/ nenergy, emin, emax

    nenergy = 200; emin = 1.d-3 ; emax = 8.d2

    inquire(file=namelist_file,exist=exis)
    if(.not.exis) then
        write(*,'(a,a)') 'WRITE_BB_H_E: Input file does not exist: ',trim(namelist_file)
        write(*,'(a)') 'Continuing with default settings...'
    else
        open(13,file=namelist_file)
        read(13,NML=H_e_cross)
        close(13)
    endif

    allocate(ebarr(nenergy))
    allocate(ioniz(n_max,nenergy))
    allocate(excit(n_max,m_max,nenergy))

    ebarr = 0.d0
    ioniz = 0.d0
    excit = 0.d0

    write(*,'(a)') "---- H-e cross sections settings ----"
    write(*,'(T2,"Emin = ",e9.2, " keV")') emin
    write(*,'(T2,"Emax = ",e9.2, " keV")') emax
    write(*,'(T2,"Nenergy = ", i4)') nenergy
    write(*,*) ''

    cnt = 0
    dlogE = (log10(emax) - log10(emin))/(nenergy - 1)
    !$OMP PARALLEL DO private(i, eb)
    do i=1, nenergy
        eb = 10.d0**(log10(emin) + (i-1)*dlogE)
        ebarr(i) = eb

        excit(:,:,i) = e_excit(eb, n_max, m_max)
        ioniz(:,i) = e_ioniz(eb, n_max)

        cnt = cnt + 1
        WRITE(*,'(f7.2,"%",a,$)') 100*cnt/real(nenergy),char(13)
    enddo
    !$OMP END PARALLEL DO

    call h5gcreate_f(id, "H_e", gid, error)

    dim1 = [1]
    dim2 = [n_max, nenergy]
    dim3 = [n_max, m_max, nenergy]

    call h5ltmake_dataset_int_f(gid, "nenergy", 0, dim1, [nenergy], error)
    call h5ltmake_dataset_int_f(gid, "n_max", 0, dim1, [n_max], error)
    call h5ltmake_dataset_int_f(gid, "m_max", 0, dim1, [m_max], error)
    call h5ltmake_dataset_double_f(gid, "dlogE", 0, dim1, [dlogE], error)
    call h5ltmake_dataset_double_f(gid, "emin", 0, dim1, [emin], error)
    call h5ltmake_dataset_double_f(gid, "emax", 0, dim1, [emax], error)


    dim1 = [nenergy]
    call h5ltmake_compressed_dataset_double_f(gid, "energy", 1, dim1, ebarr, error)
    call h5ltmake_compressed_dataset_double_f(gid, "ionization", 2, dim2, ioniz, error)
    call h5ltmake_compressed_dataset_double_f(gid, "excitation", 3, dim3, excit, error)

    call h5ltset_attribute_string_f(id, "H_e", "description", &
         "Cross sections for Hydrogen-Electron interactions", error)
    call h5ltset_attribute_string_f(gid, "nenergy", "description", &
         "Number of nucleon energy values", error)
    call h5ltset_attribute_string_f(gid, "n_max", "description", &
         "Number of initial energy levels", error)
    call h5ltset_attribute_string_f(gid, "m_max", "description", &
         "Number of final energy levels", error)
    call h5ltset_attribute_string_f(gid, "energy", "description", &
         "Nucleon energy values", error)
    call h5ltset_attribute_string_f(gid, "energy", "units", "keV/amu", error)
    call h5ltset_attribute_string_f(gid, "dlogE", "description", &
         "Energy spacing in log-10", error)
    call h5ltset_attribute_string_f(gid, "dlogE", "units", "log10(keV/amu)", error)
    call h5ltset_attribute_string_f(gid, "emin","description", &
         "Minimum Energy", error)
    call h5ltset_attribute_string_f(gid, "emin", "units", "keV/amu", error)
    call h5ltset_attribute_string_f(gid, "emax","description", &
         "Maximum Energy", error)
    call h5ltset_attribute_string_f(gid, "emax", "units", "keV/amu", error)

    call h5ltset_attribute_string_f(gid, "excitation", "description", &
         "n/m resolved excitation cross sections: excitation(n,m,energy)", error)
    call h5ltset_attribute_string_f(gid, "excitation", "units", "cm^2", error)
    call h5ltset_attribute_string_f(gid, "excitation", "reaction", &
         "e- + H(n) -> e- + H(m), m > n", error)

    call h5ltset_attribute_string_f(gid, "ionization", "description", &
         "n resolved ionization cross sections: ionization(n,energy)", error)
    call h5ltset_attribute_string_f(gid, "ionization", "units", "cm^2", error)
    call h5ltset_attribute_string_f(gid, "ionization", "reaction", &
         "e- + H(n) -> e- + H(+) + e-", error)

    call h5gclose_f(gid, error)

    deallocate(ebarr, ioniz, excit)

end subroutine write_bb_H_e