advance_FP3Dinterp_vars Subroutine

public subroutine advance_FP3Dinterp_vars(params, X_X, X_Y, X_Z, V_X, V_Y, V_Z, g, m_cache, B_X, B_Y, B_Z, E_X, E_Y, E_Z, flagCon, flagCol, P, F, PSIp)

Arguments

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

Core KORC simulation parameters.

real(kind=rp), intent(in), DIMENSION(params%pchunk):: X_X
real(kind=rp), intent(in), DIMENSION(params%pchunk):: X_Y
real(kind=rp), intent(in), DIMENSION(params%pchunk):: X_Z
real(kind=rp), intent(inout), DIMENSION(params%pchunk):: V_X
real(kind=rp), intent(inout), DIMENSION(params%pchunk):: V_Y
real(kind=rp), intent(inout), DIMENSION(params%pchunk):: V_Z
real(kind=rp), intent(inout), DIMENSION(params%pchunk):: g
real(kind=rp), intent(in) :: m_cache
real(kind=rp), intent(in), DIMENSION(params%pchunk):: B_X
real(kind=rp), intent(in), DIMENSION(params%pchunk):: B_Y
real(kind=rp), intent(in), DIMENSION(params%pchunk):: B_Z
real(kind=rp), intent(in), DIMENSION(params%pchunk):: E_X
real(kind=rp), intent(in), DIMENSION(params%pchunk):: E_Y
real(kind=rp), intent(in), DIMENSION(params%pchunk):: E_Z
integer(kind=is), intent(inout), DIMENSION(params%pchunk):: flagCon
integer(kind=is), intent(inout), DIMENSION(params%pchunk):: flagCol
type(PROFILES), intent(in) :: P
type(FIELDS), intent(in) :: F
real(kind=rp), intent(in), DIMENSION(params%pchunk):: PSIp

Calls

proc~~advance_fp3dinterp_vars~~CallsGraph proc~advance_fp3dinterp_vars advance_FP3Dinterp_vars proc~include_coulombcollisions_fo_p include_CoulombCollisions_FO_p proc~advance_fp3dinterp_vars->proc~include_coulombcollisions_fo_p proc~unitvectors_p unitVectors_p proc~include_coulombcollisions_fo_p->proc~unitvectors_p proc~ca_sd CA_SD proc~include_coulombcollisions_fo_p->proc~ca_sd proc~dca_sd dCA_SD proc~include_coulombcollisions_fo_p->proc~dca_sd proc~interp_focollision_p interp_FOcollision_p proc~include_coulombcollisions_fo_p->proc~interp_focollision_p proc~get_random get_random proc~include_coulombcollisions_fo_p->proc~get_random proc~cf_sd CF_SD proc~include_coulombcollisions_fo_p->proc~cf_sd proc~cb_ei_sd CB_ei_SD proc~include_coulombcollisions_fo_p->proc~cb_ei_sd proc~cb_ee_sd CB_ee_SD proc~include_coulombcollisions_fo_p->proc~cb_ee_sd proc~analytical_profiles_p analytical_profiles_p proc~include_coulombcollisions_fo_p->proc~analytical_profiles_p cart_to_cyl_p cart_to_cyl_p proc~include_coulombcollisions_fo_p->cart_to_cyl_p proc~psi psi proc~ca_sd->proc~psi proc~vte VTe proc~ca_sd->proc~vte proc~gammacee Gammacee proc~ca_sd->proc~gammacee proc~dca_sd->proc~psi proc~dca_sd->proc~vte proc~dca_sd->proc~gammacee ezspline_error ezspline_error proc~interp_focollision_p->ezspline_error ezspline_interp ezspline_interp proc~interp_focollision_p->ezspline_interp proc~check_if_in_profiles_domain_p check_if_in_profiles_domain_p proc~interp_focollision_p->proc~check_if_in_profiles_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~cf_sd->proc~psi proc~cf_sd->proc~vte proc~cf_sd->proc~gammacee proc~clogee CLogee proc~cf_sd->proc~clogee proc~h_j h_j proc~cf_sd->proc~h_j 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~cb_ei_sd->proc~clogee 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~~advance_fp3dinterp_vars~~CalledByGraph proc~advance_fp3dinterp_vars advance_FP3Dinterp_vars proc~adv_fointerp_top adv_FOinterp_top proc~adv_fointerp_top->proc~advance_fp3dinterp_vars proc~adv_fofio_top adv_FOfio_top proc~adv_fofio_top->proc~advance_fp3dinterp_vars program~main main program~main->proc~adv_fointerp_top program~main->proc~adv_fofio_top

Contents


Source Code

  subroutine advance_FP3Dinterp_vars(params,X_X,X_Y,X_Z,V_X,V_Y,V_Z,g, &
       m_cache,B_X,B_Y,B_Z,E_X,E_Y,E_Z,flagCon,flagCol,P,F,PSIp)    
    TYPE(KORC_PARAMS), INTENT(INOUT)                              :: params
    !! Core KORC simulation parameters.
    TYPE(PROFILES), INTENT(IN)                                 :: P
    TYPE(FIELDS), INTENT(IN)      :: F
    INTEGER                                                    :: cc,pchunk
    !! Chunk iterator.
    INTEGER(ip)                                                    :: tt
    !! time iterator.
    REAL(rp),DIMENSION(params%pchunk), INTENT(IN)  :: X_X,X_Y,X_Z
    REAL(rp),DIMENSION(params%pchunk)  :: Y_R,Y_PHI,Y_Z
    REAL(rp),DIMENSION(params%pchunk), INTENT(IN)  :: E_X,E_Y,E_Z,PSIp
    REAL(rp),DIMENSION(params%pchunk), INTENT(IN)  :: B_X,B_Y,B_Z
    REAL(rp),DIMENSION(params%pchunk) :: U_X,U_Y,U_Z
    REAL(rp),DIMENSION(params%pchunk), INTENT(INOUT)  :: V_X,V_Y,V_Z
    REAL(rp),DIMENSION(params%pchunk) :: ne,Te,Zeff
    REAL(rp),DIMENSION(params%pchunk),INTENT(INOUT) :: g
    INTEGER(is),DIMENSION(params%pchunk),INTENT(INOUT) :: flagCon,flagCol
    REAL(rp),intent(in) :: m_cache

    pchunk=params%pchunk  

    !$OMP SIMD
    !    !$OMP& aligned(U_X,U_Y,U_Z,V_X,V_Y,V_Z,g)
    do cc=1_idef,pchunk
       U_X(cc)=V_X(cc)*g(cc)
       U_Y(cc)=V_Y(cc)*g(cc)
       U_Z(cc)=V_Z(cc)*g(cc)
    end do
    !$OMP END SIMD

    do tt=1_ip,params%t_skip

       call include_CoulombCollisions_FO_p(tt,params,X_X,X_Y,X_Z, &
            U_X,U_Y,U_Z,B_X,B_Y,B_Z,m_cache,P,F,flagCon,flagCol,PSIp)

    end do

    !$OMP SIMD
    !    !$OMP& aligned(U_X,U_Y,U_Z,V_X,V_Y,V_Z,g)
    do cc=1_idef,pchunk

       g(cc)=sqrt(1._rp+U_X(cc)*U_X(cc)+U_Y(cc)*U_Y(cc)+U_Z(cc)*U_Z(cc))

       V_X(cc)=U_X(cc)/g(cc)
       V_Y(cc)=U_Y(cc)/g(cc)
       V_Z(cc)=U_Z(cc)/g(cc)
    end do
    !$OMP END SIMD

  end subroutine advance_FP3Dinterp_vars