write_bb_H_Aq Subroutine

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

Write Hydrogen-Impurity 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_aq~~CallsGraph proc~write_bb_h_aq write_bb_H_Aq proc~aq_cx Aq_cx proc~write_bb_h_aq->proc~aq_cx proc~aq_excit Aq_excit proc~write_bb_h_aq->proc~aq_excit h5ltmake_dataset_int_f h5ltmake_dataset_int_f proc~write_bb_h_aq->h5ltmake_dataset_int_f h5ltmake_dataset_double_f h5ltmake_dataset_double_f proc~write_bb_h_aq->h5ltmake_dataset_double_f h5gclose_f h5gclose_f proc~write_bb_h_aq->h5gclose_f h5ltset_attribute_string_f h5ltset_attribute_string_f proc~write_bb_h_aq->h5ltset_attribute_string_f interface~parallel_sum parallel_sum proc~write_bb_h_aq->interface~parallel_sum proc~aq_ioniz Aq_ioniz proc~write_bb_h_aq->proc~aq_ioniz interface~h5ltmake_compressed_dataset_double_f h5ltmake_compressed_dataset_double_f proc~write_bb_h_aq->interface~h5ltmake_compressed_dataset_double_f proc~aq_cx_n Aq_cx_n proc~aq_cx->proc~aq_cx_n proc~aq_excit_n Aq_excit_n proc~aq_excit->proc~aq_excit_n proc~parallel_sum_i1 parallel_sum_i1 interface~parallel_sum->proc~parallel_sum_i1 proc~parallel_sum_i0 parallel_sum_i0 interface~parallel_sum->proc~parallel_sum_i0 proc~parallel_sum_i2 parallel_sum_i2 interface~parallel_sum->proc~parallel_sum_i2 proc~parallel_sum_d4 parallel_sum_d4 interface~parallel_sum->proc~parallel_sum_d4 proc~parallel_sum_d5 parallel_sum_d5 interface~parallel_sum->proc~parallel_sum_d5 proc~parallel_sum_d2 parallel_sum_d2 interface~parallel_sum->proc~parallel_sum_d2 proc~parallel_sum_d3 parallel_sum_d3 interface~parallel_sum->proc~parallel_sum_d3 proc~parallel_sum_d0 parallel_sum_d0 interface~parallel_sum->proc~parallel_sum_d0 proc~parallel_sum_d1 parallel_sum_d1 interface~parallel_sum->proc~parallel_sum_d1 proc~aq_ioniz_n Aq_ioniz_n proc~aq_ioniz->proc~aq_ioniz_n 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_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_1 h5ltmake_compressed_dataset_double_f_1 interface~h5ltmake_compressed_dataset_double_f->proc~h5ltmake_compressed_dataset_double_f_1 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_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_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_3 h5ltmake_compressed_dataset_double_f_3 interface~h5ltmake_compressed_dataset_double_f->proc~h5ltmake_compressed_dataset_double_f_3 proc~aq_cx_n_adas Aq_cx_n_adas proc~aq_cx_n->proc~aq_cx_n_adas proc~aq_cx_n_janev Aq_cx_n_janev proc~aq_cx_n->proc~aq_cx_n_janev mpi_allreduce mpi_allreduce proc~parallel_sum_i1->mpi_allreduce proc~parallel_sum_i0->mpi_allreduce proc~parallel_sum_i2->mpi_allreduce proc~h5ltmake_compressed_dataset_double_f_4->h5ltmake_dataset_double_f h5dwrite_f h5dwrite_f proc~h5ltmake_compressed_dataset_double_f_4->h5dwrite_f h5dclose_f h5dclose_f proc~h5ltmake_compressed_dataset_double_f_4->h5dclose_f h5sclose_f h5sclose_f proc~h5ltmake_compressed_dataset_double_f_4->h5sclose_f h5pset_chunk_f h5pset_chunk_f proc~h5ltmake_compressed_dataset_double_f_4->h5pset_chunk_f h5pcreate_f h5pcreate_f proc~h5ltmake_compressed_dataset_double_f_4->h5pcreate_f proc~chunk_size chunk_size proc~h5ltmake_compressed_dataset_double_f_4->proc~chunk_size h5pset_deflate_f h5pset_deflate_f proc~h5ltmake_compressed_dataset_double_f_4->h5pset_deflate_f h5pset_shuffle_f h5pset_shuffle_f proc~h5ltmake_compressed_dataset_double_f_4->h5pset_shuffle_f h5pclose_f h5pclose_f proc~h5ltmake_compressed_dataset_double_f_4->h5pclose_f h5screate_simple_f h5screate_simple_f proc~h5ltmake_compressed_dataset_double_f_4->h5screate_simple_f h5dcreate_f h5dcreate_f proc~h5ltmake_compressed_dataset_double_f_4->h5dcreate_f proc~h5ltmake_compressed_dataset_double_f_5->h5ltmake_dataset_double_f proc~h5ltmake_compressed_dataset_double_f_5->h5dwrite_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->h5pclose_f proc~h5ltmake_compressed_dataset_double_f_5->h5screate_simple_f proc~h5ltmake_compressed_dataset_double_f_5->h5dcreate_f proc~h5ltmake_compressed_dataset_double_f_1->h5ltmake_dataset_double_f proc~h5ltmake_compressed_dataset_double_f_1->h5dwrite_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->h5pclose_f proc~h5ltmake_compressed_dataset_double_f_1->h5screate_simple_f proc~h5ltmake_compressed_dataset_double_f_1->h5dcreate_f proc~aq_excit_2_janev Aq_excit_2_janev proc~aq_excit_n->proc~aq_excit_2_janev proc~aq_excit_3_janev Aq_excit_3_janev proc~aq_excit_n->proc~aq_excit_3_janev proc~aq_excit_n_janev Aq_excit_n_janev proc~aq_excit_n->proc~aq_excit_n_janev proc~aq_excit_1_janev Aq_excit_1_janev proc~aq_excit_n->proc~aq_excit_1_janev proc~b5_ioniz_1_janev B5_ioniz_1_janev proc~aq_ioniz_n->proc~b5_ioniz_1_janev proc~aq_ioniz_n_janev Aq_ioniz_n_janev proc~aq_ioniz_n->proc~aq_ioniz_n_janev proc~c6_ioniz_1_janev C6_ioniz_1_janev proc~aq_ioniz_n->proc~c6_ioniz_1_janev proc~parallel_sum_d4->mpi_allreduce proc~parallel_sum_d5->mpi_allreduce proc~parallel_sum_d2->mpi_allreduce proc~parallel_sum_d3->mpi_allreduce proc~parallel_sum_d0->mpi_allreduce proc~parallel_sum_d1->mpi_allreduce proc~h5ltmake_compressed_dataset_double_f_6->h5ltmake_dataset_double_f proc~h5ltmake_compressed_dataset_double_f_6->h5dwrite_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->h5pclose_f proc~h5ltmake_compressed_dataset_double_f_6->h5screate_simple_f proc~h5ltmake_compressed_dataset_double_f_6->h5dcreate_f proc~h5ltmake_compressed_dataset_double_f_7->h5ltmake_dataset_double_f proc~h5ltmake_compressed_dataset_double_f_7->h5dwrite_f proc~h5ltmake_compressed_dataset_double_f_7->h5dclose_f proc~h5ltmake_compressed_dataset_double_f_7->h5sclose_f proc~h5ltmake_compressed_dataset_double_f_7->h5pset_chunk_f proc~h5ltmake_compressed_dataset_double_f_7->h5pcreate_f proc~h5ltmake_compressed_dataset_double_f_7->proc~chunk_size proc~h5ltmake_compressed_dataset_double_f_7->h5pset_deflate_f proc~h5ltmake_compressed_dataset_double_f_7->h5pset_shuffle_f proc~h5ltmake_compressed_dataset_double_f_7->h5pclose_f proc~h5ltmake_compressed_dataset_double_f_7->h5screate_simple_f proc~h5ltmake_compressed_dataset_double_f_7->h5dcreate_f proc~h5ltmake_compressed_dataset_double_f_2->h5ltmake_dataset_double_f proc~h5ltmake_compressed_dataset_double_f_2->h5dwrite_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->h5pclose_f proc~h5ltmake_compressed_dataset_double_f_2->h5screate_simple_f proc~h5ltmake_compressed_dataset_double_f_2->h5dcreate_f proc~h5ltmake_compressed_dataset_double_f_3->h5ltmake_dataset_double_f proc~h5ltmake_compressed_dataset_double_f_3->h5dwrite_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->h5pclose_f proc~h5ltmake_compressed_dataset_double_f_3->h5screate_simple_f proc~h5ltmake_compressed_dataset_double_f_3->h5dcreate_f proc~aq_excit_2_4_janev Aq_excit_2_4_janev proc~aq_excit_2_janev->proc~aq_excit_2_4_janev proc~aq_excit_2_10_janev Aq_excit_2_10_janev proc~aq_excit_2_janev->proc~aq_excit_2_10_janev proc~aq_excit_2_6_janev Aq_excit_2_6_janev proc~aq_excit_2_janev->proc~aq_excit_2_6_janev proc~aq_excit_2_9_janev Aq_excit_2_9_janev proc~aq_excit_2_janev->proc~aq_excit_2_9_janev proc~aq_excit_2_7_janev Aq_excit_2_7_janev proc~aq_excit_2_janev->proc~aq_excit_2_7_janev proc~aq_excit_2_5_janev Aq_excit_2_5_janev proc~aq_excit_2_janev->proc~aq_excit_2_5_janev proc~aq_excit_2_3_janev Aq_excit_2_3_janev proc~aq_excit_2_janev->proc~aq_excit_2_3_janev proc~aq_excit_2_8_janev Aq_excit_2_8_janev proc~aq_excit_2_janev->proc~aq_excit_2_8_janev proc~b5_cx_1_adas B5_cx_1_adas proc~aq_cx_n_adas->proc~b5_cx_1_adas proc~c6_cx_2_adas C6_cx_2_adas proc~aq_cx_n_adas->proc~c6_cx_2_adas proc~c6_cx_1_adas C6_cx_1_adas proc~aq_cx_n_adas->proc~c6_cx_1_adas proc~c6_cx_3_adas C6_cx_3_adas proc~aq_cx_n_adas->proc~c6_cx_3_adas proc~b5_cx_2_adas B5_cx_2_adas proc~aq_cx_n_adas->proc~b5_cx_2_adas proc~c6_cx_1_janev C6_cx_1_janev proc~aq_cx_n_janev->proc~c6_cx_1_janev proc~b5_cx_1_janev B5_cx_1_janev proc~aq_cx_n_janev->proc~b5_cx_1_janev proc~aq_excit_3_7_janev Aq_excit_3_7_janev proc~aq_excit_3_janev->proc~aq_excit_3_7_janev proc~aq_excit_3_6_janev Aq_excit_3_6_janev proc~aq_excit_3_janev->proc~aq_excit_3_6_janev proc~aq_excit_3_10_janev Aq_excit_3_10_janev proc~aq_excit_3_janev->proc~aq_excit_3_10_janev proc~aq_excit_3_9_janev Aq_excit_3_9_janev proc~aq_excit_3_janev->proc~aq_excit_3_9_janev proc~aq_excit_3_8_janev Aq_excit_3_8_janev proc~aq_excit_3_janev->proc~aq_excit_3_8_janev proc~aq_excit_3_4_janev Aq_excit_3_4_janev proc~aq_excit_3_janev->proc~aq_excit_3_4_janev proc~aq_excit_3_5_janev Aq_excit_3_5_janev proc~aq_excit_3_janev->proc~aq_excit_3_5_janev proc~aq_excit_1_2_janev Aq_excit_1_2_janev proc~aq_excit_1_janev->proc~aq_excit_1_2_janev proc~aq_excit_1_3_janev Aq_excit_1_3_janev proc~aq_excit_1_janev->proc~aq_excit_1_3_janev proc~aq_excit_1_5_janev Aq_excit_1_5_janev proc~aq_excit_1_janev->proc~aq_excit_1_5_janev proc~aq_excit_1_4_janev Aq_excit_1_4_janev proc~aq_excit_1_janev->proc~aq_excit_1_4_janev proc~aq_excit_3_7_janev->proc~aq_excit_3_6_janev proc~aq_excit_3_10_janev->proc~aq_excit_3_6_janev proc~aq_excit_2_10_janev->proc~aq_excit_2_5_janev proc~aq_excit_2_6_janev->proc~aq_excit_2_5_janev proc~aq_excit_3_9_janev->proc~aq_excit_3_6_janev proc~aq_excit_2_9_janev->proc~aq_excit_2_5_janev proc~aq_excit_3_8_janev->proc~aq_excit_3_6_janev proc~aq_excit_2_7_janev->proc~aq_excit_2_5_janev proc~aq_excit_2_8_janev->proc~aq_excit_2_5_janev

Called by

proc~~write_bb_h_aq~~CalledByGraph proc~write_bb_h_aq write_bb_H_Aq program~generate_tables generate_tables program~generate_tables->proc~write_bb_h_aq

Contents

Source Code


Source Code

subroutine write_bb_H_Aq(id, namelist_file, n_max, m_max)
    !+ Write Hydrogen-Impurity 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

    logical :: calculate
    integer :: q(10) = 0
    real(Float64) :: emin
    real(Float64) :: emax
    integer :: nenergy

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

    real(Float64), dimension(:,:), allocatable :: cx, 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

    character(len=10) :: aname
    character(len=5) :: asym
    integer :: i, iq, cnt, error
    logical :: exis

    NAMELIST /H_Aq_cross/ calculate, q, nenergy, emin, emax

    nenergy = 200; emin = 1.d-3 ; emax = 8.d2
    calculate = .True.; q(1) = 6

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

    if(.not.calculate) return

    q_loop: do iq=1, size(q)
        if(q(iq).eq.0) cycle q_loop

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

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

        select case (q(iq))
            case (5)
                aname = "Boron"
                asym = "H_B5"
            case (6)
                aname = "Carbon"
                asym = "H_C6"
            case DEFAULT
                write(aname,'("Impurity-",i1)') q(iq)
                write(asym,'("H_A",i1)') q(iq)
        end select

        if(verbose) then
            write(*,'(a)') "---- H-"//trim(adjustl(aname))//" cross sections settings ----"
            write(*,'(T2,"q = ", i2)') q(iq)
            write(*,'(T2,"Emin = ",e9.2, " keV")') emin
            write(*,'(T2,"Emax = ",e9.2, " keV")') emax
            write(*,'(T2,"Nenergy = ", i4)') nenergy
            write(*,*) ''
        endif

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

            cx(:,i) = Aq_cx(eb, q(iq), n_max)
            ioniz(:,i) = Aq_ioniz(eb, q(iq), n_max)
            excit(:,:,i) = Aq_excit(eb, q(iq), n_max, m_max)
            cnt = cnt + 1
            if(verbose) WRITE(*,'(f7.2,"%",a,$)') 100*cnt*istep/real(nenergy),char(13)
        enddo

#ifdef _MPI
        call parallel_sum(ebarr)
        call parallel_sum(cx)
        call parallel_sum(excit)
        call parallel_sum(ioniz)
#endif

        if(verbose) then
            call h5gcreate_f(id, trim(adjustl(asym)), 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, "cx", 2, dim2, cx, 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, trim(adjustl(asym)), "description", &
                 "Cross sections for Hydrogen-"//trim(adjustl(aname))//" 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, "cx", "description", &
                 "n resolved charge exchange / electron capture cross sections: cx(n,energy)", error)
            call h5ltset_attribute_string_f(gid, "cx", "units", "cm^2", error)
            call h5ltset_attribute_string_f(gid, "cx", "reaction", &
                 "A(q+) + H(n) -> A((q-1)+) + H(+)", 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", &
                 "A(q+) + H(n) -> A(q+) + 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", &
                 "A(q+) + H(n) -> A(q+) + H(+) + e-", error)

            call h5gclose_f(gid, error)
        endif

        deallocate(ebarr, ioniz, cx, excit)
    enddo q_loop

end subroutine write_bb_H_Aq