Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(KORC_PARAMS), | intent(inout) | :: | params | Core KORC simulation parameters. |
||
type(FIELDS), | intent(inout) | :: | F | An instance of the KORC derived type FIELDS. |
||
type(PROFILES), | intent(in) | :: | P | An instance of the KORC derived type PROFILES. |
||
type(SPECIES), | intent(inout), | DIMENSION(:), ALLOCATABLE | :: | spp | An instance of the derived type SPECIES containing all the parameters and simulation variables of the different species in the simulation. |
subroutine adv_GCeqn_top(params,F,P,spp)
TYPE(KORC_PARAMS), INTENT(INOUT) :: params
!! Core KORC simulation parameters.
TYPE(FIELDS), INTENT(INOUT) :: F
!! An instance of the KORC derived type FIELDS.
TYPE(PROFILES), INTENT(IN) :: P
!! An instance of the KORC derived type PROFILES.
TYPE(SPECIES), DIMENSION(:), ALLOCATABLE, INTENT(INOUT) :: spp
!! An instance of the derived type SPECIES containing all the parameters
!! and simulation variables of the different species in the simulation.
REAL(rp), DIMENSION(params%pchunk) :: Bmag
REAL(rp),DIMENSION(params%pchunk) :: Y_R,Y_PHI,Y_Z
REAL(rp),DIMENSION(params%pchunk) :: B_R,B_PHI,B_Z,E_PHI
REAL(rp),DIMENSION(params%pchunk) :: PSIp,ne,Te
REAL(rp),DIMENSION(params%pchunk) :: V_PLL,V_MU
REAL(rp) :: B0,EF0,R0,q0,lam,ar,m_cache,q_cache,ne0,Te0,Zeff0
INTEGER(is),DIMENSION(params%pchunk) :: flagCon,flagCol,flagRE
LOGICAL :: ss_collisions
!! Logical variable that indicates if collisions are included in
!! the simulation.
INTEGER :: ii
!! Species iterator.
INTEGER :: pp
!! Particles iterator.
INTEGER :: cc,pchunk,achunk
!! Chunk iterator.
INTEGER(ip) :: tt
INTEGER(ip) :: ttt,tttt
!! time iterator.
real(rp),dimension(F%dim_1D) :: Vden,Vdenave,VdenOMP
INTEGER :: newREs
do ii = 1_idef,params%num_species
pchunk=params%pchunk
q_cache=spp(ii)%q
m_cache=spp(ii)%m
do ttt=1_ip,params%t_it_SC
VdenOMP=0._rp
if(.not.params%LargeCollisions) then
!$OMP PARALLEL DO default(none) &
!$OMP& FIRSTPRIVATE(E0,q_cache,m_cache,pchunk) &
!$OMP& shared(F,P,params,ii,spp) &
!$OMP& PRIVATE(pp,tt,ttt,Bmag,cc,Y_R,Y_PHI,Y_Z,V_PLL,V_MU, &
!$OMP& flagCon,flagCol,B_R,B_PHI,B_Z,E_PHI,PSIp,ne, &
!$OMP& Vden,Vdenave) &
!$OMP& REDUCTION(+:VdenOMP)
do pp=1_idef,spp(ii)%ppp,pchunk
!$OMP SIMD
do cc=1_idef,pchunk
Y_R(cc)=spp(ii)%vars%Y(pp-1+cc,1)
Y_PHI(cc)=spp(ii)%vars%Y(pp-1+cc,2)
Y_Z(cc)=spp(ii)%vars%Y(pp-1+cc,3)
V_PLL(cc)=spp(ii)%vars%V(pp-1+cc,1)
V_MU(cc)=spp(ii)%vars%V(pp-1+cc,2)
PSIp(cc)=spp(ii)%vars%PSI_p(pp-1+cc)
ne(cc)=spp(ii)%vars%ne(pp-1+cc)
flagCon(cc)=spp(ii)%vars%flagCon(pp-1+cc)
flagCol(cc)=spp(ii)%vars%flagCol(pp-1+cc)
end do
!$OMP END SIMD
if (.not.params%FokPlan) then
Vdenave=0._rp
do tt=1_ip,params%t_skip
! write(output_unit_write,*) params%mpi_params%rank,'Y_R',Y_R
call advance_GCeqn_vars(spp(ii)%vars,pp, &
tt+params%t_skip*(ttt-1),params, &
Y_R,Y_PHI, Y_Z,V_PLL,V_MU,flagCon,flagCol,q_cache,m_cache, &
B_R,B_PHI,B_Z,F,P,PSIp,E_PHI)
! write(output_unit_write,*) params%mpi_params%rank,'Y_R',Y_R
if (params%collisions) then
call include_CoulombCollisions_GC_p(tt+params%t_skip*(ttt-1),params, &
Y_R,Y_PHI,Y_Z,V_PLL,V_MU,m_cache,flagCon,flagCol, &
F,P,E_PHI,ne,PSIp)
end if
if (params%SC_E) then
call calculate_SC_p(params,F,B_R,B_PHI,B_Z,Y_R,Y_Z, &
V_PLL,V_MU,m_cache,flagCon,flagCol,Vden)
Vdenave=(Vdenave*REAL(tt-1_ip)+Vden)/REAL(tt)
end if
end do !timestep iterator
VdenOMP=VdenOMP+Vdenave
!$OMP SIMD
do cc=1_idef,pchunk
spp(ii)%vars%Y(pp-1+cc,1)=Y_R(cc)
spp(ii)%vars%Y(pp-1+cc,2)=Y_PHI(cc)
spp(ii)%vars%Y(pp-1+cc,3)=Y_Z(cc)
spp(ii)%vars%V(pp-1+cc,1)=V_PLL(cc)
spp(ii)%vars%V(pp-1+cc,2)=V_MU(cc)
spp(ii)%vars%flagCon(pp-1+cc)=flagCon(cc)
spp(ii)%vars%flagCol(pp-1+cc)=flagCol(cc)
spp(ii)%vars%B(pp-1+cc,1) = B_R(cc)
spp(ii)%vars%B(pp-1+cc,2) = B_PHI(cc)
spp(ii)%vars%B(pp-1+cc,3) = B_Z(cc)
spp(ii)%vars%PSI_P(pp-1+cc) = PSIp(cc)
spp(ii)%vars%E(pp-1+cc,2) = E_PHI(cc)
end do
!$OMP END SIMD
else
do tt=1_ip,params%t_skip
!if (mod(tt,50).eq.0) then
!write(output_unit_write,*) 'iteration',tt
!flush(output_unit_write)
!endif
call include_CoulombCollisions_GC_p(tt,params, &
Y_R,Y_PHI,Y_Z,V_PLL,V_MU,m_cache,flagCon,flagCol, &
F,P,E_PHI,ne,PSIp)
end do
!$OMP SIMD
do cc=1_idef,pchunk
spp(ii)%vars%V(pp-1+cc,1)=V_PLL(cc)
spp(ii)%vars%V(pp-1+cc,2)=V_MU(cc)
spp(ii)%vars%flagCol(pp-1+cc)=flagCol(cc)
spp(ii)%vars%ne(pp-1+cc) = ne(cc)
end do
!$OMP END SIMD
end if
call analytical_fields_Bmag_p(pchunk,F,Y_R,Y_PHI,Y_Z, &
Bmag,E_PHI)
!$OMP SIMD
do cc=1_idef,pchunk
spp(ii)%vars%g(pp-1+cc)=sqrt(1+V_PLL(cc)**2+ &
2*V_MU(cc)*Bmag(cc)*m_cache)
spp(ii)%vars%eta(pp-1+cc) = rad2deg(atan2(sqrt(2*m_cache* &
Bmag(cc)*spp(ii)%vars%V(pp-1+cc,2)), &
spp(ii)%vars%V(pp-1+cc,1)))
end do
!$OMP END SIMD
end do !particle chunk iterator
!$OMP END PARALLEL DO
else if (params%LargeCollisions.and.sample_test) then
do tt=1_ip,params%coll_per_dump
!$OMP PARALLEL DO default(none) &
!$OMP& FIRSTPRIVATE(m_cache,pchunk) &
!$OMP& shared(F,P,params,ii,spp,tt) &
!$OMP& PRIVATE(pp,Bmag,cc,Y_R,Y_PHI,Y_Z,V_PLL,V_MU, &
!$OMP& flagCon,flagCol,B_R,B_PHI,B_Z,E_PHI,PSIp,ne, &
!$OMP& achunk,Te)
do pp=1_idef,spp(ii)%pinit,pchunk
if ((spp(ii)%pRE-pp).lt.pchunk) then
achunk=spp(ii)%pRE-pp+1
else
achunk=pchunk
end if
!$OMP SIMD
do cc=1_idef,achunk
Y_R(cc)=spp(ii)%vars%Y(pp-1+cc,1)
Y_PHI(cc)=spp(ii)%vars%Y(pp-1+cc,2)
Y_Z(cc)=spp(ii)%vars%Y(pp-1+cc,3)
V_PLL(cc)=spp(ii)%vars%V(pp-1+cc,1)
V_MU(cc)=spp(ii)%vars%V(pp-1+cc,2)
PSIp(cc)=spp(ii)%vars%PSI_p(pp-1+cc)
ne(cc)=spp(ii)%vars%ne(pp-1+cc)
flagCon(cc)=spp(ii)%vars%flagCon(pp-1+cc)
flagCol(cc)=spp(ii)%vars%flagCol(pp-1+cc)
end do
!$OMP END SIMD
if (params%FokPlan) then
!if (mod(tt,50).eq.0) then
!write(6,*) 'particle',pp,'iteration',tt, &
! ' of',params%t_skip
!flush(6)
!endif
call include_CoulombCollisionsLA_GC_p(spp(ii),achunk, &
tt,params,Y_R,Y_PHI,Y_Z,V_PLL,V_MU,m_cache, &
flagCon,flagCol,F,P,E_PHI,ne,Te,PSIp)
end if
!$OMP SIMD
do cc=1_idef,achunk
spp(ii)%vars%V(pp-1+cc,1)=V_PLL(cc)
spp(ii)%vars%V(pp-1+cc,2)=V_MU(cc)
spp(ii)%vars%ne(pp-1+cc)=ne(cc)
spp(ii)%vars%flagCol(pp-1+cc)=flagCol(cc)
end do
!$OMP END SIMD
end do !particle chunk iterator
!$OMP END PARALLEL DO
end do
!$OMP PARALLEL DO default(none) &
!$OMP& FIRSTPRIVATE(m_cache,pchunk) &
!$OMP& shared(F,ii,spp) &
!$OMP& PRIVATE(pp,Bmag,cc,Y_R,Y_PHI,Y_Z,V_PLL,V_MU,achunk,E_PHI)
do pp=1_idef,spp(ii)%pRE,pchunk
if ((spp(ii)%pRE-pp).lt.pchunk) then
achunk=spp(ii)%pRE-pp+1
else
achunk=pchunk
end if
!$OMP SIMD
do cc=1_idef,achunk
Y_R(cc)=spp(ii)%vars%Y(pp-1+cc,1)
Y_PHI(cc)=spp(ii)%vars%Y(pp-1+cc,2)
Y_Z(cc)=spp(ii)%vars%Y(pp-1+cc,3)
end do
!$OMP END SIMD
!write(6,*) 'Y_R',Y_R
call analytical_fields_Bmag_p(achunk,F,Y_R,Y_PHI,Y_Z, &
Bmag,E_PHI)
!$OMP SIMD
do cc=1_idef,achunk
V_PLL(cc)=spp(ii)%vars%V(pp-1+cc,1)
V_MU(cc)=spp(ii)%vars%V(pp-1+cc,2)
spp(ii)%vars%g(pp-1+cc)=sqrt(1+V_PLL(cc)**2+ &
2*V_MU(cc)*Bmag(cc)*m_cache)
spp(ii)%vars%eta(pp-1+cc) = rad2deg(atan2(sqrt(2*m_cache* &
Bmag(cc)*spp(ii)%vars%V(pp-1+cc,2)), &
spp(ii)%vars%V(pp-1+cc,1)))
end do
!$OMP END SIMD
end do !particle chunk iterator
!$OMP END PARALLEL DO
else
do tt=1_ip,params%coll_per_dump
!$OMP PARALLEL DO default(none) &
!$OMP& FIRSTPRIVATE(m_cache,pchunk,q_cache) &
!$OMP& shared(F,P,params,ii,spp,tt) &
!$OMP& PRIVATE(pp,Bmag,cc,Y_R,Y_PHI,Y_Z,V_PLL,V_MU, &
!$OMP& flagCon,flagCol,B_R,B_PHI,B_Z,E_PHI,PSIp,ne, &
!$OMP& achunk,tttt,Te)
do pp=1_idef,spp(ii)%pRE,pchunk
if ((spp(ii)%pRE-pp).lt.pchunk) then
achunk=spp(ii)%pRE-pp+1
else
achunk=pchunk
end if
!$OMP SIMD
do cc=1_idef,achunk
Y_R(cc)=spp(ii)%vars%Y(pp-1+cc,1)
Y_PHI(cc)=spp(ii)%vars%Y(pp-1+cc,2)
Y_Z(cc)=spp(ii)%vars%Y(pp-1+cc,3)
V_PLL(cc)=spp(ii)%vars%V(pp-1+cc,1)
V_MU(cc)=spp(ii)%vars%V(pp-1+cc,2)
PSIp(cc)=spp(ii)%vars%PSI_p(pp-1+cc)
ne(cc)=spp(ii)%vars%ne(pp-1+cc)
flagCon(cc)=spp(ii)%vars%flagCon(pp-1+cc)
flagCol(cc)=spp(ii)%vars%flagCol(pp-1+cc)
end do
!$OMP END SIMD
!if (mod(tt,50).eq.0) then
!write(6,*) 'particle',pp,'iteration',tt, &
! ' of',params%t_skip
!flush(6)
!endif
if (.not.params%FokPlan) then
do tttt=1_ip,params%orbits_per_coll
call advance_GCeqn_vars(spp(ii)%vars,pp, &
tttt,params, &
Y_R,Y_PHI,Y_Z,V_PLL,V_MU,flagCon,flagCol,q_cache,m_cache, &
B_R,B_PHI,B_Z,F,P,PSIp,E_PHI)
end do
endif
call include_CoulombCollisionsLA_GC_p(spp(ii),achunk, &
tt,params,Y_R,Y_PHI,Y_Z,V_PLL,V_MU,m_cache, &
flagCon,flagCol,F,P,E_PHI,ne,Te,PSIp)
!$OMP SIMD
do cc=1_idef,achunk
spp(ii)%vars%ne(pp-1+cc)=ne(cc)
spp(ii)%vars%Y(pp-1+cc,1)=Y_R(cc)
spp(ii)%vars%Y(pp-1+cc,2)=Y_PHI(cc)
spp(ii)%vars%Y(pp-1+cc,3)=Y_Z(cc)
spp(ii)%vars%V(pp-1+cc,1)=V_PLL(cc)
spp(ii)%vars%V(pp-1+cc,2)=V_MU(cc)
spp(ii)%vars%flagCon(pp-1+cc)=flagCon(cc)
spp(ii)%vars%flagCol(pp-1+cc)=flagCol(cc)
spp(ii)%vars%B(pp-1+cc,1) = B_R(cc)
spp(ii)%vars%B(pp-1+cc,2) = B_PHI(cc)
spp(ii)%vars%B(pp-1+cc,3) = B_Z(cc)
spp(ii)%vars%PSI_P(pp-1+cc) = PSIp(cc)
spp(ii)%vars%E(pp-1+cc,2) = E_PHI(cc)
end do
!$OMP END SIMD
end do !particle chunk iterator
!$OMP END PARALLEL DO
end do
!$OMP PARALLEL DO default(none) &
!$OMP& FIRSTPRIVATE(m_cache,pchunk) &
!$OMP& shared(F,ii,spp) &
!$OMP& PRIVATE(pp,Bmag,cc,Y_R,Y_PHI,Y_Z,V_PLL,V_MU,achunk,E_PHI)
do pp=1_idef,spp(ii)%pRE,pchunk
if ((spp(ii)%pRE-pp).lt.pchunk) then
achunk=spp(ii)%pRE-pp+1
else
achunk=pchunk
end if
!$OMP SIMD
do cc=1_idef,achunk
Y_R(cc)=spp(ii)%vars%Y(pp-1+cc,1)
Y_PHI(cc)=spp(ii)%vars%Y(pp-1+cc,2)
Y_Z(cc)=spp(ii)%vars%Y(pp-1+cc,3)
end do
!$OMP END SIMD
!write(6,*) 'Y_R',Y_R
call analytical_fields_Bmag_p(achunk,F,Y_R,Y_PHI,Y_Z, &
Bmag,E_PHI)
!$OMP SIMD
do cc=1_idef,achunk
V_PLL(cc)=spp(ii)%vars%V(pp-1+cc,1)
V_MU(cc)=spp(ii)%vars%V(pp-1+cc,2)
spp(ii)%vars%g(pp-1+cc)=sqrt(1+V_PLL(cc)**2+ &
2*V_MU(cc)*Bmag(cc)*m_cache)
spp(ii)%vars%eta(pp-1+cc) = rad2deg(atan2(sqrt(2*m_cache* &
Bmag(cc)*spp(ii)%vars%V(pp-1+cc,2)), &
spp(ii)%vars%V(pp-1+cc,1)))
end do
!$OMP END SIMD
end do !particle chunk iterator
!$OMP END PARALLEL DO
endif
end do
if (params%SC_E) then
call calculate_SC_E1D(params,F,VdenOMP)
end if
end do !species iterator
end subroutine adv_GCeqn_top