write_npa_weights Subroutine

public subroutine write_npa_weights()

Writes nweight to a HDF5 file

Arguments

None

Calls

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

Called by

proc~~write_npa_weights~~CalledByGraph proc~write_npa_weights write_npa_weights proc~npa_weights npa_weights proc~npa_weights->proc~write_npa_weights program~fidasim fidasim program~fidasim->proc~npa_weights

Contents

Source Code


Source Code

subroutine write_npa_weights
    !+ Writes [[libfida:nweight]] to a HDF5 file
    character(charlim) :: filename
    integer :: i
    real(Float64), dimension(:), allocatable :: ebarr,ptcharr

    !! HDF5 variables
    integer(HID_T) :: fid
    integer(HSIZE_T), dimension(5) :: dim5
    integer(HSIZE_T), dimension(3) :: dim3
    integer(HSIZE_T), dimension(2) :: dim2
    integer(HSIZE_T), dimension(1) :: d
    integer :: error

    !! define energy - array
    allocate(ebarr(inputs%ne_wght))
    do i=1,inputs%ne_wght
        ebarr(i)=real(i-0.5)*inputs%emax_wght/real(inputs%ne_wght)
    enddo

    !! define pitch - array
    allocate(ptcharr(inputs%np_wght))
    do i=1,inputs%np_wght
        ptcharr(i)=real(i-0.5)*2./real(inputs%np_wght)-1.
    enddo

    filename=trim(adjustl(inputs%result_dir))//"/"//trim(adjustl(inputs%runid))//"_npa_weights.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
    dim2 = [inputs%ne_wght, npa_chords%nchan]
    dim3 = [inputs%ne_wght, inputs%np_wght, npa_chords%nchan]
    dim5 = [inputs%ne_wght, beam_grid%nx, beam_grid%ny, beam_grid%nz, npa_chords%nchan]
    call h5ltmake_dataset_int_f(fid, "/nchan", 0, d, [npa_chords%nchan], error)
    call h5ltmake_dataset_int_f(fid, "/nenergy", 0, d, [inputs%ne_wght], error)
    call h5ltmake_dataset_int_f(fid, "/npitch", 0, d, [inputs%np_wght], error)
    call h5ltmake_compressed_dataset_double_f(fid, "/radius", 1, &
         dim2(2:2), npa_chords%radius, error)
    call h5ltmake_compressed_dataset_double_f(fid, "/energy", 1, &
         dim2(1:1), ebarr, error)
    call h5ltmake_compressed_dataset_double_f(fid, "/pitch", 1, &
         dim3(2:2), ptcharr, error)
    call h5ltmake_compressed_dataset_double_f(fid, "/flux", 2, &
         dim2, nweight%flux, error)
    call h5ltmake_compressed_dataset_double_f(fid, "/weight", 3, &
         dim3, nweight%weight, error)

    !Add attributes
    call h5ltset_attribute_string_f(fid, "/", "version", version, error)
    call h5ltset_attribute_string_f(fid,"/", "description", &
         "NPA E-p space sensitivity/weights and Flux calculated by FIDASIM", error)
    call h5ltset_attribute_string_f(fid,"/nchan", "description", &
         "Number of channels", error)
    call h5ltset_attribute_string_f(fid,"/nenergy", "description", &
         "Number of energy values", error)
    call h5ltset_attribute_string_f(fid,"/npitch", "description", &
         "Number of pitch value", error)
    call h5ltset_attribute_string_f(fid,"/energy","description", &
         "Energy array", error)
    call h5ltset_attribute_string_f(fid,"/energy", "units","keV", error)
    call h5ltset_attribute_string_f(fid,"/pitch", "description", &
         "Pitch array: p = v_parallel/v  w.r.t. the magnetic field", 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)
    call h5ltset_attribute_string_f(fid,"/flux", "description", &
         "Neutral flux: flux(energy,chan)", error)
    call h5ltset_attribute_string_f(fid,"/flux", "units", &
         "neutrals/(s*dE)", error)
    call h5ltset_attribute_string_f(fid,"/weight", "description", &
         "E-p space sensivity/weight of NPA diagnostics: weight(energy,pitch,chan)", error)
    call h5ltset_attribute_string_f(fid,"/weight","units", &
         "neutrals/(s*fast-ion*dE*dP)",error)

    if(inputs%calc_npa_wght.ge.2) then !Write diagnostic variables
        call write_beam_grid(fid, error)
        call h5ltmake_compressed_dataset_double_f(fid, "/emissivity", 4, &
             dim5(2:5), nweight%emissivity, error)
        call h5ltmake_compressed_dataset_double_f(fid, "/attenuation", 5, &
             dim5, nweight%attenuation, error)
        call h5ltmake_compressed_dataset_double_f(fid, "/cx", 5, &
             dim5, nweight%cx, error)
        call h5ltmake_compressed_dataset_double_f(fid, "/phit", 4, &
             dim5(2:5), npa_chords%phit%p, error)

        call h5ltset_attribute_string_f(fid,"/emissivity", "description", &
             "Neutral emissivity: emissivity(x,y,z,chan)", error)
        call h5ltset_attribute_string_f(fid,"/emissivity", "units", &
             "neutrals/(s*dV)", error)
        call h5ltset_attribute_string_f(fid,"/cx", "description", &
             "Charge-exchange rate: cx(energy,x,y,z,chan)", error)
        call h5ltset_attribute_string_f(fid,"/cx", "units", "s^(-1)", error)
        call h5ltset_attribute_string_f(fid,"/attenuation","description", &
             "Attenuation factor i.e. survival probability: attenuation(energy,x,y,z,chan)", error)
        call h5ltset_attribute_string_f(fid,"/phit","description", &
             "Probability of hitting the detector given an isotropic source: phit(x,y,z,chan)", error)
    endif

    !Close file
    call h5fclose_f(fid, error)

    !Close HDF5 interface
    call h5close_f(error)
    if(inputs%verbose.ge.1) then
        write(*,'(T4,a,a)') 'NPA weights written to: ',trim(filename)
    endif

end subroutine write_npa_weights