10 extern "C" void set_mom_module_ptrs(
double* den_local_cpp,
double* upara_local_cpp,
double* uperp_local_cpp,
double* temp_local_cpp,
double* Kperp_local_cpp,
11 double* Kpara_local_cpp,
double* temp_para_local_cpp,
double* den_all_cpp,
double* temp_all_cpp,
double* temp_perp_all_cpp,
double* temp_para_all_cpp,
double* upara_all_cpp);
22 View<double**, CLayout, HostType>
den_all;
31 :
den_local(
"den_local", f.n_nodes(), f.n_species()),
40 den_all(
NoInit(
"den_all"), (do_gather && pol_decomp.pol_decomp) ? grid.nnode : 1, f.n_species()),
84 for (
int ivr=0; ivr<vgrid.
nvr; ivr++){
87 double smu = vgrid.
vperp_n(ivr);
89 double en_perp= 0.5 * mu;
92 for(
int ivz=0; ivz<vgrid.
nvz; ivz++){
93 double vp = vgrid.
vpar_n(ivz);
94 double en_para2 = vp*vp;
96 double wt = mu_vol*f(isp, ivr, inode, ivz);
98 upara(inode,isp) += wt*vp;
99 uperp(inode,isp) += wt*smu;
100 Kperp(inode,isp) += wt*en_perp;
101 Kpara(inode,isp) += wt*en_para2;
108 const double UPERP_DEFAULT_COEFF = sqrt(0.25*
TWOPI);
112 double den0 = 1.0e-4*species.
f0.den_h(inode);
113 if (!std::isfinite(den(inode,isp))) den(inode,isp) = den0;
114 den(inode,isp) = max(den0,den(inode,isp));
117 double t_norm_ev = species.
f0.fg_temp_ev_h(inode);
118 double vth = 1./species.
f0.fg_vth_inv_h(inode);
121 double inv_density = 1.0/den(inode,isp);
122 Kperp(inode,isp) *= t_norm_ev*inv_density;
123 Kpara(inode,isp) *= t_norm_ev*inv_density;
124 upara(inode,isp) *= vth*inv_density;
125 uperp(inode,isp) *= vth*inv_density;
126 den(inode,isp) *= species.
f0.grid_vol_vonly_h(inode);
130 temp(inode,isp) = ( temp_para(inode,isp) + 2.0*Kperp(inode,isp) )/3.0;
133 if (!std::isfinite(den(inode,isp))) den(inode,isp) = species.
f0.den_h(inode);
134 if (!std::isfinite(Kperp(inode,isp))) Kperp(inode,isp) = t_norm_ev;
135 if (!std::isfinite(Kpara(inode,isp))) Kpara(inode,isp) = t_norm_ev;
136 if (!std::isfinite(temp_para(inode,isp))) temp_para(inode,isp) = t_norm_ev;
137 if (!std::isfinite(temp(inode,isp))) temp(inode,isp) = t_norm_ev;
138 if (!std::isfinite(upara(inode,isp))) upara(inode,isp) = species.
f0.flow_h(inode);
139 if (!std::isfinite(uperp(inode,isp))) uperp(inode,isp) = UPERP_DEFAULT_COEFF*vth;
142 den (inode,isp) = min(max(den (inode,isp),1.0e-2*species.
f0.den_h(inode)),1.0e2*species.
f0.den_h(inode));
143 Kperp (inode,isp) = min(max(Kperp (inode,isp),1.0e-2*t_norm_ev ),1.0e2*t_norm_ev );
144 Kpara (inode,isp) = min(max(Kpara (inode,isp),1.0e-2*t_norm_ev ),1.0e2*t_norm_ev );
145 temp (inode,isp) = min(max(temp (inode,isp),1.0e-2*t_norm_ev ),1.0e2*t_norm_ev );
146 temp_para (inode,isp) = min(max(temp_para (inode,isp),1.0e-2*t_norm_ev ),1.0e2*t_norm_ev );
147 upara (inode,isp) = min(max(upara (inode,isp),-0.75*vgrid.
vp_max*vth ),0.75*vgrid.
vp_max*vth );
148 uperp (inode,isp) = min(max(uperp (inode,isp),0.05*vgrid.
vp_max*vth ),0.75*vgrid.
vp_max*vth );
157 View<int*,CLayout,HostType> cnts(
NoInit(
"cnts"), pol_decomp.mpi.n_plane_ranks);
158 View<int*,CLayout,HostType> displs(
NoInit(
"cnts"), pol_decomp.mpi.n_plane_ranks);
162 for(
int i = 0; i<cnts.size(); i++){
167 int my_send_count = cnts(pol_decomp.mpi.my_plane_rank);
168 MPI_Allgatherv(
den_local.data(), my_send_count, MPI_DOUBLE,
den_all.data(), cnts.data(), displs.data(), MPI_DOUBLE, pol_decomp.mpi.plane_comm);
169 MPI_Allgatherv(
temp_local.data(), my_send_count, MPI_DOUBLE,
temp_all.data(), cnts.data(), displs.data(), MPI_DOUBLE, pol_decomp.mpi.plane_comm);
170 MPI_Allgatherv(
temp_para_local.data(),my_send_count, MPI_DOUBLE,
temp_para_all.data(), cnts.data(), displs.data(), MPI_DOUBLE, pol_decomp.mpi.plane_comm);
171 MPI_Allgatherv(
Kperp_local.data(), my_send_count, MPI_DOUBLE,
temp_perp_all.data(), cnts.data(), displs.data(), MPI_DOUBLE, pol_decomp.mpi.plane_comm);
172 MPI_Allgatherv(
upara_local.data(), my_send_count, MPI_DOUBLE,
upara_all.data(), cnts.data(), displs.data(), MPI_DOUBLE, pol_decomp.mpi.plane_comm);
View< double **, CLayout, HostType > Kperp_local
perpendicular kinetic energy, , local to compute node [Nsp,Nlocalnode]
Definition: moments.hpp:18
void calculate_moments_of_distribution(const VelocityGrid &vgrid, Plasma &plasma, const VGridDistribution< HostType > &f)
Definition: moments.hpp:69
Moments()
Definition: moments.hpp:28
Distribution< Device > f0
Species distribution in velocity space on local mesh nodes.
Definition: species.hpp:128
View< double **, CLayout, HostType > temp_local
temperature, , local to compute node [Nsp,Nlocalnode]
Definition: moments.hpp:17
View< double **, CLayout, HostType > temp_all
temperature, , on all nodes [Nsp,Nnode]
Definition: moments.hpp:24
Definition: velocity_grid.hpp:8
View< double **, CLayout, HostType > upara_all
parallel bulk velocity, , on all nodes [Nsp,Nnode]
Definition: moments.hpp:23
subroutine plasma(grid, itr, p, dene_out, deni_out, Te_out, Ti_out, Vparai_out)
Calculate the plasma density, temperature, and parallel velocity for a point in triangle itr using pl...
Definition: neutral_totalf.F90:1224
Definition: plasma.hpp:107
Definition: moments.hpp:13
int nvr
full grid size (including zero)
Definition: velocity_grid.hpp:19
void for_all_nonadiabatic_species(F func, DevicePtlOpt device_ptl_opt=UseDevicePtl)
Definition: plasma.hpp:126
View< double **, CLayout, HostType > den_all
density, , on all nodes [Nsp,Nnode]
Definition: moments.hpp:22
int nonadiabatic_idx
Index of species skipping adiabatic species (for compatibility with fortran arrays) ...
Definition: species.hpp:81
KOKKOS_INLINE_FUNCTION double vpar_n(int ivz) const
Definition: velocity_grid.hpp:72
View< double **, CLayout, HostType > upara_local
parallel bulk velocity, local to compute node [Nsp,Nlocalnode]
Definition: moments.hpp:15
View< double **, CLayout, HostType > Kpara_local
parllel kinetic energy, , local to compute node [Nsp,Nlocalnode]
Definition: moments.hpp:19
constexpr double J_2_EV
Conversion rate J to ev.
Definition: constants.hpp:6
double mass
Particle mass.
Definition: species.hpp:83
KOKKOS_INLINE_FUNCTION double mu_vol_fac(int ivr) const
Definition: velocity_grid.hpp:82
KOKKOS_INLINE_FUNCTION double kinetic_energy(double mass, double v)
Definition: basic_physics.hpp:8
View< double **, CLayout, HostType > temp_perp_all
perpendicular kinetic energy, , on all nodes [Nsp,Nnode]
Definition: moments.hpp:25
View< double **, CLayout, HostType > den_local
density, local to compute node [Nsp,Nlocalnode]
Definition: moments.hpp:14
double vp_max
max parallel velocity
Definition: velocity_grid.hpp:10
Definition: plasma.hpp:13
View< double **, CLayout, HostType > temp_para_local
parallel temperature, , local to compute node [Nsp,Nlocalnode]
Definition: moments.hpp:20
void gather(const DomainDecomposition< DeviceType > &pol_decomp)
Definition: moments.hpp:155
void set_mom_module_ptrs(double *den_local_cpp, double *upara_local_cpp, double *uperp_local_cpp, double *temp_local_cpp, double *Kperp_local_cpp, double *Kpara_local_cpp, double *temp_para_local_cpp, double *den_all_cpp, double *temp_all_cpp, double *temp_perp_all_cpp, double *temp_para_all_cpp, double *upara_all_cpp)
View< double **, CLayout, HostType > temp_para_all
parllel kinetic energy, , on all nodes [Nsp,Nnode]
Definition: moments.hpp:26
Definition: species.hpp:75
void parallel_for(const std::string name, int n_ptl, Function func, Option option, HostAoSoA aosoa_h, DeviceAoSoA aosoa_d)
Definition: streamed_parallel_for.hpp:252
int nvz
full grid size (including negative and zero)
Definition: velocity_grid.hpp:20
KOKKOS_INLINE_FUNCTION int n_nodes() const
Definition: vgrid_distribution.hpp:159
View< double **, CLayout, HostType > uperp_local
perpendicular bulk velocity, local to compute node [Nsp,Nlocalnode]
Definition: moments.hpp:16
bool pol_decomp
Use poloidal decomposition.
Definition: domain_decomposition.hpp:46
Moments(const Grid< DeviceType > &grid, const DomainDecomposition< DeviceType > &pol_decomp, const VelocityGrid &vgrid, Plasma &plasma, const VGridDistribution< HostType > &f, bool do_gather)
Definition: moments.hpp:30
Kokkos::ViewAllocateWithoutInitializing NoInit
Definition: space_settings.hpp:68
constexpr double TWOPI
Definition: constants.hpp:9
KOKKOS_INLINE_FUNCTION double vperp_n(int ivr) const
Definition: velocity_grid.hpp:77
Kokkos::View< int *, Kokkos::LayoutRight, HostType > gvid0_pid_h
Which processors get which vertices (host)
Definition: domain_decomposition.hpp:63