adv_GCinterp_2DBdB_top Subroutine

public subroutine adv_GCinterp_2DBdB_top(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(in) :: F

Calls

proc~~adv_gcinterp_2dbdb_top~~CallsGraph proc~adv_gcinterp_2dbdb_top adv_GCinterp_2DBdB_top proc~advance_fpinterp_vars advance_FPinterp_vars proc~adv_gcinterp_2dbdb_top->proc~advance_fpinterp_vars proc~include_coulombcollisions_gc_p include_CoulombCollisions_GC_p proc~advance_fpinterp_vars->proc~include_coulombcollisions_gc_p 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~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_gcinterp_2dbdb_top~~CalledByGraph proc~adv_gcinterp_2dbdb_top adv_GCinterp_2DBdB_top program~main main program~main->proc~adv_gcinterp_2dbdb_top

Contents


Source Code

  subroutine adv_GCinterp_2DBdB_top(params,spp,P,F)

    TYPE(KORC_PARAMS), INTENT(INOUT)                           :: params
    !! Core KORC simulation parameters.
    TYPE(PROFILES), INTENT(IN)                                 :: P
    TYPE(FIELDS), INTENT(IN)                                   :: 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(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
    REAL(rp),DIMENSION(params%pchunk) :: E_R,E_PHI,E_Z
    REAL(rp),DIMENSION(params%pchunk) :: ne,Te,Zeff    
    REAL(rp),DIMENSION(params%pchunk) :: V_PLL,V_MU
    REAL(rp),DIMENSION(params%pchunk) :: PSIp
    REAL(rp),DIMENSION(params%pchunk) :: curlb_R,curlb_PHI,curlb_Z
    REAL(rp),DIMENSION(params%pchunk) :: gradB_R,gradB_PHI,gradB_Z
    INTEGER(is),DIMENSION(params%pchunk) :: flagCon,flagCol
    REAL(rp) :: m_cache,q_cache,B0,EF0,R0,q0,lam,ar


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


    do ii = 1_idef,params%num_species      

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

       !$OMP PARALLEL DO default(none) &
       !$OMP& FIRSTPRIVATE(q_cache,m_cache,pchunk) &
       !$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& flagCon,flagCol,E_PHI,PSIp,curlb_R,curlb_PHI,curlb_Z, &
       !$OMP& gradB_R,gradB_PHI,gradB_Z,ne)
       do pp=1_idef,spp(ii)%ppp,pchunk

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

          !$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)

             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       
             do tt=1_ip,params%t_skip
                call advance_GCinterp_2DBdB_vars(spp(ii)%vars,pp,tt,params, &
                     Y_R,Y_PHI,Y_Z,V_PLL,V_MU,q_cache,m_cache,flagCon,flagCol, &
                     F,P,B_R,B_PHI,B_Z,E_PHI,curlb_R,curlb_PHI,curlb_Z, &
                     gradB_R,gradB_PHI,gradB_Z,PSIp)
             end do !timestep iterator


             !$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%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,flagCon,flagCol,F,P,E_PHI,ne,PSIp)

             !$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%flagCon(pp-1+cc)=flagCon(cc)
                spp(ii)%vars%flagCol(pp-1+cc)=flagCol(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,pchunk
             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

    end do !species iterator

  end subroutine adv_GCinterp_2DBdB_top