interp_2D_bfields Subroutine

private subroutine interp_2D_bfields(params, Y, B, flag)

Arguments

Type IntentOptional AttributesName
type(KORC_PARAMS), intent(in) :: params

Core KORC simulation parameters.

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

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

real(kind=rp), intent(inout), DIMENSION(:,:), ALLOCATABLE:: B

Cartesian components of interpolated magnetic field components. B(1,:)=, B(2,:)=, and B(3,:)=.

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

Flag that indicates whether particles are followed in the simulation (flag=1), or not (flag=0).


Calls

proc~~interp_2d_bfields~~CallsGraph proc~interp_2d_bfields interp_2D_bfields ezspline_error ezspline_error proc~interp_2d_bfields->ezspline_error

Called by

proc~~interp_2d_bfields~~CalledByGraph proc~interp_2d_bfields interp_2D_bfields proc~interp_fields interp_fields proc~interp_fields->proc~interp_2d_bfields proc~get_fields get_fields proc~get_fields->proc~interp_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 interp_2D_bfields(params,Y,B,flag)
  !! @note Subroutine for interpolating the pre-computed, axisymmetric magnetic
  !! field to the particles' position. @endnote
  TYPE(KORC_PARAMS), INTENT(IN)      :: params
  !! Core KORC simulation parameters.
  REAL(rp), DIMENSION(:,:), ALLOCATABLE, INTENT(IN)      :: Y
  !! Particles' position in cylindrical coordinates, Y(1,:) = \(R\),
  !! Y(2,:) = \(\phi\), and Y(3,:) = \(Z\).
  REAL(rp), DIMENSION(:,:), ALLOCATABLE, INTENT(INOUT)   :: B
  !! Cartesian components of interpolated magnetic field components.
  !! B(1,:)=\(B_x\), B(2,:)=\(B_y\), and B(3,:)=\(B_z\).
  REAL(rp), DIMENSION(:,:), ALLOCATABLE                  :: F
  !! Cylindrical components of interpolated magnetic field components.
  !! F(1,:)=\(B_R\), F(2,:)=\(B_\phi\), and F(3,:)=\(B_Z\).
  INTEGER(is), DIMENSION(:), ALLOCATABLE, INTENT(INOUT)  :: flag
  !! Flag that indicates whether particles are followed in the simulation
  !! (flag=1), or not (flag=0).
  INTEGER                                                :: pp
  !! Particle iterator.
  INTEGER                                                :: ss
  !! Species iterator.
  
  ss = size(Y,1)

  ALLOCATE(F(ss,3))
  !$OMP PARALLEL DO FIRSTPRIVATE(ss) PRIVATE(pp,ezerr) &
  !$OMP& SHARED(F,Y,B,flag,bfield_2d)
  do pp=1_idef,ss
     if ( flag(pp) .EQ. 1_is ) then
        call EZspline_interp(bfield_2d%R, Y(pp,1), Y(pp,3), F(pp,1), ezerr)
        call EZspline_error(ezerr)

        if (ezerr .NE. 0) then ! We flag the particle as lost
           flag(pp) = 0_is
        end if

        call EZspline_interp(bfield_2d%PHI, Y(pp,1), Y(pp,3), F(pp,2), ezerr)
        call EZspline_error(ezerr)

        call EZspline_interp(bfield_2d%Z, Y(pp,1), Y(pp,3), F(pp,3), ezerr)
        call EZspline_error(ezerr)

        if (.not.params%GC_coords) then
           B(pp,1) = F(pp,1)*COS(Y(pp,2)) - F(pp,2)*SIN(Y(pp,2))
           B(pp,2) = F(pp,1)*SIN(Y(pp,2)) + F(pp,2)*COS(Y(pp,2))
           B(pp,3) = F(pp,3)
        else
           B(pp,1) = F(pp,1)
           B(pp,2) = F(pp,2)
           B(pp,3) = F(pp,3)
        end if
     end if
  end do
  !$OMP END PARALLEL DO
  DEALLOCATE(F)
end subroutine interp_2D_bfields