XGC1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
varying_background.hpp
Go to the documentation of this file.
1 #ifndef VARYING_BACKGROUND_HPP
2 #define VARYING_BACKGROUND_HPP
3 
4 #include "magnetic_field.hpp"
5 #include "plasma.hpp"
6 #include "get_monte_num.hpp"
7 #include "get_volume.hpp"
8 
9 template<class Device>
11 
12  enum{
19  };
20 
21  enum{
25  N
26  };
27 
28  double pin;
29  double pout;
30  int m;
31  int mtheta;
32  double inv_dp;
33  double inv_dtheta;
34 
35  View<double**, CLayout, Device> vol;
36  View<double****, CLayout, Device> bg;
37 
38  public:
39 
40  int period;
41 
43 
45  nlr.use_namelist("col_param");
46  period = nlr.get<int>("col_vb_period", 1);
47  pin = nlr.get<double>("col_vb_pin", magnetic_field.inpsi);
48  pout = nlr.get<double>("col_vb_pout", std::min(magnetic_field.outpsi, 1.03));
49  m = nlr.get<int>("col_vb_m", 50);
50  mtheta = nlr.get<int>("col_vb_mtheta", 8);
51 
52  // Normalize
53  if(nlr.present("col_vb_pin")) pin *= magnetic_field.equil.xpt_psi;
54  if(nlr.present("col_vb_pout")) pout *= magnetic_field.equil.xpt_psi;
55 
56  double dp = (pout - pin)/m;
57  double dtheta = TWOPI/mtheta;
58  inv_dp = 1.0/dp;
59  inv_dtheta = 1.0/dtheta;
60 
61  // Allocate views
62  bg = View<double****, CLayout, Device>(NoInit("vb_bg"), n_nonadiabatic_species, m+1, mtheta, N);
63  vol = View<double**, CLayout, Device>("vb_vol", m+1, mtheta);
64 
65  int n_monte_carlo = get_monte_num(nlr, grid.nnode, grid.nplanes); // Note: n_ptl is needed for sml_monte_num
66  bool use_nonrandom_sampling = false;
67  GPTLstart("COL_VB_VOL");
68  monte_carlo_col_vb_vol(grid, magnetic_field, n_monte_carlo, use_nonrandom_sampling,
69  pin, pout, m, mtheta, inv_dp, inv_dtheta, vol);
70  GPTLstart("COL_VB_VOL");
71  }
72 
73  KOKKOS_INLINE_FUNCTION bool is_in_range(const MagneticField<DeviceType>& magnetic_field, double r, double z, double psi) const {
74  return pin<psi && psi<pout && z>magnetic_field.equil.xpt_z && magnetic_field.equil.is_in_region_1_or_2(r,z,psi);
75  }
76 
77  KOKKOS_INLINE_FUNCTION void interp_bg_profile(double psi, double theta, int isp, double& den, double& temp, double& up) const {
78  double tmp = (psi-pin)*inv_dp;
79  int j = int(tmp);
80  j = min(m-1,max(0,j));
81  int l = ((int)(floor(theta*inv_dtheta + 0.5)))%mtheta;
82  double bb = tmp - j;
83  double aa = 1.0 - bb;
84 
85  den = bg(isp,j,l,Density)*aa + bg(isp,j+1,l,Density)*bb;
86  temp = bg(isp,j,l,Temp)*aa + bg(isp,j+1,l,Temp)*bb;
87 
88  j = floor((psi-pin)*inv_dp + 0.5);
89  up = bg(isp,j,l,Flow);
90  }
91 
92  // get ion/electron density and temperature for collision routine
94  View<double***, CLayout, DeviceType> dum(NoInit("dum"),m+1,mtheta,N);
95  Kokkos::deep_copy(bg, 1.0e-99); // to avoid divide by zero
96  Kokkos::deep_copy(dum, 1.0e-99); // to avoid divide by zero
97 
98  auto vb = *this; // Use a local shallow copy so that the lambda function knows to copy it
99 
100  species.for_all_particles("col_snapshot", KOKKOS_LAMBDA( const int idx ){
101  // Get index info
102  AoSoAIndices<Device> inds(idx);
103 
104  // Copy from AoSoA to SoA
105  SimdParticles part;
106  AoSoA_2_simd(part, species.ptl(), inds.a, inds.s);
107 
108  Simd<double> psi;
109  magnetic_field.get_psi(part.ph.x(),psi);
110 
111  Simd<double> Bmag;
112  magnetic_field.bmag_interpol(part.ph.v(), Bmag);
113 
114  Simd<double> theta;
115  magnetic_field.equil.get_theta(part.ph.x(), theta);
116 
117  FORCE_SIMD
118  for (int i_simd = 0; i_simd<SIMD_SIZE; i_simd++){
119  if(part.is_active(i_simd)){
120  if(vb.pin<psi[i_simd] && psi[i_simd]<vb.pout && magnetic_field.equil.is_in_region_1_or_2(part.ph.r[i_simd],part.ph.z[i_simd],psi[i_simd])){
121  //energy calculation
122  double weight = part.ct.w0[i_simd];
123  double upar = part.ph.rho[i_simd]*Bmag[i_simd]*species.c_m;
124  double energy = 0.5*species.mass*upar*upar + part.ct.mu[i_simd]*Bmag[i_simd];
125 
126  double tmp = (psi[i_simd] - vb.pin)*vb.inv_dp;
127  int j = (int)(tmp);
128  j = min(vb.m-1,max(0,j));
129  int l = ((int)(floor(theta[i_simd]*vb.inv_dtheta + 0.5)))%vb.mtheta;
130  double bb = tmp-j;
131  double aa = 1.0-bb;
132 
133  Kokkos::atomic_add(&(dum(j,l,Density1)), weight*aa);
134  Kokkos::atomic_add(&(dum(j,l,Energy)), weight*energy*aa);
135  Kokkos::atomic_add(&(dum(j,l,VPara)), weight*upar*aa);
136  Kokkos::atomic_add(&(dum(j+1,l,Density1)), weight*bb);
137  Kokkos::atomic_add(&(dum(j+1,l,Energy)), weight*energy*bb);
138  Kokkos::atomic_add(&(dum(j+1,l,VPara)), weight*upar*bb);
139 
140 
141  j = floor((psi[i_simd]-vb.pin)*vb.inv_dp + 0.5);
142  Kokkos::atomic_add(&(dum(j,l,VPara2)), weight);
143  Kokkos::atomic_add(&(dum(j,l,Density2)), weight*upar);
144  }
145  }
146  }
147  });
148 
149 #ifdef USE_MPI
150  TIMER("COL_SNAP_RED",
151  MPI_Allreduce(MPI_IN_PLACE, dum.data(), dum.size(), MPI_DOUBLE, MPI_SUM, SML_COMM_WORLD) );
152 #endif
153 
154  int isp = species.nonadiabatic_idx;
155  int n = (m+1)*mtheta;
156  Kokkos::parallel_for("snapshot_convert", Kokkos::RangePolicy<ExSpace>(0, n), KOKKOS_LAMBDA( const int idx){
157  int i = idx%vb.mtheta; // Faster index
158  int j = idx/vb.mtheta;
159  constexpr double prefac = 2.0/3.0*J_2_EV; // Convert to eV
160  const double vpara_normed = dum(j,i,VPara)/dum(j,i,Density1);
161  const double epara = 0.5*species.mass*(vpara_normed*vpara_normed);
162  vb.bg(isp,j,i,Temp) = prefac*(dum(j,i,Energy) - epara)/dum(j,i,Density1);
163  vb.bg(isp,j,i,Density) = dum(j,i,Density1)/vb.vol(j,i);
164  vb.bg(isp,j,i,Flow) = dum(j,i,VPara2)/dum(j,i,Density2); // mean flow for moving frame collision
165  });
166  }
167 };
168 
169 #endif
double inpsi
Boundary condition used in a few spots.
Definition: magnetic_field.hpp:25
static int GPTLstart(const char *name)
Definition: timer_macro.hpp:9
MPI_Comm SML_COMM_WORLD
Definition: my_mpi.cpp:4
T get(const string &param, const T default_val, int val_ind=0)
Definition: NamelistReader.hpp:373
double c_m
c/m
Definition: species.hpp:86
VaryingBackground()
Definition: varying_background.hpp:42
View< double ****, CLayout, Device > bg
Density, Temperature and velocity of each shell for background update.
Definition: varying_background.hpp:36
void monte_carlo_col_vb_vol(const Grid< DeviceType > &grid, const MagneticField< DeviceType > &magnetic_field, long long int n_monte_carlo, bool use_nonrandom_sampling, double psi_min, double psi_max, int m, int mtheta, double inv_dp, double inv_dtheta, View< double **, CLayout, DeviceType > &vol)
Definition: get_volume.cpp:384
bool present(const string &param)
Definition: NamelistReader.hpp:363
KOKKOS_INLINE_FUNCTION VecParticles * ptl() const
Definition: species.hpp:586
int m
Definition: varying_background.hpp:30
View< double **, CLayout, Device > vol
Volume.
Definition: varying_background.hpp:35
Definition: NamelistReader.hpp:193
Definition: magnetic_field.hpp:12
KOKKOS_INLINE_FUNCTION bool is_in_range(const MagneticField< DeviceType > &magnetic_field, double r, double z, double psi) const
Definition: varying_background.hpp:73
int a
The index in the inner array of the AoSoA.
Definition: particles.hpp:150
Equilibrium equil
The object containing information about the magnetic equilibrium.
Definition: magnetic_field.hpp:32
KOKKOS_INLINE_FUNCTION void bmag_interpol(const SimdVector &v, Simd< double > &bmag) const
Definition: magnetic_field.tpp:239
int nonadiabatic_idx
Index of species skipping adiabatic species (for compatibility with fortran arrays) ...
Definition: species.hpp:81
int nplanes
Number of planes.
Definition: grid.hpp:191
int mtheta
Definition: varying_background.hpp:31
KOKKOS_INLINE_FUNCTION bool is_active(int i_simd) const
Definition: particles.hpp:66
Definition: varying_background.hpp:16
KOKKOS_INLINE_FUNCTION void AoSoA_2_simd(SimdParticles &part_one, const VecParticles *part, int a_vec, int s_vec)
Definition: particles.tpp:80
Simd< double > rho
Definition: particles.hpp:21
KOKKOS_INLINE_FUNCTION void get_psi(const SimdVector2D &x, Simd< double > &psi_out) const
Definition: magnetic_field.tpp:57
constexpr double J_2_EV
Conversion rate J to ev.
Definition: constants.hpp:6
double mass
Particle mass.
Definition: species.hpp:83
void for_all_particles(const std::string label, F lambda_func) const
Definition: species.hpp:418
#define TIMER(N, F)
Definition: timer_macro.hpp:24
double outpsi
Boundary condition used in a few spots.
Definition: magnetic_field.hpp:26
void use_namelist(const string &namelist)
Definition: NamelistReader.hpp:355
Definition: varying_background.hpp:14
Simd< double > r
Definition: particles.hpp:18
KOKKOS_INLINE_FUNCTION SimdVector2D & x()
Definition: particles.hpp:27
KOKKOS_INLINE_FUNCTION SimdVector & v()
Definition: particles.hpp:39
KOKKOS_INLINE_FUNCTION void interp_bg_profile(double psi, double theta, int isp, double &den, double &temp, double &up) const
Definition: varying_background.hpp:77
Definition: varying_background.hpp:18
SimdPhase ph
Definition: particles.hpp:59
double pin
Definition: varying_background.hpp:28
void update(const MagneticField< DeviceType > &magnetic_field, Species< DeviceType > &species)
Definition: varying_background.hpp:93
Definition: particles.hpp:58
int period
Period of background updating.
Definition: varying_background.hpp:40
Definition: varying_background.hpp:24
int s
The index in the outer array of the AoSoA.
Definition: particles.hpp:149
#define FORCE_SIMD
Definition: simd.hpp:9
Definition: varying_background.hpp:23
Simd< double > z
Definition: particles.hpp:19
KOKKOS_INLINE_FUNCTION void get_theta(const SimdVector2D &x, Simd< double > &theta) const
Definition: equil.tpp:129
Definition: varying_background.hpp:13
Simd< double > w0
Definition: particles.hpp:53
Definition: varying_background.hpp:15
Definition: varying_background.hpp:25
Definition: magnetic_field.F90:1
VaryingBackground(NLReader::NamelistReader &nlr, const MagneticField< DeviceType > &magnetic_field, const Grid< DeviceType > &grid, int n_nonadiabatic_species)
Definition: varying_background.hpp:44
SimdConstants ct
Definition: particles.hpp:60
double xpt_z
z coordinate of 1st X-point
Definition: equil.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
Definition: varying_background.hpp:17
int nnode
Number of grid nodes.
Definition: grid.hpp:190
double pout
Definition: varying_background.hpp:29
Simd< double > mu
Definition: particles.hpp:52
Definition: varying_background.hpp:10
double inv_dtheta
Definition: varying_background.hpp:33
Kokkos::ViewAllocateWithoutInitializing NoInit
Definition: space_settings.hpp:68
constexpr double TWOPI
Definition: constants.hpp:9
Definition: varying_background.hpp:22
Definition: particles.hpp:148
KOKKOS_INLINE_FUNCTION bool is_in_region_1_or_2(double r, double z, double psi) const
Definition: equil.tpp:82
double inv_dp
Definition: varying_background.hpp:32
double xpt_psi
Psi coordinate of 1st X-point.
Definition: equil.hpp:79
int get_monte_num(NLReader::NamelistReader &nlr, int nnode, int nplanes)
Definition: get_monte_num.cpp:6