adv_GCinterp_psi_top_FS Subroutine

public subroutine adv_GCinterp_psi_top_FS(params, spp, P, F)

Arguments

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

Core KORC simulation parameters.

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.

type(PROFILES), intent(in) :: P
type(FIELDS), intent(inout) :: F

Contents


Source Code

  subroutine adv_GCinterp_psi_top_FS(params,spp,P,F)
    
    TYPE(KORC_PARAMS), INTENT(INOUT)                           :: params
    !! Core KORC simulation parameters.
    TYPE(PROFILES), INTENT(IN)                                 :: P
    TYPE(FIELDS), INTENT(INOUT)                                   :: F
    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(8)               :: Bmag
    REAL(rp),DIMENSION(8) :: Y_R,Y_PHI,Y_Z
    REAL(rp),DIMENSION(8) :: B_R,B_PHI,B_Z
    REAL(rp),DIMENSION(8) :: E_R,E_PHI,E_Z
    REAL(rp),DIMENSION(8) :: ne,Te,Zeff    
    REAL(rp),DIMENSION(8) :: V_PLL,V_MU
    REAL(rp),DIMENSION(8) :: PSIp
    REAL(rp),DIMENSION(8) :: curlb_R,curlb_PHI,curlb_Z
    REAL(rp),DIMENSION(8) :: gradB_R,gradB_PHI,gradB_Z
    INTEGER(is),DIMENSION(8) :: flag_cache
    REAL(rp) :: m_cache,q_cache,B0,EF0,R0,q0,lam,ar

    
    INTEGER                                                    :: ii
    !! Species iterator.
    INTEGER                                                    :: pp
    !! Particles iterator.
    INTEGER                                                    :: cc
    !! Chunk iterator.
    INTEGER(ip)                                                    :: tt
    INTEGER(ip)                                                    :: ttt
    !! time iterator.

    real(rp),dimension(F%dim_1D) :: Vden,Vdenave,VdenOMP
    !! time iterator.
 

    do ii = 1_idef,params%num_species      

       q_cache=spp(ii)%q
       m_cache=spp(ii)%m

       do ttt=1_ip,params%t_it_SC

          VdenOMP=0._rp

          !$OMP PARALLEL DO default(none) &
          !$OMP& FIRSTPRIVATE(q_cache,m_cache) &
          !$OMP& SHARED(params,ii,spp,P,F) &
          !$OMP& PRIVATE(pp,tt,Bmag,cc,Y_R,Y_PHI,Y_Z,V_PLL,V_MU,B_R,B_PHI,B_Z, &
          !$OMP& flag_cache,E_PHI,PSIp,curlb_R,curlb_PHI,curlb_Z, &
          !$OMP& gradB_R,gradB_PHI,gradB_Z,ne, &
          !$OMP& Vden,Vdenave) &
          !$OMP& REDUCTION(+:VdenOMP)
          do pp=1_idef,spp(ii)%ppp,8

             !          write(6,'("pp: ",I16)') pp

             !$OMP SIMD
             do cc=1_idef,8_idef
                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)

                flag_cache(cc)=spp(ii)%vars%flag(pp-1+cc)           
             end do
             !$OMP END SIMD

             if (.not.params%FokPlan) then
                Vdenave=0._rp
                do tt=1_ip,params%t_skip
                   
                   call advance_GCinterp_psi_vars_FS(spp(ii)%vars,pp,tt, &
                        params, &
                        Y_R,Y_PHI,Y_Z,V_PLL,V_MU,q_cache,m_cache,flag_cache, &
                        F,P,B_R,B_PHI,B_Z,E_PHI,PSIp,curlb_R,curlb_PHI, &
                        curlb_Z,gradB_R,gradB_PHI,gradB_Z)

                   call calculate_SC_p_FS(params,F,B_R,B_PHI,B_Z,PSIp, &
                        V_PLL,V_MU,m_cache,flag_cache,Vden)

!                   write(6,*) 'pre-Vdenave',Vdenave(F%dim_1D)
                   Vdenave=(Vdenave*REAL(tt-1_ip)+Vden)/REAL(tt)

!                   write(6,*) 'Vden',Vden(F%dim_1D)
!                   write(6,*) 'post-Vdenave',Vdenave(F%dim_1D)
!                   if (pp.eq.9_idef) write(6,*) 'Vdenave',Vdenave(F%dim_1D)

                end do !timestep iterator

!                write(6,*) 'Vdenave',Vdenave(F%dim_1D)

                VdenOMP=VdenOMP+Vdenave

!                write(6,*) 'VdenOMP',VdenOMP(F%dim_1D)
                
                !$OMP SIMD
                do cc=1_idef,8_idef
                   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%flag(pp-1+cc)=flag_cache(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%gradB(pp-1+cc,1) = gradB_R(cc)
                   spp(ii)%vars%gradB(pp-1+cc,2) = gradB_PHI(cc)
                   spp(ii)%vars%gradB(pp-1+cc,3) = gradB_Z(cc)

                   spp(ii)%vars%curlb(pp-1+cc,1) = curlb_R(cc)
                   spp(ii)%vars%curlb(pp-1+cc,2) = curlb_PHI(cc)
                   spp(ii)%vars%curlb(pp-1+cc,3) = curlb_Z(cc)

                   spp(ii)%vars%E(pp-1+cc,2) = E_PHI(cc)
                   spp(ii)%vars%PSI_P(pp-1+cc) = PSIp(cc)                
                end do
                !$OMP END SIMD

             else

                call advance_FPinterp_vars(params,Y_R,Y_PHI, &
                     Y_Z,V_PLL,V_MU,m_cache,flag_cache,F,P,E_PHI,ne,PSIp)

                !$OMP SIMD
                do cc=1_idef,8_idef
                   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%flag(pp-1+cc)=flag_cache(cc)

                   spp(ii)%vars%E(pp-1+cc,2) = E_PHI(cc)
                end do
                !$OMP END SIMD

             end if


             !$OMP SIMD
             do cc=1_idef,8_idef
                B_R(cc)=spp(ii)%vars%B(pp-1+cc,1)
                B_PHI(cc)=spp(ii)%vars%B(pp-1+cc,2)
                B_Z(cc)=spp(ii)%vars%B(pp-1+cc,3)

                Bmag(cc)=sqrt(B_R(cc)*B_R(cc)+B_PHI(cc)*B_PHI(cc)+B_Z(cc)*B_Z(cc))

                spp(ii)%vars%g(pp-1+cc)=sqrt(1+V_PLL(cc)**2+ &
                     2*V_MU(cc)*Bmag(cc))

                spp(ii)%vars%eta(pp-1+cc) = atan2(sqrt(2*m_cache*Bmag(cc)* &
                     spp(ii)%vars%V(pp-1+cc,2)),spp(ii)%vars%V(pp-1+cc,1))* &
                     180.0_rp/C_PI
             end do
             !$OMP END SIMD

          end do !particle chunk iterator
          !$OMP END PARALLEL DO

          !write(6,*) 'VdenOMP',VdenOMP(F%dim_1D)


          call calculate_SC_E1D_FS(params,F,VdenOMP)             


       end do

    end do !species iterator
    
  end subroutine adv_GCinterp_psi_top_FS