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