write_bb_D_T Subroutine

public subroutine write_bb_D_T(id, namelist_file)

Write Deuterium-Tritium 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


Calls

proc~~write_bb_d_t~~CallsGraph proc~write_bb_d_t write_bb_D_T interface~h5ltmake_compressed_dataset_double_f h5ltmake_compressed_dataset_double_f proc~write_bb_d_t->interface~h5ltmake_compressed_dataset_double_f h5ltmake_dataset_int_f h5ltmake_dataset_int_f proc~write_bb_d_t->h5ltmake_dataset_int_f h5ltmake_dataset_double_f h5ltmake_dataset_double_f proc~write_bb_d_t->h5ltmake_dataset_double_f h5gclose_f h5gclose_f proc~write_bb_d_t->h5gclose_f h5ltset_attribute_string_f h5ltset_attribute_string_f proc~write_bb_d_t->h5ltset_attribute_string_f proc~d_t_fusion d_t_fusion proc~write_bb_d_t->proc~d_t_fusion interface~parallel_sum parallel_sum proc~write_bb_d_t->interface~parallel_sum h5gcreate_f h5gcreate_f proc~write_bb_d_t->h5gcreate_f 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~parallel_sum_d0 parallel_sum_d0 interface~parallel_sum->proc~parallel_sum_d0 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_d5 parallel_sum_d5 interface~parallel_sum->proc~parallel_sum_d5 proc~parallel_sum_d4 parallel_sum_d4 interface~parallel_sum->proc~parallel_sum_d4 proc~parallel_sum_d1 parallel_sum_d1 interface~parallel_sum->proc~parallel_sum_d1 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 mpi_allreduce mpi_allreduce proc~parallel_sum_d0->mpi_allreduce proc~parallel_sum_i1->mpi_allreduce proc~parallel_sum_i0->mpi_allreduce proc~parallel_sum_i2->mpi_allreduce proc~parallel_sum_d5->mpi_allreduce proc~parallel_sum_d4->mpi_allreduce proc~h5ltmake_compressed_dataset_double_f_7->h5ltmake_dataset_double_f h5sclose_f h5sclose_f proc~h5ltmake_compressed_dataset_double_f_7->h5sclose_f h5dclose_f h5dclose_f proc~h5ltmake_compressed_dataset_double_f_7->h5dclose_f h5dwrite_f h5dwrite_f proc~h5ltmake_compressed_dataset_double_f_7->h5dwrite_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 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->h5sclose_f proc~h5ltmake_compressed_dataset_double_f_6->h5dclose_f proc~h5ltmake_compressed_dataset_double_f_6->h5dwrite_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->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->h5sclose_f proc~h5ltmake_compressed_dataset_double_f_5->h5dclose_f proc~h5ltmake_compressed_dataset_double_f_5->h5dwrite_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->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->h5sclose_f proc~h5ltmake_compressed_dataset_double_f_4->h5dclose_f proc~h5ltmake_compressed_dataset_double_f_4->h5dwrite_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->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->h5sclose_f proc~h5ltmake_compressed_dataset_double_f_3->h5dclose_f proc~h5ltmake_compressed_dataset_double_f_3->h5dwrite_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->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->h5sclose_f proc~h5ltmake_compressed_dataset_double_f_2->h5dclose_f proc~h5ltmake_compressed_dataset_double_f_2->h5dwrite_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->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->h5sclose_f proc~h5ltmake_compressed_dataset_double_f_1->h5dclose_f proc~h5ltmake_compressed_dataset_double_f_1->h5dwrite_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->h5pclose_f proc~h5ltmake_compressed_dataset_double_f_1->h5dcreate_f proc~parallel_sum_d1->mpi_allreduce proc~parallel_sum_d2->mpi_allreduce proc~parallel_sum_d3->mpi_allreduce

Contents

Source Code


Source Code

subroutine write_bb_D_T(id, namelist_file)
    !+ Write Deuterium-Tritium 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

    logical :: calculate
    integer :: nbranch = 1
    real(Float64) :: emin
    real(Float64) :: emax
    integer :: nenergy

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

    real(Float64), dimension(:,:), allocatable :: fusion

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

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

    NAMELIST /D_T_cross/ calculate, nenergy, emin, emax

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

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

    if(.not.calculate) return

    allocate(ebarr(nenergy))
    allocate(fusion(nenergy,nbranch))

    ebarr = 0.d0
    fusion = 0.d0

    if(verbose) then
        write(*,'(a)') "---- D-T cross sections settings ----"
        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

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

#ifdef _MPI
    call parallel_sum(ebarr)
    call parallel_sum(fusion)
#endif

    if(verbose) then
        call h5gcreate_f(id, "D_T", gid, error)

        dim1 = [1]

        call h5ltmake_dataset_int_f(gid, "nbranch", 0, dim1, [nbranch], error)
        call h5ltmake_dataset_int_f(gid, "nenergy", 0, dim1, [nenergy], 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]
        dim2 = [nenergy, nbranch]
        call h5ltmake_compressed_dataset_double_f(gid, "energy", 1, dim1, ebarr, error)
        call h5ltmake_compressed_dataset_double_f(gid, "fusion", 2, dim2, fusion, error)

        call h5ltset_attribute_string_f(id, "D_T", "description", &
             "Cross sections for Deuterium-Tritium interactions", error)
        call h5ltset_attribute_string_f(gid, "nbranch", "description", &
             "Number of reaction branches", error)
        call h5ltset_attribute_string_f(gid, "nenergy", "description", &
             "Number of energy values", error)
        call h5ltset_attribute_string_f(gid, "energy", "description", &
             "Deuterium energy values", error)
        call h5ltset_attribute_string_f(gid, "energy", "units", "keV", 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)", error)
        call h5ltset_attribute_string_f(gid, "emin","description", &
             "Minimum energy", error)
        call h5ltset_attribute_string_f(gid, "emin", "units", "keV", error)
        call h5ltset_attribute_string_f(gid, "emax","description", &
             "Maximum energy", error)
        call h5ltset_attribute_string_f(gid, "emax", "units", "keV", error)

        call h5ltset_attribute_string_f(gid, "fusion", "description", &
             "Total cross sections for D-T nuclear reactions: fusion(deuterium energy, branch)", error)
        call h5ltset_attribute_string_f(gid, "fusion", "units", "cm^2", error)
        call h5ltset_attribute_string_f(gid, "fusion", "reaction", &
             "D + T -> He4(3.5 MeV) + n(14.1 MeV)", error)

        call h5gclose_f(gid, error)
    endif

    deallocate(ebarr, fusion)

end subroutine write_bb_D_T