WarpX
Interpolate_K.H
Go to the documentation of this file.
1 #ifndef WARPX_INTERP_K_H_
2 #define WARPX_INTERP_K_H_
3 
4 #include <AMReX_FArrayBox.H>
5 
6 namespace Interpolate {
7 
9 void interp (int j, int k, int l,
10  amrex::Array4<amrex::Real > const& fine,
12  int r_ratio, IntVect const& type) noexcept
13 {
14  using amrex::Real;
15  Real const rr = 1.0_rt/static_cast<Real>(r_ratio);
16 
17  int const jg = amrex::coarsen(j,r_ratio);
18  Real const wx = static_cast<Real>(type[0]) * static_cast<amrex::Real>(j-jg*r_ratio) * rr;
19  Real const owx = 1.0_rt-wx;
20 
21 #if (AMREX_SPACEDIM >= 2)
22  int const kg = amrex::coarsen(k,r_ratio);
23  Real const wy = static_cast<Real>(type[1]) * static_cast<amrex::Real>(k-kg*r_ratio) * rr;
24  Real const owy = 1.0_rt-wy;
25 #endif
26 
27 #if defined(WARPX_DIM_1D_Z)
28  fine(j,k,l) = owx * crse(jg ,0,0)
29  + wx * crse(jg+1,0,0);
30 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
31  fine(j,k,l) = owx * owy * crse(jg ,kg ,0)
32  + owx * wy * crse(jg ,kg+1,0)
33  + wx * owy * crse(jg+1,kg ,0)
34  + wx * wy * crse(jg+1,kg+1,0);
35 #else
36  int const lg = amrex::coarsen(l,r_ratio);
37  Real const wz = static_cast<Real>(type[2]) * static_cast<amrex::Real>(l-lg*r_ratio) * rr;
38  Real const owz = 1.0_rt-wz;
39  fine(j,k,l) = owx * owy * owz * crse(jg ,kg ,lg )
40  + wx * owy * owz * crse(jg+1,kg ,lg )
41  + owx * wy * owz * crse(jg ,kg+1,lg )
42  + wx * wy * owz * crse(jg+1,kg+1,lg )
43  + owx * owy * wz * crse(jg ,kg ,lg+1)
44  + wx * owy * wz * crse(jg+1,kg ,lg+1)
45  + owx * wy * wz * crse(jg ,kg+1,lg+1)
46  + wx * wy * wz * crse(jg+1,kg+1,lg+1);
47 #endif
48 }
49 
50 }
51 
52 #endif
Definition: Interpolate.cpp:21
#define AMREX_FORCE_INLINE
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Box coarsen(const Box &b, int ref_ratio) noexcept
#define AMREX_GPU_DEVICE
Array4< Real > fine
type
Definition: run_alltests_1node.py:72
Array4< Real const > crse
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void interp(int j, int k, int l, amrex::Array4< amrex::Real > const &fine, amrex::Array4< amrex::Real const > const &crse, int r_ratio, IntVect const &type) noexcept
Definition: Interpolate_K.H:9