XGC1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
col_grid.hpp
Go to the documentation of this file.
1 #ifndef COL_GRID_HPP
2 #define COL_GRID_HPP
3 
4 #include <vector>
5 #include <stdio.h>
6 #include <Kokkos_Core.hpp>
7 #include <Kokkos_DualView.hpp>
8 
10 #include "velocity_grid.hpp"
11 #include "vgrid_distribution.hpp"
12 #include "magnetic_field.hpp"
13 
14 #include "col_grid_matrix.hpp"
15 #include "col_species.hpp"
16 #include "col_vgrids.hpp"
17 
18 //> XGC's multi-species Fokker-Planck-Landau collision operator
19 // see R. Hager et al., Journal of Computational Physics 315 (2016), pp. 644-660
20 // E. S. Yoon et al., Physics of Plasmas 21 (2014), 032503
21 //
22 // The headers of the subroutines in this module refer to the corresponding parts of Hager et al.
23 // where applicable in order to make the source code easier to understand.
24 
25 
26 // Enum for derivatives of ED array
27 namespace ED{
28 enum{
29  RR=0,
30  RZ,
31  ZZ,
32  ER,
33  EZ,
34  N
35 };
36 }
37 
38 // <SOLVER-related>
39 const int vpic_inner_iter_max=20;
40 
41 namespace Convergence{
42 
43 // Ideally this would be an enum class,
44 // but the are MPI and I/O operations that expect ints
45 enum Status{
54 };
55 
56 inline bool is_okay(Status convergence_status){
57  return (convergence_status!=MaxedIter &&
58  convergence_status!=NaNOrInfFound &&
59  convergence_status!=FirstStepFail);
60 }
61 
62 class Moments{
63  // Values
64  double w_sum; //> sum of energy over species
65  double p_sum; //> sum of momentum over species
66  double min_p_thres; //> momentum error criterion
67 
68  // Deltas
69  double dw_sum; //> sum of energy change over species (i.e. error in energy)
70  double dp_sum; //> sum of momentum change over species (i.e. error in momentum)
71  double dn_n_max; //> maximum fraction of density changes among species (i.e. maximum relative error in density in all species )
73  bool en_exit_ok;
75 
76  //> Checks the convergence criterion (i.e. conservation of mass, energy and momentum)
77  template<class Device>
78  inline void eval(int isp, int mesh_ind, const CollisionVelocityGrids<Device>& cvg, const CollisionSpecies<Device>& cs,double& dw,double& dp,double& dn_n) const{
79  double vpic_dn = 0.0;
80  double vpic_dw = 0.0;
81  double vpic_dp = 0.0;
82 
83  int gri = cs.s.h_view(isp, mesh_ind).grid_index;
84 
85  for (int index_i = 0; index_i<cvg.mesh_z.extent(2); index_i++){
86  double mesh_z2 = cvg.mesh_z(mesh_ind,gri,index_i)*cvg.mesh_z(mesh_ind,gri,index_i);
87  for (int index_j = 0; index_j<cvg.mesh_r.extent(2); index_j++){
88  double mesh_r2 = cvg.mesh_r(mesh_ind,gri,index_j)*cvg.mesh_r(mesh_ind,gri,index_j);
89  double vpic_dfc = (cs.pdf_np1.h_view(mesh_ind,isp,index_i,index_j) - cs.pdf_n.h_view(mesh_ind,isp,index_i,index_j)) * cvg.vol_h(mesh_ind,gri,index_j);
90 
91  vpic_dn += vpic_dfc;
92  vpic_dw += vpic_dfc*(mesh_z2 + mesh_r2);
93  vpic_dp += vpic_dfc*cvg.mesh_z(mesh_ind,gri,index_i);
94  }
95  }
96 
97  dw = vpic_dw*cs.s.h_view(isp,mesh_ind).mass;
98  dp = vpic_dp*cs.s.h_view(isp,mesh_ind).mass;
99  dn_n = std::abs(vpic_dn/cs.s.h_view(isp,mesh_ind).dens);
100  }
101 
102  public:
103 
104  template<class Device>
105  inline void set(int mesh_ind, const CollisionSpecies<Device>& col_spall){
106  w_sum = 0;
107  p_sum = 0;
108  min_p_thres = 1e99; //just large value to be reinitialized with a small value
109  for (int spi = 0; spi<col_spall.n(); spi++){
110  w_sum += col_spall.s.h_view(spi,mesh_ind).ens;
111  p_sum += col_spall.s.h_view(spi,mesh_ind).mom;
112  min_p_thres = std::min( col_spall.s.h_view(spi,mesh_ind).mass*col_spall.s.h_view(spi,mesh_ind).dens*sqrt(col_spall.s.h_view(spi,mesh_ind).numeric_vth2),
113  min_p_thres );
114  }
115  }
116 
117  template<class Device>
118  inline bool evaluate_conservation(int mesh_ind, const CollisionVelocityGrids<Device>& col_vgrids, const CollisionSpecies<Device>& col_spall, bool& values_are_all_finite){
119  dw_sum = 0.0;
120  dp_sum = 0.0;
121  dn_n_max = 0.0;
122  dens_exit_ok = true;
123 
124  for (int spi = 0; spi<col_spall.n(); spi++){
125  // Re-compute conservation errors after correction
126  double dw,dp,dn_n;
127  eval(spi, mesh_ind, col_vgrids, col_spall,dw,dp,dn_n);// moments evalulation from difference betw. pdf_n and pdf_np1
128  dw_sum += dw;
129  dp_sum += dp;
130  if(dn_n_max<dn_n) dn_n_max = dn_n;
131  dens_exit_ok = dens_exit_ok && (dn_n<1e-10);
132  }
133 
134  // Check for NaNs and Infs
135  values_are_all_finite = (std::isfinite(dw_sum) && std::isfinite(dp_sum) && std::isfinite(dn_n_max));
136 
137  en_exit_ok = std::abs(dw_sum/w_sum) <= 1e-7;
138  // bool mom_exit_ok = std::abs(col_dp_sum)/(std::max(std::abs(col_p_sum),min_p_thres)) <= (col_f_vp_max/col_f_nvz);// v_par = max_v_par of simulation domain / vth
139  // ALS : momentum tolerance changed to the tolerance in XGC's former collisions module:
140  // For momentum tolerance use main ion mass and electron density as normalization (when there are more than 2 species).
141  double mass_use = (col_spall.s.h_view(0,mesh_ind).is_electron ? col_spall.s.h_view(1,mesh_ind).mass : col_spall.s.h_view(0,mesh_ind).mass);
142  mom_exit_ok = std::abs(dp_sum)/(std::max(std::abs(p_sum),1e-3*mass_use*col_spall.s.h_view(0,mesh_ind).dens)) <= 1.0e-7;
143 
144  return (dens_exit_ok && en_exit_ok && mom_exit_ok);
145  }
146 };
147 
148 }
149 
150 
151 template<class Device>
152 struct TmpColData{
153 
154  // These views are both on CPU and GPU - dual views
155  Kokkos::View<double***,Device> gammac_spall;
156  Kokkos::DualView<int[4][4],Device> mat_pos_rel_indx; // Shouldn't need to be a dual view
157 
158  Kokkos::View<double***,Device> fhalf;
159  Kokkos::View<double***,Device> dfdr;
160  Kokkos::View<double***,Device> dfdz;
161  Kokkos::View<double****,Device> EDs;
162  Kokkos::View<double******,Device> M_ab;
163  Kokkos::View<double*****,Device> M_s;
164 
165  TmpColData<Device>(int nvr, int nvz, int sp_num, int n_vgrids, int mb_n_nodes)
166  : // Note: layout left
167  gammac_spall("gammac_spall",mb_n_nodes,sp_num,sp_num),
168  fhalf("fhalf",mb_n_nodes,nvz-1,nvr-1),
169  dfdr("dfdr",mb_n_nodes,nvz-1,nvr-1),
170  dfdz("dfdz",mb_n_nodes,nvz-1,nvr-1),
171  EDs("EDs",mb_n_nodes,nvz-1,nvr-1,ED::N),
172  M_s("M_s",mb_n_nodes,n_vgrids, (nvr-1)*(nvz-1), ED::N, nvr-1),
173  M_ab("M_ab",mb_n_nodes,n_vgrids, n_vgrids-1, (nvr-1)*(nvz-1), 3, (nvr-1)*(nvz-1)),
174  mat_pos_rel_indx("mat_pos_rel_indx")
175  {
176  // Set up mat_pos_rel_indx and send to GPU
177  int mat_loc[16] = {4,5,7,8,
178  1,2,4,5,
179  3,4,6,7,
180  0,1,3,4};
181 
182  for (int i=0;i<16;i++)
183  mat_pos_rel_indx.h_view(i/4,i%4)=mat_loc[i];
184 
185 #ifdef USE_GPU
186  Kokkos::deep_copy(mat_pos_rel_indx.d_view, mat_pos_rel_indx.h_view);
187 #endif
188  }
189 };
190 
191 
192 template<class Device>
194 
195  public:
196 
198 
200 
203 
206 
208 
210 
211  //> Set up global variables of the collision operator
213  nlr.use_namelist("col_param");
214  inner_psi_bound = nlr.get<double>("col_pin",magnetic_field.inpsi);
215  outer_psi_bound = nlr.get<double>("col_pout",magnetic_field.outpsi);
216 
217  // Normalize by xpt
218  // Default values are already normalized, so only normalize if not default
219  if(nlr.present("col_pin")) inner_psi_bound *= magnetic_field.equil.xpt_psi;
220  if(nlr.present("col_pout")) outer_psi_bound *= magnetic_field.equil.xpt_psi;
221 
222  start_step = nlr.get<int>("col_f_start", 1);
223 
226  // These defaults have not been systematically optimized
227 #ifdef USE_GPU
228  int default_batch_size = 64;
229 #else
230  int default_batch_size = 1;//omp_get_max_threads();
231 #endif
232  std::string linalg_backend_str;
233  if(nlr.namelist_present("performance_param")){
234  nlr.use_namelist("performance_param");
235  batch_size = nlr.get<int>("mesh_batch_size", default_batch_size);
236  linalg_backend_str = nlr.get<std::string>("collisions_solver", "lapack");
237  ginkgo_residual_reduction = nlr.get<double>("ginkgo_residual_reduction", 1.0e-16);
238  ginkgo_max_iterations = nlr.get<int>("ginkgo_max_iterations", 300);
239  }else{
240  batch_size = default_batch_size;
241  linalg_backend_str = "lapack";
242  ginkgo_residual_reduction = 1.0e-16;
243  ginkgo_max_iterations = 300;
244  }
245 
247  }
248 
250  nlr.use_namelist("ptl_param");
251  int plasma_nspecies = nlr.get<int>("ptl_nsp",1);
252 
253  nlr.use_namelist("sml_param");
254  int sml_special = nlr.get<int>("sml_special", 0);
255  bool sml_electron_on = nlr.get<bool>("sml_electron_on", false);
256  if(sml_special==4) sml_electron_on = false;
257  if(sml_electron_on) plasma_nspecies += 1; // Add electrons
258 
259  VelocityGrid vgrid(nlr);
260 #ifdef USE_GPU
261  int default_batch_size = 64;
262 #else
263  int default_batch_size = 1;//omp_get_max_threads();
264 #endif
265  int col_grid_batch_size;
266  if(nlr.namelist_present("performance_param")){
267  nlr.use_namelist("performance_param");
268  col_grid_batch_size = nlr.get<int>("mesh_batch_size", default_batch_size);
269  }else{
270  col_grid_batch_size = default_batch_size;
271  }
272 
273  double device_col_GB_per_vertex = plasma_nspecies*std::max(plasma_nspecies-1,1)*(vgrid.nvp*vgrid.nmu)*(vgrid.nvp*vgrid.nmu)*8*BYTES_TO_GB;
274  double gpu_memory_usage = (col_grid_batch_size*device_col_GB_per_vertex);
275 
276  MemoryPrediction memory_prediction{"Fields", 0.0, gpu_memory_usage};
277 
278  return memory_prediction;
279  }
280 
281  // Core init
282  void core_init(int isp, int mesh_ind, const CollisionVelocityGrids<Device>& col_vgrids, const CollisionSpecies<Device>& col_spall ) const;
283 
284  // Get maxw_fac
285  static KOKKOS_INLINE_FUNCTION double get_maxw_fac(double mesh_dr, double mesh_r, double numeric_vth2) ;
286 
287  // Get numeric v_thermal equilibrium for a given species against a given grid
288  View<double*,Device> get_numeric_v_thermal_equil(int mb_n_nodes, int spi, int grj, const CollisionVelocityGrids<Device>& col_vgrids, const CollisionSpecies<Device>& col_spall) const;
289 
290  // Delta init
291  void core_delta_init(int mb_n_nodes, int gri, int grj, int spi, CollisionVelocityGrids<Device>& col_vgrids, const CollisionSpecies<Device>& col_spall) const;
292 
293  // LU_matrix_ftn
294  static KOKKOS_INLINE_FUNCTION void LU_matrix_ftn(int mesh_ind, int gri, int grj, int spi, int cell_i,int cell_j,const CollisionVelocityGrids<Device>& col_vgrids,int mprl_col, int mat_pos,double coeff1,double coeff2, const TmpColData<Device>& tcd, const Kokkos::View<int**,Device>& index_map_LU_d,typename Collisions::GridMatrix<Device>::values_array_t LU_values);
295 
296  // LU_matrix
297  void LU_matrix(int mb_n_nodes, int gri, int grj, int spi, const CollisionVelocityGrids<Device>& col_vgrids,const TmpColData<Device>& tcd, Collisions::GridMatrix<Device>* const mtx) const;
298 
299  // the core loop including the picard loop is now inside this function
300  void picard_loop(int vpic_inner_iter_max, const CollisionVelocityGrids<Device>& col_vgrids, CollisionSpecies<Device>& col_spall, TmpColData<Device>& tcd, int mb_n_nodes, Kokkos::View<Convergence::Status*,HostType>& convergence_status, Kokkos::View<Convergence::Moments*,HostType>& moments) const;
301 
302  // dispatch for E_and_D calc
303  void E_and_D(int mb_n_nodes, int gri, int grj, const CollisionVelocityGrids<Device>& col_vgrids, TmpColData<Device>& tcd) const;
304 
305  // E_and_D_s
306  static KOKKOS_INLINE_FUNCTION void E_and_D_s(int idx, int mb_n_nodes, int nvrm1, int nvzm1, const TmpColData<Device>& tcd, int gri) ;
307 
308  // dispatch for angle_avg calc
309  void angle_avg(int mb_n_nodes, int gri, int grj, CollisionVelocityGrids<Device>& col_vgrids, TmpColData<Device>& tcd) const;
310 
311  // angle_avg_s
312  static KOKKOS_INLINE_FUNCTION void angle_avg_s(int idx, int mb_n_nodes, int nvrm1, int nvzm1, const CollisionVelocityGrids<Device>& col_vgrids, const TmpColData<Device>& tcd, int gri);
313 
314  // E_and_D_ab
315  static KOKKOS_INLINE_FUNCTION void E_and_D_ab(int idx, int mb_n_nodes, int nvrm1, int nvzm1, const CollisionVelocityGrids<Device>& col_vgrids, const TmpColData<Device>& tcd, int gri, int grj);
316 
317  // angle_avg_ab
318  static KOKKOS_INLINE_FUNCTION void angle_avg_ab(int idx, int mb_n_nodes, int nvrm1, int nvzm1, const CollisionVelocityGrids<Device>& col_vgrids, const TmpColData<Device>& tcd, int gri, int grj);
319 
320  void f_df(int mb_n_nodes, const CollisionVelocityGrids<Device>& col_vgrids, const CollisionSpecies<Device>& col_spall,int spi,int grj,TmpColData<Device>& tcd) const;
321 
322  // Main collision algorithm
323  Kokkos::View<Convergence::Status*,HostType> core(CollisionVelocityGrids<Device>& col_vgrids, CollisionSpecies<Device>& col_spall, TmpColData<Device>& tcd, int mb_n_nodes) const;
324 
325 
326  // Carries out collisions
327  void collision(std::vector<Species<DeviceType>> &all_species, const Moments& moments, const DomainDecomposition<DeviceType>& pol_decomp, const View<bool*, HostType>& in_range, bool symmetric_f, const VGridDistribution<HostType>& f0_f, VGridDistribution<HostType>& df_out, double dt, View<int*,CLayout,HostType>& converged_all, View<double*,CLayout,HostType>& node_cost) const;
328 };
329 
330 #include "col_grid.tpp"
331 
332 #endif
double min_p_thres
Definition: col_grid.hpp:66
double inpsi
Boundary condition used in a few spots.
Definition: magnetic_field.hpp:25
Definition: col_grid_matrix.hpp:23
double w_sum
Definition: col_grid.hpp:64
int nmu
n points in mu (not including zero)
Definition: velocity_grid.hpp:13
int nvp
n points in parallel velocity (not including zero)
Definition: velocity_grid.hpp:9
Kokkos::View< double ******, Device > M_ab
Definition: col_grid.hpp:162
static KOKKOS_INLINE_FUNCTION void LU_matrix_ftn(int mesh_ind, int gri, int grj, int spi, int cell_i, int cell_j, const CollisionVelocityGrids< Device > &col_vgrids, int mprl_col, int mat_pos, double coeff1, double coeff2, const TmpColData< Device > &tcd, const Kokkos::View< int **, Device > &index_map_LU_d, typename Collisions::GridMatrix< Device >::values_array_t LU_values)
Definition: col_grid.tpp:838
int batch_size
Definition: col_grid.hpp:197
Kokkos::View< double ***, Kokkos::LayoutRight, HostType > mesh_z
Definition: col_vgrids.hpp:21
T get(const string &param, const T default_val, int val_ind=0)
Definition: NamelistReader.hpp:373
constexpr double BYTES_TO_GB
Definition: constants.hpp:12
Kokkos::View< Convergence::Status *, HostType > core(CollisionVelocityGrids< Device > &col_vgrids, CollisionSpecies< Device > &col_spall, TmpColData< Device > &tcd, int mb_n_nodes) const
Definition: col_grid.tpp:138
Kokkos::View< double *****, Device > M_s
Definition: col_grid.hpp:163
Definition: col_grid.hpp:193
Kokkos::DualView< double ****, Kokkos::LayoutRight, Device > pdf_np1
Definition: col_species.hpp:88
Definition: velocity_grid.hpp:8
void set(int mesh_ind, const CollisionSpecies< Device > &col_spall)
Definition: col_grid.hpp:105
bool present(const string &param)
Definition: NamelistReader.hpp:363
Definition: col_grid.hpp:46
Kokkos::View< double ***, Device > fhalf
Definition: col_grid.hpp:158
Status
Definition: col_grid.hpp:45
Definition: col_grid.hpp:51
void collision(std::vector< Species< DeviceType >> &all_species, const Moments &moments, const DomainDecomposition< DeviceType > &pol_decomp, const View< bool *, HostType > &in_range, bool symmetric_f, const VGridDistribution< HostType > &f0_f, VGridDistribution< HostType > &df_out, double dt, View< int *, CLayout, HostType > &converged_all, View< double *, CLayout, HostType > &node_cost) const
Definition: col_grid.tpp:22
Collisions::LinAlgBackend labackend
Definition: col_grid.hpp:199
double p_sum
Definition: col_grid.hpp:65
Kokkos::View< value_type ***, Kokkos::LayoutRight, Device, Kokkos::MemoryTraits< Kokkos::Unmanaged >> values_array_t
Definition: col_grid_matrix.hpp:31
Kokkos::View< double ***, Device > gammac_spall
Definition: col_grid.hpp:155
Definition: NamelistReader.hpp:193
Definition: magnetic_field.hpp:12
Definition: moments.hpp:13
Equilibrium equil
The object containing information about the magnetic equilibrium.
Definition: magnetic_field.hpp:32
Kokkos::View< double ***, Kokkos::LayoutRight, HostType > vol_h
Definition: col_vgrids.hpp:23
bool en_exit_ok
Definition: col_grid.hpp:73
Definition: col_grid.hpp:29
CollisionGrid(NLReader::NamelistReader &nlr, const MagneticField< DeviceType > &magnetic_field)
Definition: col_grid.hpp:212
Kokkos::DualView< CollisionSpeciesScalars **, Device > s
Definition: col_species.hpp:85
double dn_n_max
Definition: col_grid.hpp:71
CollisionGrid()
Definition: col_grid.hpp:209
int ginkgo_max_iterations
Definition: col_grid.hpp:202
double inner_psi_bound
Inner psi bound; collisions skipped below this.
Definition: col_grid.hpp:204
Definition: col_vgrids.hpp:13
void core_delta_init(int mb_n_nodes, int gri, int grj, int spi, CollisionVelocityGrids< Device > &col_vgrids, const CollisionSpecies< Device > &col_spall) const
Definition: col_grid.tpp:777
View< double *, Device > get_numeric_v_thermal_equil(int mb_n_nodes, int spi, int grj, const CollisionVelocityGrids< Device > &col_vgrids, const CollisionSpecies< Device > &col_spall) const
Definition: col_grid.tpp:745
static KOKKOS_INLINE_FUNCTION void E_and_D_ab(int idx, int mb_n_nodes, int nvrm1, int nvzm1, const CollisionVelocityGrids< Device > &col_vgrids, const TmpColData< Device > &tcd, int gri, int grj)
Definition: col_grid.tpp:515
Definition: col_grid.hpp:47
Definition: col_grid.hpp:152
void LU_matrix(int mb_n_nodes, int gri, int grj, int spi, const CollisionVelocityGrids< Device > &col_vgrids, const TmpColData< Device > &tcd, Collisions::GridMatrix< Device > *const mtx) const
Definition: col_grid.tpp:883
const int vpic_inner_iter_max
Definition: col_grid.hpp:39
static KOKKOS_INLINE_FUNCTION void E_and_D_s(int idx, int mb_n_nodes, int nvrm1, int nvzm1, const TmpColData< Device > &tcd, int gri)
Definition: col_grid.tpp:684
void core_init(int isp, int mesh_ind, const CollisionVelocityGrids< Device > &col_vgrids, const CollisionSpecies< Device > &col_spall) const
Definition: col_grid.tpp:356
double dw_sum
Definition: col_grid.hpp:69
double outpsi
Boundary condition used in a few spots.
Definition: magnetic_field.hpp:26
idx
Definition: diag_f0_df_port1.hpp:32
Kokkos::DualView< double ****, Kokkos::LayoutRight, Device > pdf_n
Definition: col_species.hpp:87
void use_namelist(const string &namelist)
Definition: NamelistReader.hpp:355
Definition: col_grid.hpp:62
Definition: col_grid.hpp:31
Definition: memory_prediction.hpp:4
void eval(int isp, int mesh_ind, const CollisionVelocityGrids< Device > &cvg, const CollisionSpecies< Device > &cs, double &dw, double &dp, double &dn_n) const
Definition: col_grid.hpp:78
void E_and_D(int mb_n_nodes, int gri, int grj, const CollisionVelocityGrids< Device > &col_vgrids, TmpColData< Device > &tcd) const
Definition: col_grid.tpp:491
Definition: col_grid.hpp:48
Kokkos::View< double ***, Device > dfdz
Definition: col_grid.hpp:160
double outer_psi_bound
Outer psi bound; collisions skipped above this.
Definition: col_grid.hpp:205
void picard_loop(int vpic_inner_iter_max, const CollisionVelocityGrids< Device > &col_vgrids, CollisionSpecies< Device > &col_spall, TmpColData< Device > &tcd, int mb_n_nodes, Kokkos::View< Convergence::Status *, HostType > &convergence_status, Kokkos::View< Convergence::Moments *, HostType > &moments) const
Definition: col_grid.tpp:239
Definition: col_grid.hpp:34
double ginkgo_residual_reduction
Definition: col_grid.hpp:201
static KOKKOS_INLINE_FUNCTION void angle_avg_s(int idx, int mb_n_nodes, int nvrm1, int nvzm1, const CollisionVelocityGrids< Device > &col_vgrids, const TmpColData< Device > &tcd, int gri)
Definition: col_grid.tpp:568
bool dens_exit_ok
Definition: col_grid.hpp:72
bool namelist_present(const string &namelist)
Definition: NamelistReader.hpp:351
Definition: col_grid.hpp:30
Kokkos::View< double ****, Device > EDs
Definition: col_grid.hpp:161
Definition: col_grid.hpp:32
Definition: col_grid.hpp:53
Definition: magnetic_field.F90:1
LinAlgBackend
Definition: col_grid_matrix.hpp:221
int n() const
Definition: col_species.hpp:181
Definition: col_grid.hpp:50
Kokkos::View< double ***, Device > dfdr
Definition: col_grid.hpp:159
static KOKKOS_INLINE_FUNCTION void angle_avg_ab(int idx, int mb_n_nodes, int nvrm1, int nvzm1, const CollisionVelocityGrids< Device > &col_vgrids, const TmpColData< Device > &tcd, int gri, int grj)
Definition: col_grid.tpp:421
int start_step
starting step for collisions; should be part of step trigger
Definition: col_grid.hpp:207
bool mom_exit_ok
Definition: col_grid.hpp:74
Definition: species.hpp:75
static MemoryPrediction estimate_memory_usage(NLReader::NamelistReader &nlr)
Definition: col_grid.hpp:249
Definition: col_species.hpp:81
Kokkos::DualView< int[4][4], Device > mat_pos_rel_indx
Definition: col_grid.hpp:156
void angle_avg(int mb_n_nodes, int gri, int grj, CollisionVelocityGrids< Device > &col_vgrids, TmpColData< Device > &tcd) const
Definition: col_grid.tpp:403
bool is_okay(Status convergence_status)
Definition: col_grid.hpp:56
Definition: col_grid.hpp:52
double dp_sum
Definition: col_grid.hpp:70
void f_df(int mb_n_nodes, const CollisionVelocityGrids< Device > &col_vgrids, const CollisionSpecies< Device > &col_spall, int spi, int grj, TmpColData< Device > &tcd) const
Definition: col_grid.tpp:633
Kokkos::View< double ***, Kokkos::LayoutRight, HostType > mesh_r
Definition: col_vgrids.hpp:20
static KOKKOS_INLINE_FUNCTION double get_maxw_fac(double mesh_dr, double mesh_r, double numeric_vth2)
Definition: col_grid.tpp:729
Definition: col_grid.hpp:49
bool evaluate_conservation(int mesh_ind, const CollisionVelocityGrids< Device > &col_vgrids, const CollisionSpecies< Device > &col_spall, bool &values_are_all_finite)
Definition: col_grid.hpp:118
double xpt_psi
Psi coordinate of 1st X-point.
Definition: equil.hpp:84
Definition: col_grid.hpp:33