read_chords Subroutine

public subroutine read_chords()

Reads the spectral geometry and stores the quantities in spec_chords

Arguments

None

Calls

proc~~read_chords~~CallsGraph proc~read_chords read_chords h5ltread_dataset_string_f h5ltread_dataset_string_f proc~read_chords->h5ltread_dataset_string_f h5close_f h5close_f proc~read_chords->h5close_f h5gopen_f h5gopen_f proc~read_chords->h5gopen_f h5ltpath_valid_f h5ltpath_valid_f proc~read_chords->h5ltpath_valid_f h5gclose_f h5gclose_f proc~read_chords->h5gclose_f proc~h5ltread_dataset_int_scalar_f h5ltread_dataset_int_scalar_f proc~read_chords->proc~h5ltread_dataset_int_scalar_f h5open_f h5open_f proc~read_chords->h5open_f h5ltread_dataset_double_f h5ltread_dataset_double_f proc~read_chords->h5ltread_dataset_double_f h5fclose_f h5fclose_f proc~read_chords->h5fclose_f h5fopen_f h5fopen_f proc~read_chords->h5fopen_f h5ltread_dataset_int_f h5ltread_dataset_int_f proc~h5ltread_dataset_int_scalar_f->h5ltread_dataset_int_f

Called by

proc~~read_chords~~CalledByGraph proc~read_chords read_chords program~fidasim fidasim program~fidasim->proc~read_chords

Contents

Source Code


Source Code

subroutine read_chords
    !+ Reads the spectral geometry and stores the quantities in [[libfida:spec_chords]]
    integer(HID_T) :: fid, gid
    integer(HSIZE_T), dimension(2) :: dims
    logical :: path_valid

    real(Float64), dimension(:,:), allocatable :: lenses
    real(Float64), dimension(:,:), allocatable :: axes
    real(Float64) :: xyz_lens(3), xyz_axis(3)
    character(len=20) :: system = ''

    integer :: i
    integer :: error

    if(inputs%verbose.ge.1) then
        write(*,'(a)') '---- FIDA/BES settings ----'
    endif

    !!Initialize HDF5 interface
    call h5open_f(error)

    !!Open HDF5 file
    call h5fopen_f(inputs%geometry_file, H5F_ACC_RDONLY_F, fid, error)

    !!Check if SPEC group exists
    call h5ltpath_valid_f(fid, "/spec", .True., path_valid, error)
    if(.not.path_valid) then
        if(inputs%verbose.ge.1) then
            write(*,'(a)') 'FIDA/BES geometry is not in the geometry file'
            write(*,'(a)') 'Continuing without spectral diagnostics'
        endif
        inputs%calc_spec = 0
        inputs%tot_spectra=0
        inputs%calc_fida = 0
        inputs%calc_pfida = 0
        inputs%calc_bes = 0
        inputs%calc_dcx = 0
        inputs%calc_halo = 0
        inputs%calc_cold = 0
        inputs%calc_brems = 0
        inputs%calc_fida_wght = 0
        call h5fclose_f(fid, error)
        call h5close_f(error)
        return
    endif

    !!Open SPEC group
    call h5gopen_f(fid, "/spec", gid, error)

    call h5ltread_dataset_string_f(gid, "/spec/system", system, error)
    call h5ltread_dataset_int_scalar_f(gid, "/spec/nchan", spec_chords%nchan, error)

    allocate(lenses(3, spec_chords%nchan))
    allocate(axes(3, spec_chords%nchan))
    allocate(spec_chords%los(spec_chords%nchan))
    allocate(spec_chords%radius(spec_chords%nchan))

    dims = [3,spec_chords%nchan]
    call h5ltread_dataset_double_f(gid, "/spec/lens", lenses, dims, error)
    call h5ltread_dataset_double_f(gid, "/spec/axis", axes, dims, error)
    call h5ltread_dataset_double_f(gid, "/spec/spot_size", spec_chords%los%spot_size, dims(2:2), error)
    call h5ltread_dataset_double_f(gid, "/spec/sigma_pi", spec_chords%los%sigma_pi, dims(2:2), error)
    call h5ltread_dataset_double_f(gid, "/spec/radius", spec_chords%radius, dims(2:2), error)

    !!Close SPEC group
    call h5gclose_f(gid, error)

    !!Close file id
    call h5fclose_f(fid, error)

    !!Close HDF5 interface
    call h5close_f(error)

    chan_loop: do i=1,spec_chords%nchan
        call uvw_to_xyz(lenses(:,i), xyz_lens)
        xyz_axis = matmul(beam_grid%inv_basis, axes(:,i))
        spec_chords%los(i)%lens = xyz_lens
        spec_chords%los(i)%axis = xyz_axis
        spec_chords%los(i)%lens_uvw = lenses(:,i)
        spec_chords%los(i)%axis_uvw = axes(:,i)
    enddo chan_loop

    if(inputs%verbose.ge.1) then
        write(*,'(T2,"FIDA/BES System: ",a)') trim(adjustl(system))
        write(*,'(T2,"Number of channels: ",i5)') spec_chords%nchan
        write(*,*) ''
    endif

    deallocate(lenses,axes)

end subroutine read_chords