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