Subroutine for interpolating the pre-computed, axisymmetric gradient of the magnitude of themagnetic field to the particles' position. Stored as cylindrical components of field.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
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 | :: | BR | Cylindrical components of gradient of R-component of magnetic field. |
|
real(kind=rp), | intent(inout), | DIMENSION(:,:), ALLOCATABLE | :: | BPHI | Cylindrical components of gradient of R-component of magnetic field. |
|
real(kind=rp), | intent(inout), | DIMENSION(:,:), ALLOCATABLE | :: | BZ | Cylindrical components of gradient of R-component of magnetic field. |
|
integer(kind=is), | intent(inout), | DIMENSION(:), ALLOCATABLE | :: | flag | Flag that indicates whether particles are followed in the simulation (flag=1), or not (flag=0). |
subroutine gradient_2D_Bfields(Y,BR,BPHI,BZ,flag)
!! @note Subroutine for interpolating the pre-computed, axisymmetric
!! gradient of the magnitude of themagnetic field to the particles'
!! position. Stored as cylindrical components of field. @endnote
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) :: BR
!! Cylindrical components of gradient of R-component of magnetic field.
REAL(rp), DIMENSION(:,:), ALLOCATABLE, INTENT(INOUT) :: BPHI
!! Cylindrical components of gradient of R-component of magnetic field.
REAL(rp), DIMENSION(:,:), ALLOCATABLE, INTENT(INOUT) :: BZ
!! Cylindrical components of gradient of R-component of magnetic field.
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(2,ss))
!$OMP PARALLEL DO FIRSTPRIVATE(ss) PRIVATE(pp,ezerr) &
!$OMP& SHARED(Y,BR,BPHI,BZ,flag,bfield_2d)
do pp=1_idef,ss
if ( flag(pp) .EQ. 1_is ) then
call EZspline_gradient(bfield_2d%R, Y(pp,1), Y(pp,3), F(:,pp), ezerr)
call EZspline_error(ezerr)
BR(pp,1) = F(pp,1)
BR(pp,2) = 0._rp
BR(pp,3) = F(pp,2)
if (ezerr .NE. 0) then ! We flag the particle as lost
flag(pp) = 0_is
end if
call EZspline_gradient(bfield_2d%PHI, Y(pp,1), Y(pp,3), F(:,pp), &
ezerr)
call EZspline_error(ezerr)
BPHI(pp,1) = F(pp,1)
BPHI(pp,2) = 0._rp
BPHI(pp,3) = F(pp,2)
call EZspline_gradient(bfield_2d%Z, Y(pp,1), Y(pp,3), F(:,pp), ezerr)
call EZspline_error(ezerr)
BZ(pp,1) = F(pp,1)
BZ(pp,2) = 0._rp
BZ(pp,3) = F(pp,2)
end if
end do
!$OMP END PARALLEL DO
DEALLOCATE(F)
end subroutine gradient_2D_Bfields