gradient_2D_Bfields Subroutine

private subroutine gradient_2D_Bfields(Y, BR, BPHI, BZ, flag)

Arguments

Type IntentOptional AttributesName
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).


Calls

proc~~gradient_2d_bfields~~CallsGraph proc~gradient_2d_bfields gradient_2D_Bfields ezspline_error ezspline_error proc~gradient_2d_bfields->ezspline_error

Called by

proc~~gradient_2d_bfields~~CalledByGraph proc~gradient_2d_bfields gradient_2D_Bfields proc~interp_fields interp_fields proc~interp_fields->proc~gradient_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 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