XGCa
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
col_vgrids.hpp
Go to the documentation of this file.
1 #ifndef COL_VGRIDS_HPP
2 #define COL_VGRIDS_HPP
3 
4 #include <vector>
5 #include <stdio.h>
6 #include <Kokkos_Core.hpp>
7 #include <Kokkos_DualView.hpp>
8 
9 #include "col_species.hpp"
10 
11 
12 template<class Device>
14  Kokkos::DualView<int*,Kokkos::LayoutRight,Device> nspecies_in_grid;
15  const int n;
16  const int mb_n_nodes;
17  const int n_species;
18  const int nvr;
19  const int nvz;
20  Kokkos::View<double***,Kokkos::LayoutRight,HostType> mesh_r;
21  Kokkos::View<double***,Kokkos::LayoutRight,HostType> mesh_z;
22  Kokkos::View<double***,Kokkos::LayoutRight,Device> vol;
23  Kokkos::View<double***,Kokkos::LayoutRight,HostType> vol_h;
24  Kokkos::View<double****,Device> delta_r;
25  Kokkos::View<double****,Device> delta_z;
26 
27  // These three are sent to GPU, so use a dualview
28  Kokkos::DualView<double***,Device> mesh_r_half;
29  Kokkos::DualView<double***,Device> mesh_z_half;
30  Kokkos::DualView<double***,Device> local_center_volume;
31  Kokkos::DualView<double**,Device> vpar_beg;
32  Kokkos::DualView<double**,Device> mesh_dz;
33  Kokkos::DualView<double**,Device> mesh_dr;
34 
35  Kokkos::DualView<int**,Kokkos::LayoutRight,Device> map_grid_to_species;
36 
37  // Count how many species are using each grid
38  inline Kokkos::DualView<int*,Kokkos::LayoutRight,Device> count_species_in_grid(const CollisionSpecies<Device>& col_spall) const{
39  Kokkos::DualView<int*,Kokkos::LayoutRight,Device> nspecies_in_grid_tmp("nspecies_in_grid", col_spall.s.h_view.extent(0));
40  // Loop over species and set up mapping
41  for (int isp = 0; isp < col_spall.s.h_view.extent(0); isp++){
42  int gri = col_spall.s.h_view(isp,0).grid_index; // Which grid this species is on
43  nspecies_in_grid_tmp.h_view(gri)++;
44  }
45  return nspecies_in_grid_tmp;
46  }
47 
48  // Count how many grids are used
49  inline int count_grids() const{
50  int ngrids=0;
51  for(int i=0; i<nspecies_in_grid.h_view.size(); i++){
52  if(nspecies_in_grid.h_view(i)>0) ngrids++;
53  }
54  return ngrids;
55  }
56 
58 
61  n(count_grids()),
62  mb_n_nodes(col_spall.pdf_n.h_view.extent(0)), // Use pdf_n dimensions: mb_n_nodes, n_species, nvz, nvr
63  n_species(col_spall.pdf_n.h_view.extent(1)),
64  nvz(col_spall.pdf_n.h_view.extent(2)),
65  nvr(col_spall.pdf_n.h_view.extent(3)),
66  map_grid_to_species("map_grid_to_species",n,n_species), // n_species is the max possible
67  mesh_dr("mesh_dr", mb_n_nodes, n),
68  mesh_dz("mesh_dz", mb_n_nodes, n),
69  vpar_beg("vpar_beg", mb_n_nodes, n),
70  mesh_r("mesh_r", mb_n_nodes, n, nvr),
71  mesh_z("mesh_z", mb_n_nodes, n, nvz),
72  mesh_r_half("mesh_r_half",mb_n_nodes,nvr-1,n),
73  mesh_z_half("mesh_z_half",mb_n_nodes,nvz-1,n),
74  local_center_volume("local_center_volume",mb_n_nodes,nvr-1,n),
75  vol("vol",nvr,n, mb_n_nodes),
76  vol_h("vol_h",mb_n_nodes,n, nvr), // Must have nvr fastest index for passing into fortran negative_f_correction routine
77  delta_r("delta_r",n,nvr-1,n_species, mb_n_nodes),
78  delta_z("delta_z",n,nvz-1,n_species, mb_n_nodes)
79  {
80  // Loop over species and set up mapping
81  int isp_in_grid=0;
82  int gri=0;
83  for (int isp = 0; isp < n_species; isp++){
84  // Move to next species or to next grid
85  if(gri!=col_spall.s.h_view(isp,0).grid_index) isp_in_grid=0;
86 
87  gri = col_spall.s.h_view(isp,0).grid_index; // Which grid this species is on
88  map_grid_to_species.h_view(gri,isp_in_grid) = isp;
89  isp_in_grid++;
90  }
91 
92 #ifdef USE_GPU
93  // Copy mapping views to device
94  Kokkos::deep_copy(map_grid_to_species.d_view, map_grid_to_species.h_view);
95  Kokkos::deep_copy(nspecies_in_grid.d_view, nspecies_in_grid.h_view);
96 #endif
97  }
98 
99  //> Convert the distribution function to the units used
100  // internally by the collision operator
101  // The distribution function set up with f0_module is in normalized velocity space
102  // whereas the collision operator uses SI units.
103  void setup_one(int isp, int mesh_ind, double vth, double dsmu, double dvp, double vp_max) {
104 
105  // vol
106  double pi4bb0vth3_dd = TWOPI*vth*vth*vth*dsmu*dvp;
107 
108  // --- 0th component
109  double smu_n = dsmu/3.0;
110  vol_h(mesh_ind,isp,0) = 0.5*pi4bb0vth3_dd*smu_n;
111  // --- 1st to (last-1)
112  for (int imu = 1; imu < nvr-1; imu++){
113  smu_n = dsmu*imu;
114  vol_h(mesh_ind,isp,imu) = pi4bb0vth3_dd*smu_n;
115  }
116  // --- last component
117  smu_n = dsmu*((nvr-1)-1.0/3.0);//rh imu or f0_nmu???
118  vol_h(mesh_ind,isp,nvr-1) = 0.5*pi4bb0vth3_dd*smu_n;
119 
120  //vth_par = sqrt(sp_t_par(i) * UNIT_CHARGE / mass)
121  //vth_perp = sqrt(sp_t_perp(i) * UNIT_CHARGE / mass)
122 
123  // < Mesh-related quantities>
124  double vpar_beg_l = -vp_max*vth;// ( = dlx)
125  double mesh_dz_l = dvp*vth;// ( = mesh_dz)
126  double mesh_dr_l = dsmu*vth;// ( = mesh_dr)
127 
128  vpar_beg.h_view(mesh_ind, isp) = vpar_beg_l;
129  mesh_dz.h_view(mesh_ind, isp) = mesh_dz_l;
130  mesh_dr.h_view(mesh_ind, isp) = mesh_dr_l;
131 
132  // mesh_r, mesh_r_half, mesh_z, mesh_z_half
133  for (int j = 0; j < nvr-1; j++){
134  mesh_r(mesh_ind,isp,j) = mesh_dr_l*j;
135  mesh_r_half.h_view(mesh_ind,j,isp) = mesh_dr_l*(j+0.5);
136  }
137  mesh_r(mesh_ind,isp,nvr-1) = mesh_dr_l*(nvr-1);
138  for (int j = 0; j < nvz-1; j++){
139  mesh_z(mesh_ind,isp,j) = vpar_beg_l + mesh_dz_l*j;
140  mesh_z_half.h_view(mesh_ind,j,isp) = vpar_beg_l + mesh_dz_l*(j+0.5);
141  }
142  mesh_z(mesh_ind,isp,nvz-1) = vpar_beg_l + mesh_dz_l*(nvz-1);
143  // local_center_volume (mesh_r_half, mesh_dr, mesh_dz)
144  // : volume centered at a cell
145  for (int j = 0; j < (nvr-1); j++){
146  local_center_volume.h_view(mesh_ind,j,isp) = mesh_r_half.h_view(mesh_ind,j,isp)*mesh_dr_l*mesh_dz_l;
147  }
148 
149  }
150 
151  void setup_all(const std::vector<Species<Device>>& all_species, const CollisionSpecies<Device>& col_species,
152  const Kokkos::View<int**,HostType>& mesh_nodes, int mesh_batch_ind){
153 
154  for (int mesh_ind = 0; mesh_ind<col_species.s.h_view.extent(1); mesh_ind++){
155  for (int gri = 0; gri < n; gri++){
156  const int FIRST_SPECIES_OF_GRID = 0; // Initialize grid using the first species assigned to the grid
157  int isp = map_grid_to_species.h_view(gri,FIRST_SPECIES_OF_GRID);
158 
159  // Use equilibrium thermal velocity to construct velocity grid
160  int inode = mesh_nodes(mesh_batch_ind,mesh_ind);
161  double vth = all_species[col_species.s.h_view(isp,mesh_ind).gid].get_f0_eq_thermal_velocity_lnode_h(inode);
162  setup_one(isp, mesh_ind, vth, col_species.f.dsmu, col_species.f.dvp, col_species.f.vp_max);
163  }
164  }
165 
166 #ifdef USE_GPU
167  // Copy dual views to device
168  Kokkos::deep_copy(mesh_dr.d_view, mesh_dr.h_view);
169  Kokkos::deep_copy(mesh_dz.d_view, mesh_dz.h_view);
170  Kokkos::deep_copy(vpar_beg.d_view, vpar_beg.h_view);
171  Kokkos::deep_copy(mesh_r_half.d_view, mesh_r_half.h_view);
172  Kokkos::deep_copy(mesh_z_half.d_view, mesh_z_half.h_view);
173  Kokkos::deep_copy(local_center_volume.d_view, local_center_volume.h_view);
174 #endif
175  // Transpose vol
176  // First create a mirror on device of vol_h
177  Kokkos::View<double***,Kokkos::LayoutRight,Device> vol_tmp("vol_tmp",vol_h.extent(0),vol_h.extent(1),vol_h.extent(2));
178  // Deep copy that
179  Kokkos::deep_copy(vol_tmp, vol_h);
180 
181  // Make a local copy of the View itself for the lambda function to copy to device
182  Kokkos::View<double***,Kokkos::LayoutRight,Device> lvol = vol;
183 
184  // Naive transpose copy on device - unoptimized, should be small though
185  int n2 = vol_h.extent(1)*nvr;
186  int nvr_loc = nvr;
187  Kokkos::parallel_for("transpose_vol", Kokkos::RangePolicy<ExSpace>( 0, n2*vol_h.extent(0)), KOKKOS_LAMBDA (const int idx){
188  int i0 = idx/n2; // slowest index
189  int rem = idx%n2;
190  int i1 = rem/nvr_loc;
191  int i2 = rem%nvr_loc; // fastest index
192  lvol(i2,i1,i0) = vol_tmp(i0,i1,i2);
193  });
194  }
195 };
196 
197 #endif
Kokkos::View< double ****, Device > delta_r
Definition: col_vgrids.hpp:24
Kokkos::View< double ***, Kokkos::LayoutRight, Device > vol
Definition: col_vgrids.hpp:22
const int n
Number of velocity grids.
Definition: col_vgrids.hpp:15
Kokkos::DualView< double **, Device > mesh_dr
Definition: col_vgrids.hpp:33
Kokkos::DualView< double **, Device > mesh_dz
Definition: col_vgrids.hpp:32
Kokkos::View< double ***, Kokkos::LayoutRight, HostType > mesh_z
Definition: col_vgrids.hpp:21
Kokkos::View< double ****, Device > delta_z
Definition: col_vgrids.hpp:25
void setup_all(const std::vector< Species< Device >> &all_species, const CollisionSpecies< Device > &col_species, const Kokkos::View< int **, HostType > &mesh_nodes, int mesh_batch_ind)
Definition: col_vgrids.hpp:151
Kokkos::DualView< double ***, Device > mesh_z_half
Definition: col_vgrids.hpp:29
Kokkos::View< double ***, Kokkos::LayoutRight, HostType > vol_h
Definition: col_vgrids.hpp:23
double dvp
grid spacing in parallel velocity
Definition: vgrid_distribution.hpp:25
const int nvz
grid points in z
Definition: col_vgrids.hpp:19
Kokkos::DualView< int **, Kokkos::LayoutRight, Device > map_grid_to_species
Definition: col_vgrids.hpp:35
Kokkos::DualView< CollisionSpeciesScalars **, Device > s
Definition: col_species.hpp:85
Definition: col_vgrids.hpp:13
const int nvr
grid points in r
Definition: col_vgrids.hpp:18
Kokkos::DualView< int *, Kokkos::LayoutRight, Device > count_species_in_grid(const CollisionSpecies< Device > &col_spall) const
Definition: col_vgrids.hpp:38
idx
Definition: diag_f0_df_port1.hpp:32
const int mb_n_nodes
Batch size.
Definition: col_vgrids.hpp:16
void setup_one(int isp, int mesh_ind, double vth, double dsmu, double dvp, double vp_max)
Definition: col_vgrids.hpp:103
Kokkos::DualView< int *, Kokkos::LayoutRight, Device > nspecies_in_grid
Definition: col_vgrids.hpp:14
Kokkos::DualView< double ***, Device > mesh_r_half
Definition: col_vgrids.hpp:28
Kokkos::DualView< double **, Device > vpar_beg
Definition: col_vgrids.hpp:31
const int n_species
Number of species total.
Definition: col_vgrids.hpp:17
Kokkos::DualView< double ***, Device > local_center_volume
Definition: col_vgrids.hpp:30
double dsmu
grid spacing in mu
Definition: vgrid_distribution.hpp:27
const VGridDistribution< HostType > f
Definition: col_species.hpp:83
Definition: species.hpp:75
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
int count_grids() const
Definition: col_vgrids.hpp:49
Definition: col_species.hpp:81
constexpr double TWOPI
Definition: constants.hpp:9
Kokkos::View< double ***, Kokkos::LayoutRight, HostType > mesh_r
Definition: col_vgrids.hpp:20
double vp_max
max parallel velocity
Definition: vgrid_distribution.hpp:24