write_dcx Subroutine

public subroutine write_dcx()

Writes the direct charge exchange (DCX) neutrals and spectra to a HDF5 file

Arguments

None

Calls

proc~~write_dcx~~CallsGraph proc~write_dcx write_dcx h5fcreate_f h5fcreate_f proc~write_dcx->h5fcreate_f proc~write_beam_grid write_beam_grid proc~write_dcx->proc~write_beam_grid h5open_f h5open_f proc~write_dcx->h5open_f h5ltset_attribute_string_f h5ltset_attribute_string_f proc~write_dcx->h5ltset_attribute_string_f h5close_f h5close_f proc~write_dcx->h5close_f h5ltmake_dataset_int_f h5ltmake_dataset_int_f proc~write_dcx->h5ltmake_dataset_int_f interface~h5ltmake_compressed_dataset_double_f h5ltmake_compressed_dataset_double_f proc~write_dcx->interface~h5ltmake_compressed_dataset_double_f h5fclose_f h5fclose_f proc~write_dcx->h5fclose_f proc~write_beam_grid->h5ltset_attribute_string_f proc~write_beam_grid->interface~h5ltmake_compressed_dataset_double_f h5gclose_f h5gclose_f proc~write_beam_grid->h5gclose_f h5gcreate_f h5gcreate_f proc~write_beam_grid->h5gcreate_f proc~xyz_to_uvw xyz_to_uvw proc~write_beam_grid->proc~xyz_to_uvw 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

Contents

Source Code


Source Code

subroutine write_dcx
    !+ Writes the direct charge exchange (DCX) neutrals and spectra to a HDF5 file
    integer(HID_T) :: fid
    integer(HSIZE_T), dimension(4) :: dims
    integer(HSIZE_T), dimension(1) :: d
    integer :: error

    character(charlim) :: filename
    character(15) :: spec_str
    integer :: i
    real(Float64), dimension(:),   allocatable :: lambda_arr
    real(Float64), dimension(:,:), allocatable :: dcx_spec

    filename=trim(adjustl(inputs%result_dir))//"/"//trim(adjustl(inputs%runid))//"_dcx.h5"
    spec_str = ""
    if(inputs%calc_spec.ge.1) then
        spec_str = " spectra and"
        allocate(lambda_arr(inputs%nlambda))
        do i=1,inputs%nlambda
            lambda_arr(i) = (i-0.5)*inputs%dlambda + inputs%lambdamin ! [nm]
        enddo

        allocate(dcx_spec(inputs%nlambda,spec_chords%nchan))
        !! convert [Ph/(s*wavel_bin*cm^2*all_directions)] to [Ph/(s*nm*sr*m^2)]!
        dcx_spec = spec%bes(:,:,halo_type)/(inputs%dlambda)/(4.d0*pi)*1.d4
    endif

    !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
    call write_beam_grid(fid, error)

    d(1) =1
    call h5ltmake_dataset_int_f(fid,"/nlevel", 0, d, [nlevs], error)
    dims = [nlevs, beam_grid%nx, beam_grid%ny, beam_grid%nz ]
    call h5ltmake_compressed_dataset_double_f(fid, "/dens", 4, dims, &
         neut%dens(:,halo_type,:,:,:), error)

    if(inputs%calc_spec.ge.1) then
        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
        call h5ltmake_compressed_dataset_double_f(fid, "/spec", 2, dims(1:2), &
             dcx_spec, error)
        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)
    endif

    !Add attributes
    call h5ltset_attribute_string_f(fid,"/nlevel","description", &
         "Number of atomic energy levels", error)
    call h5ltset_attribute_string_f(fid,"/dens", "description", &
         "Direct Charge Exchange (DCX) neutral density: dcx(level,x,y,z)", error)
    call h5ltset_attribute_string_f(fid,"/dens","units","neutrals*cm^-3", error)

    if(inputs%calc_spec.ge.1) then
        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,"/spec","description", &
             "Direct Charge Exchange (DCX) beam emission: spec(lambda, chan)", error)
        call h5ltset_attribute_string_f(fid,"/spec","units","Ph/(s*nm*sr*m^2)",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)
    endif

    call h5ltset_attribute_string_f(fid, "/", "version", version, error)
    call h5ltset_attribute_string_f(fid,"/","description", &
         "Direct Charge Exchange (DCX)"//trim(spec_str)//" neutral density calculated by FIDASIM", error)

    !Close file
    call h5fclose_f(fid, error)

    !Close HDF5 interface
    call h5close_f(error)

    if(inputs%calc_spec.ge.1) then
        deallocate(dcx_spec,lambda_arr)
    endif

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

end subroutine write_dcx