interp_FOfields1_p Subroutine

public subroutine interp_FOfields1_p(pchunk, F, Y_R, Y_PHI, Y_Z, B_X, B_Y, B_Z, E_X, E_Y, E_Z, PSIp, flag_cache)

Arguments

Type IntentOptional AttributesName
integer, intent(in) :: pchunk
type(FIELDS), intent(in) :: F
real(kind=rp), intent(in), DIMENSION(pchunk):: Y_R
real(kind=rp), intent(in), DIMENSION(pchunk):: Y_PHI
real(kind=rp), intent(in), DIMENSION(pchunk):: Y_Z
real(kind=rp), intent(out), DIMENSION(pchunk):: B_X
real(kind=rp), intent(out), DIMENSION(pchunk):: B_Y
real(kind=rp), intent(out), DIMENSION(pchunk):: B_Z
real(kind=rp), intent(out), DIMENSION(pchunk):: E_X
real(kind=rp), intent(out), DIMENSION(pchunk):: E_Y
real(kind=rp), intent(out), DIMENSION(pchunk):: E_Z
real(kind=rp), intent(out), DIMENSION(pchunk):: PSIp
integer(kind=is), intent(inout), DIMENSION(pchunk):: flag_cache

Calls

proc~~interp_fofields1_p~~CallsGraph proc~interp_fofields1_p interp_FOfields1_p ezspline_interp ezspline_interp proc~interp_fofields1_p->ezspline_interp proc~check_if_in_fields_domain_p check_if_in_fields_domain_p proc~interp_fofields1_p->proc~check_if_in_fields_domain_p ezspline_error ezspline_error proc~interp_fofields1_p->ezspline_error proc~calculate_magnetic_field_p calculate_magnetic_field_p proc~interp_fofields1_p->proc~calculate_magnetic_field_p proc~calculate_magnetic_field_p->ezspline_interp proc~calculate_magnetic_field_p->ezspline_error ezspline_gradient ezspline_gradient proc~calculate_magnetic_field_p->ezspline_gradient

Called by

proc~~interp_fofields1_p~~CalledByGraph proc~interp_fofields1_p interp_FOfields1_p proc~adv_fointerp_top adv_FOinterp_top proc~adv_fointerp_top->proc~interp_fofields1_p proc~fo_init FO_init proc~fo_init->proc~interp_fofields1_p program~main main program~main->proc~adv_fointerp_top program~main->proc~fo_init

Contents

Source Code


Source Code

subroutine interp_FOfields1_p(pchunk,F,Y_R,Y_PHI,Y_Z,B_X,B_Y,B_Z, &
     E_X,E_Y,E_Z,PSIp,flag_cache)
  INTEGER, INTENT(IN)  :: pchunk
  TYPE(FIELDS), INTENT(IN)                               :: F
  REAL(rp),DIMENSION(pchunk),INTENT(IN)   :: Y_R,Y_PHI,Y_Z
  REAL(rp),DIMENSION(pchunk),INTENT(OUT)   :: B_X,B_Y,B_Z
  REAL(rp),DIMENSION(pchunk)   :: B_R,B_PHI  
  REAL(rp),DIMENSION(pchunk),INTENT(OUT)   :: E_X,E_Y,E_Z
  REAL(rp),DIMENSION(pchunk),INTENT(OUT)   :: PSIp
  REAL(rp),DIMENSION(pchunk)   :: E_R,E_PHI
  REAL(rp),DIMENSION(pchunk)   :: cP,sP  
  !  INTEGER(ip) :: ezerr
  INTEGER                                      :: cc
  !! Particle chunk iterator.
  INTEGER(is),DIMENSION(pchunk),INTENT(INOUT)   :: flag_cache

  call check_if_in_fields_domain_p(pchunk,F,Y_R,Y_PHI,Y_Z,flag_cache)

  call EZspline_interp(bfield_2d%A,pchunk,Y_R, Y_Z,PSIp, ezerr)
  call EZspline_error(ezerr)
  
  call calculate_magnetic_field_p(pchunk,F,Y_R,Y_Z,B_R,B_PHI,B_Z)
  
  call EZspline_interp(efield_2d%R,efield_2d%PHI,efield_2d%Z,pchunk,Y_R,Y_Z, &
       E_R,E_PHI,E_Z,ezerr)
  call EZspline_error(ezerr)
 
  !$OMP SIMD
!  !$OMP& aligned (cP,sP,B_X,B_Y,E_X,E_Y,Y_PHI,B_R,B_PHI,E_R,E_PHI)
  do cc=1_idef,pchunk
     cP(cc)=cos(Y_PHI(cc))
     sP(cc)=sin(Y_PHI(cc))
     
     B_X(cc) = B_R(cc)*cP(cc) - B_PHI(cc)*sP(cc)
     B_Y(cc) = B_R(cc)*sP(cc) + B_PHI(cc)*cP(cc)

     E_X(cc) = E_R(cc)*cP(cc) - E_PHI(cc)*sP(cc)
     E_Y(cc) = E_R(cc)*sP(cc) + E_PHI(cc)*cP(cc)

  end do
  !$OMP END SIMD
  
end subroutine interp_FOfields1_p