WarpX
Loading...
Searching...
No Matches
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
6namespace Interpolate {
7
9void interp (int j, int k, int l,
12 const amrex::IntVect r_ratio, IntVect const& type) noexcept
13{
14 using amrex::Real;
15 Real const rrx = 1.0_rt/static_cast<Real>(r_ratio[0]);
16
17 int const jg = amrex::coarsen(j,r_ratio[0]);
18 Real const wx = static_cast<Real>(type[0]) * static_cast<amrex::Real>(j-jg*r_ratio[0]) * rrx;
19 Real const owx = 1.0_rt-wx;
20
21#if (AMREX_SPACEDIM >= 2)
22 Real const rry = 1.0_rt/static_cast<Real>(r_ratio[1]);
23 int const kg = amrex::coarsen(k,r_ratio[1]);
24 Real const wy = static_cast<Real>(type[1]) * static_cast<amrex::Real>(k-kg*r_ratio[1]) * rry;
25 Real const owy = 1.0_rt-wy;
26#endif
27
28#if AMREX_SPACEDIM == 1
29 fine(j,k,l) = owx * crse(jg ,0,0)
30 + wx * crse(jg+1,0,0);
31#elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
32 fine(j,k,l) = owx * owy * crse(jg ,kg ,0)
33 + owx * wy * crse(jg ,kg+1,0)
34 + wx * owy * crse(jg+1,kg ,0)
35 + wx * wy * crse(jg+1,kg+1,0);
36#else
37 Real const rrz = 1.0_rt/static_cast<Real>(r_ratio[2]);
38 int const lg = amrex::coarsen(l,r_ratio[2]);
39 Real const wz = static_cast<Real>(type[2]) * static_cast<amrex::Real>(l-lg*r_ratio[2]) * rrz;
40 Real const owz = 1.0_rt-wz;
41 fine(j,k,l) = owx * owy * owz * crse(jg ,kg ,lg )
42 + wx * owy * owz * crse(jg+1,kg ,lg )
43 + owx * wy * owz * crse(jg ,kg+1,lg )
44 + wx * wy * owz * crse(jg+1,kg+1,lg )
45 + owx * owy * wz * crse(jg ,kg ,lg+1)
46 + wx * owy * wz * crse(jg+1,kg ,lg+1)
47 + owx * wy * wz * crse(jg ,kg+1,lg+1)
48 + wx * wy * wz * crse(jg+1,kg+1,lg+1);
49#endif
50}
51
52}
53
54#endif
#define AMREX_FORCE_INLINE
#define AMREX_GPU_DEVICE
Array4< Real > fine
Array4< Real const > crse
Definition Interpolate.cpp:24
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, const amrex::IntVect r_ratio, IntVect const &type) noexcept
Definition Interpolate_K.H:9
__host__ __device__ BoxND< dim > coarsen(const BoxND< dim > &b, int ref_ratio) noexcept
IntVectND< 3 > IntVect