XGC1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Public Attributes | List of all members
mom_module Module Reference

Moment generating function from f0_f data. More...

Public Member Functions

subroutine set_mom_module_ptrs (den_local_cpp, upara_local_cpp, uperp_local_cpp, temp_local_cpp, Kperp_local_cpp, Kpara_local_cpp, temp_para_local_cpp, den_all_cpp, temp_all_cpp, temp_perp_all_cpp, temp_para_all_cpp, upara_all_cpp)
 

Public Attributes

real(kind=8), dimension(:,:),
pointer 
den_local
 density, \( \int f d\mathbf{v} \) local to compute node [Nsp,Nlocalnode] More...
 
real(kind=8), dimension(:,:),
pointer 
upara_local
 parallel bulk velocity, \( \int v_\parallel f d\mathbf{v} \) local to compute node [Nsp,Nlocalnode] More...
 
real(kind=8), dimension(:,:),
pointer 
uperp_local
 perpendicular bulk velocity, \( \int v_\perp f d\mathbf{v} \) local to compute node [Nsp,Nlocalnode] More...
 
real(kind=8), dimension(:,:),
pointer 
temp_local
 temperature, \( 2/3 \int 1/2 (\mathbf{v} - u_\parallel)^2 f d\mathbf{v} \), local to compute node [Nsp,Nlocalnode] More...
 
real(kind=8), dimension(:,:),
pointer 
kperp_local
 perpendicular kinetic energy, \( \int 1/2 v_\perp^2 f d\mathbf{v} \), local to compute node [Nsp,Nlocalnode] More...
 
real(kind=8), dimension(:,:),
pointer 
kpara_local
 parllel kinetic energy, \( \int v_\parallel^2 d\mathbf{v} \), local to compute node [Nsp,Nlocalnode] More...
 
real(kind=8), dimension(:,:),
pointer 
temp_para_local
 parallel temperature, \( \int (v_\parallel-u_\parallel)^2 d\mathbf{v} \), local to compute node [Nsp,Nlocalnode] More...
 
real(kind=8), dimension(:,:),
pointer 
den_all
 density, \( \int f d\mathbf{v} \), on all nodes [Nsp,Nnode] More...
 
real(kind=8), dimension(:,:),
pointer 
upara_all
 parallel bulk velocity, \( \int v_\parallel f d\mathbf{v} \), on all nodes [Nsp,Nnode] More...
 
real(kind=8), dimension(:,:),
pointer 
temp_all
 temperature, \( 2/3 \int 1/2 (\mathbf{v} - u_\parallel)^2 f d\mathbf{v} \), on all nodes [Nsp,Nnode] More...
 
real(kind=8), dimension(:,:),
pointer 
temp_perp_all
 perpendicular kinetic energy, \( \int 1/2 v_\perp^2 f d\mathbf{v} \), on all nodes [Nsp,Nnode] More...
 
real(kind=8), dimension(:,:),
pointer 
temp_para_all
 parllel kinetic energy, \( \int v_\parallel^2 d\mathbf{v} \), on all nodes [Nsp,Nnode] More...
 

Detailed Description

Moment generating function from f0_f data.

Member Function/Subroutine Documentation

subroutine mom_module::set_mom_module_ptrs ( real (kind=8), dimension(ptl_isp:ptl_nsp,f0_inode1:f0_inode2), target  den_local_cpp,
real (kind=8), dimension(ptl_isp:ptl_nsp,f0_inode1:f0_inode2), target  upara_local_cpp,
real (kind=8), dimension(ptl_isp:ptl_nsp,f0_inode1:f0_inode2), target  uperp_local_cpp,
real (kind=8), dimension(ptl_isp:ptl_nsp,f0_inode1:f0_inode2), target  temp_local_cpp,
real (kind=8), dimension(ptl_isp:ptl_nsp,f0_inode1:f0_inode2), target  Kperp_local_cpp,
real (kind=8), dimension(ptl_isp:ptl_nsp,f0_inode1:f0_inode2), target  Kpara_local_cpp,
real (kind=8), dimension(ptl_isp:ptl_nsp,f0_inode1:f0_inode2), target  temp_para_local_cpp,
real (kind=8), dimension(ptl_isp:ptl_nsp, grid_global%nnode), target  den_all_cpp,
real (kind=8), dimension(ptl_isp:ptl_nsp, grid_global%nnode), target  temp_all_cpp,
real (kind=8), dimension(ptl_isp:ptl_nsp, grid_global%nnode), target  temp_perp_all_cpp,
real (kind=8), dimension(ptl_isp:ptl_nsp, grid_global%nnode), target  temp_para_all_cpp,
real (kind=8), dimension(ptl_isp:ptl_nsp, grid_global%nnode), target  upara_all_cpp 
)

Member Data Documentation

real (kind=8), dimension(:,:), pointer mom_module::den_all

density, \( \int f d\mathbf{v} \), on all nodes [Nsp,Nnode]

real (kind=8), dimension(:,:), pointer mom_module::den_local

density, \( \int f d\mathbf{v} \) local to compute node [Nsp,Nlocalnode]

real (kind=8), dimension(:,:), pointer mom_module::kpara_local

parllel kinetic energy, \( \int v_\parallel^2 d\mathbf{v} \), local to compute node [Nsp,Nlocalnode]

real (kind=8), dimension(:,:), pointer mom_module::kperp_local

perpendicular kinetic energy, \( \int 1/2 v_\perp^2 f d\mathbf{v} \), local to compute node [Nsp,Nlocalnode]

real (kind=8), dimension(:,:), pointer mom_module::temp_all

temperature, \( 2/3 \int 1/2 (\mathbf{v} - u_\parallel)^2 f d\mathbf{v} \), on all nodes [Nsp,Nnode]

real (kind=8), dimension(:,:), pointer mom_module::temp_local

temperature, \( 2/3 \int 1/2 (\mathbf{v} - u_\parallel)^2 f d\mathbf{v} \), local to compute node [Nsp,Nlocalnode]

real (kind=8), dimension(:,:), pointer mom_module::temp_para_all

parllel kinetic energy, \( \int v_\parallel^2 d\mathbf{v} \), on all nodes [Nsp,Nnode]

real (kind=8), dimension(:,:), pointer mom_module::temp_para_local

parallel temperature, \( \int (v_\parallel-u_\parallel)^2 d\mathbf{v} \), local to compute node [Nsp,Nlocalnode]

real (kind=8), dimension(:,:), pointer mom_module::temp_perp_all

perpendicular kinetic energy, \( \int 1/2 v_\perp^2 f d\mathbf{v} \), on all nodes [Nsp,Nnode]

real (kind=8), dimension(:,:), pointer mom_module::upara_all

parallel bulk velocity, \( \int v_\parallel f d\mathbf{v} \), on all nodes [Nsp,Nnode]

real (kind=8), dimension(:,:), pointer mom_module::upara_local

parallel bulk velocity, \( \int v_\parallel f d\mathbf{v} \) local to compute node [Nsp,Nlocalnode]

real (kind=8), dimension(:,:), pointer mom_module::uperp_local

perpendicular bulk velocity, \( \int v_\perp f d\mathbf{v} \) local to compute node [Nsp,Nlocalnode]


The documentation for this module was generated from the following file: