1 #ifndef FIELD_FOLLOWING_COORDINATES_HPP
2 #define FIELD_FOLLOWING_COORDINATES_HPP
10 View<int**,CLayout,DeviceType>
tr;
11 View<double***,CLayout,DeviceType>
p;
14 View<double**,CLayout,DeviceType>
result;
22 p(
NoInit(
"p"),projection.
p.layout()),
26 Kokkos::deep_copy(
tr, projection.
tr);
27 Kokkos::deep_copy(
p, projection.
p);
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");
41 Kokkos::deep_copy(output, 0.0);
44 View<double**,CLayout,DeviceType> vol(
"vol",input.layout());
47 View<double**,CLayout,DeviceType> node_vol_ff(
NoInit(
"node_vol_ff"),input.layout());
57 for (
int dir=0;dir<2;dir++){
58 for (
int j=0;j<3;j++){
61 double wt=ff.
p(dir,i,j)*node_vol_ff(dir,nd);
64 Kokkos::atomic_add(&(output(dir,nd)), wt*input(dir,i));
65 Kokkos::atomic_add(&(vol(dir,nd)), wt);
72 for (
int dir=0;dir<2;dir++){
73 if(vol(dir, i) == 0.0){
76 output(dir, i) /= vol(dir, i);
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++){
87 for (
int j=0;j<3;j++){
89 double wt=ff.
p(cdir,i,j);
92 output(dir,i) += wt*input(dir,nd);
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{
109 Kokkos::deep_copy(tmp, input);
115 Kokkos::deep_copy(output,
result);
127 Kokkos::deep_copy(output, input);
133 Kokkos::deep_copy(output,
result);
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");
148 Kokkos::deep_copy(output, 0.0);
155 for (
int dir=0;dir<2;dir++){
156 for (
int j=0;j<3;j++){
159 double wt=ff.
p(dir,i,j);
162 output(dir,i) += wt*input(dir,nd);
178 Kokkos::deep_copy(input,
result);
187 for(
int i = 0; i<input.extent(0); i++){
201 for(
int i = 0; i<input.extent(0); i++){
215 Kokkos::View<double**,Kokkos::LayoutRight,DeviceType> dest_view(
NoInit(
"dest_view"), nphi, field_in.nnode());
218 field_in.copy_to_double_view(
my_subview(dest_view, 0));
219 field_in.copy_to_double_view(
my_subview(dest_view, 1));
225 field_out.transpose_copy_from_double_view(
result);
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:18
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:19
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:68