pfida_f Subroutine

public subroutine pfida_f()

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

Arguments

None

Calls

proc~~pfida_f~~CallsGraph proc~pfida_f pfida_f proc~get_nlaunch_pass_grid get_nlaunch_pass_grid proc~pfida_f->proc~get_nlaunch_pass_grid proc~store_fida_photons store_fida_photons proc~pfida_f->proc~store_fida_photons proc~bt_cx_rates bt_cx_rates proc~pfida_f->proc~bt_cx_rates proc~mc_fastion_pass_grid mc_fastion_pass_grid proc~pfida_f->proc~mc_fastion_pass_grid proc~ind2sub ind2sub proc~pfida_f->proc~ind2sub proc~get_plasma get_plasma proc~pfida_f->proc~get_plasma proc~parallel_merge_reservoirs parallel_merge_reservoirs proc~pfida_f->proc~parallel_merge_reservoirs proc~track_cylindrical track_cylindrical proc~pfida_f->proc~track_cylindrical proc~store_photon_birth store_photon_birth proc~pfida_f->proc~store_photon_birth proc~colrad colrad proc~pfida_f->proc~colrad interface~parallel_sum parallel_sum proc~pfida_f->interface~parallel_sum proc~gyro_correction gyro_correction proc~pfida_f->proc~gyro_correction proc~get_nlaunch_pass_grid->proc~ind2sub proc~rng_init rng_init proc~get_nlaunch_pass_grid->proc~rng_init interface~randind_cdf randind_cdf proc~get_nlaunch_pass_grid->interface~randind_cdf proc~cumsum cumsum proc~get_nlaunch_pass_grid->proc~cumsum proc~store_photons store_photons proc~store_fida_photons->proc~store_photons interface~interpol_coeff interpol_coeff proc~bt_cx_rates->interface~interpol_coeff proc~get_distribution get_distribution proc~mc_fastion_pass_grid->proc~get_distribution interface~randind randind proc~mc_fastion_pass_grid->interface~randind proc~mc_passive_grid mc_passive_grid proc~mc_fastion_pass_grid->proc~mc_passive_grid interface~randu randu proc~mc_fastion_pass_grid->interface~randu proc~get_fields get_fields proc~mc_fastion_pass_grid->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~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~track_cylindrical->proc~get_plasma proc~get_passive_grid_indices get_passive_grid_indices proc~track_cylindrical->proc~get_passive_grid_indices proc~cyl_to_uvw cyl_to_uvw proc~track_cylindrical->proc~cyl_to_uvw proc~line_cylinder_intersect line_cylinder_intersect proc~track_cylindrical->proc~line_cylinder_intersect proc~doppler_stark doppler_stark proc~track_cylindrical->proc~doppler_stark proc~track_cylindrical->proc~in_plasma proc~line_plane_intersect line_plane_intersect proc~track_cylindrical->proc~line_plane_intersect proc~plane_basis plane_basis proc~track_cylindrical->proc~plane_basis proc~track_cylindrical->proc~get_fields 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~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~parallel_sum_d0->mpi_allreduce 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 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~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 proc~plane_basis->proc~cross_product 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~mc_passive_grid->proc~cyl_to_uvw proc~mc_passive_grid->interface~randu 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~balback->proc~rswap proc~cyl_to_xyz->proc~cyl_to_uvw proc~cyl_to_xyz->proc~uvw_to_xyz proc~comdiv Comdiv proc~hqrvec->proc~comdiv proc~outerprod outerprod proc~ludcmp->proc~outerprod proc~swap swap proc~ludcmp->proc~swap

Called by

proc~~pfida_f~~CalledByGraph proc~pfida_f pfida_f program~fidasim fidasim program~fidasim->proc~pfida_f

Contents

Source Code


Source Code

subroutine pfida_f
    !+ Calculate Passive FIDA emission using a Fast-ion distribution function F(E,p,r,z)
    integer :: i,j,k,ic,ncell,is
    integer(Int64) :: iion
    real(Float64), dimension(3) :: ri      !! start position
    real(Float64), dimension(3) :: vi      !! velocity of fast ions
    real(Float64), dimension(3) :: xyz_vi
    real(Float64) :: denf !! fast-ion density
    integer, dimension(3) :: ind      !! new actual cell
    logical :: los_intersect
    !! Determination of the CX probability
    type(LocalEMFields) :: fields
    type(LocalProfiles) :: plasma
    real(Float64), dimension(nlevs) :: rates, rates_is !! CX rates

    !! Collisiional radiative model along track
    integer :: ntrack
    integer :: jj      !! counter along track
    type(ParticleTrack),dimension(pass_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) :: max_papprox, eb, ptch
    integer, dimension(pass_grid%ngrid) :: cell_ind
    real(Float64), dimension(pass_grid%nr,pass_grid%nz,pass_grid%nphi) :: papprox
    integer(Int32), dimension(pass_grid%nr,pass_grid%nz,pass_grid%nphi) :: nlaunch

    !! Estimate how many particles to launch in each cell
    papprox=0.d0
    do ic=1,pass_grid%ngrid
        call ind2sub(pass_grid%dims,ic,ind)
        i = ind(1) ; j = ind(2) ; k = ind(3)
        call get_plasma(plasma,ind=ind,input_coords=2)
        if(.not.plasma%in_plasma) cycle
        papprox(i,j,k) = sum(plasma%denn)*plasma%denf
    enddo
    max_papprox = maxval(papprox)
    where (papprox.lt.(max_papprox*1.d-6))
        papprox = 0.0
    endwhere

    ncell = 0
    do ic=1,pass_grid%ngrid
        call ind2sub(pass_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_pass_grid(inputs%n_pfida, 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,xyz_vi,ri,fields, &
    !$OMP tracks,ntrack,jj,plasma,rates,rates_is,is,denn,states,photons,denf,eb,ptch,los_intersect)
    loop_over_cells: do ic = istart, ncell, istep
        call ind2sub(pass_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_pass_grid(ind, fields, eb, ptch, denf, output_coords=1)
            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)
            xyz_vi = matmul(beam_grid%inv_basis,vi)

            !! Find the particles path through the interpolation grid
            call track_cylindrical(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 edge neutrals
            call get_plasma(plasma, pos=ri, input_coords=1)
            rates = 0.0
            do is=1,n_thermal
                call bt_cx_rates(plasma, plasma%denn(:,is), thermal_mass(is), xyz_vi, rates_is)
                rates = rates + rates_is
            enddo
            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,input_coords=1)

                call colrad(plasma, fbm%A, xyz_vi, tracks(jj)%time, states, denn, photons)

                call store_fida_photons(tracks(jj)%pos, xyz_vi, beam_lambda0, photons/nlaunch(i,j,k), passive=.True.)
                if(inputs%calc_res.ge.1) then 
                    call store_photon_birth(tracks(1)%pos, photons/nlaunch(i,j,k), spatres%pfida, passive=.True.)
                endif
            enddo loop_along_track
        enddo loop_over_fast_ions
    enddo loop_over_cells
    !$OMP END PARALLEL DO

#ifdef _MPI
    call parallel_sum(spec%pfida)
    if(inputs%calc_res.ge.1) then
        do jj=1,spec_chords%nchan
            call parallel_merge_reservoirs(spatres%pfida(jj))
        enddo
    endif
#endif

end subroutine pfida_f