analytical_fields_p Subroutine

public subroutine analytical_fields_p(B0, E0, R0, q0, lam, ar, X_X, X_Y, X_Z, B_X, B_Y, B_Z, E_X, E_Y, E_Z, flag_cache)

Arguments

Type IntentOptional AttributesName
real(kind=rp), intent(in) :: B0
real(kind=rp), intent(in) :: E0
real(kind=rp), intent(in) :: R0
real(kind=rp), intent(in) :: q0
real(kind=rp), intent(in) :: lam
real(kind=rp), intent(in) :: ar
real(kind=rp), intent(in), DIMENSION(8):: X_X
real(kind=rp), intent(in), DIMENSION(8):: X_Y
real(kind=rp), intent(in), DIMENSION(8):: X_Z
real(kind=rp), intent(out), DIMENSION(8):: B_X
real(kind=rp), intent(out), DIMENSION(8):: B_Y
real(kind=rp), intent(out), DIMENSION(8):: B_Z
real(kind=rp), intent(out), DIMENSION(8):: E_X
real(kind=rp), intent(out), DIMENSION(8):: E_Y
real(kind=rp), intent(out), DIMENSION(8):: E_Z
integer(kind=is), intent(inout), DIMENSION(8):: flag_cache

Contents

Source Code


Source Code

  subroutine analytical_fields_p(B0,E0,R0,q0,lam,ar,X_X,X_Y,X_Z, &
       B_X,B_Y,B_Z,E_X,E_Y,E_Z,flag_cache)
    REAL(rp),  INTENT(IN)      :: R0,B0,lam,q0,E0,ar
    REAL(rp),  INTENT(IN),DIMENSION(8)      :: X_X,X_Y,X_Z
    REAL(rp),  INTENT(OUT),DIMENSION(8)     :: B_X,B_Y,B_Z
    REAL(rp),  INTENT(OUT),DIMENSION(8)     :: E_X,E_Y,E_Z
    INTEGER(is),  INTENT(INOUT),DIMENSION(8)     :: flag_cache
    REAL(rp),DIMENSION(8)     :: T_R,T_T,T_Z
    REAL(rp),DIMENSION(8)                               :: Ezeta
    !! Toroidal electric field \(E_\zeta\).
    REAL(rp),DIMENSION(8)                               :: Bzeta
    !! Toroidal magnetic field \(B_\zeta\).
    REAL(rp),DIMENSION(8)                              :: Bp
    !! Poloidal magnetic field \(B_\theta(r)\).
    REAL(rp),DIMENSION(8)                               :: eta
    !! Aspect ratio \(\eta\).
    REAL(rp),DIMENSION(8)                                :: q
    !! Safety profile \(q(r)\).
    REAL(rp),DIMENSION(8)                             :: cT,sT,cZ,sZ
    INTEGER                                      :: cc
    !! Particle chunk iterator.

    call cart_to_tor_check_if_confined_p(ar,R0,X_X,X_Y,X_Z, &
         T_R,T_T,T_Z,flag_cache)

    !$OMP SIMD
    !    !$OMP& aligned(cT,sT,cZ,sZ,eta,q,Bp,Bzeta,B_X,B_Y,B_Z, &
    !    !$OMP& Ezeta,E_X,E_Y,E_Z,T_T,T_Z,T_R)
    do cc=1_idef,8
       cT(cc)=cos(T_T(cc))
       sT(cc)=sin(T_T(cc))
       cZ(cc)=cos(T_Z(cc))
       sZ(cc)=sin(T_Z(cc))

       eta(cc) = T_R(cc)/R0
       q(cc) = q0*(1.0_rp + (T_R(cc)*T_R(cc)/(lam*lam)))
       Bp(cc) = -eta(cc)*B0/(q(cc)*(1.0_rp + eta(cc)*cT(cc)))
       Bzeta(cc) = B0/( 1.0_rp + eta(cc)*cT(cc))


       B_X(cc) =  Bzeta(cc)*cZ(cc) - Bp(cc)*sT(cc)*sZ(cc)
       B_Y(cc) = -Bzeta(cc)*sZ(cc) - Bp(cc)*sT(cc)*cZ(cc)
       B_Z(cc) = Bp(cc)*cT(cc)

       Ezeta(cc) = -E0/( 1.0_rp + eta(cc)*cT(cc))

       E_X(cc) = Ezeta(cc)*cZ(cc)
       E_Y(cc) = -Ezeta(cc)*sZ(cc)
       E_Z(cc) = 0.0_rp
    end do
    !$OMP END SIMD

  end subroutine analytical_fields_p