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