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  Kokkos::deep_copy(node_vol_ff, grid.volumes_and_areas.node_vol_ff_h);
49 
50  // Make a (SHALLOW) copy so the lambda doesn't have to implicitly capture *this
51  FieldFollowingCoordinates ff = *this;
52 
53  // for all node point of ff
54  Kokkos::parallel_for("cnvt_grid_ff2real", Kokkos::RangePolicy<ExSpace>(0,grid.nnode), KOKKOS_LAMBDA( const int i ){
55  for (int dir=0;dir<2;dir++){
56  for (int j=0;j<3;j++){
57  // node and weight
58  int nd=grid.get_node_index(ff.tr(dir,i), j);
59  double wt=ff.p(dir,i,j)*node_vol_ff(dir,nd);
60 
61  // accumulate values
62  Kokkos::atomic_add(&(output(dir,nd)), wt*input(dir,i));
63  Kokkos::atomic_add(&(vol(dir,nd)), wt);
64  }
65  }
66  });
67 
68  // Normalize to vol
69  Kokkos::parallel_for("cnvt_grid_ff2real_norm", Kokkos::RangePolicy<ExSpace>(0,grid.nnode), KOKKOS_LAMBDA( const int i ){
70  for (int dir=0;dir<2;dir++){
71  if(vol(dir, i) == 0.0){
72  output(dir, i) = 0.0;
73  }else{
74  output(dir, i) /= vol(dir, i);
75  }
76  }
77  });
78 
79  // fill empty area --- replace nan values
80  Kokkos::parallel_for("cnvt_grid_ff2real_fillempty", Kokkos::RangePolicy<ExSpace>(0,grid.nnode), KOKKOS_LAMBDA( const int i ){
81  for (int dir=0; dir<2; dir++){
82  if(vol(dir,i)==0.0){
83  int cdir = 1-dir;
84  output(dir,i)=0.0;
85  for (int j=0;j<3;j++){
86  int nd=grid.get_node_index(ff.tr(cdir,i), j);
87  double wt=ff.p(cdir,i,j);
88 
89  // obtain interpolated values
90  output(dir,i) += wt*input(dir,nd);
91  }
92  }
93  }
94  });
95  Kokkos::fence();
96  }
97 
98  /* Conversion from field-following coordinates to real coordinates when views are on host
99  * User also provides a temporary device allocation
100  * @param[in] grid is the Grid object needed here to map elements to nodes
101  * @param[in] input is the array one wants to convert
102  * @param[out] output is the resulting converted array
103  * @param[in] tmp is an array with the same size as input/output used for temporary storage
104  * */
105  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{
106  // Copy input to tmp
107  Kokkos::deep_copy(tmp, input);
108 
109  // Convert
110  cnvt_grid_ff2real(grid, tmp, result);
111 
112  // Copy to output
113  Kokkos::deep_copy(output, result);
114  }
115 
116  /* Conversion from real coordinates to field-following coordinates
117  * @param[in] grid is the Grid object needed here to map elements to nodes
118  * @param[in] input is the array one wants to convert
119  * @param[out] output is the resulting converted array
120  * */
121  void cnvt_grid_real2ff(const Grid<DeviceType>& grid, const View<double**,CLayout,DeviceType>& input, const View<double**,Kokkos::LayoutRight,DeviceType>& output) const{
122  // Check the input and output views are the correct size
123  if(input.extent(0) != 2 || input.extent(1) != grid.nnode ||
124  output.extent(0)!= 2 || output.extent(1)!= grid.nnode ) exit_XGC("\nError in cnvt_grid_real2ff: Unexpected View sizes. \n");
125 
126  // Initialize output to zero
127  Kokkos::deep_copy(output, 0.0);
128 
129 
130  // Make a (SHALLOW) copy so the lambda doesn't have to implicitly capture *this
131  FieldFollowingCoordinates ff = *this;
132 
133  Kokkos::parallel_for("cnvt_grid_real2ff", Kokkos::RangePolicy<ExSpace>(0,grid.nnode), KOKKOS_LAMBDA( const int i ){
134  for (int dir=0;dir<2;dir++){
135  for (int j=0;j<3;j++){
136  int nd=grid.get_node_index(ff.tr(dir,i), j);
137 
138  double wt=ff.p(dir,i,j);
139 
140  // obtain interpolated values
141  output(dir,i) += wt*input(dir,nd);
142  }
143  }
144  });
145  Kokkos::fence();
146  }
147 
148  /* In-place conversion from real coordinates to field-following coordinates, using a buffer array to store the result temporarily
149  * before copying back to the input array
150  * @param[in] grid is the Grid object needed here to map elements to nodes
151  * @param[inout] input is the array one wants to convert, and the resulting converted array
152  * */
153  void cnvt_grid_real2ff(const Grid<DeviceType>& grid, const View<double**,CLayout,DeviceType>& input) const{
154  cnvt_grid_real2ff(grid, input, result);
155 
156  // Copy back to original array
157  Kokkos::deep_copy(input, result);
158  }
159 
160  /* Conversion from real coordinates to field-following coordinates for vectors
161  * @param[in] grid is the Grid object needed here to map elements to nodes
162  * @param[in] input is the array one wants to convert
163  * @param[out] output is the resulting converted array
164  * */
165  void cnvt_grid_real2ff(const Grid<DeviceType>& grid, const View<double***,CLayout,DeviceType>& input, const View<double***,Kokkos::LayoutRight,DeviceType>& output) const{
166  for(int i = 0; i<input.extent(0); i++){
167  auto input_i = my_subview(input, i);
168  auto output_i = my_subview(output, i);
169 
170  cnvt_grid_real2ff(grid, input_i, output_i);
171  }
172  }
173 
174  /* In-place conversion from real coordinates to field-following coordinates for vectors, using a buffer array to store the result temporarily
175  * before copying back to the input array
176  * @param[in] grid is the Grid object needed here to map elements to nodes
177  * @param[inout] input is the array one wants to convert, and the resulting converted array
178  * */
179  void cnvt_grid_real2ff(const Grid<DeviceType>& grid, const View<double***,CLayout,DeviceType>& input)const{
180  for(int i = 0; i<input.extent(0); i++){
181  auto input_i = my_subview(input, i);
182  cnvt_grid_real2ff(grid, input_i);
183  }
184  }
185 
186 
187  /* Conversion from real coordinates to field-following coordinates from a GridField with PhiInterpType::None to a GridField with PhiInterpType::Planes.
188  * 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. */
189  void cnvt_grid_real2ff(const Grid<DeviceType>& grid, const ScalarGridField& field_in,
191 
192  // Field-following algorithm requires argument to be nphi=2, and inode as the contiguous index
193  int nphi = 2;
194  Kokkos::View<double**,Kokkos::LayoutRight,DeviceType> dest_view(NoInit("dest_view"), nphi, field_in.nnode());
195 
196  // Since the input is a ScalarGridField (nphi=1), copy the field into both the iphi=0 and iphi=1 indices
197  field_in.copy_to_double_view(my_subview(dest_view, 0)); // iphi=0;
198  field_in.copy_to_double_view(my_subview(dest_view, 1)); // iphi=1;
199 
200  // Convert to field-following
201  cnvt_grid_real2ff(grid, dest_view, result);
202 
203  // Now transpose and copy into the output argument:
204  field_out.transpose_copy_from_double_view(result);
205  }
206 
207 };
208 
209 #endif
VolumesAndAreas volumes_and_areas
Definition: grid.hpp:313
Definition: field_following_coordinates.hpp:9
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:121
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:25
View< double ***, CLayout, Device > p
Definition: grid.hpp:55
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:165
void cnvt_grid_real2ff(const Grid< DeviceType > &grid, const View< double ***, CLayout, DeviceType > &input) const
Definition: field_following_coordinates.hpp:179
KOKKOS_INLINE_FUNCTION int get_node_index(int triangle_index, int tri_vertex_index) const
Definition: grid.tpp:891
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:105
void cnvt_grid_real2ff(const Grid< DeviceType > &grid, const View< double **, CLayout, DeviceType > &input) const
Definition: field_following_coordinates.hpp:153
void exit_XGC(std::string msg)
Definition: globals.hpp:37
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
Definition: grid.hpp:54
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:271
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:189
Kokkos::ViewAllocateWithoutInitializing NoInit
Definition: space_settings.hpp:68
View< double **, CLayout, HostType > node_vol_ff_h
Definition: grid.hpp:23