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