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