Calculate Active FIDA emission using a Fast-ion distribution function F(E,p,r,z)
subroutine fida_f
!+ Calculate Active FIDA emission using a Fast-ion distribution function F(E,p,r,z)
integer :: i,j,k,ic,ncell
integer(Int64) :: iion
real(Float64), dimension(3) :: ri !! start position
real(Float64), dimension(3) :: vi !! velocity of fast ions
real(Float64) :: denf !! fast-ion density
integer, dimension(3) :: ind !! new actual cell
integer, dimension(5) :: neut_types=[1,2,3,4,5]
logical :: los_intersect
!! Determination of the CX probability
type(LocalEMFields) :: fields
type(LocalProfiles) :: plasma
real(Float64), dimension(nlevs) :: rates !! CX rates
!! Collisiional radiative model along track
integer :: ntrack
integer :: jj !! counter along track
type(ParticleTrack),dimension(beam_grid%ntrack) :: tracks
real(Float64) :: photons !! photon flux
real(Float64), dimension(nlevs) :: states !! Density of n-states
real(Float64), dimension(nlevs) :: denn
!! Number of particles to launch
real(Float64) :: eb, ptch
integer, dimension(beam_grid%ngrid) :: cell_ind
real(Float64), dimension(beam_grid%nx,beam_grid%ny,beam_grid%nz) :: papprox
integer(Int32), dimension(beam_grid%nx,beam_grid%ny,beam_grid%nz) :: nlaunch
!! Estimate how many particles to launch in each cell
papprox=0.d0
do ic=1,beam_grid%ngrid
call ind2sub(beam_grid%dims,ic,ind)
i = ind(1) ; j = ind(2) ; k = ind(3)
call get_plasma(plasma,ind=ind)
if(.not.plasma%in_plasma) cycle
papprox(i,j,k) = (sum(neut%full%dens(:,i,j,k)) + &
sum(neut%half%dens(:,i,j,k)) + &
sum(neut%third%dens(:,i,j,k)) + &
sum(neut%dcx%dens(:,i,j,k)) + &
sum(neut%halo%dens(:,i,j,k)))* &
plasma%denf
enddo
ncell = 0
do ic=1,beam_grid%ngrid
call ind2sub(beam_grid%dims,ic,ind)
i = ind(1) ; j = ind(2) ; k = ind(3)
if(papprox(i,j,k).gt.0.0) then
ncell = ncell + 1
cell_ind(ncell) = ic
endif
enddo
call get_nlaunch(inputs%n_fida, papprox, nlaunch)
if(inputs%verbose.ge.1) then
write(*,'(T6,"# of markers: ",i10)') sum(nlaunch)
endif
!! Loop over all cells that have neutrals
!$OMP PARALLEL DO schedule(dynamic,1) private(ic,i,j,k,ind,iion,vi,ri,fields, &
!$OMP tracks,ntrack,jj,plasma,rates,denn,states,photons,denf,eb,ptch,los_intersect)
loop_over_cells: do ic = istart, ncell, istep
call ind2sub(beam_grid%dims,cell_ind(ic),ind)
i = ind(1) ; j = ind(2) ; k = ind(3)
loop_over_fast_ions: do iion=1, nlaunch(i, j, k)
!! Sample fast ion distribution for velocity and position
call mc_fastion(ind, fields, eb, ptch, denf)
if(denf.le.0.0) cycle loop_over_fast_ions
!! Correct for gyro motion and get particle position and velocity
call gyro_correction(fields, eb, ptch, fbm%A, ri, vi)
!! Find the particles path through the beam grid
call track(ri, vi, tracks, ntrack,los_intersect)
if(.not.los_intersect) cycle loop_over_fast_ions
if(ntrack.eq.0) cycle loop_over_fast_ions
!! Calculate CX probability with beam and halo neutrals
call get_total_cx_rate(tracks(1)%ind, ri, vi, neut_types, rates)
if(sum(rates).le.0.) cycle loop_over_fast_ions
!! Weight CX rates by ion source density
states=rates*denf
!! Calculate the spectra produced in each cell along the path
loop_along_track: do jj=1,ntrack
call get_plasma(plasma,pos=tracks(jj)%pos)
call colrad(plasma, fbm%A, vi, tracks(jj)%time, states, denn, photons)
call store_fida_photons(tracks(jj)%pos, vi, beam_lambda0, photons/nlaunch(i,j,k))
if(inputs%calc_res.ge.1) call store_photon_birth(tracks(1)%pos, photons/nlaunch(i,j,k), spatres%fida)
enddo loop_along_track
enddo loop_over_fast_ions
enddo loop_over_cells
!$OMP END PARALLEL DO
#ifdef _MPI
call parallel_sum(spec%fida)
if(inputs%calc_res.ge.1) then
do jj=1,spec_chords%nchan
call parallel_merge_reservoirs(spatres%fida(jj))
enddo
endif
#endif
end subroutine fida_f