XGCa
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
species.hpp
Go to the documentation of this file.
1 #ifndef SPECIES_HPP
2 #define SPECIES_HPP
3 #include <Cabana_AoSoA.hpp>
4 #include <Cabana_DeepCopy.hpp>
5 #include <Kokkos_Core.hpp>
6 #include "NamelistReader.hpp"
7 #include "timer_macro.hpp"
8 #include "magnetic_field.hpp"
9 #include "grid.hpp"
10 #include "domain_decomposition.hpp"
11 #include "particles.hpp"
12 #include "space_settings.hpp"
13 #include "distribution.hpp"
14 #include "profile.hpp"
15 #include "gyro_avg_mat.hpp"
17 #include "basic_physics.hpp"
18 #include "memory_prediction.hpp"
19 #include "xgc_io.hpp"
20 
21 extern "C" void set_spall_num_and_ptr(int idx, int n_ptl, int n_vecs, VecParticles* ptl);
22 extern "C" void set_min_max_num(int isp, int n_ptl);
23 extern "C" void adjust_n_ptl_for_core_ptl(int* n_ptl);
24 
29  return false;
30 }
31 
36  return false;
37 }
38 
39 // Used for Cabana slices (getting one particle property at a time)
40 namespace PtlSlice{
41 #ifdef ESC_PTL
42 enum{Ph=0,Ct,Gid,Flag};
43 #else
44 enum{Ph=0,Ct,Gid};
45 #endif
46 }
47 
48 struct PtlMvmt{
49  // Options for particle location management when looping over particles
50  enum SendOpt{
51  NoSend=0,
54  };
55 
56  enum ReturnOpt{
60  };
61 
64 
65  PtlMvmt(SendOpt send_opt, ReturnOpt return_opt) : send_opt(send_opt), return_opt(return_opt){}
66 };
67 
71 };
72 
73 // Species class
74 template<class Device>
75 class Species{
76  public:
77 
78  int idx;
79  bool is_electron;
80  bool is_adiabatic;
83  double mass;
84  double charge;
85  double charge_eu;
86  double c_m;
87  double c2_2m;
88 
90  bool is_deltaf;
91 
92  int ncycles;
94 
96  int n_ptl;
97  Cabana::AoSoA<ParticleDataTypes,HostType,VEC_LEN> particles;
98 
99  // Device particles
100  Cabana::AoSoA<ParticleDataTypes,Device,VEC_LEN> particles_d;
101 
103 
106 
107  /*** Could be its own class inside species? ***/
111 
112  // phase0 (for ion restoration)
113  Cabana::AoSoA<PhaseDataTypes,HostType,VEC_LEN> phase0;
114  Cabana::AoSoA<PhaseDataTypes,Device,VEC_LEN> phase0_d;
115 
116  // For electron restoration
117  Cabana::AoSoA<ParticleDataTypes,HostType,VEC_LEN> backup_particles;
118  int n_backup_particles; // Number of particles stored in backup_particles (can't be deduced from its size due to buffer)
119  /****/
120 
122 
124 
125  Eq::Profile<Device> eq_temp; // Equilibrium temperature
126  Eq::Profile<Device> eq_den; // Equilibrium density
127  Eq::Profile<Device> eq_flow; // Equilibrium flow
128  int eq_flow_type; // Type of Equilibirum flow
129 
130  Eq::Profile<Device> eq_fg_temp; // Equilibrium temperature
131  Eq::Profile<Device> eq_fg_den; // Equilibrium density
132  Eq::Profile<Device> eq_fg_flow; // Equilibrium flow
133  int eq_fg_flow_type; // Type of Equilibirum flow
134 
135  Eq::Profile<Device> eq_mk_temp; // Equilibrium temperature
136  Eq::Profile<Device> eq_mk_den; // Equilibrium density
137  Eq::Profile<Device> eq_mk_flow; // Equilibrium flow
138  int eq_mk_flow_type; // Type of Equilibirum flow
139 
140 
142 
143  /*** Constructors ***/
144 
145  Species(int idx_in, int nonadiabatic_idx_in, bool is_electron_in, bool is_adiabatic_in, KinType kintype_in, double mass_in, double charge_in, double charge_eu_in, bool is_deltaf_in,
146  int ncycles_in);
147 
148  Species(NLReader::NamelistReader& nlr, const Grid<DeviceType> &grid, const MagneticField<DeviceType> &magnetic_field, const DomainDecomposition<DeviceType>& pol_decomp, int idx_in, int nonadiabatic_idx_in);
149 
150  // Electron or ion default constructor
151  Species(SpeciesType sp_type, int n_ptl)
152  : idx(sp_type==ELECTRON ? 0 : 1),
153  is_electron(sp_type==ELECTRON),
154  mass(is_electron ? 3.344e-30 : PROTON_MASS),
156  charge_eu(is_electron ? -1.0 : 1.0),
157  is_deltaf(true),
158  is_adiabatic(false),
159  nonadiabatic_idx(idx), // Since is_adiabatic is false above
161  ncycles(is_electron ? 70 : 1),
162  c_m(charge/mass),
163  c2_2m(0.5*charge*charge/mass),
165  n_ptl(n_ptl),
166  backup_particles("backup_particles", 0),
167  particles("particles", add_vec_buffer(n_ptl)),
170  eq_temp(1.0e3,-0.1),
171  eq_den(1.0e19,-0.1),
172  eq_flow_type(2),
173  owns_particles_d(false),
177 
178  // Special constructor for tests that involve tracking particles in memory
179  // The idea is to use these particles to test functions that reorder particles,
180  // e.g. sort, shift, and cleaning
181  Species(int n_ptl_in)
182  : n_ptl(n_ptl_in),
183  is_electron(true),
184  is_adiabatic(false),
185  particles("particles", add_vec_buffer(n_ptl_in)),
187  owns_particles_d(false),
192  {
193 
194  // Slice particle properties
195  auto ph = Cabana::slice<PtlSlice::Ph>(particles);
196  auto ct = Cabana::slice<PtlSlice::Ct>(particles);
197  auto gid = Cabana::slice<PtlSlice::Gid>(particles);
198 #ifdef ESC_PTL
199  auto flag = Cabana::slice<PtlSlice::Flag>(particles);
200 #endif
201 
202  // Offset gid if using MPI
203 #ifdef USE_MPI
204  long long int gid_offset = n_ptl*SML_COMM_RANK;
205 #else
206  long long int gid_offset = 0;
207 #endif
208 
209  // Assign trackable values
210  for (int i=0;i<n_ptl;i++){
211  // Set GID in order
212  gid(i) = gid_offset + i+1; // 1-indexed
213 
214  // Value of properties is gid + 0.1*(property index)
215  // First particle: (1.0, 1.1, ... 1.8)
216  // Second particle: (2.0, 2.1, ... 2.8)
217  for (int j=0;j<6;j++) ph(i, j) = gid(i) + (j)*0.1;
218  for (int j=0;j<3;j++) ct(i, j) = gid(i) + (j+6)*0.1;
219  }
220 
221  // Buffer particles: same but with gid = -1
222  if(n_ptl>0){
223  for (int i=n_ptl;i<add_vec_buffer(n_ptl);i++){
224  gid(i) = -1;
225  for (int j=0;j<6;j++) ph(i, j) = gid(i) + (j)*0.1;
226  for (int j=0;j<3;j++) ct(i, j) = gid(i) + (j+6)*0.1;
227  }
228  }
229  }
230 
231  static std::vector<MemoryPrediction> estimate_memory_usage(NLReader::NamelistReader& nlr, const Grid<DeviceType> &grid, const DomainDecomposition<DeviceType>& pol_decomp, int species_idx);
232 
233  static int get_initial_n_ptl(NLReader::NamelistReader& nlr, const Grid<DeviceType> &grid, const DomainDecomposition<DeviceType>& pol_decomp, int sml_special, int species_idx, bool verbose);
234 
235  void resize_particles(int new_n_ptl){
236  n_ptl = new_n_ptl;
237 
238 #ifndef USE_GPU
239  // If CPU-only, particles_d points to the same location as particles. If particles is resized, then Cabana will not deallocate the first allocation
240  // since it is still used by particles_d. So, reset particles_d before resize, and point it back to particles only afterwards
241  particles_d = Cabana::AoSoA<ParticleDataTypes,Device,VEC_LEN>();
242 #endif
243 
244  particles.reserve(minimum_ptl_reservation); // Can only raise reservation (no-op if AoSoA is already larger)
245  particles.resize(add_vec_buffer(n_ptl));
246 
247 #ifndef USE_GPU
248  // Point particles_d back to particles after resize
250 #endif
251 
253  }
254 
256 #ifdef USE_GPU
257  particles.reserve(minimum_ptl_reservation); // Can only raise reservation (no-op if AoSoA is already larger)
258  particles.resize(particles_d.size());
259 #else
260  // Point particles back to particles_d after resize
262 #endif
263 
265  }
266 
267  /* If using CPU-only, then "device" particles are a shallow copy of host particles so that
268  * there is no unnecessary duplication. When "device" particles are resized, Cabana will keep the
269  * original allocation if there is a second reference (i.e. host particles).
270  * To resolve this, we free the host particles here so that there is no second reference
271  * */
273  particles = Cabana::AoSoA<ParticleDataTypes,HostType,VEC_LEN>();
274  }
275 
276  /* Resizes device particles if on GPU, or just creates a shallow copy if CPU only
277  * */
279  if(!owns_particles_d) exit_XGC("\nSpecies tried to resize device particles, but doesn't own the device array.");
280 
281 #ifdef USE_GPU
282  particles_d.reserve(minimum_ptl_reservation); // Can only raise reservation (no-op if AoSoA is already larger)
283  // Resize device particles to match n particles
285 #else
286  // If kernels are on CPU, do shallow copy
288 #endif
289  }
290 
291  /* Resizes device particles
292  * */
293  void resize_device_particles(int new_n_ptl){
294  if(!owns_particles_d) exit_XGC("\nSpecies tried to resize device particles, but doesn't own the device array.");
295 
296  n_ptl = new_n_ptl;
297 
298  particles_d.reserve(minimum_ptl_reservation); // Can only raise reservation (no-op if AoSoA is already larger)
299 
300  // Resize device particles to match host particles
302  }
303 
304  /* Copies particles to device - deep copy if using GPU, otherwise shallow copy
305  * Also takes the opportunity to set the buffer particles to realistic values
306  * */
308  if(!owns_particles_d) exit_XGC("\nSpecies tried to copy particles to device, but doesn't own the device array.");
309 
310 #ifdef USE_GPU
311  // Copy to device
312  Cabana::deep_copy(particles_d, particles);
313 #else
314  // No operation required if CPU-only
315 #endif
316 
317  // Copy last particle to fill remainder of trailing vector in AoSoA
319  }
320 
321  /* Copies particles from device - deep copy if using GPU, otherwise no copy is necessary
322  * */
324  if(!owns_particles_d) exit_XGC("\nSpecies tried to copy particles from device, but doesn't own the device array.");
325 
326 #ifdef USE_GPU
327  // Copy particles to host
328  Cabana::deep_copy(particles, particles_d);
329 #else
330  // No operation required if CPU-only
331 #endif
332  }
333 
334  /* Copies particles to device if they are resident on the device
335  * */
338  }
339 
340  /* Copies particles from device if they are resident on the device
341  * */
344  }
345 
346  /* Copies particles to device if they are NOT resident on the device
347  * */
350  }
351 
352  /* Copies particles from device if they are NOT resident on the device
353  * */
356  }
357 
366  if (n_ptl>0){
367  int last_ptl_index = n_ptl - 1;
368  auto ph = Cabana::slice<PtlSlice::Ph>(particles_d);
369  auto ct = Cabana::slice<PtlSlice::Ct>(particles_d);
370  auto gid = Cabana::slice<PtlSlice::Gid>(particles_d);
371 #ifdef ESC_PTL
372  auto flag = Cabana::slice<PtlSlice::Flag>(particles_d);
373 #endif
374 
375  Kokkos::parallel_for("set_buffer_particles_d", Kokkos::RangePolicy<ExSpace>( n_ptl, add_vec_buffer(n_ptl) ), KOKKOS_LAMBDA( const int i ){
376  // Buffer particles: same as last particle, gid = -1
377  for (int j=0;j<6;j++) ph(i, j) = ph(last_ptl_index, j);
378  for (int j=0;j<3;j++) ct(i, j) = ct(last_ptl_index, j);
379  gid(i) = -1;
380 #ifdef ESC_PTL
381  flag(i) = flag(last_ptl_index);
382 #endif
383  });
384  }
385  }
386 
395  if (n_ptl>0){
396  int last_ptl_index = n_ptl - 1;
397  auto ph = Cabana::slice<PtlSlice::Ph>(phase0_d);
398 
399  Kokkos::parallel_for("set_buffer_phase0", Kokkos::RangePolicy<ExSpace>( n_ptl, add_vec_buffer(n_ptl) ), KOKKOS_LAMBDA( const int i ){
400  // copy final real particle
401  for (int j=0;j<6;j++) ph(i, j) = ph(last_ptl_index, j);
402  });
403  }
404  }
405 
406  // Options for custom launch bounds since kokkos defaults are suboptimal for electron push kernel
407  enum class LaunchBounds{
408  Default,
409  Custom
410  };
411 
417  template<typename F>
418  inline void for_all_particles(const std::string label, F lambda_func) const {
419  Kokkos::RangePolicy<ExSpace> particle_range_policy( 0, p_range<DeviceType>(n_ptl) );
420  Kokkos::parallel_for(label, Opt::require(particle_range_policy, Async), lambda_func);
421  }
422 
423  inline void back_up_SoA(Cabana::AoSoA<ParticleDataTypes,Device,VEC_LEN>& backup_SoA, int offset, int n) const{
424  auto ph_b = Cabana::slice<PtlSlice::Ph>(backup_SoA);
425  auto ct_b = Cabana::slice<PtlSlice::Ct>(backup_SoA);
426  auto gid_b = Cabana::slice<PtlSlice::Gid>(backup_SoA);
427 #ifdef ESC_PTL
428  auto flag_b = Cabana::slice<PtlSlice::Flag>(backup_SoA);
429 #endif
430 
431  auto ph = Cabana::slice<PtlSlice::Ph>(particles_d);
432  auto ct = Cabana::slice<PtlSlice::Ct>(particles_d);
433  auto gid = Cabana::slice<PtlSlice::Gid>(particles_d);
434 #ifdef ESC_PTL
435  auto flag = Cabana::slice<PtlSlice::Flag>(particles_d);
436 #endif
437 
438  Kokkos::parallel_for("backup_first_soa", Kokkos::RangePolicy<ExSpace>( 0, n ), KOKKOS_LAMBDA( const int i ){
439  int i_offset = i + offset;
440  // Make backup copy
441  for (int j=0;j<6;j++) ph_b(i, j) = ph(i_offset, j);
442  for (int j=0;j<3;j++) ct_b(i, j) = ct(i_offset, j);
443  gid_b(i) = gid(i_offset);
444 #ifdef ESC_PTL
445  flag_b(i) = flag(i_offset);
446 #endif
447  // Deactivate
448  gid(i_offset) = -1;
449  });
450  }
451 
452  inline void restore_backup_SoA(Cabana::AoSoA<ParticleDataTypes,Device,VEC_LEN>& backup_SoA, int offset, int n) const{
453  auto ph_b = Cabana::slice<PtlSlice::Ph>(backup_SoA);
454  auto ct_b = Cabana::slice<PtlSlice::Ct>(backup_SoA);
455  auto gid_b = Cabana::slice<PtlSlice::Gid>(backup_SoA);
456 #ifdef ESC_PTL
457  auto flag_b = Cabana::slice<PtlSlice::Flag>(backup_SoA);
458 #endif
459 
460  auto ph = Cabana::slice<PtlSlice::Ph>(particles_d);
461  auto ct = Cabana::slice<PtlSlice::Ct>(particles_d);
462  auto gid = Cabana::slice<PtlSlice::Gid>(particles_d);
463 #ifdef ESC_PTL
464  auto flag = Cabana::slice<PtlSlice::Flag>(particles_d);
465 #endif
466 
467  Kokkos::parallel_for("backup_first_soa", Kokkos::RangePolicy<ExSpace>( 0, n ), KOKKOS_LAMBDA( const int i ){
468  int i_offset = i + offset;
469  // Restore from backup copy
470  for (int j=0;j<6;j++) ph(i_offset, j) = ph_b(i, j);
471  for (int j=0;j<3;j++) ct(i_offset, j) = ct_b(i, j);
472  gid(i_offset) = gid_b(i);
473 #ifdef ESC_PTL
474  flag(i_offset) = flag_b(i);
475 #endif
476  });
477  }
478 
484  template<typename F>
485  inline void for_particle_range(int begin_idx, int end_idx, const std::string label, F lambda_func) const {
486  if(end_idx <= begin_idx) return; // Return if range is 0 or less
487 
488  // Still need the subset to line up with the AoSoA vector length
489  int first_soa = begin_idx/VEC_LEN;
490  int n_other_ptl_in_first_soa = begin_idx - first_soa*VEC_LEN;
491  bool first_soa_is_partial = (n_other_ptl_in_first_soa>0);
492  int last_soa = (end_idx-1)/VEC_LEN;
493  int n_other_ptl_in_last_soa = (last_soa+1)*VEC_LEN - end_idx;
494  bool last_soa_is_partial = (n_other_ptl_in_last_soa>0);
495 
496 #ifdef USE_GPU
497  int first_item_in_shifted_range = first_soa*VEC_LEN;
498 #else
499  int first_item_in_shifted_range = first_soa;
500 #endif
501 
502  Cabana::AoSoA<ParticleDataTypes,Device,VEC_LEN> ptl_first_soa;
503  Cabana::AoSoA<ParticleDataTypes,Device,VEC_LEN> ptl_last_soa;
504  if(first_soa_is_partial){
505  // Make a backup of the first SoA and set the particles_d GIDs to -1
506  ptl_first_soa = Cabana::AoSoA<ParticleDataTypes,Device,VEC_LEN>("ptl_first_soa", n_other_ptl_in_first_soa);
507  back_up_SoA(ptl_first_soa, first_soa*VEC_LEN, n_other_ptl_in_first_soa);
508  }
509  if(last_soa_is_partial){
510  // Make a backup of the last SoA and set the particles_d GIDs to -1
511  ptl_last_soa = Cabana::AoSoA<ParticleDataTypes,Device,VEC_LEN>("ptl_last_soa", n_other_ptl_in_last_soa);
512  back_up_SoA(ptl_last_soa, end_idx, n_other_ptl_in_last_soa);
513  }
514 
515  // Finally, do the parallel_for
516  Kokkos::RangePolicy<ExSpace> particle_range_policy( first_item_in_shifted_range, p_range<DeviceType>(end_idx) );
517  Kokkos::parallel_for(label, Opt::require(particle_range_policy, Async), lambda_func);
518 
519  // Restore from backup
520  if(first_soa_is_partial){
521  restore_backup_SoA(ptl_first_soa, first_soa*VEC_LEN, n_other_ptl_in_first_soa);
522  }
523  if(last_soa_is_partial){
524  restore_backup_SoA(ptl_last_soa, end_idx, n_other_ptl_in_last_soa);
525  }
526  }
527 
535  template<typename F>
536  inline void for_all_particles(const std::string label, F lambda_func,
537  const PtlMvmt mvmt, LaunchBounds launch_bounds=LaunchBounds::Default) {
538  if(!owns_particles_d) exit_XGC("\nSpecies tried to loop over particles on device, but doesn't own the device array.");
539 
540  bool use_streaming = stream_particles;
541 #ifndef USE_STREAMS
542  use_streaming = false; // Just to be safe, turn streams off here
543 #endif
544 
545  bool send_ptl = ( mvmt.send_opt==PtlMvmt::Send ||
547  bool return_ptl = ( mvmt.return_opt==PtlMvmt::Return ||
549 
550  // Don't need to stream if particles are already present and don't need to be returned
551  if((!send_ptl) && (!return_ptl)) use_streaming = false;
552 
553  if(use_streaming){
554 #ifdef USE_STREAMS
555  Streamed::Option stream_option = Streamed::Normal; // Send to device and back
556  if(!send_ptl) stream_option = Streamed::NoSend;
557  if(!return_ptl) stream_option = Streamed::NoReturn;
558 
559  // Execute streaming parallel_for
560  Streamed::parallel_for(label, n_ptl, lambda_func, stream_option, particles, particles_d);
561 #endif
562  }else{
563  if(send_ptl) TIMER("copy_ptl_to_device",copy_particles_to_device() );
564 
565  if (launch_bounds==LaunchBounds::Custom) {
566 #ifdef USE_EPUSH_LAUNCH_BOUNDS
567 # if !defined(PUSH_MAX_THREADS_PER_BLOCK) || !defined(PUSH_MIN_WARPS_PER_EU)
568 # error "USE_EPUSH_LAUNCH_BOUNDS requires PUSH_MAX_THREADS_PER_BLOCK and PUSH_MIN_WARPS_PER_EU to be defined"
569 # endif
570  Kokkos::RangePolicy<ExSpace, Kokkos::LaunchBounds<PUSH_MAX_THREADS_PER_BLOCK, PUSH_MIN_WARPS_PER_EU>>
571  particle_range_policy( 0, p_range<DeviceType>(n_ptl) );
572  Kokkos::parallel_for(label, Opt::require(particle_range_policy, Async), lambda_func);
573 #else
574  exit_XGC("\nERROR: LaunchBounds::Custom specified, but USE_EPUSH_LAUNCH_BOUNDS is not defined\n");
575 #endif
576  } else {
577  Kokkos::RangePolicy<ExSpace>
578  particle_range_policy( 0, p_range<DeviceType>(n_ptl) );
579  Kokkos::parallel_for(label, Opt::require(particle_range_policy, Async), lambda_func);
580  }
581 
582  if(return_ptl) TIMER("copy_ptl_from_device", copy_particles_from_device() );
583  }
584  }
585 
586  KOKKOS_INLINE_FUNCTION VecParticles* ptl() const{
587  return (VecParticles*)(&particles_d.access(0));
588  }
589 
590  KOKKOS_INLINE_FUNCTION VecPhase* ph0() const{
591  return (VecPhase*)(&phase0_d.access(0));
592  }
593 
595  for_all_particles("copy_to_phase0", KOKKOS_LAMBDA( const int idx ){
596  AoSoAIndices<DeviceType> inds(idx);
597  VecParticles* ptl_loc = species.ptl();
598  VecPhase* ph0_loc = species.ph0();
599  for (int i_simd = 0; i_simd<SIMD_SIZE; i_simd++){
600  int p_vec = inds.a + i_simd;
601  ph0_loc[inds.s].r[p_vec] = ptl_loc[inds.s].ph.r[p_vec];
602  ph0_loc[inds.s].z[p_vec] = ptl_loc[inds.s].ph.z[p_vec];
603  ph0_loc[inds.s].phi[p_vec] = ptl_loc[inds.s].ph.phi[p_vec];
604  ph0_loc[inds.s].rho[p_vec] = ptl_loc[inds.s].ph.rho[p_vec];
605  ph0_loc[inds.s].w1[p_vec] = ptl_loc[inds.s].ph.w1[p_vec];
606  ph0_loc[inds.s].w2[p_vec] = ptl_loc[inds.s].ph.w2[p_vec];
607  }
608  });
609  }
610 
613  // If particles are resident on the device, then use the host particle allocation as the backup
614  // Otherwise, resize the backup_particle array
617  }else{
618  backup_particles.resize(particles_d.size());
619  }
620 
621  // Copy particle data to backup
622  Cabana::deep_copy(backup_particles, particles_d);
624  } else {
625  // For ions, save phase to phase0
626  phase0_d.resize(particles_d.size());
627  copy_to_phase0(*this); // Separate function to avoid implicit copy into lambda
628  }
630  }
631 
634  // If particles are resident on device, don't need to resize host particles since they are already used as the backup.
635  // If not, resize host particle array to fit backup particles
639  }else{
640  // If particles are not resident on device, resize host particle array to fit backup particles
642  }
643 
644  // Resize device particle array to fit backed up particles
646 
647  // Restore particle data from backup
648  Cabana::deep_copy(particles_d, backup_particles);
649 
651  // If the backup particles are simply pointing to the host particles, then
652  // point them to their own 0-sized allocation when finished using them
653  // so that they don't make a copy when host particles get resized
654  backup_particles = Cabana::AoSoA<ParticleDataTypes,HostType,VEC_LEN>("backup_particles", 0);
655  }
656  } else {
657  // When ipc==2, copy phase0 to device
658  // First, resize device phase0 and reset pointer
659  phase0_d.resize(phase0.size());
660 
661  // Next copy data to phase0 on device
662  Cabana::deep_copy(phase0_d, phase0);
663 
664  // Fill buffer with realistic ptl data
666  }
667  particles_are_backed_up = false;
668  }
669 
670  KOKKOS_INLINE_FUNCTION void restore_phase_from_phase0(const AoSoAIndices<Device>& inds, SimdParticles& part_one) const {
671  VecPhase* ph0_loc = ph0();
672  for (int i_simd = 0; i_simd<SIMD_SIZE; i_simd++){
673  int p_vec = inds.a + i_simd;
674  part_one.ph.r[i_simd] = ph0_loc[inds.s].r[p_vec];
675  part_one.ph.z[i_simd] = ph0_loc[inds.s].z[p_vec];
676  part_one.ph.phi[i_simd] = ph0_loc[inds.s].phi[p_vec];
677  part_one.ph.rho[i_simd] = ph0_loc[inds.s].rho[p_vec];
678  part_one.ph.w1[i_simd] = ph0_loc[inds.s].w1[p_vec];
679  part_one.ph.w2[i_simd] = ph0_loc[inds.s].w2[p_vec];
680  }
681  }
682 
683  long long int get_total_n_ptl(){
684 #ifdef USE_MPI
685  long long int tmp_n_ptl = n_ptl;
686  long long int out_n_ptl = 0;
687  MPI_Allreduce(&tmp_n_ptl, &out_n_ptl, 1, MPI_LONG_LONG_INT, MPI_SUM, SML_COMM_WORLD);
688  return out_n_ptl;
689 #else
690  return (long long int)(n_ptl);
691 #endif
692  }
693 
695 #ifdef USE_MPI
696  int tmp_n_ptl = n_ptl;
697  int out_n_ptl = 0;
698  MPI_Allreduce(&tmp_n_ptl, &out_n_ptl, 1, MPI_INT, MPI_MAX, SML_COMM_WORLD);
699  return out_n_ptl;
700 #else
701  return n_ptl;
702 #endif
703  }
704 
705  // Gets the gyro_radius of a species based on equilibrium temperature
706  // inode is the LOCAL (poloidally decomposed) grid node index to get temperature
707  // smu_n is the normalized sqrt(mu)
708  // bfield is the magnetic field at inode
709  KOKKOS_INLINE_FUNCTION double get_f0_eq_gyro_radius(int inode, double smu_n, double bfield) const{
710  // Should replace UNIT_CHARGE*charge_eu with charge(?)
711  return smu_n*sqrt(mass*f0.temp_ev(inode)*EV_2_J) / (UNIT_CHARGE*charge_eu*bfield);
712  }
713 
714  // Gets the equilibrium thermal velocity of a species based on f0 temperature
715  // inode is the GLOBAL node index to get temperature
716  KOKKOS_INLINE_FUNCTION double get_f0_eq_thermal_velocity(int inode) const{
717  return thermal_velocity(mass, f0.temp_global(inode));
718  }
719 
720  // Gets the equilibrium thermal velocity of a species based on f0 temperature, on device
721  // inode is the local node index to get temperature
722  KOKKOS_INLINE_FUNCTION double get_f0_eq_thermal_velocity_lnode(int inode) const{
723  return thermal_velocity(mass, f0.temp_ev(inode));
724  }
725 
726  // Gets the equilibrium thermal velocity of a species based on f0 temperature, on host
727  // inode is the local node index to get temperature
728  KOKKOS_INLINE_FUNCTION double get_f0_eq_thermal_velocity_lnode_h(int inode) const{
729  return thermal_velocity(mass, f0.temp_ev_h(inode));
730  }
731 
732  // Get species velocity
733  KOKKOS_INLINE_FUNCTION void get_particle_velocity_and_nearest_node(const Grid<DeviceType>& grid, const MagneticField<DeviceType>& magnetic_field, const DomainDecomposition<DeviceType>& pol_decomp, SimdParticles& part, Simd<double>& smu, Simd<double>& vp, Simd<int>& nearest_node, Simd<bool>& not_in_triangle, Simd<bool>& not_in_poloidal_domain) const{
734 
735  // This modulo surely doesnt need to be here (at least, should be elsewhere).
736  // Modulo phi coordinate
737  grid.wedge_modulo_phi(part.ph.phi);
738 
740  grid.get_grid_weights(magnetic_field, part.ph.v(), grid_wts0);
741 
742  // Output argument
743  for (int i_simd = 0; i_simd<SIMD_SIZE; i_simd++){
744  not_in_triangle[i_simd] = !grid_wts0.is_valid(i_simd);
745  }
746 
747  Simd<double> bmag;
748  magnetic_field.bmag_interpol(part.ph.v(), bmag);
749 
750  for (int i_simd = 0; i_simd<SIMD_SIZE; i_simd++){
751  if(!grid_wts0.is_valid(i_simd)) continue;
752 
753  nearest_node[i_simd]=grid_wts0.node[i_simd] - pol_decomp.node_offset;
754  not_in_poloidal_domain[i_simd] = (nearest_node[i_simd]<0 || nearest_node[i_simd]>=pol_decomp.nnodes);
755 
756  double temp_ev_norm = not_in_poloidal_domain[i_simd] ? f0.temp_ev(0) : f0.temp_ev(nearest_node[i_simd]);
757 
758  // get vp and smu
759  const double& B = bmag[i_simd];
760  vp[i_simd] = normalized_v_para(c_m, mass, B, temp_ev_norm, part.ph.rho[i_simd]);
761  smu[i_simd] = normalized_sqrt_mu(B, temp_ev_norm, part.ct.mu[i_simd]);
762  }
763  }
764 
766  const Grid<DeviceType>& grid,
768  const VelocityGrid& vgrid);
769 
772 
773  void write_ptl_checkpoint_files(const DomainDecomposition<DeviceType>& pol_decomp, const XGC_IO_Stream& stream, std::string sp_name);
774  void write_f0_checkpoint_files(const DomainDecomposition<DeviceType>& pol_decomp, const XGC_IO_Stream& stream, std::string sp_name);
775  void read_f0_checkpoint_files(const DomainDecomposition<DeviceType>& pol_decomp, const XGC_IO_Stream& stream, std::string sp_name);
776  void read_ptl_checkpoint_files(const DomainDecomposition<DeviceType>& pol_decomp, const XGC_IO_Stream& stream, std::string sp_name, bool n_ranks_is_same, int version);
777 
778  long long int get_max_gid() const;
779  void get_ptl_write_total_and_offsets(const DomainDecomposition<DeviceType>& pol_decomp, long long int& inum_total, long long int& ioff) const;
780 };
781 
782 #endif
void calculate_global_f0_arrays(const Grid< DeviceType > &grid, const MagneticField< DeviceType > &magnetic_field)
Definition: species.cpp:453
Cabana::AoSoA< PhaseDataTypes, HostType, VEC_LEN > phase0
Definition: species.hpp:113
bool stream_particles
Whether to stream particles between host and device if possible.
Definition: species.hpp:105
Definition: globals.hpp:84
KOKKOS_INLINE_FUNCTION VecPhase * ph0() const
Definition: species.hpp:590
KOKKOS_INLINE_FUNCTION int divide_and_round_up(int a, int b)
Definition: globals.hpp:193
bool owns_particles_d
Whether the species owns the device particle allocation right now.
Definition: species.hpp:102
void back_up_SoA(Cabana::AoSoA< ParticleDataTypes, Device, VEC_LEN > &backup_SoA, int offset, int n) const
Definition: species.hpp:423
KOKKOS_INLINE_FUNCTION double normalized_v_para(double c_m, double mass, double B, double temp_ev, double rho)
Definition: basic_physics.hpp:80
KOKKOS_INLINE_FUNCTION double get_f0_eq_gyro_radius(int inode, double smu_n, double bfield) const
Definition: species.hpp:709
subroutine adjust_n_ptl_for_core_ptl(n_ptl)
Definition: load.F90:472
void set_spall_num_and_ptr(int idx, int n_ptl, int n_vecs, VecParticles *ptl)
void for_particle_range(int begin_idx, int end_idx, const std::string label, F lambda_func) const
Definition: species.hpp:485
constexpr double PROTON_MASS
Definition: constants.hpp:7
MarkerType
Definition: globals.hpp:110
Distribution< Device > f0
Species distribution in velocity space on local mesh nodes.
Definition: species.hpp:121
MPI_Comm SML_COMM_WORLD
Definition: my_mpi.cpp:4
Cabana::AoSoA< ParticleDataTypes, HostType, VEC_LEN > backup_particles
Copy of particles to be restored for RK2.
Definition: species.hpp:117
bool is_electron
Whether this species is the electrons.
Definition: species.hpp:79
Eq::Profile< Device > eq_mk_temp
Definition: species.hpp:135
void for_all_particles(const std::string label, F lambda_func, const PtlMvmt mvmt, LaunchBounds launch_bounds=LaunchBounds::Default)
Definition: species.hpp:536
void save_backup_particles()
Definition: species.hpp:611
void update_decomposed_f0_calculations(const DomainDecomposition< DeviceType > &pol_decomp, const Grid< DeviceType > &grid, const MagneticField< DeviceType > &magnetic_field, const VelocityGrid &vgrid)
Definition: species.cpp:406
double c2_2m
c2/2m
Definition: species.hpp:87
double rho[VEC_LEN]
Definition: particles.hpp:96
void copy_to_phase0(Species< Device > &species)
Definition: species.hpp:594
Definition: species.hpp:58
Simd< double > w1
Definition: particles.hpp:22
double c_m
c/m
Definition: species.hpp:86
Definition: species.hpp:57
constexpr double EV_2_J
Conversion rate ev to J.
Definition: constants.hpp:5
Definition: velocity_grid.hpp:8
bool default_streaming_option()
Definition: species.hpp:28
Eq::Profile< Device > eq_den
Definition: species.hpp:126
Definition: globals.hpp:89
KOKKOS_INLINE_FUNCTION VecParticles * ptl() const
Definition: species.hpp:586
Definition: grid_weights.hpp:47
KOKKOS_INLINE_FUNCTION double thermal_velocity(double mass, double temp_ev)
Definition: basic_physics.hpp:58
Definition: NamelistReader.hpp:193
KinType kintype
Whether the species is gyrokinetic or drift kinetic.
Definition: species.hpp:82
Eq::Profile< Device > eq_fg_temp
Definition: species.hpp:130
Definition: magnetic_field.hpp:12
void get_ptl_write_total_and_offsets(const DomainDecomposition< DeviceType > &pol_decomp, long long int &inum_total, long long int &ioff) const
Definition: species.cpp:538
int add_vec_buffer(int n_ptl)
Definition: particles.hpp:194
int idx
Index in all_species.
Definition: species.hpp:78
Definition: particles.hpp:92
int a
The index in the inner array of the AoSoA.
Definition: particles.hpp:150
Definition: particles.hpp:109
KOKKOS_INLINE_FUNCTION void bmag_interpol(const SimdVector &v, Simd< double > &bmag) const
Definition: magnetic_field.tpp:239
bool particles_are_backed_up
Whether particles are currently backed up.
Definition: species.hpp:110
int nonadiabatic_idx
Index of species skipping adiabatic species (for compatibility with fortran arrays) ...
Definition: species.hpp:81
bool default_residence_option()
Definition: species.hpp:35
int n_ptl
Number of particles.
Definition: species.hpp:96
Eq::Profile< Device > eq_mk_flow
Definition: species.hpp:137
KOKKOS_INLINE_FUNCTION void get_particle_velocity_and_nearest_node(const Grid< DeviceType > &grid, const MagneticField< DeviceType > &magnetic_field, const DomainDecomposition< DeviceType > &pol_decomp, SimdParticles &part, Simd< double > &smu, Simd< double > &vp, Simd< int > &nearest_node, Simd< bool > &not_in_triangle, Simd< bool > &not_in_poloidal_domain) const
Definition: species.hpp:733
Definition: streamed_parallel_for.hpp:16
int node_offset
Offset of first mesh node belonging to this MPI rank.
Definition: domain_decomposition.hpp:55
Definition: streamed_parallel_for.hpp:14
Eq::Profile< Device > eq_fg_flow
Definition: species.hpp:132
void set_buffer_phase0_d()
Definition: species.hpp:394
long long int get_total_n_ptl()
Definition: species.hpp:683
void set_buffer_particles_d()
Definition: species.hpp:365
Simd< double > rho
Definition: particles.hpp:21
Definition: species.hpp:59
int p_range< DeviceType >(int num_particle)
Definition: particles.hpp:187
int eq_flow_type
Definition: species.hpp:128
double charge_eu
Particle charge in eu.
Definition: species.hpp:85
Definition: species.hpp:51
int nnodes
Number of nodes belonging to this MPI rank.
Definition: domain_decomposition.hpp:56
void resize_particles(int new_n_ptl)
Definition: species.hpp:235
double mass
Particle mass.
Definition: species.hpp:83
KOKKOS_INLINE_FUNCTION double get_f0_eq_thermal_velocity(int inode) const
Definition: species.hpp:716
Species(int idx_in, int nonadiabatic_idx_in, bool is_electron_in, bool is_adiabatic_in, KinType kintype_in, double mass_in, double charge_in, double charge_eu_in, bool is_deltaf_in, int ncycles_in)
Definition: species.cpp:25
void write_ptl_checkpoint_files(const DomainDecomposition< DeviceType > &pol_decomp, const XGC_IO_Stream &stream, std::string sp_name)
Definition: species.cpp:563
void for_all_particles(const std::string label, F lambda_func) const
Definition: species.hpp:418
Cabana::AoSoA< ParticleDataTypes, Device, VEC_LEN > particles_d
Particles on device.
Definition: species.hpp:100
Definition: species.hpp:48
double w2[VEC_LEN]
Definition: particles.hpp:98
#define TIMER(N, F)
Definition: timer_macro.hpp:24
Eq::Profile< Device > eq_mk_den
Definition: species.hpp:136
RKRestorationMethod
Definition: species.hpp:68
void copy_particles_to_device_if_not_resident()
Definition: species.hpp:348
RKRestorationMethod RK_restoration_method
Currently, electrons must use first method and ions must use second.
Definition: species.hpp:108
Simd< double > r
Definition: particles.hpp:18
void resize_host_particles_to_match_device()
Definition: species.hpp:255
Definition: species.hpp:69
KOKKOS_INLINE_FUNCTION SimdVector & v()
Definition: particles.hpp:39
ReturnOpt return_opt
Definition: species.hpp:63
ReturnOpt
Definition: species.hpp:56
Option
Definition: streamed_parallel_for.hpp:13
void restore_particles_from_backup()
Definition: species.hpp:632
MarkerType marker_type
Marker type: reduced delta-f, total-f, full-f, or none (placeholder for adiabatic species) ...
Definition: species.hpp:89
Definition: globals.hpp:90
static std::vector< MemoryPrediction > estimate_memory_usage(NLReader::NamelistReader &nlr, const Grid< DeviceType > &grid, const DomainDecomposition< DeviceType > &pol_decomp, int species_idx)
Definition: species.cpp:56
SendOpt send_opt
Definition: species.hpp:62
double charge
Particle charge.
Definition: species.hpp:84
SimdPhase ph
Definition: particles.hpp:59
void copy_particles_from_device()
Definition: species.hpp:323
KOKKOS_INLINE_FUNCTION double get_f0_eq_thermal_velocity_lnode(int inode) const
Definition: species.hpp:722
void copy_particles_from_device_if_not_resident()
Definition: species.hpp:354
KOKKOS_INLINE_FUNCTION void wedge_modulo_phi(Simd< double > &phi_mod) const
Definition: grid.tpp:86
void unassign_host_particles()
Definition: species.hpp:272
void read_ptl_checkpoint_files(const DomainDecomposition< DeviceType > &pol_decomp, const XGC_IO_Stream &stream, std::string sp_name, bool n_ranks_is_same, int version)
Definition: species.cpp:613
int ncycles_between_sorts
Number of subcycles between sorts.
Definition: species.hpp:93
Definition: particles.hpp:58
Cabana::AoSoA< PhaseDataTypes, Device, VEC_LEN > phase0_d
Definition: species.hpp:114
int SML_COMM_RANK
Definition: my_mpi.cpp:5
KinType
Definition: globals.hpp:88
Species(SpeciesType sp_type, int n_ptl)
Definition: species.hpp:151
bool is_deltaf
Whether this species is deltaf.
Definition: species.hpp:90
Definition: xgc_io.hpp:24
VecPhase ph
Definition: particles.hpp:110
Definition: species.hpp:44
Definition: species.hpp:70
void set_min_max_num(int isp, int n_ptl)
int minimum_ptl_reservation
The minimum reservation size for particles.
Definition: species.hpp:95
int s
The index in the outer array of the AoSoA.
Definition: particles.hpp:149
KOKKOS_INLINE_FUNCTION double normalized_sqrt_mu(double B, double temp_ev, double mu)
Definition: basic_physics.hpp:90
Simd< double > z
Definition: particles.hpp:19
void copy_particles_to_device_if_resident()
Definition: species.hpp:336
Definition: species.hpp:52
void resize_device_particles(int new_n_ptl)
Definition: species.hpp:293
Definition: species.hpp:53
constexpr double UNIT_CHARGE
Charge of an electron (C)
Definition: constants.hpp:4
long long int get_max_gid() const
Definition: species.cpp:515
KOKKOS_INLINE_FUNCTION double get_f0_eq_thermal_velocity_lnode_h(int inode) const
Definition: species.hpp:728
void exit_XGC(std::string msg)
Definition: globals.hpp:37
void copy_particles_from_device_if_resident()
Definition: species.hpp:342
bool is_adiabatic
Whether this species is adiabatic.
Definition: species.hpp:80
Simd< double > phi
Definition: particles.hpp:20
Definition: magnetic_field.F90:1
static constexpr const Kokkos::Experimental::WorkItemProperty::HintLightWeight_t Async
Definition: space_settings.hpp:82
int n_backup_particles
Definition: species.hpp:118
Eq::Profile< Device > eq_flow
Definition: species.hpp:127
Definition: streamed_parallel_for.hpp:15
SendOpt
Definition: species.hpp:50
void read_f0_checkpoint_files(const DomainDecomposition< DeviceType > &pol_decomp, const XGC_IO_Stream &stream, std::string sp_name)
Definition: species.cpp:708
SimdConstants ct
Definition: particles.hpp:60
int eq_fg_flow_type
Definition: species.hpp:133
void write_f0_checkpoint_files(const DomainDecomposition< DeviceType > &pol_decomp, const XGC_IO_Stream &stream, std::string sp_name)
Definition: species.cpp:672
void copy_particles_to_device()
Definition: species.hpp:307
KOKKOS_INLINE_FUNCTION void restore_phase_from_phase0(const AoSoAIndices< Device > &inds, SimdParticles &part_one) const
Definition: species.hpp:670
Species(int n_ptl_in)
Definition: species.hpp:181
double phi[VEC_LEN]
Definition: particles.hpp:95
Simd< double > w2
Definition: particles.hpp:23
double r[VEC_LEN]
Definition: particles.hpp:93
GyroAverageMatrices< HostType > gyro_avg_matrices
Definition: species.hpp:141
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
KOKKOS_INLINE_FUNCTION void get_grid_weights(const MagneticField< Device > &magnetic_field, const SimdVector &v, const Simd< double > &psi, SimdVector2D &xff, SimdGridWeights< Order::One, PIT > &grid_wts) const
Definition: grid.tpp:32
Definition: species.hpp:44
Eq::Profile< Device > eq_temp
Definition: species.hpp:125
bool particles_resident_on_device
Whether the particles can reside on device.
Definition: species.hpp:104
Simd< double > mu
Definition: particles.hpp:52
Eq::Profile< Device > eq_fg_den
Definition: species.hpp:131
PtlMvmt(SendOpt send_opt, ReturnOpt return_opt)
Definition: species.hpp:65
static int get_initial_n_ptl(NLReader::NamelistReader &nlr, const Grid< DeviceType > &grid, const DomainDecomposition< DeviceType > &pol_decomp, int sml_special, int species_idx, bool verbose)
Definition: species.cpp:113
int ncycles
Number of subcycles.
Definition: species.hpp:92
int eq_mk_flow_type
Definition: species.hpp:138
Definition: profile.hpp:171
int collision_grid_index
Which collision grid to use.
Definition: species.hpp:123
Definition: particles.hpp:148
SpeciesType
Definition: globals.hpp:83
double z[VEC_LEN]
Definition: particles.hpp:94
void restore_backup_SoA(Cabana::AoSoA< ParticleDataTypes, Device, VEC_LEN > &backup_SoA, int offset, int n) const
Definition: species.hpp:452
Definition: distribution.hpp:10
void resize_device_particles()
Definition: species.hpp:278
int get_max_n_ptl()
Definition: species.hpp:694
LaunchBounds
Definition: species.hpp:407
Cabana::AoSoA< ParticleDataTypes, HostType, VEC_LEN > particles
Particles.
Definition: species.hpp:97
double w1[VEC_LEN]
Definition: particles.hpp:97
Definition: species.hpp:44