1 #ifndef WARPX_COARSEN_IO_H_ 2 #define WARPX_COARSEN_IO_H_ 4 #include <AMReX_Array.H> 5 #include <AMReX_Array4.H> 6 #include <AMReX_Extension.H> 7 #include <AMReX_GpuQualifiers.H> 8 #include <AMReX_IntVect.H> 9 #include <AMReX_REAL.H> 11 #include <AMReX_BaseFwd.H> 17 using namespace amrex;
38 Real
Interp ( Array4<Real const>
const& arr_src,
39 GpuArray<int,3>
const& sf,
40 GpuArray<int,3>
const& sc,
41 GpuArray<int,3>
const&
cr,
48 const int ic[3] = {
i, j, k };
51 int np[3], idx_min[3];
54 for (
int l = 0; l < 3; ++l ) {
55 if ( cr[l] == 1 ) np[l] = 1+amrex::Math::abs(sf[l]-sc[l]);
60 for (
int l = 0; l < 3; ++l ) {
61 if ( cr[l] == 1 ) idx_min[l] = ic[l]-sc[l]*(1-sf[l]);
62 else idx_min[l] = ic[l]*cr[l]+
static_cast<int>(cr[l]/2)*(1-sc[l])-(1-sf[l]);
66 const int numx = np[0];
67 const int numy = np[1];
68 const int numz = np[2];
69 const int imin = idx_min[0];
70 const int jmin = idx_min[1];
71 const int kmin = idx_min[2];
77 for (
int kref = 0; kref < numz; ++kref) {
78 for (
int jref = 0; jref < numy; ++jref) {
79 for (
int iref = 0; iref < numx; ++iref) {
83 wx = 1.0_rt/
static_cast<Real
>(numx);
84 wy = 1.0_rt/
static_cast<Real
>(numy);
85 wz = 1.0_rt/
static_cast<Real
>(numz);
86 c += wx*wy*wz*arr_src(ii,jj,kk,comp);
112 void Loop ( MultiFab& mf_dst,
113 const MultiFab& mf_src,
118 const IntVect crse_ratio=IntVect(1) );
139 void Coarsen ( MultiFab& mf_dst,
140 const MultiFab& mf_src,
145 const IntVect crse_ratio=IntVect(1) );
146 void Coarsen ( MultiFab& mf_dst,
147 const MultiFab& mf_src,
151 const IntVect ngrowvect,
152 const IntVect crse_ratio=IntVect(1) );
155 #endif // WARPX_COARSEN_IO_H_ imin
Definition: check_interp_points_and_weights.py:121
void Coarsen(MultiFab &mf_dst, const MultiFab &mf_src, const int dcomp, const int scomp, const int ncomp, const int ngrow, const IntVect crse_ratio=IntVect(1))
Stores in the coarsened MultiFab mf_dst the values obtained by interpolating the data contained in th...
Definition: CoarsenIO.cpp:87
Definition: CoarsenIO.H:15
ii
Definition: check_interp_points_and_weights.py:145
i
Definition: check_interp_points_and_weights.py:171
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: CoarsenIO.H:38
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
void Loop(MultiFab &mf_dst, const MultiFab &mf_src, const int dcomp, const int scomp, const int ncomp, const IntVect ngrow, const IntVect crse_ratio=IntVect(1))
Loops over the boxes of the coarsened MultiFab mf_dst and fills them by interpolating the data contai...
Definition: CoarsenIO.cpp:20