calculate_3DBdBfields1_p Subroutine

public subroutine calculate_3DBdBfields1_p(pchunk, F, Y_R, Y_PHI, Y_Z, B_R, B_PHI, B_Z, E_R, E_PHI, E_Z, curlb_R, curlb_PHI, curlb_Z, gradB_R, gradB_PHI, gradB_Z, flag_cache, PSIp)

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_R
real(kind=rp), intent(out), DIMENSION(pchunk):: B_PHI
real(kind=rp), intent(out), DIMENSION(pchunk):: B_Z
real(kind=rp), intent(out), DIMENSION(pchunk):: E_R
real(kind=rp), intent(out), DIMENSION(pchunk):: E_PHI
real(kind=rp), intent(out), DIMENSION(pchunk):: E_Z
real(kind=rp), intent(out), DIMENSION(pchunk):: curlb_R
real(kind=rp), intent(out), DIMENSION(pchunk):: curlb_PHI
real(kind=rp), intent(out), DIMENSION(pchunk):: curlb_Z
real(kind=rp), intent(out), DIMENSION(pchunk):: gradB_R
real(kind=rp), intent(out), DIMENSION(pchunk):: gradB_PHI
real(kind=rp), intent(out), DIMENSION(pchunk):: gradB_Z
integer(kind=is), intent(inout), DIMENSION(pchunk):: flag_cache
real(kind=rp), intent(out), DIMENSION(pchunk):: PSIp

Calls

proc~~calculate_3dbdbfields1_p~~CallsGraph proc~calculate_3dbdbfields1_p calculate_3DBdBfields1_p ezspline_interp ezspline_interp proc~calculate_3dbdbfields1_p->ezspline_interp proc~check_if_in_fields_domain_p check_if_in_fields_domain_p proc~calculate_3dbdbfields1_p->proc~check_if_in_fields_domain_p ezspline_error ezspline_error proc~calculate_3dbdbfields1_p->ezspline_error

Called by

proc~~calculate_3dbdbfields1_p~~CalledByGraph proc~calculate_3dbdbfields1_p calculate_3DBdBfields1_p proc~advance_gcinterp_3dbdb1_vars advance_GCinterp_3DBdB1_vars proc~advance_gcinterp_3dbdb1_vars->proc~calculate_3dbdbfields1_p

Contents


Source Code

subroutine calculate_3DBdBfields1_p(pchunk,F,Y_R,Y_PHI,Y_Z,B_R,B_PHI,B_Z, &
     E_R,E_PHI,E_Z,curlb_R,curlb_PHI,curlb_Z,gradB_R,gradB_PHI,gradB_Z, &
     flag_cache,PSIp)
  INTEGER, INTENT(IN)  :: pchunk
  REAL(rp), DIMENSION(pchunk), INTENT(IN)      :: Y_R,Y_PHI,Y_Z
  real(rp), DIMENSION(pchunk) :: Y_PHI_mod
  TYPE(FIELDS), INTENT(IN)                               :: F
  REAL(rp), DIMENSION(pchunk),  INTENT(OUT)   :: B_R,B_PHI,B_Z
  REAL(rp), DIMENSION(pchunk,4)   :: BR,BPHI,BZ
  REAL(rp), DIMENSION(pchunk)   :: dBRdR,dBPHIdR,dBZdR
  REAL(rp), DIMENSION(pchunk)   :: dBRdPHI,dBPHIdPHI,dBZdPHI
  REAL(rp), DIMENSION(pchunk)   :: dBRdZ,dBPHIdZ,dBZdZ
  REAL(rp), DIMENSION(pchunk),  INTENT(OUT)   :: gradB_R,gradB_PHI,gradB_Z
  REAL(rp), DIMENSION(pchunk),  INTENT(OUT)   :: curlb_R,curlb_PHI,curlb_Z
  REAL(rp), DIMENSION(pchunk),  INTENT(OUT)   :: E_R,E_PHI,E_Z
  REAL(rp), DIMENSION(pchunk),  INTENT(OUT)   :: PSIp
  REAL(rp), DIMENSION(pchunk)   :: Bmag
  INTEGER                                                :: cc
  REAL(rp), DIMENSION(pchunk,6)  :: A
  INTEGER(is),DIMENSION(pchunk),INTENT(INOUT)   :: flag_cache

  Y_PHI_mod=modulo(Y_PHI,2._rp*C_PI)
  
  call check_if_in_fields_domain_p(pchunk,F,Y_R,Y_PHI_mod,Y_Z,flag_cache)

  call EZspline_interp(bfield_3d%R,bfield_3d%PHI,bfield_3d%Z, &
       efield_3d%R,efield_3d%PHI,efield_3d%Z,bfield_3d%A, &
       pchunk,Y_R,Y_PHI_mod,Y_Z,BR,BPHI,BZ,E_R,E_PHI,E_Z,PSIp,ezerr)
  call EZspline_error(ezerr)

  
!  call EZspline_interp(bfield_2d%R,bfield_2d%PHI,bfield_2d%Z, &
!       dbdR_2d%R,dbdR_2d%PHI,dBdR_2d%Z, &
!       dbdPHI_2d%R,dbdPHI_2d%PHI,dbdPHI_2d%Z, &
!       dbdZ_2d%R,dbdZ_2d%PHI,dbdZ_2d%Z, &
!       efield_2d%R,efield_2d%PHI,efield_2d%Z,8,Y_R,Y_Z,B_R,B_PHI,B_Z, &
!       dBRdR,dBPHIdR,dBZdR,dBRdPHI,dBPHIdPHI,dBZdPHI,dBRdZ,dBPHIdZ,dBZdZ, &
!       E_R,E_PHI,E_Z,ezerr)
!  call EZspline_error(ezerr)

  !$OMP SIMD
!    !$OMP& aligned(PSIp,A,B_R,Y_R,B_PHI,B_Z,Bmag,gradB_R,gradB_PHI,gradB_Z, &
!    !$OMP& curlb_R,curlb_PHI,curlb_Z,E_R,E_PHI,E_Z)
  do cc=1_idef,pchunk

     B_R(cc)=BR(cc,1)
     B_PHI(cc)=BPHI(cc,1)
     B_Z(cc)=BZ(cc,1)

     dBRdR(cc)=BR(cc,2)
     dBRdPHI(cc)=BR(cc,3)
     dBRdZ(cc)=BR(cc,4)

     dBPHIdR(cc)=BPHI(cc,2)
     dBPHIdPHI(cc)=BPHI(cc,3)
     dBPHIdZ(cc)=BPHI(cc,4)

     dBZdR(cc)=BZ(cc,2)
     dBZdPHI(cc)=BZ(cc,3)
     dBZdZ(cc)=BZ(cc,4)
     
     Bmag(cc)=sqrt(B_R(cc)*B_R(cc)+B_PHI(cc)*B_PHI(cc)+B_Z(cc)*B_Z(cc))

     gradB_R(cc)=(B_R(cc)*dBRdR(cc)+B_PHI(cc)*dBPHIdR(cc)+ &
          B_Z(cc)*dBZdR(cc))/Bmag(cc)
     gradB_PHI(cc)=(B_R(cc)*dBRdPHI(cc)+B_PHI(cc)*dBPHIdPHI(cc)+ &
          B_Z(cc)*dBZdPHI(cc))/(Y_R(cc)*Bmag(cc))
     gradB_Z(cc)=(B_R(cc)*dBRdZ(cc)+B_PHI(cc)*dBPHIdZ(cc)+ &
          B_Z(cc)*dBZdZ(cc))/Bmag(cc)

     curlb_R(cc)=(Bmag(cc)*dBZdPHI(cc)/Y_R(cc)-B_Z(cc)*gradB_PHI(cc)- &
          Bmag(cc)*dBPHIdZ(cc)+B_PHI(cc)*gradB_Z(cc))/(Bmag(cc)*Bmag(cc))
     curlb_PHI(cc)=(Bmag(cc)*dBRdZ(cc)-B_R(cc)*gradB_Z(cc)- &
          Bmag(cc)*dBZdR(cc)+B_Z(cc)*gradB_R(cc))/(Bmag(cc)*Bmag(cc))
     curlb_Z(cc)=(Bmag(cc)*B_PHI(cc)/Y_R(cc)+Bmag(cc)*dBPHIdR(cc)- &
          B_PHI(cc)*gradB_R(cc)-Bmag(cc)*dBRdPHI(cc)/Y_R(cc)+ &
          B_R(cc)*gradB_PHI(cc))/(Bmag(cc)*Bmag(cc))
     
  end do
  !$OMP END SIMD
  
   
!  write(output_unit_write,'("PSIp: ",E17.10)') PSIp
  
!  write(output_unit_write,'("Y_R: ",E17.10)') Y_R
!  write(output_unit_write,'("Y_Z: ",E17.10)') Y_Z
  
!  write(output_unit_write,'("B_R: ",E17.10)') B_R
!  write(output_unit_write,'("B_PHIinterp: ",E17.10)') B_PHI
!  write(output_unit_write,'("B_Z: ",E17.10)') B_Z

end subroutine calculate_3DBdBfields1_p