fida_f Subroutine

public subroutine fida_f()

Calculate Active FIDA emission using a Fast-ion distribution function F(E,p,r,z)

Arguments

None

Calls

proc~~fida_f~~CallsGraph proc~fida_f fida_f proc~get_nlaunch get_nlaunch proc~fida_f->proc~get_nlaunch proc~ind2sub ind2sub proc~fida_f->proc~ind2sub proc~get_plasma get_plasma proc~fida_f->proc~get_plasma proc~track track proc~fida_f->proc~track proc~gyro_correction gyro_correction proc~fida_f->proc~gyro_correction interface~parallel_sum parallel_sum proc~fida_f->interface~parallel_sum proc~mc_fastion mc_fastion proc~fida_f->proc~mc_fastion interface~randind_cdf randind_cdf proc~get_nlaunch->interface~randind_cdf proc~cumsum cumsum proc~get_nlaunch->proc~cumsum proc~rng_init rng_init proc~get_nlaunch->proc~rng_init proc~get_position get_position proc~get_plasma->proc~get_position proc~uvw_to_xyz uvw_to_xyz proc~get_plasma->proc~uvw_to_xyz proc~xyz_to_uvw xyz_to_uvw proc~get_plasma->proc~xyz_to_uvw proc~in_plasma in_plasma proc~get_plasma->proc~in_plasma proc~get_indices get_indices proc~track->proc~get_indices proc~track->proc~in_plasma proc~get_fields get_fields proc~track->proc~get_fields proc~doppler_stark doppler_stark proc~track->proc~doppler_stark proc~pitch_to_vec pitch_to_vec proc~gyro_correction->proc~pitch_to_vec proc~gyro_step gyro_step proc~gyro_correction->proc~gyro_step interface~randu randu proc~gyro_correction->interface~randu proc~parallel_sum_d4 parallel_sum_d4 interface~parallel_sum->proc~parallel_sum_d4 proc~parallel_sum_d5 parallel_sum_d5 interface~parallel_sum->proc~parallel_sum_d5 proc~parallel_sum_d2 parallel_sum_d2 interface~parallel_sum->proc~parallel_sum_d2 proc~parallel_sum_d3 parallel_sum_d3 interface~parallel_sum->proc~parallel_sum_d3 proc~parallel_sum_d0 parallel_sum_d0 interface~parallel_sum->proc~parallel_sum_d0 proc~parallel_sum_d1 parallel_sum_d1 interface~parallel_sum->proc~parallel_sum_d1 proc~parallel_sum_i1 parallel_sum_i1 interface~parallel_sum->proc~parallel_sum_i1 proc~parallel_sum_i0 parallel_sum_i0 interface~parallel_sum->proc~parallel_sum_i0 proc~parallel_sum_i2 parallel_sum_i2 interface~parallel_sum->proc~parallel_sum_i2 proc~get_distribution get_distribution proc~mc_fastion->proc~get_distribution proc~mc_fastion->interface~randu proc~mc_fastion->proc~get_fields interface~randind randind proc~mc_fastion->interface~randind proc~get_distribution->proc~xyz_to_uvw interface~interpol interpol proc~get_distribution->interface~interpol proc~cyl_to_xyz cyl_to_xyz proc~get_position->proc~cyl_to_xyz mpi_allreduce mpi_allreduce proc~parallel_sum_d4->mpi_allreduce proc~cross_product cross_product proc~gyro_step->proc~cross_product proc~parallel_sum_d5->mpi_allreduce proc~parallel_sum_d2->mpi_allreduce proc~parallel_sum_d3->mpi_allreduce proc~parallel_sum_d0->mpi_allreduce proc~parallel_sum_d1->mpi_allreduce proc~parallel_sum_i1->mpi_allreduce proc~parallel_sum_i0->mpi_allreduce proc~parallel_sum_i2->mpi_allreduce proc~in_plasma->proc~xyz_to_uvw interface~interpol_coeff interpol_coeff proc~in_plasma->interface~interpol_coeff proc~cyl_to_uvw cyl_to_uvw proc~in_plasma->proc~cyl_to_uvw proc~get_fields->proc~uvw_to_xyz proc~get_fields->proc~xyz_to_uvw proc~get_fields->proc~in_plasma proc~calc_perp_vectors calc_perp_vectors proc~get_fields->proc~calc_perp_vectors proc~rng_seed rng_seed proc~rng_init->proc~rng_seed proc~my_rank my_rank proc~rng_init->proc~my_rank proc~cyl_interpol3d_coeff cyl_interpol3D_coeff interface~interpol_coeff->proc~cyl_interpol3d_coeff proc~interpol1d_coeff interpol1D_coeff interface~interpol_coeff->proc~interpol1d_coeff proc~interpol2d_coeff_arr interpol2D_coeff_arr interface~interpol_coeff->proc~interpol2d_coeff_arr proc~interpol2d_coeff interpol2D_coeff interface~interpol_coeff->proc~interpol2d_coeff proc~cyl_interpol3d_coeff_arr cyl_interpol3D_coeff_arr interface~interpol_coeff->proc~cyl_interpol3d_coeff_arr proc~interpol1d_coeff_arr interpol1D_coeff_arr interface~interpol_coeff->proc~interpol1d_coeff_arr proc~cyl_to_xyz->proc~uvw_to_xyz proc~cyl_to_xyz->proc~cyl_to_uvw proc~interpol2d_arr interpol2D_arr interface~interpol->proc~interpol2d_arr proc~interpol3d_arr interpol3D_arr interface~interpol->proc~interpol3d_arr proc~interpol2d_2d_arr interpol2D_2D_arr interface~interpol->proc~interpol2d_2d_arr proc~interpol1d_arr interpol1D_arr interface~interpol->proc~interpol1d_arr proc~interpol3d_2d_arr interpol3D_2D_arr interface~interpol->proc~interpol3d_2d_arr proc~interpol2d_arr->interface~interpol_coeff proc~cyl_interpol3d_coeff->proc~interpol2d_coeff proc~interpol3d_arr->interface~interpol_coeff proc~interpol2d_2d_arr->interface~interpol_coeff proc~interpol1d_arr->interface~interpol_coeff proc~interpol2d_coeff_arr->proc~interpol2d_coeff proc~cyl_interpol3d_coeff_arr->proc~cyl_interpol3d_coeff proc~cyl_interpol3d_coeff_arr->proc~interpol2d_coeff proc~interpol1d_coeff_arr->proc~interpol1d_coeff proc~interpol3d_2d_arr->interface~interpol_coeff

Called by

proc~~fida_f~~CalledByGraph proc~fida_f fida_f program~fidasim fidasim program~fidasim->proc~fida_f

Contents

Source Code


Source Code

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(:,i,j,k)) + &
                          sum(neut%half(:,i,j,k)) + &
                          sum(neut%third(:,i,j,k)) + &
                          sum(neut%dcx(:,i,j,k)) + &
                          sum(neut%halo(:,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, 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_beam_cx_rate(tracks(1)%ind, ri, vi, beam_ion, 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,beam_ion, vi, tracks(jj)%time, states, denn, photons)

                call store_fida_photons(tracks(jj)%pos, vi, photons/nlaunch(i,j,k))
            enddo loop_along_track
        enddo loop_over_fast_ions
    enddo loop_over_cells
    !$OMP END PARALLEL DO

#ifdef _MPI
    call parallel_sum(spec%fida)
#endif

end subroutine fida_f