XGCa
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
get_potential_grad.hpp
Go to the documentation of this file.
1 #ifndef GET_POTENTIAL_GRAD
2 #define GET_POTENTIAL_GRAD
3 
4 #include "sml.hpp"
5 #include "grid.hpp"
7 #include "magnetic_field.hpp"
8 #include "my_mirror_view.hpp"
9 #include "gradparx2.hpp"
11 #include "gradient_matrices.hpp"
12 #include "gyro_avg_mat.hpp"
13 #include "task_group.hpp"
14 #include "grid_field.hpp"
15 #include "grid_deriv.hpp"
16 
17 // Fortran pointer retrievals
18 
19 extern "C" int* get_psn_ff_1dp_tr_loc();
20 extern "C" double* get_psn_ff_1dp_p_loc();
21 extern "C" double* get_psn_ff_1dp_dx_loc();
22 extern "C" rtype* get_psn_pbd0_2_iseg_loc();
23 
24 extern "C" int get_psn_pbd0_2_nseg();
25 
26 // Fortran routines
27 extern "C" void get_field_Ah_cv_ff(Field<VarType::Scalar,PIT_GLOBAL>* psn_Ah_cv_ff);
28 extern "C" void get_pot_epar_em_filter(double* tmp, double* E_para_em);
29 
30 
31 template<class Device>
33  View<double**,CLayout,Device> field;
34 
35  public:
36 
37  AlternatingStorage(const std::string& name, int nnode)
38  : field(NoInit(name), 2, nnode)
39  {}
40 
41  // Functions to access the "alternating" Views
42  View<double*,CLayout, Device> left(int i_plane){
43  int i_left_alt = i_plane%2;
44  return my_subview(field,i_left_alt);
45  }
46 
47  View<double*,CLayout, Device> right(int i_plane){
48  int i_right_alt = (i_plane+1)%2;
49  return my_subview(field,i_right_alt);
50  }
51 };
52 
53 template<class Device, class DeviceOut>
55 
56  int nnode;
57  int nrho;
58  int nphi;
59  int ndim;
60 
61  View<double**,CLayout,Device> input_potential;
65  View<double**,CLayout,Device> potential;
66  View<double***,CLayout, Device> gradient;
67  View<double****,CLayout, Device> gradient_rho;
68  View<double***,CLayout, Device> potential_rho;
69 
70  View<double*,CLayout,Device> scratch;
71 
72  // Additional objects that should go elsewhere eventually, but are here
73  // because they are needed for get_pot_grad and can be reused for the different fields
78 
79  GetPotentialGradTemp(const Simulation<DeviceType>& sml, const Grid<DeviceType>& grid, const DomainDecomposition<DeviceType>& pol_decomp, const MagneticField<DeviceType>& magnetic_field, int n_input_potential_planes, bool gyroavg_requested, const GyroAverageMatrices<HostType>& gyro_avg_matrices_in = GyroAverageMatrices<HostType>())
80  : nnode(grid.nnode),
81  nrho(gyroavg_requested ? gyro_avg_matrices_in.nrho : 0),
82  nphi((PIT_GLOBAL==PhiInterpType::Planes) ? 2 : 1),
84  potential_alt("potential_alt", nnode),
85  gradient_r_alt("gradient_r_alt", nnode),
86  gradient_z_alt("gradient_r_alt", nnode),
87  potential(NoInit("potential"), nphi, nnode),
88  gradient(NoInit("gradient"), ndim, nphi, nnode),
89  potential_rho(NoInit("potential_rho"), nrho+1, nphi, nnode),
90  gradient_rho(NoInit("gradient_rho"), nrho+1, nphi, nnode, ndim),
91  gyro_avg_matrices(gyro_avg_matrices_in.template mirror<DeviceType>()) // Send gyromatrices to GPU
92  {
93  // Use scratch for transpose on device if the fields can't be written directly
94  if(!std::is_same<DeviceOut,Device>()){
95  scratch = View<double*,CLayout,Device>(NoInit("scratch"), gradient_rho.size());
96  input_potential = View<double**,CLayout,Device>(NoInit("input_potential"), n_input_potential_planes, nnode);
97  }
98 
99  GPTLstart("GET_POT_GRAD_MAT_SETUP");
100  grad_matrices = GradientMatrices<DeviceType>(true); // Copy in fortran data
101  GPTLstop("GET_POT_GRAD_MAT_SETUP");
102 
103  // This could be set up only once (if it doesnt take up too much device memory)
104  GPTLstart("GET_POT_GRAD_FF_SETUP");
105  if (!sml.is_XGCa) ff = FieldFollowingCoordinates(grid);
106  GPTLstop("GET_POT_GRAD_FF_SETUP");
107 
108  // Set up for parallel gradient
109  GPTLstart("GET_POT_GRAD_GPTX_SETUP");
110  if(!sml.is_XGCa) gptx = GradParXTmp(grid, get_psn_pbd0_2_nseg(), magnetic_field.bt_sign, get_psn_ff_1dp_tr_loc(),
112  GPTLstop("GET_POT_GRAD_GPTX_SETUP");
113  }
114 };
115 
116 template<class Device>
117 struct Field00{
120 
121  View<Field<VarType::Scalar,PhiInterpType::None>*,CLayout, Device> pot_managed; // Use if device != host
122  View<double*,CLayout, Device, Kokkos::MemoryTraits<Kokkos::Unmanaged>> pot;
123  View<double*,CLayout, Device> r;
124  View<double*,CLayout, Device> z;
125 
127  : is_provided(false)
128  {}
129 
130  template<class DeviceIn>
132  : r(NoInit("r"), field_in.f.extent(0)),
133  z(NoInit("z"), field_in.f.extent(0)),
134  is_provided(true),
135  discard_when_basis_is_one(discard_when_basis_is_one_in)
136  {
137  // Create a mirror if Device is different from the input DeviceIn
138  pot_managed = my_mirror_view(field_in.f, Device());
139  mirror_copy(pot_managed, field_in.f);
140  pot = View<double*,CLayout,Device, Kokkos::MemoryTraits<Kokkos::Unmanaged>>((double*)(pot_managed.data()), pot_managed.layout());
141  }
142 
143  void calculate_gradient(const Grid<DeviceType>& grid, const GradientMatrices<DeviceType>& grad_matrices){
144  grid_deriv(grad_matrices, grid,pot,r,z, discard_when_basis_is_one);
145 
146  auto r_tmp = r; // Use these to avoid using member function in kernel
147  auto z_tmp = z; // Use these to avoid using member function in kernel
148 
149  // rh This is general now --->
150  // We compute the electric field --> flip sign
151  // E=-E
152  Kokkos::parallel_for("get_pot_grad_flip_00", Kokkos::RangePolicy<ExSpace>(0,grid.nnode), KOKKOS_LAMBDA( const int i ){
153  r_tmp(i)=-r_tmp(i);
154  z_tmp(i)=-z_tmp(i);
155  });
156  }
157 
158  template<class DeviceOut, PhiInterpType PIT>
159  void set_field(const Grid<DeviceType>& grid, FieldFollowingCoordinates& ff, View<Field<VarType::Vector2D,PIT>*,CLayout,DeviceOut>& field00_ff_h){
160  if (is_provided){
161  // Copy to two planes if field following is needed
162  int nphi = (PIT==PhiInterpType::Planes) ? 2 : 1;
163  View<double**,CLayout, Device> r_tmp(NoInit("r_tmp"), nphi, r.size());
164  View<double**,CLayout, Device> z_tmp(NoInit("z_tmp"), nphi, z.size());
165 
166  for(int iphi=0; iphi<nphi; iphi++){
167  Kokkos::deep_copy(my_subview(r_tmp, iphi), r);
168  Kokkos::deep_copy(my_subview(z_tmp, iphi), z);
169  }
170 
171  if(PIT==PhiInterpType::Planes){
172  ff.cnvt_grid_real2ff(grid,r_tmp);
173  ff.cnvt_grid_real2ff(grid,z_tmp);
174  }
175 
176  auto field00_ff = my_mirror_view(field00_ff_h, Device());
177  for(int iphi=0; iphi<nphi; iphi++){
178  Kokkos::parallel_for("set_field", Kokkos::RangePolicy<ExSpace>(0,field00_ff_h.extent(0)), KOKKOS_LAMBDA( const int i ){
179 #ifdef XGC1
180  field00_ff(i).V[iphi][0]=r_tmp(iphi,i);
181  field00_ff(i).V[iphi][1]=z_tmp(iphi,i);
182 #else
183  field00_ff(i).E[0]=r_tmp(iphi,i);
184  field00_ff(i).E[1]=z_tmp(iphi,i);
185 #endif
186  });
187  }
188  Kokkos::fence();
189  // Copy if not a mirror
190  mirror_copy(field00_ff_h, field00_ff);
191  }else{
192  Kokkos::parallel_for("set_field", Kokkos::RangePolicy<typename DeviceOut::execution_space>(0,field00_ff_h.extent(0)), KOKKOS_LAMBDA( const int i ){
193 #ifdef XGC1
194  field00_ff_h(i).V[0][0]=0.0;
195  field00_ff_h(i).V[0][1]=0.0;
196  field00_ff_h(i).V[1][0]=0.0;
197  field00_ff_h(i).V[1][1]=0.0;
198 #else
199  field00_ff_h(i).E[0]=0.0;
200  field00_ff_h(i).E[1]=0.0;
201 #endif
202  });
203  }
204  }
205 };
206 
207 template<class Device, class DeviceOut>
208 struct EMParField{
209  bool requested;
210 
211  View<double**,CLayout, Device> field;
212  View<double***,CLayout, Device> field_rho;
213 
215 
217  : requested(false)
218  {}
219 
221  requested = true;
222  field_out = output_field;
223  field_rho = View<double***,CLayout, Device>(NoInit("E_para_em_rho"), field_out.f.extent(1), 2, field_out.f.extent(0));
224  field = View<double**,CLayout, Device>(NoInit("E_para_em"), 2, field_out.f.extent(0));
225  }
226 
227  void calculate(const Grid<DeviceType>& grid, const View<double**,CLayout, Device>& field_para){
228  // Obtain E_para with the same filters that are applied in push_As
229  // This filtered E_para must be used in the equation of motion
230  // for dA_s/dt in case of pullback-mode 4.
231  GPTLstart("GET_POT_EPAR_EM");
232  // tmp_copy gets copied in and modified internally
233  View<double**,CLayout, HostType> tmp_copy(NoInit("tmp"),2,grid.nnode);
234  Kokkos::deep_copy(tmp_copy, field_para);
235 
236  // Output:
237  auto field_h = my_mirror_view(field, HostType());
238 
239  get_pot_epar_em_filter(tmp_copy.data(), field_h.data());
240 
241  // Copy result of get_pot_epar_em_filter back to device
242  mirror_copy(field, field_h);
243  GPTLstop("GET_POT_EPAR_EM");
244  }
245 };
246 
247 template<class DeviceIn, class DeviceOut, VarType VT, PhiInterpType PIT, TorType TT, KinType KT>
249  View<double**,CLayout,DeviceIn, Kokkos::MemoryTraits<Kokkos::Unmanaged>> input_potential;
257 
258  // Constructor if passed an unmanaged view
259  GetPotGradFieldArgs(const View<double**,CLayout,DeviceIn, Kokkos::MemoryTraits<Kokkos::Unmanaged>>& input_potential_in, bool ignore_poloidal_dpot_in=false)
260  : input_potential(input_potential_in),
261  ignore_poloidal_dpot(ignore_poloidal_dpot_in),
262  potential_is_requested(false),
263  gradient_is_requested(false)
264  {}
265 
266  // Constructor if passed a 1D managed view
267  GetPotGradFieldArgs(const View<double*,CLayout,DeviceIn>& input_potential_in, bool ignore_poloidal_dpot_in=false)
268  : input_potential(View<double**,CLayout,DeviceIn,Kokkos::MemoryTraits<Kokkos::Unmanaged>>((double*)(input_potential_in.data()), 1, input_potential_in.extent(0))),
269  ignore_poloidal_dpot(ignore_poloidal_dpot_in),
270  potential_is_requested(false),
271  gradient_is_requested(false)
272  {}
273 
275  potential = potential_in;
276  potential_is_requested = true;
277  }
278 
280  gradient = gradient_in;
281  gradient_is_requested = true;
282  }
283 
285  E_para_em.request(output_field);
286  }
287 };
288 
289 template<class DeviceIn, class DeviceOut, VarType VT, PhiInterpType PIT, TorType TT, KinType KT>
290 void get_field_grad(const Grid<DeviceType>& grid,
291  const DomainDecomposition<DeviceType>& pol_decomp,
295 
296 #endif
VarType
Definition: field.hpp:16
AlternatingStorage(const std::string &name, int nnode)
Definition: get_potential_grad.hpp:37
View< double **, CLayout, Device > input_potential
Definition: get_potential_grad.hpp:61
Definition: field_following_coordinates.hpp:18
static int GPTLstart(const char *name)
Definition: timer_macro.hpp:9
int nrho
Definition: get_potential_grad.hpp:57
View< double *, CLayout, Device > left(int i_plane)
Definition: get_potential_grad.hpp:42
Definition: gradparx2.hpp:7
constexpr VarType vec2d_if_axisym()
View< double *, CLayout, Device > scratch
Definition: get_potential_grad.hpp:70
int ndim
Definition: get_potential_grad.hpp:59
View< Field< VarType::Scalar, PhiInterpType::None > *, CLayout, Device > pot_managed
Definition: get_potential_grad.hpp:121
void mirror_copy(T1 &view_dest, const T2 &view_src)
Definition: my_mirror_view.hpp:122
bool ignore_poloidal_dpot
Definition: get_potential_grad.hpp:252
void request(const GridField< DeviceOut, VarType::Scalar, PIT_GLOBAL, TorType::OnePlane, KinType::GyroKin > &output_field)
Definition: get_potential_grad.hpp:220
void get_pot_epar_em_filter(double *tmp, double *E_para_em)
bool is_provided
Definition: get_potential_grad.hpp:118
Field00()
Definition: get_potential_grad.hpp:126
Definition: sml.hpp:8
Kokkos::Device< HostExSpace, HostMemSpace > HostType
Definition: space_settings.hpp:56
FieldFollowingCoordinates ff
Definition: get_potential_grad.hpp:76
View< double **, CLayout, Device > potential
Definition: get_potential_grad.hpp:65
int * get_psn_ff_1dp_tr_loc()
int nphi
Definition: get_potential_grad.hpp:58
Definition: magnetic_field.hpp:12
Definition: get_potential_grad.hpp:248
int get_psn_pbd0_2_nseg()
bool discard_when_basis_is_one
Definition: get_potential_grad.hpp:119
Definition: gradparx2.hpp:43
Definition: get_potential_grad.hpp:117
void cnvt_grid_real2ff(const Grid< DeviceType > &grid, const View< double **, CLayout, DeviceType > &input, const View< double **, Kokkos::LayoutRight, DeviceType > &output)
Definition: field_following_coordinates.hpp:373
Kokkos::LayoutRight CLayout
Definition: space_settings.hpp:67
EMParField< DeviceType, DeviceOut > E_para_em
Definition: get_potential_grad.hpp:256
GradientMatrices< DeviceType > grad_matrices
Definition: get_potential_grad.hpp:75
GridField< DeviceOut, VT, PIT, TT, KT > gradient
Definition: get_potential_grad.hpp:250
double * get_psn_ff_1dp_dx_loc()
bool potential_is_requested
Definition: get_potential_grad.hpp:253
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:428
Definition: grid_field.hpp:25
View< double *, CLayout, Device, Kokkos::MemoryTraits< Kokkos::Unmanaged > > pot
Definition: get_potential_grad.hpp:122
AlternatingStorage< Device > gradient_z_alt
Definition: get_potential_grad.hpp:64
bool gradient_is_requested
Definition: get_potential_grad.hpp:254
View< double **, CLayout, DeviceIn, Kokkos::MemoryTraits< Kokkos::Unmanaged > > input_potential
Definition: get_potential_grad.hpp:249
void get_field_Ah_cv_ff(Field< VarType::Scalar, PIT_GLOBAL > *psn_Ah_cv_ff)
View< double ***, CLayout, Device > gradient
Definition: get_potential_grad.hpp:66
PhiInterpType
Definition: globals.hpp:88
View< double ***, CLayout, Device > potential_rho
Definition: get_potential_grad.hpp:68
View< double *, CLayout, Device > right(int i_plane)
Definition: get_potential_grad.hpp:47
subroutine grid_deriv(grid, qty, qty_deriv_x, qty_deriv_y, psi_only)
Definition: search.F90:3994
constexpr PhiInterpType PIT_GLOBAL
Definition: globals.hpp:98
AlternatingStorage< Device > gradient_r_alt
Definition: get_potential_grad.hpp:63
double bt_sign
Whether toroidal field is reversed?
Definition: magnetic_field.hpp:19
View< double ***, CLayout, Device > field_rho
Definition: get_potential_grad.hpp:212
Field00< DeviceType > field00
Definition: get_potential_grad.hpp:255
GetPotGradFieldArgs(const View< double *, CLayout, DeviceIn > &input_potential_in, bool ignore_poloidal_dpot_in=false)
Definition: get_potential_grad.hpp:267
EMParField()
Definition: get_potential_grad.hpp:216
void request_potential(const GridField< DeviceOut, VarType::Scalar, PIT, TT, KT > &potential_in)
Definition: get_potential_grad.hpp:274
rtype * get_psn_pbd0_2_iseg_loc()
void request_para_em(const GridField< DeviceOut, VarType::Scalar, PIT_GLOBAL, TorType::OnePlane, KinType::GyroKin > &output_field)
Definition: get_potential_grad.hpp:284
GyroAverageMatrices< DeviceType > gyro_avg_matrices
Definition: get_potential_grad.hpp:74
View< double *, CLayout, Device > z
Definition: get_potential_grad.hpp:124
Field00(const GridField< DeviceIn, VarType::Scalar, PhiInterpType::None, TorType::OnePlane, KinType::DriftKin > &field_in, bool discard_when_basis_is_one_in)
Definition: get_potential_grad.hpp:131
int nnode
Definition: get_potential_grad.hpp:56
Definition: field.hpp:50
Definition: get_potential_grad.hpp:208
View< T *, CLayout, Device > my_mirror_view(const View< T *, CLayout, Device > &view, Device nd)
Definition: my_mirror_view.hpp:14
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
Kokkos::Device< ExSpace, MemSpace > DeviceType
Definition: space_settings.hpp:47
const bool is_XGCa
Equivalent to the preprocessor flag for now.
Definition: sml.hpp:15
GridField< DeviceOut, VarType::Scalar, PIT, TT, KT > potential
Definition: get_potential_grad.hpp:251
void calculate(const Grid< DeviceType > &grid, const View< double **, CLayout, Device > &field_para)
Definition: get_potential_grad.hpp:227
Definition: get_potential_grad.hpp:32
AlternatingStorage< Device > potential_alt
Definition: get_potential_grad.hpp:62
void request_gradient(const GridField< DeviceOut, VT, PIT, TT, KT > &gradient_in)
Definition: get_potential_grad.hpp:279
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
void calculate_gradient(const Grid< DeviceType > &grid, const GradientMatrices< DeviceType > &grad_matrices)
Definition: get_potential_grad.hpp:143
View< double **, CLayout, Device > field
Definition: get_potential_grad.hpp:211
GridField< DeviceOut, VarType::Scalar, PIT_GLOBAL, TorType::OnePlane, KinType::GyroKin > field_out
Definition: get_potential_grad.hpp:214
int nnode
Number of grid nodes.
Definition: grid.hpp:216
GradParXTmp gptx
Definition: get_potential_grad.hpp:77
View< double *, CLayout, Device > r
Definition: get_potential_grad.hpp:123
double * get_psn_ff_1dp_p_loc()
View< double **, CLayout, Device > field
Definition: get_potential_grad.hpp:33
Kokkos::ViewAllocateWithoutInitializing NoInit
Definition: space_settings.hpp:68
bool requested
Definition: get_potential_grad.hpp:209
GetPotGradFieldArgs(const View< double **, CLayout, DeviceIn, Kokkos::MemoryTraits< Kokkos::Unmanaged >> &input_potential_in, bool ignore_poloidal_dpot_in=false)
Definition: get_potential_grad.hpp:259
View< double ****, CLayout, Device > gradient_rho
Definition: get_potential_grad.hpp:67
static int GPTLstop(const char *name)
Definition: timer_macro.hpp:10
void set_field(const Grid< DeviceType > &grid, FieldFollowingCoordinates &ff, View< Field< VarType::Vector2D, PIT > *, CLayout, DeviceOut > &field00_ff_h)
Definition: get_potential_grad.hpp:159
Definition: get_potential_grad.hpp:54