1 #ifndef WARPX_COARSEN_MR_H_ 2 #define WARPX_COARSEN_MR_H_ 5 #include <AMReX_Math.H> 30 Real
Interp ( Array4<Real const>
const& arr_src,
31 GpuArray<int,3>
const& sf,
32 GpuArray<int,3>
const& sc,
33 GpuArray<int,3>
const&
cr,
40 const int ic[3] = {
i, j, k };
43 int np[3], idx_min[3];
46 for (
int l = 0; l < 3; ++l ) {
47 if ( cr[l] == 1 ) np[l] = 1;
48 else np[l] = cr[l]*(1-sf[l])*(1-sc[l])
49 +(2*(cr[l]-1)+1)*sf[l]*sc[l];
53 for (
int l = 0; l < 3; ++l ) {
54 if ( cr[l] == 1 ) idx_min[l] = ic[l];
55 else idx_min[l] = ic[l]*cr[l]*(1-sf[l])*(1-sc[l])
56 +(ic[l]*cr[l]-cr[l]+1)*sf[l]*sc[l];
60 const int numx = np[0];
61 const int numy = np[1];
62 const int numz = np[2];
63 const int imin = idx_min[0];
64 const int jmin = idx_min[1];
65 const int kmin = idx_min[2];
66 const int sfx = sf[0];
67 const int sfy = sf[1];
68 const int sfz = sf[2];
69 const int scx = sc[0];
70 const int scy = sc[1];
71 const int scz = sc[2];
72 const int crx = cr[0];
73 const int cry = cr[1];
74 const int crz = cr[2];
79 auto const arr_src_safe = [arr_src]
80 AMREX_GPU_DEVICE (
int const ix,
int const iy,
int const iz,
int const n) noexcept
82 return arr_src.contains( ix, iy, iz ) ? arr_src(ix,iy,iz,n) : 0.0_rt;
94 for (
int kref = 0; kref < numz; ++kref) {
95 for (
int jref = 0; jref < numy; ++jref) {
96 for (
int iref = 0; iref < numx; ++iref) {
100 wx = (1.0_rt/
static_cast<Real
>(numx))*(1-sfx)*(1-scx)
101 +((amrex::Math::abs(crx-amrex::Math::abs(ii-i*crx)))/static_cast<Real>(crx*crx))*sfx*scx;
102 wy = (1.0_rt/
static_cast<Real
>(numy))*(1-sfy)*(1-scy)
103 +((amrex::Math::abs(cry-amrex::Math::abs(jj-j*cry)))/static_cast<Real>(cry*cry))*sfy*scy;
104 wz = (1.0_rt/
static_cast<Real
>(numz))*(1-sfz)*(1-scz)
105 +((amrex::Math::abs(crz-amrex::Math::abs(kk-k*crz)))/static_cast<Real>(crz*crz))*sfz*scz;
106 c += wx*wy*wz*arr_src_safe(ii,jj,kk,comp);
126 void Loop ( MultiFab& mf_dst,
127 const MultiFab& mf_src,
130 const IntVect crse_ratio );
142 void Coarsen ( MultiFab& mf_dst,
143 const MultiFab& mf_src,
144 const IntVect crse_ratio );
147 #endif // WARPX_COARSEN_MR_H_ void Loop(MultiFab &mf_dst, const MultiFab &mf_src, const int ncomp, const IntVect ngrow, const IntVect crse_ratio)
Loops over the boxes of the coarsened MultiFab mf_dst and fills them by interpolating the data contai...
Definition: CoarsenMR.cpp:6
imin
Definition: check_interp_points_and_weights.py:121
ii
Definition: check_interp_points_and_weights.py:145
int n
Definition: run_libensemble_on_warpx.py:68
void Coarsen(MultiFab &mf_dst, const MultiFab &mf_src, const IntVect crse_ratio)
Stores in the coarsened MultiFab mf_dst the values obtained by interpolating the data contained in th...
Definition: CoarsenMR.cpp:65
i
Definition: check_interp_points_and_weights.py:171
jj
Definition: check_interp_points_and_weights.py:157
cr
TODO Coarsening for IO: interpolation points and weights def coarsening_points_and_weights_for_IO( i...
Definition: check_interp_points_and_weights.py:98
Definition: CoarsenMR.H:7
AMREX_GPU_DEVICE AMREX_FORCE_INLINE Real Interp(Array4< Real const > const &arr_src, GpuArray< int, 3 > const &sf, GpuArray< int, 3 > const &sc, GpuArray< int, 3 > const &cr, const int i, const int j, const int k, const int comp)
Interpolates the floating point data contained in the source Array4 arr_src, extracted from a fine Mu...
Definition: CoarsenMR.H:30