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