adv_GCeqn_top Subroutine

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

Arguments

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


Calls

proc~~adv_gceqn_top~~CallsGraph proc~adv_gceqn_top adv_GCeqn_top proc~calculate_sc_e1d calculate_SC_E1D proc~adv_gceqn_top->proc~calculate_sc_e1d proc~analytical_fields_bmag_p analytical_fields_Bmag_p proc~adv_gceqn_top->proc~analytical_fields_bmag_p proc~include_coulombcollisions_gc_p include_CoulombCollisions_GC_p proc~adv_gceqn_top->proc~include_coulombcollisions_gc_p proc~rad2deg~2 rad2deg proc~adv_gceqn_top->proc~rad2deg~2 proc~calculate_sc_p calculate_SC_p proc~adv_gceqn_top->proc~calculate_sc_p mpi_allreduce mpi_allreduce proc~calculate_sc_e1d->mpi_allreduce proc~initialize_sc1d_field_interpolant initialize_SC1D_field_interpolant proc~calculate_sc_e1d->proc~initialize_sc1d_field_interpolant proc~calculate_gcfieldswe_p calculate_GCfieldswE_p proc~include_coulombcollisions_gc_p->proc~calculate_gcfieldswe_p proc~dca_sd dCA_SD proc~include_coulombcollisions_gc_p->proc~dca_sd proc~interp_fields_3d_p interp_fields_3D_p proc~include_coulombcollisions_gc_p->proc~interp_fields_3d_p proc~get_random get_random proc~include_coulombcollisions_gc_p->proc~get_random proc~calculate_gcfields_p calculate_GCfields_p proc~include_coulombcollisions_gc_p->proc~calculate_gcfields_p proc~analytical_profiles_p analytical_profiles_p proc~include_coulombcollisions_gc_p->proc~analytical_profiles_p proc~interp_hcollision_p interp_Hcollision_p proc~include_coulombcollisions_gc_p->proc~interp_hcollision_p proc~cb_ei_sd CB_ei_SD proc~include_coulombcollisions_gc_p->proc~cb_ei_sd proc~interp_fields_p interp_fields_p proc~include_coulombcollisions_gc_p->proc~interp_fields_p proc~calculate_gcfields_p_fs calculate_GCfields_p_FS proc~include_coulombcollisions_gc_p->proc~calculate_gcfields_p_fs proc~interp_focollision_p interp_FOcollision_p proc~include_coulombcollisions_gc_p->proc~interp_focollision_p proc~analytical_fields_gc_p analytical_fields_GC_p proc~include_coulombcollisions_gc_p->proc~analytical_fields_gc_p proc~calculate_2dbdbfields_p calculate_2DBdBfields_p proc~include_coulombcollisions_gc_p->proc~calculate_2dbdbfields_p proc~cb_ei_sd_fio CB_ei_SD_FIO proc~include_coulombcollisions_gc_p->proc~cb_ei_sd_fio proc~cf_sd_fio CF_SD_FIO proc~include_coulombcollisions_gc_p->proc~cf_sd_fio proc~ca_sd CA_SD proc~include_coulombcollisions_gc_p->proc~ca_sd proc~calculate_gcfields_2x1t_p calculate_GCfields_2x1t_p proc~include_coulombcollisions_gc_p->proc~calculate_gcfields_2x1t_p proc~add_analytical_e_p add_analytical_E_p proc~include_coulombcollisions_gc_p->proc~add_analytical_e_p proc~cf_sd CF_SD proc~include_coulombcollisions_gc_p->proc~cf_sd proc~cb_ee_sd CB_ee_SD proc~include_coulombcollisions_gc_p->proc~cb_ee_sd ezspline_error ezspline_error proc~calculate_gcfieldswe_p->ezspline_error proc~check_if_in_fields_domain_p check_if_in_fields_domain_p proc~calculate_gcfieldswe_p->proc~check_if_in_fields_domain_p ezspline_derivative ezspline_derivative proc~calculate_gcfieldswe_p->ezspline_derivative proc~psi psi proc~dca_sd->proc~psi proc~vte VTe proc~dca_sd->proc~vte proc~gammacee Gammacee proc~dca_sd->proc~gammacee proc~interp_fields_3d_p->ezspline_error ezspline_interp ezspline_interp proc~interp_fields_3d_p->ezspline_interp proc~interp_fields_3d_p->proc~check_if_in_fields_domain_p omp_get_thread_num omp_get_thread_num proc~get_random->omp_get_thread_num interface~random_get_number random_get_number proc~get_random->interface~random_get_number proc~calculate_gcfields_p->ezspline_error proc~calculate_gcfields_p->proc~check_if_in_fields_domain_p proc~calculate_gcfields_p->ezspline_derivative proc~interp_hcollision_p->ezspline_error proc~interp_hcollision_p->ezspline_interp proc~check_if_in_profiles_domain_p check_if_in_profiles_domain_p proc~interp_hcollision_p->proc~check_if_in_profiles_domain_p proc~cb_ei_sd->proc~vte proc~clogei CLogei proc~cb_ei_sd->proc~clogei proc~g_j g_j proc~cb_ei_sd->proc~g_j proc~cb_ei_sd->proc~gammacee proc~clogee CLogee proc~cb_ei_sd->proc~clogee proc~interp_fields_p->ezspline_error proc~interp_fields_p->ezspline_interp proc~interp_fields_p->proc~check_if_in_fields_domain_p proc~calculate_gcfields_p_fs->ezspline_error proc~calculate_gcfields_p_fs->proc~check_if_in_fields_domain_p proc~calculate_gcfields_p_fs->ezspline_derivative proc~initialize_sc1d_field_interpolant->ezspline_error ezspline_allocated ezspline_allocated proc~initialize_sc1d_field_interpolant->ezspline_allocated ezspline_init ezspline_init proc~initialize_sc1d_field_interpolant->ezspline_init ezspline_setup ezspline_setup proc~initialize_sc1d_field_interpolant->ezspline_setup proc~interp_focollision_p->ezspline_error proc~interp_focollision_p->ezspline_interp proc~interp_focollision_p->proc~check_if_in_profiles_domain_p proc~calculate_2dbdbfields_p->ezspline_error proc~calculate_2dbdbfields_p->ezspline_interp proc~calculate_2dbdbfields_p->proc~check_if_in_fields_domain_p proc~cb_ei_sd_fio->proc~vte proc~cb_ei_sd_fio->proc~clogei proc~cb_ei_sd_fio->proc~g_j proc~cb_ei_sd_fio->proc~gammacee proc~cb_ei_sd_fio->proc~clogee proc~cf_sd_fio->proc~psi proc~cf_sd_fio->proc~vte proc~cf_sd_fio->proc~gammacee proc~cf_sd_fio->proc~clogee proc~h_j h_j proc~cf_sd_fio->proc~h_j proc~ca_sd->proc~psi proc~ca_sd->proc~vte proc~ca_sd->proc~gammacee proc~calculate_gcfields_2x1t_p->ezspline_error proc~calculate_gcfields_2x1t_p->proc~check_if_in_fields_domain_p proc~calculate_gcfields_2x1t_p->ezspline_derivative proc~cf_sd->proc~psi proc~cf_sd->proc~vte proc~cf_sd->proc~gammacee proc~cf_sd->proc~clogee proc~cf_sd->proc~h_j proc~cb_ee_sd->proc~psi proc~cb_ee_sd->proc~vte proc~cb_ee_sd->proc~gammacee proc~delta delta proc~cb_ee_sd->proc~delta proc~clogei->proc~vte proc~clog0 CLog0 proc~clogei->proc~clog0 proc~gammacee->proc~clogee proc~clogee->proc~vte proc~clogee->proc~clog0 proc~delta->proc~vte

Called by

proc~~adv_gceqn_top~~CalledByGraph proc~adv_gceqn_top adv_GCeqn_top program~main main program~main->proc~adv_gceqn_top

Contents

Source Code


Source Code

  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