write_spectra Subroutine

public subroutine write_spectra()

Writes Spectra to a HDF5 file

Arguments

None

Calls

proc~~write_spectra~~CallsGraph proc~write_spectra write_spectra h5fcreate_f h5fcreate_f proc~write_spectra->h5fcreate_f h5open_f h5open_f proc~write_spectra->h5open_f h5ltset_attribute_string_f h5ltset_attribute_string_f proc~write_spectra->h5ltset_attribute_string_f h5close_f h5close_f proc~write_spectra->h5close_f h5ltmake_dataset_int_f h5ltmake_dataset_int_f proc~write_spectra->h5ltmake_dataset_int_f interface~h5ltmake_compressed_dataset_double_f h5ltmake_compressed_dataset_double_f proc~write_spectra->interface~h5ltmake_compressed_dataset_double_f h5fclose_f h5fclose_f proc~write_spectra->h5fclose_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 h5dwrite_f h5dwrite_f proc~h5ltmake_compressed_dataset_double_f_7->h5dwrite_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 h5ltmake_dataset_double_f h5ltmake_dataset_double_f proc~h5ltmake_compressed_dataset_double_f_7->h5ltmake_dataset_double_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->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->h5ltmake_dataset_double_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->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->h5ltmake_dataset_double_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->h5dwrite_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->h5ltmake_dataset_double_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->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->h5ltmake_dataset_double_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->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->h5ltmake_dataset_double_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->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->h5ltmake_dataset_double_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

Called by

proc~~write_spectra~~CalledByGraph proc~write_spectra write_spectra program~fidasim fidasim program~fidasim->proc~write_spectra

Contents

Source Code


Source Code

subroutine write_spectra
    !+ Writes [[libfida:spectra]] to a HDF5 file
    integer(HID_T) :: fid
    integer(HSIZE_T), dimension(3) :: dims
    integer(HSIZE_T), dimension(1) :: d
    integer :: error

    character(charlim) :: filename
    integer :: i
    real(Float64) :: factor
    real(Float64), dimension(:), allocatable :: lambda_arr

    allocate(lambda_arr(inputs%nlambda))
    do i=1,inputs%nlambda
        lambda_arr(i) = (i-0.5)*inputs%dlambda + inputs%lambdamin
    enddo

    !! convert [Ph/(s*wavel_bin*cm^2*all_directions)] to [Ph/(s*nm*sr*m^2)]!
    factor = 1.d0/(inputs%dlambda)/(4.d0*pi)*1.d4

    !! write to file
    filename=trim(adjustl(inputs%result_dir))//"/"//trim(adjustl(inputs%runid))//"_spectra.h5"

    !Open HDF5 interface
    call h5open_f(error)

    !Create file overwriting any existing file
    call h5fcreate_f(filename, H5F_ACC_TRUNC_F, fid, error)

    !Write variables
    d(1) = 1
    call h5ltmake_dataset_int_f(fid, "/nchan", 0, d, [spec_chords%nchan], error)
    call h5ltmake_dataset_int_f(fid, "/nlambda", 0, d, [inputs%nlambda], error)
    dims(1) = inputs%nlambda
    dims(2) = spec_chords%nchan
    dims(3) = particles%nclass
    call h5ltmake_compressed_dataset_double_f(fid, "/lambda", 1, dims(1:1), &
         lambda_arr, error)
    call h5ltmake_compressed_dataset_double_f(fid, "/radius", 1, dims(2:2), &
         spec_chords%radius, error)

    !Add attributes
    call h5ltset_attribute_string_f(fid,"/nchan", "description", &
         "Number of channels", error)
    call h5ltset_attribute_string_f(fid,"/nlambda", "description", &
         "Number of wavelengths", error)
    call h5ltset_attribute_string_f(fid,"/lambda","description", &
         "Wavelength array", error)
    call h5ltset_attribute_string_f(fid,"/lambda","units","nm", error)
    call h5ltset_attribute_string_f(fid,"/radius", "description", &
         "Line of sight radius at midplane or tangency point", error)
    call h5ltset_attribute_string_f(fid,"/radius","units","cm", error)

    if(inputs%calc_brems.ge.1) then
        spec%brems = factor*spec%brems
        !Write variables
        call h5ltmake_compressed_dataset_double_f(fid, "/brems", 2, &
             dims(1:2), spec%brems, error)
        !Add attributes
        call h5ltset_attribute_string_f(fid,"/brems","description", &
             "Visible Bremsstrahlung: brems(lambda,chan)", error)
        call h5ltset_attribute_string_f(fid,"/brems","units",&
             "Ph/(s*nm*sr*m^2)",error )
    endif

    if(inputs%calc_bes.ge.1) then
        spec%full  = factor*spec%full
        spec%half  = factor*spec%half
        spec%third = factor*spec%third
        !Write variables
        call h5ltmake_compressed_dataset_double_f(fid, "/full", 2, dims(1:2), &
             spec%full, error)
        call h5ltmake_compressed_dataset_double_f(fid, "/half", 2, dims(1:2), &
             spec%half, error)
        call h5ltmake_compressed_dataset_double_f(fid, "/third", 2, dims(1:2),&
             spec%third, error)
        !Add attributes
        call h5ltset_attribute_string_f(fid,"/full","description", &
             "Full energy component of the beam emmision: full(lambda,chan)", error)
        call h5ltset_attribute_string_f(fid,"/full","units","Ph/(s*nm*sr*m^2)",error )
        call h5ltset_attribute_string_f(fid,"/half","description", &
             "Half energy component of the beam emmision: half(lambda,chan)", error)
        call h5ltset_attribute_string_f(fid,"/half","units","Ph/(s*nm*sr*m^2)",error )
        call h5ltset_attribute_string_f(fid,"/third","description", &
             "Third energy component of the beam emmision: third(lambda,chan)", error)
        call h5ltset_attribute_string_f(fid,"/third","units","Ph/(s*nm*sr*m^2)",error )
    endif

    if(inputs%calc_dcx.ge.1) then
        spec%dcx   = factor*spec%dcx
        call h5ltmake_compressed_dataset_double_f(fid, "/dcx", 2, dims(1:2), &
             spec%dcx, error)
        call h5ltset_attribute_string_f(fid,"/dcx","description", &
             "Direct Charge Exchange (DCX) emission: dcx(lambda,chan)", error)
        call h5ltset_attribute_string_f(fid,"/dcx","units","Ph/(s*nm*sr*m^2)",error )
    endif

    if(inputs%calc_halo.ge.1) then
        spec%halo  = factor*spec%halo
        call h5ltmake_compressed_dataset_double_f(fid, "/halo", 2, dims(1:2), &
             spec%halo, error)
        call h5ltset_attribute_string_f(fid,"/halo","description", &
             "Halo component of the beam emmision: halo(lambda,chan)", error)
        call h5ltset_attribute_string_f(fid,"/halo","units","Ph/(s*nm*sr*m^2)",error )
    endif

    if(inputs%calc_cold.ge.1) then
        spec%cold  = factor*spec%cold
        call h5ltmake_compressed_dataset_double_f(fid, "/cold", 2, dims(1:2), &
             spec%cold, error)
        call h5ltset_attribute_string_f(fid,"/cold","description", &
             "Cold D-alpha emission: cold(lambda,chan)", error)
        call h5ltset_attribute_string_f(fid,"/cold","units","Ph/(s*nm*sr*m^2)",error )
    endif

    if(inputs%calc_fida.ge.1) then
        spec%fida  = factor*spec%fida
        !Write variables
        if(particles%nclass.le.1) then
            call h5ltmake_compressed_dataset_double_f(fid, "/fida", 2, &
                 dims(1:2), spec%fida(:,:,1), error)
            !Add attributes
            call h5ltset_attribute_string_f(fid,"/fida","description", &
                 "Active Fast-ion D-alpha (FIDA) emmision: fida(lambda,chan)", error)
        else
            call h5ltmake_dataset_int_f(fid,"/nclass", 0, d, [particles%nclass], error)
            call h5ltmake_compressed_dataset_double_f(fid, "/fida", 3, &
                 dims, spec%fida, error)
            !Add attributes
            call h5ltset_attribute_string_f(fid,"/fida","description", &
                 "Active Fast-ion D-alpha (FIDA) emmision: fida(lambda,chan,class)", error)
       endif
        call h5ltset_attribute_string_f(fid,"/fida","units","Ph/(s*nm*sr*m^2)",error )
    endif

    if(inputs%calc_pfida.ge.1) then
        spec%pfida = factor*spec%pfida
        !Write variables
        if(particles%nclass.le.1) then
            call h5ltmake_compressed_dataset_double_f(fid, "/pfida", 2, &
                 dims(1:2), spec%pfida(:,:,1), error)
            !Add attributes
            call h5ltset_attribute_string_f(fid,"/pfida","description", &
                 "Passive Fast-ion D-alpha (p-FIDA) emmision: pfida(lambda,chan)", error)
        else
            if(inputs%calc_fida.le.0) then
                call h5ltmake_dataset_int_f(fid,"/nclass", 0, d, [particles%nclass], error)
            endif
            call h5ltmake_compressed_dataset_double_f(fid, "/pfida", 3, &
                 dims, spec%pfida, error)
            !Add attributes
            call h5ltset_attribute_string_f(fid,"/pfida","description", &
                 "Passive Fast-ion D-alpha (p-FIDA) emmision: pfida(lambda,chan,class)", error)
       endif
        call h5ltset_attribute_string_f(fid,"/pfida","units","Ph/(s*nm*sr*m^2)",error )
    endif

    call h5ltset_attribute_string_f(fid, "/", "version", version, error)
    call h5ltset_attribute_string_f(fid,"/","description",&
         "Spectra calculated by FIDASIM", error)

    !Close file
    call h5fclose_f(fid, error)

    !Close HDF5 interface
    call h5close_f(error)

    if(inputs%verbose.ge.1) then
        write(*,'(T4,a,a)') 'Spectra written to: ', trim(filename)
    endif

end subroutine write_spectra