1 #ifndef WARPX_COARSEN_MR_H_ 2 #define WARPX_COARSEN_MR_H_ 4 #include <AMReX_Array.H> 5 #include <AMReX_Array4.H> 6 #include <AMReX_Extension.H> 7 #include <AMReX_GpuQualifiers.H> 8 #include <AMReX_REAL.H> 10 #include <AMReX_BaseFwd.H> 16 using namespace amrex;
37 Real
Interp ( Array4<Real const>
const& arr_src,
38 GpuArray<int,3>
const& sf,
39 GpuArray<int,3>
const& sc,
40 GpuArray<int,3>
const&
cr,
47 const int ic[3] = {
i, j, k };
50 int np[3], idx_min[3];
53 for (
int l = 0; l < 3; ++l ) {
54 if ( cr[l] == 1 ) np[l] = 1;
55 else np[l] = cr[l]*(1-sf[l])*(1-sc[l])
56 +(2*(cr[l]-1)+1)*sf[l]*sc[l];
60 for (
int l = 0; l < 3; ++l ) {
61 if ( cr[l] == 1 ) idx_min[l] = ic[l];
62 else idx_min[l] = ic[l]*cr[l]*(1-sf[l])*(1-sc[l])
63 +(ic[l]*cr[l]-cr[l]+1)*sf[l]*sc[l];
67 const int numx = np[0];
68 const int numy = np[1];
69 const int numz = np[2];
70 const int imin = idx_min[0];
71 const int jmin = idx_min[1];
72 const int kmin = idx_min[2];
73 const int sfx = sf[0];
74 const int sfy = sf[1];
75 const int sfz = sf[2];
76 const int scx = sc[0];
77 const int scy = sc[1];
78 const int scz = sc[2];
79 const int crx = cr[0];
80 const int cry = cr[1];
81 const int crz = cr[2];
86 auto const arr_src_safe = [arr_src]
87 AMREX_GPU_DEVICE (
int const ix,
int const iy,
int const iz,
int const n) noexcept
89 return arr_src.contains( ix, iy, iz ) ? arr_src(ix,iy,iz,n) : 0.0_rt;
101 for (
int kref = 0; kref < numz; ++kref) {
102 for (
int jref = 0; jref < numy; ++jref) {
103 for (
int iref = 0; iref < numx; ++iref) {
107 wx = (1.0_rt/
static_cast<Real
>(numx))*(1-sfx)*(1-scx)
108 +((amrex::Math::abs(crx-amrex::Math::abs(ii-i*crx)))/static_cast<Real>(crx*crx))*sfx*scx;
109 wy = (1.0_rt/
static_cast<Real
>(numy))*(1-sfy)*(1-scy)
110 +((amrex::Math::abs(cry-amrex::Math::abs(jj-j*cry)))/static_cast<Real>(cry*cry))*sfy*scy;
111 wz = (1.0_rt/
static_cast<Real
>(numz))*(1-sfz)*(1-scz)
112 +((amrex::Math::abs(crz-amrex::Math::abs(kk-k*crz)))/static_cast<Real>(crz*crz))*sfz*scz;
113 c += wx*wy*wz*arr_src_safe(ii,jj,kk,comp);
133 void Loop ( MultiFab& mf_dst,
134 const MultiFab& mf_src,
137 const IntVect crse_ratio );
149 void Coarsen ( MultiFab& mf_dst,
150 const MultiFab& mf_src,
151 const IntVect crse_ratio );
154 #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:17
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:76
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: BreitWheelerEngineWrapper.H:35
Definition: CoarsenMR.H:14
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:37