npa_mc Subroutine

public subroutine npa_mc()

Calculate NPA flux using a Monte Carlo fast-ion distribution

Arguments

None

Calls

proc~~npa_mc~~CallsGraph proc~npa_mc npa_mc proc~get_fields get_fields proc~npa_mc->proc~get_fields proc~npa_gyro_range npa_gyro_range proc~npa_mc->proc~npa_gyro_range proc~store_npa store_npa proc~npa_mc->proc~store_npa proc~randu randu proc~npa_mc->proc~randu proc~uvw_to_xyz uvw_to_xyz proc~npa_mc->proc~uvw_to_xyz proc~get_indices get_indices proc~npa_mc->proc~get_indices proc~hit_npa_detector hit_npa_detector proc~npa_mc->proc~hit_npa_detector proc~gyro_trajectory gyro_trajectory proc~npa_mc->proc~gyro_trajectory proc~gyro_surface gyro_surface proc~npa_mc->proc~gyro_surface proc~attenuate attenuate proc~npa_mc->proc~attenuate proc~get_beam_cx_rate get_beam_cx_rate proc~npa_mc->proc~get_beam_cx_rate proc~in_plasma in_plasma proc~get_fields->proc~in_plasma proc~calc_perp_vectors calc_perp_vectors proc~get_fields->proc~calc_perp_vectors proc~approx_eq approx_eq proc~npa_gyro_range->proc~approx_eq proc~gyro_range gyro_range proc~npa_gyro_range->proc~gyro_range proc~store_npa->proc~get_fields proc~xyz_to_uvw xyz_to_uvw proc~store_npa->proc~xyz_to_uvw omp_get_thread_num omp_get_thread_num proc~randu->omp_get_thread_num proc~rng_uniform rng_uniform proc~randu->proc~rng_uniform proc~line_plane_intersect line_plane_intersect proc~hit_npa_detector->proc~line_plane_intersect proc~in_boundary in_boundary proc~hit_npa_detector->proc~in_boundary proc~get_plasma get_plasma proc~attenuate->proc~get_plasma proc~colrad colrad proc~attenuate->proc~colrad proc~bt_cx_rates bt_cx_rates proc~get_beam_cx_rate->proc~bt_cx_rates proc~get_beam_cx_rate->proc~get_plasma proc~bb_cx_rates bb_cx_rates proc~get_beam_cx_rate->proc~bb_cx_rates interface~interpol_coeff interpol_coeff proc~bt_cx_rates->interface~interpol_coeff proc~get_plasma->proc~in_plasma proc~gyro_range->proc~line_plane_intersect proc~gyro_range->proc~in_boundary proc~boundary_edge boundary_edge proc~gyro_range->proc~boundary_edge proc~line_gyro_surface_intersect line_gyro_surface_intersect proc~gyro_range->proc~line_gyro_surface_intersect proc~gyro_surface_coordinates gyro_surface_coordinates proc~gyro_range->proc~gyro_surface_coordinates proc~in_gyro_surface in_gyro_surface proc~gyro_range->proc~in_gyro_surface proc~get_rate_matrix get_rate_matrix proc~colrad->proc~get_rate_matrix proc~eigen eigen proc~colrad->proc~eigen proc~matinv matinv proc~colrad->proc~matinv proc~in_plasma->proc~xyz_to_uvw proc~in_plasma->interface~interpol_coeff proc~bb_cx_rates->interface~interpol_coeff proc~get_rate_matrix->interface~interpol_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~interpol1d_coeff_arr interpol1D_coeff_arr interface~interpol_coeff->proc~interpol1d_coeff_arr proc~hqr2 hqr2 proc~eigen->proc~hqr2 proc~balback balback proc~eigen->proc~balback proc~elmtrans elmtrans proc~eigen->proc~elmtrans proc~elmhes elmhes proc~eigen->proc~elmhes proc~balance balance proc~eigen->proc~balance proc~ludcmp ludcmp proc~matinv->proc~ludcmp proc~outerprod outerprod proc~ludcmp->proc~outerprod proc~hqrvec hqrvec proc~hqr2->proc~hqrvec proc~interpol2d_coeff_arr->proc~interpol2d_coeff proc~interpol1d_coeff_arr->proc~interpol1d_coeff

Called by

proc~~npa_mc~~CalledByGraph proc~npa_mc npa_mc program~fidasim fidasim program~fidasim->proc~npa_mc

Contents

Source Code


Source Code

subroutine npa_mc
    !+ Calculate NPA flux using a Monte Carlo fast-ion distribution
    integer :: iion,iphi,nphi
    type(FastIon) :: fast_ion
    real(Float64) :: phi,theta,dtheta
    real(Float64), dimension(3) :: ri, rf, rg, vi
    integer :: det,j,ichan,ir,nrange
    type(LocalEMFields) :: fields
    type(GyroSurface) :: gs
    real(Float64), dimension(nlevs) :: rates
    real(Float64), dimension(nlevs) :: states
    real(Float64) :: flux
    integer, dimension(4) :: neut_types=[1,2,3,4]
    integer, dimension(3) :: ind
    real(Float64), dimension(3) :: uvw, uvw_vi
    real(Float64), dimension(2,4) :: gyrange
    real(Float64) :: s,c
    real(Float64) :: maxcnt, inv_maxcnt, cnt
    real(Float64), dimension(1) :: randomu

    maxcnt=particles%nparticle
    inv_maxcnt = 100.d0/maxcnt
    nphi = ceiling(dble(inputs%n_npa)/particles%nparticle)
    if(inputs%verbose.ge.1) then
        write(*,'(T6,"# of markers: ",i9)') int(particles%nparticle*nphi,Int64)
    endif

    cnt=0.0
    !$OMP PARALLEL DO schedule(guided) private(iion,iphi,ind,fast_ion,vi,ri,rf,phi,s,c,ir, &
    !$OMP& randomu,rg,fields,uvw,uvw_vi,rates,states,flux,det,ichan,gs,nrange,gyrange,theta,dtheta)
    loop_over_fast_ions: do iion=1,particles%nparticle
        cnt=cnt+1
        fast_ion = particles%fast_ion(iion)
        if(fast_ion%vabs.eq.0)cycle loop_over_fast_ions
        if(.not.fast_ion%cross_grid) cycle loop_over_fast_ions
        phi_loop: do iphi=1,nphi
            !! Pick random torodial angle
            call randu(randomu)
            phi = fast_ion%phi_enter + fast_ion%delta_phi*randomu(1)
            s = sin(phi) ; c = cos(phi)

            !! Calculate position in machine coordinates
            uvw(1) = fast_ion%r*c
            uvw(2) = fast_ion%r*s
            uvw(3) = fast_ion%z

            if(inputs%dist_type.eq.2) then
                !! Convert to beam grid coordinates
                call uvw_to_xyz(uvw, rg)

                !! Get electomagnetic fields
                call get_fields(fields, pos=rg)

                !! Correct for gyro motion and get position and velocity
                call gyro_surface(fields, fast_ion%energy, fast_ion%pitch, gs)

                detector_loop: do ichan=1,npa_chords%nchan
                    call npa_gyro_range(ichan, gs, gyrange, nrange)
                    if(nrange.eq.0) cycle detector_loop

                    gyro_range_loop: do ir=1,nrange
                        dtheta = gyrange(2,ir)
                        theta = gyrange(1,ir) + 0.5*dtheta
                        call gyro_trajectory(gs, theta, ri, vi)

                        !! Check if particle hits a NPA detector
                        call hit_npa_detector(ri, vi ,det, rf, det=ichan)
                        if(det.ne.ichan) then
                            if (inputs%verbose.ge.0)then
                                write(*,*) "NPA_MC: Missed Detector ",ichan
                            endif
                            cycle gyro_range_loop
                        endif

                        !! Get beam grid indices at ri
                        call get_indices(ri,ind)

                        !! Calculate CX probability with beam and halo neutrals
                        call get_beam_cx_rate(ind,ri,vi,beam_ion,neut_types,rates)
                        if(sum(rates).le.0.) cycle gyro_range_loop

                        !! Weight CX rates by ion source density
                        states=rates*fast_ion%weight/nphi

                        !! Attenuate states as the particle move through plasma
                        call attenuate(ri,rf,vi,states)

                        !! Store NPA Flux
                        flux = (dtheta/(2*pi))*sum(states)*beam_grid%dv
                        call store_npa(det,ri,rf,vi,flux,fast_ion%class)
                    enddo gyro_range_loop
                enddo detector_loop
            else !! Full Orbit
                !! Convert to beam grid coordinates
                call uvw_to_xyz(uvw, ri)

                !! Calculate velocity vector
                uvw_vi(1) = c*fast_ion%vr - s*fast_ion%vt
                uvw_vi(2) = s*fast_ion%vr + c*fast_ion%vt
                uvw_vi(3) = fast_ion%vz
                vi = matmul(beam_grid%inv_basis,uvw_vi)

                !! Check if particle hits a NPA detector
                call hit_npa_detector(ri, vi ,det, rf)
                if(det.eq.0) cycle phi_loop

                !! Get beam grid indices at ri
                call get_indices(ri,ind)

                !! Calculate CX probability with beam and halo neutrals
                call get_beam_cx_rate(ind,ri,vi,beam_ion,neut_types,rates)
                if(sum(rates).le.0.) cycle phi_loop

                !! Weight CX rates by ion source density
                states=rates*fast_ion%weight/nphi

                !! Attenuate states as the particle moves though plasma
                call attenuate(ri,rf,vi,states)

                !! Store NPA Flux
                flux = sum(states)*beam_grid%dv
                call store_npa(det,ri,rf,vi,flux,fast_ion%class)
            endif
        enddo phi_loop
        if (inputs%verbose.ge.2)then
          WRITE(*,'(f7.2,"% completed",a,$)') cnt*inv_maxcnt,char(13)
        endif
    enddo loop_over_fast_ions
    !$OMP END PARALLEL DO
    if(inputs%verbose.ge.1) then
        write(*,'(T4,"Number of NPA particles that hit a detector: ",i8)') npa%npart
    endif

end subroutine npa_mc