1 #ifndef ERF_PlaneAverage_H
2 #define ERF_PlaneAverage_H
5 #include "AMReX_iMultiFab.H"
6 #include "AMReX_MultiFab.H"
7 #include "AMReX_GpuContainers.H"
18 amrex::Geometry geom_in,
20 bool inc_ghost=
false);
28 [[nodiscard]] AMREX_FORCE_INLINE
34 [[nodiscard]] amrex::Real
dx ()
const {
return m_dx; }
35 [[nodiscard]] amrex::Real
xlo ()
const {
return m_xlo; }
43 [[nodiscard]]
const amrex::Vector<amrex::Real>&
line_average ()
const
49 void line_average (
int comp, amrex::Gpu::HostVector<amrex::Real>& l_vec);
56 [[nodiscard]]
const amrex::MultiFab&
field ()
const {
return *
m_field; }
80 amrex::IntVect
m_ng = amrex::IntVect(0);
84 template <
typename IndexSelector>
91 amrex::Geometry geom_in,
94 : m_field(field_in), m_geom(geom_in), m_axis(axis_in), m_inc_ghost(inc_ghost)
96 AMREX_ALWAYS_ASSERT(
m_axis >= 0 &&
m_axis < AMREX_SPACEDIM);
103 amrex::Box domain =
m_geom.Domain();
105 m_ng = field_in->nGrowVect();
106 domain.grow(axis_in,
m_ng[axis_in]);
108 amrex::IntVect dom_lo(domain.loVect());
109 amrex::IntVect dom_hi(domain.hiVect());
114 auto period =
m_geom.periodicity();
115 for (
int i = 0; i < AMREX_SPACEDIM; ++i) {
116 int p_fac = (!period.isPeriodic(i)) ? 1 : 0;
131 AMREX_ALWAYS_ASSERT(comp >= 0 && comp <
m_ncomp);
137 ind =
static_cast<int>(floor((
x -
m_xlo) /
m_dx - 0.5));
138 const amrex::Real x1 =
m_xlo + (ind + 0.5) *
m_dx;
147 AMREX_ALWAYS_ASSERT(ind >= 0 && ind + 1 <
m_ncell_line);
156 AMREX_ALWAYS_ASSERT(comp >= 0 && comp <
m_ncomp);
177 amrex::Abort(
"axis must be equal to 0, 1, or 2");
182 template <
typename IndexSelector>
188 amrex::Real* line_avg = lavg.data();
193 amrex::IntVect ng = amrex::IntVect(0);
197 std::unique_ptr<amrex::iMultiFab> mask = OwnerMask(*
m_field,
m_geom.periodicity(), ng);
200 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
202 for (amrex::MFIter mfi(mfab, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
203 amrex::Box tbx = mfi.tilebox();
206 amrex::Box pbx = PerpendicularBox<IndexSelector>(tbx, amrex::IntVect(0));
208 const amrex::Array4<const amrex::Real>& fab_arr = mfab.const_array(mfi);
209 const amrex::Array4<const int >& mask_arr = mask->const_array(mfi);
211 amrex::ParallelFor(amrex::Gpu::KernelInfo().setReduction(
true), pbx, [=]
212 AMREX_GPU_DEVICE(
int p_i,
int p_j,
int p_k,
213 amrex::Gpu::Handler
const& handler) noexcept
217 amrex::Box lbx = ParallelBox<IndexSelector>(tbx, amrex::IntVect{p_i, p_j, p_k});
219 for (
int k = lbx.smallEnd(2); k <= lbx.bigEnd(2); ++k) {
220 for (
int j = lbx.smallEnd(1); j <= lbx.bigEnd(1); ++j) {
221 for (
int i = lbx.smallEnd(0); i <= lbx.bigEnd(0); ++i) {
222 int ind = idxOp.getIndx(i, j, k) +
offset;
223 for (
int n = 0; n <
ncomp; ++n) {
227 amrex::Real fac = (mask_arr(i,j,k)) ? 1.0 : 0.0;
228 amrex::Gpu::deviceReduceSum(&line_avg[
ncomp * ind + n],
229 fab_arr(i, j, k, n) * denom * fac, handler);
DirectionSelector< 0 > XDir
Definition: ERF_DirectionSelector.H:36
DirectionSelector< 1 > YDir
Definition: ERF_DirectionSelector.H:37
DirectionSelector< 2 > ZDir
Definition: ERF_DirectionSelector.H:38
AMREX_FORCE_INLINE IntVect offset(const int face_dir, const int normal)
Definition: ERF_ReadBndryPlanes.cpp:28
Definition: ERF_PlaneAverage.H:14
const amrex::MultiFab & field() const
Definition: ERF_PlaneAverage.H:56
int m_ncell_line
Definition: ERF_PlaneAverage.H:70
AMREX_FORCE_INLINE amrex::Real line_average_interpolated(amrex::Real x, int comp) const
Definition: ERF_PlaneAverage.H:129
int m_precision
Definition: ERF_PlaneAverage.H:72
int ncomp() const
Definition: ERF_PlaneAverage.H:39
amrex::Real m_xlo
Definition: ERF_PlaneAverage.H:67
amrex::IntVect m_ixtype
Definition: ERF_PlaneAverage.H:79
amrex::IntVect m_ng
Definition: ERF_PlaneAverage.H:80
amrex::Vector< amrex::Real > m_line_xcentroid
Definition: ERF_PlaneAverage.H:64
const amrex::Vector< amrex::Real > & line_centroids() const
Definition: ERF_PlaneAverage.H:51
AMREX_FORCE_INLINE void compute_averages(const IndexSelector &idxOp, const amrex::MultiFab &mfab)
void set_precision(int p)
Definition: ERF_PlaneAverage.H:32
const int m_level
Definition: ERF_PlaneAverage.H:73
const amrex::MultiFab * m_field
Definition: ERF_PlaneAverage.H:75
int level() const
Definition: ERF_PlaneAverage.H:38
int m_ncell_plane
Definition: ERF_PlaneAverage.H:69
amrex::Vector< amrex::Real > m_line_average
Definition: ERF_PlaneAverage.H:62
amrex::Real xlo() const
Definition: ERF_PlaneAverage.H:35
AMREX_FORCE_INLINE void operator()()
Definition: ERF_PlaneAverage.H:163
int ncell_plane() const
Definition: ERF_PlaneAverage.H:40
const int m_axis
Definition: ERF_PlaneAverage.H:77
int ncell_line() const
Definition: ERF_PlaneAverage.H:41
bool m_inc_ghost
Definition: ERF_PlaneAverage.H:78
amrex::Real dx() const
Definition: ERF_PlaneAverage.H:34
amrex::Real m_dx
Definition: ERF_PlaneAverage.H:66
amrex::Geometry m_geom
Definition: ERF_PlaneAverage.H:76
int axis() const
Definition: ERF_PlaneAverage.H:37
int m_ncomp
Definition: ERF_PlaneAverage.H:59
const amrex::Vector< amrex::Real > & line_average() const
Definition: ERF_PlaneAverage.H:43