cyl_check_if_confined Subroutine

public subroutine cyl_check_if_confined(F, Xcyl, flag)

Arguments

Type IntentOptional AttributesName
type(FIELDS), intent(in) :: F

An instance of the KORC derived type FIELDS.

real(kind=rp), intent(in), DIMENSION(:,:), ALLOCATABLE:: Xcyl

Particles' position in cylindrical coordinates. Xcyl(1,:) = , Xcyl(2,:) = , Xcyl(3,:) =

integer(kind=is), intent(inout), DIMENSION(:), ALLOCATABLE:: flag

Called by

proc~~cyl_check_if_confined~~CalledByGraph proc~cyl_check_if_confined cyl_check_if_confined proc~get_analytical_fields get_analytical_fields proc~get_analytical_fields->proc~cyl_check_if_confined proc~get_fields get_fields proc~get_fields->proc~get_analytical_fields proc~mh_psi MH_psi proc~mh_psi->proc~get_fields proc~sample_hollmann_distribution_3d_psi sample_Hollmann_distribution_3D_psi proc~sample_hollmann_distribution_3d_psi->proc~get_fields proc~unitvectors unitVectors proc~unitvectors->proc~get_fields proc~get_hollmann_distribution_3d_psi get_Hollmann_distribution_3D_psi proc~get_hollmann_distribution_3d_psi->proc~sample_hollmann_distribution_3d_psi proc~gyro_distribution gyro_distribution proc~gyro_distribution->proc~unitvectors

Contents

Source Code


Source Code

  subroutine cyl_check_if_confined(F,Xcyl,flag)
    implicit none
    TYPE(FIELDS), INTENT(IN)                               :: F
    !! An instance of the KORC derived type FIELDS.
    REAL(rp), DIMENSION(:,:), ALLOCATABLE, INTENT(IN)      :: Xcyl
    !! Particles' position in cylindrical coordinates. Xcyl(1,:) = \(R\),
    !! Xcyl(2,:) = \(\phi\), Xcyl(3,:) = \(Z\)
    INTEGER(is), DIMENSION(:), ALLOCATABLE, INTENT(INOUT)   :: flag
    REAL(rp)                                               :: a
    !! Distance to plasma edge as measured from the magnetic axis.
    REAL(rp)                                               :: Ro
    !! Radial position of the magnetic axis.
    INTEGER                                                :: pp
    !! Iterator.
    INTEGER                                                :: ss
    !! Iterator.
 

    if (size(Xcyl,1).eq.1) then
       ss = size(Xcyl,1)
    else
       if (Xcyl(2,1).eq.0) then
          ss=1_idef
       else
          ss = size(Xcyl,1)
       end if
    endif

    
    a = F%AB%a
    Ro = F%AB%Ro
    
    !$OMP PARALLEL DO FIRSTPRIVATE(ss,a,Ro) PRIVATE(pp) SHARED(Xcyl,flag)
    do pp=1_idef,ss
       if (sqrt((Xcyl(pp,1)-Ro)**2+Xcyl(pp,3)**2) .gt. a) then
          flag(pp)=0_is
       endif       
    end do
    !$OMP END PARALLEL DO
  end subroutine cyl_check_if_confined