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~store_fida_photons store_fida_photons proc~fida_f->proc~store_fida_photons proc~mc_fastion mc_fastion proc~fida_f->proc~mc_fastion proc~get_plasma get_plasma proc~fida_f->proc~get_plasma proc~get_nlaunch get_nlaunch proc~fida_f->proc~get_nlaunch proc~parallel_merge_reservoirs parallel_merge_reservoirs proc~fida_f->proc~parallel_merge_reservoirs proc~ind2sub ind2sub proc~fida_f->proc~ind2sub proc~store_photon_birth store_photon_birth proc~fida_f->proc~store_photon_birth proc~track track proc~fida_f->proc~track proc~get_total_cx_rate get_total_cx_rate proc~fida_f->proc~get_total_cx_rate proc~colrad colrad proc~fida_f->proc~colrad interface~parallel_sum parallel_sum proc~fida_f->interface~parallel_sum proc~gyro_correction gyro_correction proc~fida_f->proc~gyro_correction proc~store_photons store_photons proc~store_fida_photons->proc~store_photons proc~get_distribution get_distribution proc~mc_fastion->proc~get_distribution proc~mc_beam_grid mc_beam_grid proc~mc_fastion->proc~mc_beam_grid interface~randind randind proc~mc_fastion->interface~randind interface~randu randu proc~mc_fastion->interface~randu proc~get_fields get_fields proc~mc_fastion->proc~get_fields proc~xyz_to_uvw xyz_to_uvw proc~get_plasma->proc~xyz_to_uvw proc~uvw_to_xyz uvw_to_xyz proc~get_plasma->proc~uvw_to_xyz proc~in_plasma in_plasma proc~get_plasma->proc~in_plasma proc~get_position get_position proc~get_plasma->proc~get_position proc~get_nlaunch->proc~ind2sub proc~rng_init rng_init proc~get_nlaunch->proc~rng_init interface~randind_cdf randind_cdf proc~get_nlaunch->interface~randind_cdf proc~cumsum cumsum proc~get_nlaunch->proc~cumsum proc~parallel_merge_reservoirs->interface~parallel_sum proc~init_reservoir init_reservoir proc~parallel_merge_reservoirs->proc~init_reservoir proc~my_rank my_rank proc~parallel_merge_reservoirs->proc~my_rank proc~merge_reservoirs merge_reservoirs proc~parallel_merge_reservoirs->proc~merge_reservoirs proc~num_ranks num_ranks proc~parallel_merge_reservoirs->proc~num_ranks proc~get_passive_grid_indices get_passive_grid_indices proc~store_photon_birth->proc~get_passive_grid_indices proc~update_reservoir update_reservoir proc~store_photon_birth->proc~update_reservoir proc~get_indices get_indices proc~store_photon_birth->proc~get_indices proc~track->proc~get_plasma proc~doppler_stark doppler_stark proc~track->proc~doppler_stark proc~track->proc~in_plasma proc~track->proc~get_fields proc~track->proc~get_indices proc~neutral_cx_rate neutral_cx_rate proc~get_total_cx_rate->proc~neutral_cx_rate proc~get_rate_matrix get_rate_matrix proc~colrad->proc~get_rate_matrix proc~eigen eigen proc~colrad->proc~eigen proc~linsolve linsolve proc~colrad->proc~linsolve proc~parallel_sum_i3 parallel_sum_i3 interface~parallel_sum->proc~parallel_sum_i3 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_d5 parallel_sum_d5 interface~parallel_sum->proc~parallel_sum_d5 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~parallel_sum_d4 parallel_sum_d4 interface~parallel_sum->proc~parallel_sum_d4 proc~parallel_sum_d3 parallel_sum_d3 interface~parallel_sum->proc~parallel_sum_d3 proc~parallel_sum_d2 parallel_sum_d2 interface~parallel_sum->proc~parallel_sum_d2 proc~gyro_step gyro_step proc~gyro_correction->proc~gyro_step proc~pitch_to_vec pitch_to_vec proc~gyro_correction->proc~pitch_to_vec proc~gyro_correction->interface~randu proc~xyz_to_cyl xyz_to_cyl proc~get_passive_grid_indices->proc~xyz_to_cyl proc~parallel_sum_i3->proc~parallel_sum_i1 mpi_allreduce mpi_allreduce proc~parallel_sum_i3->mpi_allreduce proc~get_distribution->proc~xyz_to_uvw proc~get_distribution->proc~get_position interface~interpol interpol proc~get_distribution->interface~interpol proc~rng_init->proc~my_rank proc~rng_seed rng_seed proc~rng_init->proc~rng_seed proc~init_reservoir->interface~randu proc~bb_cx_rates bb_cx_rates proc~neutral_cx_rate->proc~bb_cx_rates proc~parallel_sum_d0->mpi_allreduce proc~mc_beam_grid->interface~randu proc~parallel_sum_d1->proc~parallel_sum_d1 proc~parallel_sum_d1->mpi_allreduce proc~parallel_sum_d5->proc~parallel_sum_d1 proc~parallel_sum_d5->mpi_allreduce proc~parallel_sum_i1->proc~parallel_sum_i1 proc~parallel_sum_i1->mpi_allreduce proc~parallel_sum_i0->mpi_allreduce proc~merge_reservoirs->proc~init_reservoir proc~merge_reservoirs->proc~update_reservoir proc~merge_reservoirs->interface~randu float float proc~merge_reservoirs->float proc~parallel_sum_i2->proc~parallel_sum_i1 proc~parallel_sum_i2->mpi_allreduce proc~update_reservoir->proc~init_reservoir proc~update_reservoir->interface~randind proc~update_reservoir->interface~randu interface~interpol_coeff interpol_coeff proc~get_rate_matrix->interface~interpol_coeff proc~elmhes elmhes proc~eigen->proc~elmhes proc~balance balance proc~eigen->proc~balance proc~hqr2 hqr2 proc~eigen->proc~hqr2 proc~balback balback proc~eigen->proc~balback proc~elmtrans elmtrans proc~eigen->proc~elmtrans proc~cross_product cross_product proc~gyro_step->proc~cross_product proc~in_plasma->proc~xyz_to_uvw proc~cyl_to_uvw cyl_to_uvw proc~in_plasma->proc~cyl_to_uvw proc~in_plasma->interface~interpol_coeff proc~parallel_sum_d4->proc~parallel_sum_d1 proc~parallel_sum_d4->mpi_allreduce dgetrs dgetrs proc~linsolve->dgetrs proc~matinv matinv proc~linsolve->proc~matinv dgetrf dgetrf proc~linsolve->dgetrf proc~parallel_sum_d3->proc~parallel_sum_d1 proc~parallel_sum_d3->mpi_allreduce proc~cyl_to_xyz cyl_to_xyz proc~get_position->proc~cyl_to_xyz proc~store_photons->proc~get_passive_grid_indices proc~store_photons->proc~uvw_to_xyz proc~store_photons->proc~get_fields proc~store_photons->proc~get_indices proc~spectrum spectrum proc~store_photons->proc~spectrum proc~get_fields->proc~xyz_to_uvw proc~get_fields->proc~uvw_to_xyz proc~get_fields->proc~in_plasma proc~get_fields->proc~get_position proc~calc_perp_vectors calc_perp_vectors proc~get_fields->proc~calc_perp_vectors proc~parallel_sum_d2->proc~parallel_sum_d1 proc~parallel_sum_d2->mpi_allreduce proc~rswap RSWAP proc~elmhes->proc~rswap proc~xyz_to_cyl->proc~xyz_to_uvw proc~uvw_to_cyl uvw_to_cyl proc~xyz_to_cyl->proc~uvw_to_cyl proc~balance->proc~rswap proc~hqrvec hqrvec proc~hqr2->proc~hqrvec proc~lubksb lubksb proc~matinv->proc~lubksb proc~ludcmp ludcmp proc~matinv->proc~ludcmp proc~bb_cx_rates->interface~interpol_coeff proc~balback->proc~rswap proc~cyl_to_xyz->proc~uvw_to_xyz proc~cyl_to_xyz->proc~cyl_to_uvw proc~comdiv Comdiv proc~hqrvec->proc~comdiv proc~outerprod outerprod proc~ludcmp->proc~outerprod proc~swap swap proc~ludcmp->proc~swap

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%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