XGCa
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
field_following_coordinates.hpp
Go to the documentation of this file.
1 #ifndef FIELD_FOLLOWING_COORDINATES_HPP
2 #define FIELD_FOLLOWING_COORDINATES_HPP
3 
4 #include "my_subview.hpp"
5 #include "grid.hpp"
6 #include "grid_field.hpp"
7 
8 
10  View<int**,CLayout,DeviceType> tr;
11  View<double***,CLayout,DeviceType> p;
12 
13  // View for temporarily storing result data for in-place conversions
14  View<double**,CLayout,DeviceType> result;
15 
16  public:
17 
19 
21  : tr(NoInit("tr"),projection.tr.layout()),
22  p(NoInit("p"),projection.p.layout()),
23  result(NoInit("result"),tr.layout())
24  {
25  // Copy grid conversion info to device
26  Kokkos::deep_copy(tr, projection.tr);
27  Kokkos::deep_copy(p, projection.p);
28  }
29 
30  /* Conversion from field-following coordinates to toroidal coordinates
31  * @param[in] grid is the Grid object needed here to map elements to nodes
32  * @param[in] input is the array one wants to convert
33  * @param[out] output is the resulting converted array
34  * */
35  void cnvt_grid_ff2real(const Grid<DeviceType>& grid, const View<double**,CLayout,DeviceType>& input, const View<double**,Kokkos::LayoutRight,DeviceType>& output) const{
36  // Check the input and output views are the correct size
37  if(input.extent(0) != 2 || input.extent(1) != grid.nnode ||
38  output.extent(0)!= 2 || output.extent(1)!= grid.nnode ) exit_XGC("\nError in cnvt_grid_ff2real: Unexpected View sizes. \n");
39 
40  // Initialize output to zero
41  Kokkos::deep_copy(output, 0.0);
42 
43  // Initialize vol to zero
44  View<double**,CLayout,DeviceType> vol("vol",input.layout());
45 
46  // Copy grid%node_vol_ff from fortran
47  View<double**,CLayout,DeviceType> node_vol_ff(NoInit("node_vol_ff"),input.layout());
48  // For stellarators, Dir 0 uses left plane node volume, Dir 1 uses right plane node volume
49  Kokkos::deep_copy(my_subview(node_vol_ff,0), my_subview(grid.volumes_and_areas.node_vol_ff_l_h,0));
50  Kokkos::deep_copy(my_subview(node_vol_ff,1), my_subview(grid.volumes_and_areas.node_vol_ff_r_h,1));
51 
52  // Make a (SHALLOW) copy so the lambda doesn't have to implicitly capture *this
53  FieldFollowingCoordinates ff = *this;
54 
55  // for all node point of ff
56  Kokkos::parallel_for("cnvt_grid_ff2real", Kokkos::RangePolicy<ExSpace>(0,grid.nnode), KOKKOS_LAMBDA( const int i ){
57  for (int dir=0;dir<2;dir++){
58  for (int j=0;j<3;j++){
59  // node and weight
60  int nd=grid.get_node_index(ff.tr(dir,i), j);
61  double wt=ff.p(dir,i,j)*node_vol_ff(dir,nd);
62 
63  // accumulate values
64  Kokkos::atomic_add(&(output(dir,nd)), wt*input(dir,i));
65  Kokkos::atomic_add(&(vol(dir,nd)), wt);
66  }
67  }
68  });
69 
70  // Normalize to vol
71  Kokkos::parallel_for("cnvt_grid_ff2real_norm", Kokkos::RangePolicy<ExSpace>(0,grid.nnode), KOKKOS_LAMBDA( const int i ){
72  for (int dir=0;dir<2;dir++){
73  if(vol(dir, i) == 0.0){
74  output(dir, i) = 0.0;
75  }else{
76  output(dir, i) /= vol(dir, i);
77  }
78  }
79  });
80 
81  // fill empty area --- replace nan values
82  Kokkos::parallel_for("cnvt_grid_ff2real_fillempty", Kokkos::RangePolicy<ExSpace>(0,grid.nnode), KOKKOS_LAMBDA( const int i ){
83  for (int dir=0; dir<2; dir++){
84  if(vol(dir,i)==0.0){
85  int cdir = 1-dir;
86  output(dir,i)=0.0;
87  for (int j=0;j<3;j++){
88  int nd=grid.get_node_index(ff.tr(cdir,i), j);
89  double wt=ff.p(cdir,i,j);
90 
91  // obtain interpolated values
92  output(dir,i) += wt*input(dir,nd);
93  }
94  }
95  }
96  });
97  Kokkos::fence();
98  }
99 
100  /* Conversion from field-following coordinates to real coordinates when views are on host
101  * User also provides a temporary device allocation
102  * @param[in] grid is the Grid object needed here to map elements to nodes
103  * @param[in] input is the array one wants to convert
104  * @param[out] output is the resulting converted array
105  * @param[in] tmp is an array with the same size as input/output used for temporary storage
106  * */
107  void cnvt_grid_ff2real(const Grid<DeviceType>& grid, const View<double**,CLayout,HostType>& input, const View<double**,CLayout,HostType>& output, const View<double**,CLayout,DeviceType>& tmp) const{
108  // Copy input to tmp
109  Kokkos::deep_copy(tmp, input);
110 
111  // Convert
112  cnvt_grid_ff2real(grid, tmp, result);
113 
114  // Copy to output
115  Kokkos::deep_copy(output, result);
116  }
117 
118 #ifdef USE_GPU
119  /* Conversion from field-following coordinates to real coordinates when input view is on host but output view is on device
120  * User also provides a temporary device allocation
121  * @param[in] grid is the Grid object needed here to map elements to nodes
122  * @param[in] input is the array one wants to convert
123  * @param[out] output is the resulting converted array
124  * */
125  void cnvt_grid_ff2real(const Grid<DeviceType>& grid, const View<double**,CLayout,HostType>& input, const View<double**,CLayout,DeviceType>& output) const{
126  // Copy input to output
127  Kokkos::deep_copy(output, input);
128 
129  // Convert
130  cnvt_grid_ff2real(grid, output, result);
131 
132  // Copy to output
133  Kokkos::deep_copy(output, result);
134  }
135 #endif
136 
137  /* Conversion from real coordinates to field-following coordinates
138  * @param[in] grid is the Grid object needed here to map elements to nodes
139  * @param[in] input is the array one wants to convert
140  * @param[out] output is the resulting converted array
141  * */
142  void cnvt_grid_real2ff(const Grid<DeviceType>& grid, const View<double**,CLayout,DeviceType>& input, const View<double**,Kokkos::LayoutRight,DeviceType>& output) const{
143  // Check the input and output views are the correct size
144  if(input.extent(0) != 2 || input.extent(1) != grid.nnode ||
145  output.extent(0)!= 2 || output.extent(1)!= grid.nnode ) exit_XGC("\nError in cnvt_grid_real2ff: Unexpected View sizes. \n");
146 
147  // Initialize output to zero
148  Kokkos::deep_copy(output, 0.0);
149 
150 
151  // Make a (SHALLOW) copy so the lambda doesn't have to implicitly capture *this
152  FieldFollowingCoordinates ff = *this;
153 
154  Kokkos::parallel_for("cnvt_grid_real2ff", Kokkos::RangePolicy<ExSpace>(0,grid.nnode), KOKKOS_LAMBDA( const int i ){
155  for (int dir=0;dir<2;dir++){
156  for (int j=0;j<3;j++){
157  int nd=grid.get_node_index(ff.tr(dir,i), j);
158 
159  double wt=ff.p(dir,i,j);
160 
161  // obtain interpolated values
162  output(dir,i) += wt*input(dir,nd);
163  }
164  }
165  });
166  Kokkos::fence();
167  }
168 
169  /* In-place conversion from real coordinates to field-following coordinates, using a buffer array to store the result temporarily
170  * before copying back to the input array
171  * @param[in] grid is the Grid object needed here to map elements to nodes
172  * @param[inout] input is the array one wants to convert, and the resulting converted array
173  * */
174  void cnvt_grid_real2ff(const Grid<DeviceType>& grid, const View<double**,CLayout,DeviceType>& input) const{
175  cnvt_grid_real2ff(grid, input, result);
176 
177  // Copy back to original array
178  Kokkos::deep_copy(input, result);
179  }
180 
181  /* Conversion from real coordinates to field-following coordinates for vectors
182  * @param[in] grid is the Grid object needed here to map elements to nodes
183  * @param[in] input is the array one wants to convert
184  * @param[out] output is the resulting converted array
185  * */
186  void cnvt_grid_real2ff(const Grid<DeviceType>& grid, const View<double***,CLayout,DeviceType>& input, const View<double***,Kokkos::LayoutRight,DeviceType>& output) const{
187  for(int i = 0; i<input.extent(0); i++){
188  auto input_i = my_subview(input, i);
189  auto output_i = my_subview(output, i);
190 
191  cnvt_grid_real2ff(grid, input_i, output_i);
192  }
193  }
194 
195  /* In-place conversion from real coordinates to field-following coordinates for vectors, using a buffer array to store the result temporarily
196  * before copying back to the input array
197  * @param[in] grid is the Grid object needed here to map elements to nodes
198  * @param[inout] input is the array one wants to convert, and the resulting converted array
199  * */
200  void cnvt_grid_real2ff(const Grid<DeviceType>& grid, const View<double***,CLayout,DeviceType>& input)const{
201  for(int i = 0; i<input.extent(0); i++){
202  auto input_i = my_subview(input, i);
203  cnvt_grid_real2ff(grid, input_i);
204  }
205  }
206 
207 
208  /* Conversion from real coordinates to field-following coordinates from a GridField with PhiInterpType::None to a GridField with PhiInterpType::Planes.
209  * This implementation could be improved. It uses the existing algorithm which requires the input to have nphi=2, so there is a lot of superfluous copying etc. */
210  void cnvt_grid_real2ff(const Grid<DeviceType>& grid, const ScalarGridField& field_in,
212 
213  // Field-following algorithm requires argument to be nphi=2, and inode as the contiguous index
214  int nphi = 2;
215  Kokkos::View<double**,Kokkos::LayoutRight,DeviceType> dest_view(NoInit("dest_view"), nphi, field_in.nnode());
216 
217  // Since the input is a ScalarGridField (nphi=1), copy the field into both the iphi=0 and iphi=1 indices
218  field_in.copy_to_double_view(my_subview(dest_view, 0)); // iphi=0;
219  field_in.copy_to_double_view(my_subview(dest_view, 1)); // iphi=1;
220 
221  // Convert to field-following
222  cnvt_grid_real2ff(grid, dest_view, result);
223 
224  // Now transpose and copy into the output argument:
225  field_out.transpose_copy_from_double_view(result);
226  }
227 
228 };
229 
230 #endif
VolumesAndAreas volumes_and_areas
Definition: grid.hpp:176
Definition: field_following_coordinates.hpp:9
void cnvt_grid_ff2real(const Grid< DeviceType > &grid, const View< double **, CLayout, HostType > &input, const View< double **, CLayout, DeviceType > &output) const
Definition: field_following_coordinates.hpp:125
void cnvt_grid_ff2real(const Grid< DeviceType > &grid, const View< double **, CLayout, DeviceType > &input, const View< double **, Kokkos::LayoutRight, DeviceType > &output) const
Definition: field_following_coordinates.hpp:35
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
FieldFollowingCoordinates(const Projection< HostType > &projection)
Definition: field_following_coordinates.hpp:20
FieldFollowingCoordinates()
Definition: field_following_coordinates.hpp:18
View< int **, CLayout, DeviceType > tr
Definition: field_following_coordinates.hpp:10
Definition: grid_field.hpp:22
View< double ***, CLayout, Device > p
Definition: ff_projection.hpp:24
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:186
View< double **, CLayout, HostType > node_vol_ff_l_h
Definition: volumes_and_areas.hpp:12
void cnvt_grid_real2ff(const Grid< DeviceType > &grid, const View< double ***, CLayout, DeviceType > &input) const
Definition: field_following_coordinates.hpp:200
KOKKOS_INLINE_FUNCTION int get_node_index(int triangle_index, int tri_vertex_index) const
Definition: grid.tpp:158
void cnvt_grid_ff2real(const Grid< DeviceType > &grid, const View< double **, CLayout, HostType > &input, const View< double **, CLayout, HostType > &output, const View< double **, CLayout, DeviceType > &tmp) const
Definition: field_following_coordinates.hpp:107
void cnvt_grid_real2ff(const Grid< DeviceType > &grid, const View< double **, CLayout, DeviceType > &input) const
Definition: field_following_coordinates.hpp:174
void exit_XGC(std::string msg)
Definition: globals.hpp:37
View< double **, CLayout, HostType > node_vol_ff_r_h
Definition: volumes_and_areas.hpp:13
View< double ***, CLayout, DeviceType > p
Definition: field_following_coordinates.hpp:11
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
View< int **, CLayout, Device > tr
The triangles that are mapped to by the projection. One-indexed.
Definition: ff_projection.hpp:22
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
View< double **, CLayout, DeviceType > result
Definition: field_following_coordinates.hpp:14
int nnode
Number of grid nodes.
Definition: grid.hpp:159
void cnvt_grid_real2ff(const Grid< DeviceType > &grid, const ScalarGridField &field_in, const GridField< DeviceType, VarType::Scalar, PhiInterpType::Planes, TorType::OnePlane, KinType::DriftKin > &field_out) const
Definition: field_following_coordinates.hpp:210
Kokkos::ViewAllocateWithoutInitializing NoInit
Definition: space_settings.hpp:69