1 #ifndef WARPX_COMPUTEDIAGFUNCTOR_H_
2 #define WARPX_COMPUTEDIAGFUNCTOR_H_
21 m_ncomp(ncomp), m_crse_ratio(crse_ratio) {}
38 virtual void operator() (
amrex::MultiFab& mf_dst,
int dcomp,
int i_buffer = 0)
const = 0;
41 int nComp ()
const {
return m_ncomp; }
60 amrex::Real current_z_boost,
62 const int snapshot_full) {
65 current_z_boost, buffer_box,
66 k_index_zlab, snapshot_full);
75 if (convertRZmodes2cartesian) {
80 "The RZ averaging over modes must write into one single component");
84 for (
int ic=1 ; ic < mf_src.
nComp() ; ic += 2) {
#define WARPX_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition: TextMsg.H:13
Functor to compute a diagnostic and store the result in existing MultiFab.
Definition: ComputeDiagFunctor.H:18
ComputeDiagFunctor(ComputeDiagFunctor &&)=default
ComputeDiagFunctor & operator=(ComputeDiagFunctor &&)=default
int nComp() const
Definition: ComputeDiagFunctor.H:41
virtual void PrepareFunctorData(int i_buffer, bool z_slice_in_domain, amrex::Real current_z_boost, amrex::Box buffer_box, const int k_index_zlab, const int snapshot_full)
Prepare data required to process fields in the operator() Note that this function has parameters that...
Definition: ComputeDiagFunctor.H:59
ComputeDiagFunctor & operator=(const ComputeDiagFunctor &)=default
void InterpolateMFForDiag(amrex::MultiFab &mf_dst, const amrex::MultiFab &mf_src, int dcomp, const amrex::DistributionMapping &dm, bool convertRZmodes2cartesian) const
Definition: ComputeDiagFunctor.H:70
amrex::IntVect m_crse_ratio
Definition: ComputeDiagFunctor.H:106
ComputeDiagFunctor(const ComputeDiagFunctor &)=default
virtual void InitData()
Definition: ComputeDiagFunctor.H:68
ComputeDiagFunctor(int ncomp, amrex::IntVect crse_ratio)
Definition: ComputeDiagFunctor.H:20
virtual ~ComputeDiagFunctor()=default
int m_ncomp
Definition: ComputeDiagFunctor.H:102
IntVect nGrowVect() const noexcept
const BoxArray & boxArray() const noexcept
int nComp() const noexcept
static void Add(MultiFab &dst, const MultiFab &src, int srccomp, int dstcomp, int numcomp, int nghost)
static void Copy(MultiFab &dst, const MultiFab &src, int srccomp, int dstcomp, int numcomp, int nghost)
void Coarsen(amrex::MultiFab &mf_dst, const amrex::MultiFab &mf_src, const int dcomp, const int scomp, const int ncomp, const int ngrow, const amrex::IntVect crse_ratio)
Stores in the coarsened MultiFab mf_dst the values obtained by interpolating the data contained in th...
Definition: sample.cpp:82
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void ignore_unused(const Ts &...)