XGC1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
electric_field.hpp
Go to the documentation of this file.
1 #ifndef ELECTRIC_FIELD_HPP
2 #define ELECTRIC_FIELD_HPP
3 #include "space_settings.hpp"
4 #include "NamelistReader.hpp"
5 #include "grid.hpp"
6 #include "grid_field_pack.hpp"
7 #include "timer_macro.hpp"
8 #include "smoothing.hpp"
9 #include "get_potential_grad.hpp"
10 #include "plane_field_gatherer.hpp"
12 #ifdef SONIC_GK
13 #include "get_sonic_fields.hpp"
14 #endif
15 #ifdef EXPLICIT_EM
16 #include "globals.hpp"
17 #endif
18 #include "bias_potential.hpp"
19 
20 // Fortran calls
21 #ifdef XGC1
23 # ifdef EXPLICIT_EM
25  double* dpot_es
26 # endif
27  );
28 #else
30 # ifdef SONIC_GK
34 # endif
35  );
36 #endif
37 extern "C" void set_pot_fort(int n, int npsi00, double* pot0m, double* pot00_1d);
38 
39 // Electric field class
40 template<class Device>
41 struct ElectricField {
42 
43  // Dimensions
44  int nnode;
45  int nphi;
46 
47  // Gyromatrix info (to be moved soon)
48  int nrho;
49 
50  // Optimization options
52 
53  // Field variables
54  bool turb_efield;
55  bool E00_efield;
57 
58  // Host Views
61 
62  //For EM
73 
76 //#ifdef XGC1
78 //#else
79 //# ifdef SONIC_GK
83 //# endif
84 //#endif
88  View<double***,CLayout, HostType> save_dpot0; // Used for adiabatic response
89  View<double***,CLayout, HostType> save_dpot; // Used for adiabatic response
90  View<double***,CLayout, HostType> dpotsave; // Used for split-weight scheme
91  View<double*,CLayout, HostType> dpot_es; // Used for a diagnostic
92 
93  View<double*,CLayout, HostType> pot0m;
94  View<double*,CLayout, HostType> pot00_1d;
95 
96  // Host View for the adjustable loop voltage current drive
100 
101  // Bias potential
103 
104  // Default constructor
106 
107  // Makes sense to use the Grid in the constructor here, but the Grid still copies from fortran so until that changes,
108  // it's simpler to not use it here so that ElectricField can be defined before the restart, which uses psn%As
109  ElectricField(NLReader::NamelistReader& nlr, const Grid<DeviceType>& grid, int nrho_in, double dt=1.0, double psi_norm=1.0) : //Grid<HostTypeType>& grid) :
110  nnode(grid.nnode),
111  nphi(grid.nplanes),
112  nrho(nrho_in),
113  // Host Views
114 #ifdef XGC1
115  pot("pot",nphi,nrho,nnode),
116  E("E",nphi,nrho,nnode),
117 # if defined(DELTAF_CONV) && !defined(EXPLICIT_EM)
118  E00_ff_h("E00_ff_h",nnode),
119  ddpotdt("ddpotdt",nphi,nrho,nnode),
120  dpotsave("dpotsave", 2, 2, nnode), // View. dims: RK2, ff, nnode
121 # endif
122  save_dpot("save_dpot", 2, 2, nnode), // View. dims: RK2, ff, nnode
123  save_dpot0("save_dpot0", 2, 2, nnode), // View. dims: RK2, ff, nnode
124 # ifdef EXPLICIT_EM
125  dpot_es("dpot_es", nnode),
126  dAh("dAh",nphi,nrho,nnode),
127  Ah("Ah",nphi,nrho,nnode),
128  Ah_h("Ah_h",4,nnode), // -1:2 in fortran
129  As_h("As_h",4,nnode), // -1:2 in fortran
130  As_h_backup("As_h_backup",4,nnode), // -1:2 in fortran
131 # endif
132  dpot_h("dpot_h",4,nnode), // -1:2 in fortran
133  dpot_ff_h("dpot_ff_h",nnode),
134  dpot_n0_h("dpot_n0_h",nnode),
135 #else
136  E("E",nrho, nnode),
137  dpot_h("dpot_h",1,nnode),
138 # ifdef SONIC_GK
139  dEr_B2("dEr_B2",nrho, nnode),
140  dEz_B2("dEz_B2",nrho, nnode),
141  du2_E("du2_E",nrho, nnode),
142 # endif
143  dpot_n0_h("dpot_n0_h",nnode),
144 #endif
145  pot0m("pot0m", nnode),
146  pot00_1d("pot00_1d", grid.psi00.n),
147  pot0_h("pot0_h",nnode),
148  loop_voltage_h("loop_voltage_h",nnode),
149  jpara_diff_previous_h("jpara_diff_previous_h",nnode),
150  jpara_diff_integral_h("jpara_diff_integral_h",nnode),
151  bias_potential(nlr, dt, psi_norm)
152  {
153  nlr.use_namelist("sml_param");
154  turb_efield = nlr.get<bool>("sml_turb_efield",true);
155  E00_efield = nlr.get<bool>("sml_00_efield",true);
156  n0_m_efield = nlr.get<bool>("sml_0m_efield",true);
157 
158 #ifdef EXPLICIT_EM
159  // Pass these (or sml) as arguments
160  bool em_mixed_variable = nlr.get<bool>("sml_em_mixed_variable",true);
161  bool em_control_variate = nlr.get<bool>("sml_em_control_variate",false);
162  int em_pullback_mode = nlr.get<int>("sml_em_pullback_mode",3);
163 
164  if (!em_mixed_variable)
165  exit_XGC("\nError: sml_em_mixed_variable is explicitly set to .false. which is incompatible with EXPLICIT_EM\n");
166 #endif
167 
168  if(nlr.namelist_present("performance_param")){
169  nlr.use_namelist("performance_param");
170  calculate_phi_ff_on_device = nlr.get<bool>("calculate_phi_ff_on_device", false);
171  }else{
173  }
174 
175 #ifdef EXPLICIT_EM
176  if (em_control_variate){
179  }
180 
181  if (em_mixed_variable){
184  if (em_pullback_mode==4){
186  }
187  }
188 #endif
189 
190  // Set pointers in Fortran to rho arrays
191 #ifndef NO_FORTRAN_MODULES
192  set_pot_fort(grid.nnode, grid.psi00.n, pot0m.data(), pot00_1d.data());
193 #endif
194 #ifdef XGC1
195  set_rho_ff_pointers(nrho, nnode, pot0_h.data(), dpot_n0_h.data(), dpot_h.data(), dpot_ff_h.data(), pot.data(), E.data()
196 # ifdef EXPLICIT_EM
197  , Ah_h.data(), Ah.data(), dAh.data(), As_h.data(), As_h_backup.data(), As.data(), dAs.data(), Epar_em.data(), dpot_es.data()
198 # endif
199  );
200 #else
201  set_rho_pointers(nrho, nnode, pot0_h.data(), dpot_h.data(), E.data()
202 # ifdef SONIC_GK
203  , dEr_B2.data()
204  , dEz_B2.data()
205  , du2_E.data()
206 # endif
207  );
208 #endif
209  // Initialize loop voltage
210  nlr.use_namelist("src_param");
211  const double initial_loop_vol = nlr.get<double>("src_loop_voltage",0.0);
212  //
213  for (int inode=0; inode<nnode; inode++){
214  loop_voltage_h(inode) = initial_loop_vol;
215  }
216  }
217 
218 
219  /**************************************************/
220  /************ Copy to GridFieldPack *************/
221  /**************************************************/
222 
223  // Resizes allocations for large arrays that have a phi dimension (i.e used for electron push)
224  template<PhiInterpType PIT>
226  // Re-initialize new gridpack - deallocates whatever Views were already allocated
227  GridFieldPack<Device, PIT> gfpack(kintype_in, turb_efield);
228  //
229  // Copy axisymmetric loop voltage to temporary device-view
230  //
231  ScalarGridField loop_voltage("loop_voltage",loop_voltage_h.nnode());
232  Kokkos::deep_copy(loop_voltage.f,loop_voltage_h.f);
233  //
234 #ifdef XGC1
235  // Copy field data to device (or shallow copy if no device)
236  if(kintype_in==DriftKin){
237  GPTLstart("GAT_FIELD_INFO");
238  if(pol_decomp.decompose_fields){
239  if(near_field){
240  // Gathering the "near field", defined as the field decomposition region that contains this rank's
241  // global domain decomposition
242  int local_pid = pol_decomp.field_decomp.find_domain_owner(pol_decomp.plane_index, grid.nplanes, pol_decomp.node_offset, grid.nnode);
243  gfpack.phi_offset = pol_decomp.field_decomp.all_first_plane(local_pid);
244  gfpack.node_offset = pol_decomp.field_decomp.all_first_node(local_pid);
245  }else{
246  gfpack.phi_offset = pol_decomp.field_decomp.first_plane;
247  gfpack.node_offset = pol_decomp.field_decomp.first_node;
248  }
249  }
250  PlaneFieldGatherer pfg(pol_decomp, grid, near_field);
251 
252 # if defined(DELTAF_CONV) && !defined(EXPLICIT_EM)
253  if(ddpotdt.f.size()>0) pfg.gather_phi_ff_on_device(ddpotdt.f, gfpack.template get<Label::ddpotdt>().f);
254  grid_field_copy(gfpack.template get<Label::E00>(), E00_ff_h); // No gather_phi needed due to axisymmetry
255 # endif
256 
258  // Allocate temporary views
259  GPTLstart("GET_POT_GRAD_ALLOC_VIEWS");
260  int n_input_field_planes = pot0_h.nphi();
261  bool gradient_requested = true;
262  bool gyroaverage_requested = false;
263  GetPotentialGradTemp<DeviceType, DeviceType> gpg_tmp(sml, grid, pol_decomp, magnetic_field, n_input_field_planes, gradient_requested, gyroaverage_requested);
264  GPTLstop("GET_POT_GRAD_ALLOC_VIEWS");
265 
266  // Whether to ignore dpot
267  bool ignore_poloidal_dpot = !(turb_efield);
268  bool use_field00 = E00_efield;
269  bool potential_is_requested = false; // Not used in push
271  pfg.calculate_phi_ff_on_device(sml,grid, magnetic_field, smoothing, gpg_tmp, my_subview(dpot_h.f, 2), pot0_h, gfpack.template get<Label::E>(), unused_pot, potential_is_requested, use_field00, ignore_poloidal_dpot);
272 # ifdef EXPLICIT_EM
274  pfg.calculate_phi_ff_on_device(sml,grid, magnetic_field, smoothing, gpg_tmp, my_subview(Ah_h.f, 2), pot0_em_h, gfpack.template get<Label::dAh>(), gfpack.template get<Label::Ah>());
275  if(As.f.size()>0) pfg.calculate_phi_ff_on_device(sml,grid, magnetic_field, smoothing, gpg_tmp, my_subview(As_h.f, 2), pot0_em_h, gfpack.template get<Label::dAs>(), gfpack.template get<Label::As>());
276 # endif
277  }else{
278  pfg.gather_phi_ff_on_device(E.f, gfpack.template get<Label::E>().f);
279 # ifdef EXPLICIT_EM
280  pfg.gather_phi_ff_on_device(dAh.f, gfpack.template get<Label::dAh>().f);
281  pfg.gather_phi_ff_on_device(Ah.f, gfpack.template get<Label::Ah>().f);
282  if(dAs.f.size()>0) pfg.gather_phi_ff_on_device(dAs.f, gfpack.template get<Label::dAs>().f);
283  if(As.f.size()>0) pfg.gather_phi_ff_on_device(As.f, gfpack.template get<Label::As>().f);
284 # endif
285  }
286 # ifdef EXPLICIT_EM
287  if(Epar_em.f.size()>0) pfg.gather_phi_ff_on_device(Epar_em.f, gfpack.template get<Label::Epar_em>().f);
288 # endif
289  GPTLstop("GAT_FIELD_INFO");
290  }else{
291  GPTLstart("Copy_rho_ff_fields_to_device");
292  grid_field_copy(gfpack.template get<Label::E_gyro>(), E);
293 # ifdef EXPLICIT_EM
294  grid_field_copy(gfpack.template get<Label::dAs_gyro>(), dAs);
295  grid_field_copy(gfpack.template get<Label::dAh_gyro>(), dAh);
296  grid_field_copy(gfpack.template get<Label::As_gyro>(), As);
297  grid_field_copy(gfpack.template get<Label::Ah_gyro>(), Ah);
298  grid_field_copy(gfpack.template get<Label::Epar_em_gyro>(), Epar_em);
299 # endif
300  GPTLstop("Copy_rho_ff_fields_to_device");
301  }
302  //
303  // Loop voltage is identical for gyrokinetic and drift-kinetic species,
304  // but different between XGC1 and XGCa
305  //
309 #else
310  grid_field_copy(gfpack.E_rho, E);
311 # ifdef SONIC_GK
312  grid_field_copy(gfpack.dEr_B2_rho, dEr_B2);
313  grid_field_copy(gfpack.dEz_B2_rho, dEz_B2);
314  grid_field_copy(gfpack.du2_E_rho, du2_E);
315 # endif
317 #endif
318  //grid_field_copy(gfpack.loop_voltage,loop_voltage_h);
319  return gfpack;
320  }
321 
322  /* Set of functions used to allocate/deallocate gpu memory for Ah during particle scatters */
323  // Resizes allocation of Ah_cv
324  template<PhiInterpType PIT>
326  GridFieldPack<Device, PIT> gfpack(kintype_in, turb_efield);
327  // No need to gather Ah_phi_ff, because copy_Ah_cv_to_device_as_Ah is only used in scatter, which already assumes particles are
328  // on the correct plane!
329 #ifdef EXPLICIT_EM
330  constexpr int ZERO_PHI = 0;
331  gfpack.template get<Label::Ah>().f = View<Field<VarType::Scalar,PIT_GLOBAL>**, CLayout,DeviceType>("gfpack_view",grid.nplanes, grid.nnode);
332  Kokkos::deep_copy(my_subview(gfpack.template get<Label::Ah>().f,ZERO_PHI), Ah_cv.f);
333 #endif
334  return gfpack;
335  }
336 
337  template<PhiInterpType PIT>
339  GridFieldPack<Device, PIT> gfpack(kintype_in, turb_efield);
340 #ifdef EXPLICIT_EM
341  if(kintype_in==DriftKin){
342  // No need to gather Ah_phi_ff, because copy_Ah_ff_to_device is only used in scatter and pullback, which already assume particles are
343  // on the correct plane!
344  constexpr int ZERO_GYRO = 0;
345  constexpr int ZERO_PHI = 0;
347  Kokkos::parallel_for("Ah_ff_rhozero", Kokkos::RangePolicy<HostExSpace>( 0, grid.nnode), [=] (const int inode){
348  Ah_tmp(inode) = Ah(inode,ZERO_GYRO);
349  });
350  gfpack.template get<Label::Ah>().f = View<Field<VarType::Scalar,PIT_GLOBAL>**, CLayout,DeviceType>("gfpack_view",grid.nplanes, grid.nnode);
351  Kokkos::deep_copy(my_subview(gfpack.template get<Label::Ah>().f,ZERO_PHI), Ah_tmp.f);
352  }else{
353  grid_field_copy(gfpack.template get<Label::Ah_gyro>(), Ah);
354  }
355 #endif
356  return gfpack;
357  }
358 
359 
360  // Resizes allocation of dpot_ff and copies to device
361  template<PhiInterpType PIT>
363  GridFieldPack<Device, PIT> gfpack(kintype_in, turb_efield);
364 #ifdef XGC1
365  grid_field_copy(gfpack.dpot_ff, dpot_ff_h);
366 #else
367  const int PLANE_ZERO=0;
368  grid_field_copy(gfpack.dpot, dpot_h.subfield(PLANE_ZERO));
369 #endif
370  return gfpack;
371  }
372 
383  GPTLstart("GET_AH_RHO_FF_EXCL_DESTR"); // Timer excludes the destructors which are non-negligible
384 
385  // Allocate temporary views
386  GPTLstart("GET_AH_RHO_FF_ALLOC_VIEWS");
387  int n_input_field_planes = Ah_h.nphi();
388  bool gradient_requested = false;
389  bool gyroaverage_requested = true;
391  GetPotentialGradTemp<DeviceType, HostType> tmp(sml, grid, pol_decomp, magnetic_field, n_input_field_planes, gradient_requested, gyroaverage_requested, species.gyro_avg_matrices);
392  GPTLstop("GET_AH_RHO_FF_ALLOC_VIEWS");
393 
394  /****************************************/
395  // Calculation of the Hamiltonian vector potential A_h
396  /****************************************/
398  FieldArgs gpg_field_args(Ah_h.unmanaged());
399  gpg_field_args.request_potential(Ah);
400  get_field_grad(sml, grid, pol_decomp, magnetic_field, smoothing, gpg_field_args, tmp);
401 
402  GPTLstop("GET_AH_RHO_FF_EXCL_DESTR");
403 
405  }
406 
417  void get_potential_grad(const Simulation<DeviceType>& sml, const Grid<DeviceType>& grid, const DomainDecomposition<DeviceType>& pol_decomp, const MagneticField<DeviceType>& magnetic_field, Smoothing& smoothing, const View<double*, CLayout, HostType>& spitzer_resistivity, Species<DeviceType>& species){
418  GPTLstart("GET_POT_GRAD_EXCL_DESTR"); // Timer excludes the destructors which are non-negligible
419 
420  // Allocate temporary views
421  GPTLstart("GET_POT_GRAD_ALLOC_VIEWS");
422  int n_input_field_planes = dpot_h.nphi();
423  bool gradient_requested = true;
424  bool gyroaverage_requested = true;
425 #ifdef MULTI_RATE
426  gyroaverage_requested = false; // This is a hack. This variable controls the size of the tmp arrays for get_pot_grad,
427  // but the rho_ff arrays will still exist and the rho=0 values will be correct
428 #endif
429  if(gyroaverage_requested) species.gyro_avg_matrices.copy_to_device_if_not_resident();
430  GetPotentialGradTemp<DeviceType, HostType> tmp(sml, grid, pol_decomp, magnetic_field, n_input_field_planes, gradient_requested, gyroaverage_requested, species.gyro_avg_matrices);
431  GPTLstop("GET_POT_GRAD_ALLOC_VIEWS");
432 
433  /****************************************/
434  // Calculation of the electric potential E
435  /****************************************/
436  GPTLstart("GET_POT_GRAD_PHI");
437 
439  FieldArgs gpg_field_args(dpot_h.unmanaged(), !(turb_efield));
440 
441  // Set up E00
442  if(E00_efield){
443  gpg_field_args.field00 = Field00<DeviceType>(pot0_h, sml.grad_psitheta);
444  gpg_field_args.field00.calculate_gradient(grid, tmp.grad_matrices);
445  }
446  // Write results of E00
447  if(sml.reduced_deltaf && sml.electron_on) gpg_field_args.field00.set_field(grid, tmp.ff, E00_ff_h.f);
448 
449  // Set up EMPar
450  if(sml.explicit_electromagnetic && (sml.em_mixed_variable && sml.em_pullback_mode==4)) gpg_field_args.request_para_em(Epar_em, sml.em_pullback_dampfac, spitzer_resistivity, gyroaverage_requested);
451 
452  // Request potential and/or gradient outputs
453  // XGC1 algorithm uses potential but XGCa doesnt
454  if(!sml.is_XGCa) gpg_field_args.request_potential(pot);
455  gpg_field_args.request_gradient(E);
456 
457  // Get potentials
458  get_field_grad(sml, grid, pol_decomp, magnetic_field, smoothing, gpg_field_args, tmp);
459 
460  GPTLstop("GET_POT_GRAD_PHI");
461 
462  if(sml.explicit_electromagnetic){
463  GPTLstart("GET_POT_GRAD_APAR");
464 
465  /****************************************/
466  // Calculation of the Hamiltonian vector potential A_h
467  /****************************************/
468  gpg_field_args = FieldArgs(Ah_h.unmanaged());
469  gpg_field_args.request_potential(Ah);
470 #ifndef MULTI_RATE
471  gpg_field_args.request_gradient(dAh);
472 #endif
473  get_field_grad(sml, grid, pol_decomp, magnetic_field, smoothing, gpg_field_args, tmp);
474 
475  /****************************************/
476  // Calculation of the symplectic vector potential A_s
477  /****************************************/
478  if (sml.em_mixed_variable){
479  gpg_field_args = FieldArgs(As_h.unmanaged());
480  gpg_field_args.request_potential(As);
481 #ifndef MULTI_RATE
482  gpg_field_args.request_gradient(dAs);
483 #endif
484  get_field_grad(sml, grid, pol_decomp, magnetic_field, smoothing, gpg_field_args, tmp);
485  }
486 
487  /****************************************/
488  // Calculation of the Hamiltonian vector potential A_h for the control-variate method
489  /****************************************/
490  if (sml.em_control_variate){
491  TIMER("GET_POT_APAR_CV",
492  get_field_Ah_cv_ff(grid, pol_decomp, Ah_cv_h, Ah_cv) );
493  }
494  GPTLstop("GET_POT_GRAD_APAR");
495  }
496 
497 #ifdef SONIC_GK
498  GPTLstart("GET_POT_GRAD_SONIC");
499  View<double*,CLayout,HostType> Er_B2_h("Er_B2_h",grid.nnode);
500  View<double*,CLayout,HostType> Ez_B2_h("Ez_B2_h",grid.nnode);
501  View<double*,CLayout,HostType> u2_E_h("u2_E_h",grid.nnode);
502 
503  TIMER("GET_POT_SONIC",
504  get_sonic_fields(grid, E, Er_B2_h, Ez_B2_h, u2_E_h) );
505 
506  gpg_field_args = FieldArgs(Er_B2_h);
507  gpg_field_args.request_gradient(dEr_B2);
508  get_field_grad(sml, grid, pol_decomp, magnetic_field, smoothing, gpg_field_args, tmp);
509 
510  gpg_field_args = FieldArgs(Ez_B2_h);
511  gpg_field_args.request_gradient(dEz_B2);
512  get_field_grad(sml, grid, pol_decomp, magnetic_field, smoothing, gpg_field_args, tmp);
513 
514  gpg_field_args = FieldArgs(u2_E_h);
515  gpg_field_args.request_gradient(du2_E);
516  get_field_grad(sml, grid, pol_decomp, magnetic_field, smoothing, gpg_field_args, tmp);
517 
518  GPTLstop("GET_POT_GRAD_SONIC");
519 #endif
520 
522 
523  GPTLstop("GET_POT_GRAD_EXCL_DESTR");
524  }
525 
527  nlr.use_namelist("sml_param");
528  int field_n_rho = nlr.get<int>("sml_grid_nrho", 6);
529 
530 #ifdef EXPLICIT_EM
531  int n_fields = 3 + 1; // E, dAs, dAh, + scalars As,Ah
532 #else
533  int n_fields = 1;
534 #endif
535  double host_field_GB_per_vertex = (53 + 8*grid.nplanes + 11*(field_n_rho + 1))*8*BYTES_TO_GB;
536  double cpu_memory_usage = grid.nnode*host_field_GB_per_vertex;
537 
538  double device_field_GB_per_vertex = std::max(n_fields*8*grid.nplanes, n_fields*8*(field_n_rho+1))*8*BYTES_TO_GB;
539  double gpu_memory_usage = grid.nnode*device_field_GB_per_vertex;
540 
541  MemoryPrediction memory_prediction{"Collisions", cpu_memory_usage, gpu_memory_usage};
542 
543  return memory_prediction;
544  }
545 };
546 
547 #endif
GridField< HostType, VarType::Scalar, PhiInterpType::Planes, TorType::OnePlane, KinType::DriftKin > ddpotdt
Time derivative of - Not field-following?
Definition: electric_field.hpp:60
bool decompose_fields
Whether to decompose fields.
Definition: domain_decomposition.hpp:55
GridField< HostType, VarType::Scalar, PIT_GLOBAL, TorType::OnePlane, KinType::GyroKin > pot
Definition: electric_field.hpp:74
Definition: field_following_coordinates.hpp:9
static int GPTLstart(const char *name)
Definition: timer_macro.hpp:9
BiasPotential bias_potential
Definition: electric_field.hpp:102
GridField< HostType, VarType::Scalar, PhiInterpType::None, TorType::OnePlane, KinType::DriftKin > jpara_diff_integral_h
Definition: electric_field.hpp:99
Definition: bias_potential.hpp:7
GridField< Device, VarType::Scalar, PIT, TorType::OnePlane, KinType::DriftKin > loop_voltage
Definition: grid_field_pack.hpp:91
int first_plane
First plane belonging to this rank, including ghost planes.
Definition: field_decomposition.hpp:35
void get_potential_grad(const Simulation< DeviceType > &sml, const Grid< DeviceType > &grid, const DomainDecomposition< DeviceType > &pol_decomp, const MagneticField< DeviceType > &magnetic_field, Smoothing &smoothing, const View< double *, CLayout, HostType > &spitzer_resistivity, Species< DeviceType > &species)
Definition: electric_field.hpp:417
GridField< HostType, vec2d_if_axisym< PIT_GLOBAL >), PIT_GLOBAL, TorType::OnePlane, KinType::GyroKin > dAh
Definition: electric_field.hpp:64
void set_rho_ff_pointers(int nrho, int n, Field< VarType::Scalar, PhiInterpType::None > *pot0, Field< VarType::Scalar, PhiInterpType::None > *dpot_n0, Field< VarType::Scalar, PhiInterpType::None > *dpot, Field< VarType::Scalar, PhiInterpType::Planes > *dpot_ff, Field< VarType::Scalar, PhiInterpType::Planes > *pot_rho_ff, Field< VarType::Vector, PhiInterpType::Planes > *E_rho_ff)
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
GridField< HostType, VarType::Scalar, PhiInterpType::None, TorType::OnePlane, KinType::DriftKin > loop_voltage_h
Definition: electric_field.hpp:97
GridField< HostType, VarType::Scalar, PIT_GLOBAL, TorType::OnePlane, KinType::GyroKin > Ah
Definition: electric_field.hpp:65
bool em_mixed_variable
Switch for use of mixed-variable formulation.
Definition: sml.hpp:77
GridField< HostType, VarType::Vector2D, PIT_GLOBAL, TorType::OnePlane, KinType::DriftKin > E00_ff_h
Radial electric field from in field-following format.
Definition: electric_field.hpp:59
View< double *, CLayout, HostType > dpot_es
Definition: electric_field.hpp:91
GridField< HostType, VarType::Scalar, PhiInterpType::Planes, TorType::OnePlane, KinType::DriftKin > dpot_ff_h
Definition: electric_field.hpp:77
Definition: sml.hpp:8
static constexpr bool is_XGCa
Equivalent to the preprocessor flag for now.
Definition: sml.hpp:19
static constexpr bool explicit_electromagnetic
Equivalent to the preprocessor flag for now.
Definition: sml.hpp:37
GridField< HostType, VarType::Scalar, PhiInterpType::None, TorType::MultiplePlanes, KinType::DriftKin > As_h
Definition: electric_field.hpp:66
void copy_to_device_if_not_resident()
Definition: gyro_avg_mat.hpp:50
int plane_index
Offset of local plane.
Definition: domain_decomposition.hpp:52
Definition: globals.hpp:89
Definition: plane_field_gatherer.hpp:9
FieldFollowingCoordinates ff
Definition: get_potential_grad.hpp:75
KOKKOS_INLINE_FUNCTION int find_domain_owner(int global_plane_index, int nplanes_total, int global_node_index, int nnodes_total) const
Definition: field_decomposition.hpp:136
int nrho
Definition: electric_field.hpp:48
void get_Ah_rho_ff(const Simulation< DeviceType > &sml, const Grid< DeviceType > &grid, const DomainDecomposition< DeviceType > &pol_decomp, const MagneticField< DeviceType > &magnetic_field, Smoothing &smoothing, Species< DeviceType > &species)
Definition: electric_field.hpp:382
Definition: NamelistReader.hpp:193
Definition: magnetic_field.hpp:14
static MemoryPrediction estimate_memory_usage(NLReader::NamelistReader &nlr, const Grid< DeviceType > &grid)
Definition: electric_field.hpp:526
Definition: get_potential_grad.hpp:285
int phi_offset
Offset for phi_ff field decomposition.
Definition: grid_field_pack.hpp:57
void cnvt_grid_real2ff(const Grid< DeviceType > &grid, const View< double **, CLayout, DeviceType > &input, const View< double **, Kokkos::LayoutRight, DeviceType > &output) const
Definition: field_following_coordinates.hpp:142
Definition: electric_field.hpp:41
void get_field_grad(const Simulation< DeviceType > &sml, const Grid< DeviceType > &grid, const DomainDecomposition< DeviceType > &pol_decomp, const MagneticField< DeviceType > &magnetic_field, Smoothing &smoothing, GetPotGradFieldArgs< DeviceIn, DeviceOut, VT, PIT, TT, KT > &args, GetPotentialGradTemp< DeviceType, DeviceOut > &tmp)
Definition: get_potential_grad.cpp:395
GridField< HostType, VarType::Scalar, PhiInterpType::None, TorType::OnePlane, KinType::DriftKin > pot0_h
Definition: electric_field.hpp:86
int nplanes
Number of planes.
Definition: grid.hpp:157
GridFieldPack< Device, PIT > copy_Ah_ff_to_device(KinType kintype_in, const DomainDecomposition< DeviceType > &pol_decomp, const Grid< DeviceType > &grid) const
Definition: electric_field.hpp:338
static constexpr bool reduced_deltaf
Equivalent to the preprocessor flag for now.
Definition: sml.hpp:25
Definition: grid_field_pack.hpp:21
int node_offset
Offset of first mesh node belonging to this MPI rank.
Definition: domain_decomposition.hpp:59
GridFieldPack< Device, PIT > copy_push_fields_to_device(KinType kintype_in, const Simulation< DeviceType > &sml, const DomainDecomposition< DeviceType > &pol_decomp, const MagneticField< DeviceType > &magnetic_field, const Grid< DeviceType > &grid, Smoothing &smoothing, bool near_field=false) const
Definition: electric_field.hpp:225
Kokkos::LayoutRight CLayout
Definition: space_settings.hpp:67
GradientMatrices< DeviceType > grad_matrices
Definition: get_potential_grad.hpp:74
Projection< HostType > half_plane_ff
Definition: grid.hpp:176
void deallocate_device_matrices_if_not_resident()
Definition: gyro_avg_mat.hpp:56
ElectricField(NLReader::NamelistReader &nlr, const Grid< DeviceType > &grid, int nrho_in, double dt=1.0, double psi_norm=1.0)
Definition: electric_field.hpp:109
View< double ***, CLayout, HostType > save_dpot
Definition: electric_field.hpp:89
void request_para_em(const GridField< DeviceOut, VarType::Scalar, PhiInterpType::Planes, TorType::OnePlane, KinType::GyroKin > &output_field, double em_pullback_dampfac_in, const View< double *, CLayout, HostType > &spitzer_resistivity, bool gyroaverage_requested)
Definition: get_potential_grad.hpp:321
void set_pot_fort(int n, int npsi00, double *pot0m, double *pot00_1d)
int em_pullback_mode
Definition: sml.hpp:79
bool grad_psitheta
Definition: sml.hpp:100
GridField< HostType, VarType::Scalar, PhiInterpType::Planes, TorType::OnePlane, KinType::DriftKin > Ah_cv
Definition: electric_field.hpp:71
ElectricField()
Definition: electric_field.hpp:105
#define TIMER(N, F)
Definition: timer_macro.hpp:24
GridField< HostType, VarType::Scalar, PhiInterpType::None, TorType::MultiplePlanes, KinType::DriftKin > Ah_cv_h
Definition: electric_field.hpp:70
void use_namelist(const string &namelist)
Definition: NamelistReader.hpp:355
constexpr PhiInterpType PIT_GLOBAL
Definition: globals.hpp:103
GridField< HostType, VarType::Scalar, PhiInterpType::None, TorType::MultiplePlanes, KinType::DriftKin > As_h_backup
Definition: electric_field.hpp:67
Definition: memory_prediction.hpp:4
void grid_field_copy(T &dest, const T &src)
Definition: grid_field.hpp:452
Definition: smoothing.hpp:8
GridField< HostType, VarType::Vector2D, PIT_GLOBAL, TorType::OnePlane, KinType::GyroKin > dEz_B2
Gradient of Ez/B^2.
Definition: electric_field.hpp:81
void gather_phi_ff_on_device(View< FT **, CLayout, HostType > &tmp, View< FT *, CLayout, HostType > &tmp_full, const T_h &rho_ff_h, View< FT **, CLayout, DeviceType > &phi_ff)
Definition: plane_field_gatherer.hpp:197
Field00< DeviceType > field00
Definition: get_potential_grad.hpp:292
Definition: globals.hpp:90
GridField< HostType, vec2d_if_axisym< PIT_GLOBAL >), PIT_GLOBAL, TorType::OnePlane, KinType::GyroKin > dAs
Definition: electric_field.hpp:68
void request_potential(const GridField< DeviceOut, VarType::Scalar, PIT, TT, KT > &potential_in)
Definition: get_potential_grad.hpp:311
GridField< HostType, VarType::Scalar, PhiInterpType::None, TorType::MultiplePlanes, KinType::DriftKin > Ah_h
Definition: electric_field.hpp:63
bool E00_efield
Flux-surface averaged potential not used for calculating the electric field if .false.
Definition: electric_field.hpp:55
GridField< HostType, VarType::Vector2D, PIT_GLOBAL, TorType::OnePlane, KinType::GyroKin > dEr_B2
Gradient of Er/B^2.
Definition: electric_field.hpp:80
View< double *, CLayout, HostType > pot00_1d
1D flux-surface averaged potential
Definition: electric_field.hpp:94
int node_offset
Offset for phi_ff field decomposition.
Definition: grid_field_pack.hpp:58
GridField< HostType, VarType::Vector2D, PIT_GLOBAL, TorType::OnePlane, KinType::GyroKin > du2_E
Gradient of the square (dot product) of the E x B drift.
Definition: electric_field.hpp:82
Definition: field.hpp:45
KinType
Definition: globals.hpp:88
GridField< HostType, VarType::Scalar, PhiInterpType::None, TorType::OnePlane, KinType::DriftKin > jpara_diff_previous_h
Definition: electric_field.hpp:98
View< int *, CLayout, HostType > all_first_plane
First plane of each rank.
Definition: field_decomposition.hpp:43
bool namelist_present(const string &namelist)
Definition: NamelistReader.hpp:351
bool electron_on
Use kinetic electrons.
Definition: sml.hpp:63
GridField< HostType, VarType::Scalar, PhiInterpType::Planes, TorType::OnePlane, KinType::GyroKin > Epar_em
Definition: electric_field.hpp:72
GridField< Device, VarType::Scalar, PIT, TorType::OnePlane, KinType::DriftKin > dpot_ff
Definition: grid_field_pack.hpp:79
void get_sonic_fields(const Grid< DeviceType > &grid, const GridField< HostType, vec2d_if_axisym< PIT_GLOBAL >(), PIT_GLOBAL, TorType::OnePlane, KinType::GyroKin > &E, View< double *, CLayout, HostType > &Er_B2_h, View< double *, CLayout, HostType > &Ez_B2_h, View< double *, CLayout, HostType > &u2_E_h)
bool calculate_phi_ff_on_device
Definition: electric_field.hpp:51
void exit_XGC(std::string msg)
Definition: globals.hpp:37
UniformRange psi00
Definition: grid.hpp:166
Kokkos::View< T *, Kokkos::LayoutRight, Device > my_subview(const Kokkos::View< T ****, Kokkos::LayoutRight, Device > &view, int i, int j, int k)
Definition: my_subview.hpp:8
Definition: magnetic_field.F90:1
bool turb_efield
E-field calculated only with if .false., psndpot will still contain all (n=0,|m|&gt;0) and (|n|&gt;0...
Definition: electric_field.hpp:54
Kokkos::Device< ExSpace, MemSpace > DeviceType
Definition: space_settings.hpp:47
View< double *, CLayout, HostType > pot0m
n=0 component of the electrostatic potential before smoothing
Definition: electric_field.hpp:93
View< double ***, CLayout, HostType > save_dpot0
Definition: electric_field.hpp:88
int nphi
Definition: electric_field.hpp:45
GridField< HostType, VarType::Scalar, PIT_GLOBAL, TorType::OnePlane, KinType::GyroKin > As
Definition: electric_field.hpp:69
int nnode
Definition: electric_field.hpp:44
FieldDecomposition< Device > field_decomp
Definition: domain_decomposition.hpp:56
Definition: species.hpp:75
void request_gradient(const GridField< DeviceOut, VT, PIT, TT, KT > &gradient_in)
Definition: get_potential_grad.hpp:316
void get_field_Ah_cv_ff(const Grid< DeviceType > &grid, const DomainDecomposition< DeviceType > &pol_decomp, const GridField< HostType, VarType::Scalar, PhiInterpType::None, TorType::MultiplePlanes, KinType::DriftKin > &Ah_cv_h, const GridField< HostType, VarType::Scalar, PhiInterpType::Planes, TorType::OnePlane, KinType::DriftKin > &Ah_cv_ff)
Definition: get_potential_grad.cpp:7
GridFieldPack< Device, PIT > copy_dpot_to_device(KinType kintype_in) const
Definition: electric_field.hpp:362
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
GridFieldPack< Device, PIT > copy_Ah_cv_to_device_as_Ah(KinType kintype_in, const DomainDecomposition< DeviceType > &pol_decomp, const Grid< DeviceType > &grid) const
Definition: electric_field.hpp:325
void calculate_gradient(const Grid< DeviceType > &grid, const GradientMatrices< DeviceType > &grad_matrices)
Definition: get_potential_grad.hpp:153
GyroAverageMatrices< Device > gyro_avg_matrices
Definition: species.hpp:146
GridField< HostType, VarType::Scalar, PhiInterpType::None, TorType::MultiplePlanes, KinType::DriftKin > dpot_h
Definition: electric_field.hpp:85
int nnode
Number of grid nodes.
Definition: grid.hpp:156
double em_pullback_dampfac
Damping term gamma on -b.grad(phi) in pullback mode 4.
Definition: sml.hpp:85
GridField< HostType, VarType::Scalar, PhiInterpType::None, TorType::OnePlane, KinType::DriftKin > dpot_n0_h
Definition: electric_field.hpp:87
bool n0_m_efield
Definition: electric_field.hpp:56
View< int *, CLayout, HostType > all_first_node
First node of each rank.
Definition: field_decomposition.hpp:41
bool em_control_variate
Switch for use of control variate method.
Definition: sml.hpp:78
int first_node
First mesh node belonging to this rank, including ghost nodes.
Definition: field_decomposition.hpp:32
GridField< HostType, vec2d_if_axisym< PIT_GLOBAL >), PIT_GLOBAL, TorType::OnePlane, KinType::GyroKin > E
Definition: electric_field.hpp:75
void set_field(const Grid< DeviceType > &grid, const FieldFollowingCoordinates &ff, View< Field< VarType::Vector2D, PIT > *, CLayout, DeviceOut > &field00_ff_h)
Definition: get_potential_grad.hpp:169
View< double ***, CLayout, HostType > dpotsave
Definition: electric_field.hpp:90
void calculate_phi_ff_on_device(const Simulation< DeviceType > &sml, const Grid< DeviceType > &grid, const MagneticField< DeviceType > &magnetic_field, Smoothing &smoothing, GetPotentialGradTemp< DeviceType, DeviceType > &tmp, const View< Field< VarType::Scalar, PhiInterpType::None > *, CLayout, HostType > &dpot_h, const GT0 &pot0_h, GT &phi_ff, GT2 &pot_phi_ff, bool potential_is_requested=true, bool use_field00=false, bool ignore_poloidal_dpot=false)
Definition: plane_field_gatherer.hpp:117
static int GPTLstop(const char *name)
Definition: timer_macro.hpp:10
int n
Definition: uniform_range.hpp:7
Definition: get_potential_grad.hpp:53