WarpX
ComputeDiagFunctor.H
Go to the documentation of this file.
1 #ifndef WARPX_COMPUTEDIAGFUNCTOR_H_
2 #define WARPX_COMPUTEDIAGFUNCTOR_H_
3 
5 #include "Utils/TextMsg.H"
6 
8 
9 #include <AMReX.H>
10 #include <AMReX_MultiFab.H>
11 
16 class
18 {
19 public:
20  ComputeDiagFunctor( int ncomp, amrex::IntVect crse_ratio) :
21  m_ncomp(ncomp), m_crse_ratio(crse_ratio) {}
22  //** Virtual Destructor to handle clean destruction of derived classes */
23  virtual ~ComputeDiagFunctor() = default;
24 
25  // Default move and copy operations
30 
38  virtual void operator() (amrex::MultiFab& mf_dst, int dcomp, int i_buffer = 0) const = 0;
41  int nComp () const { return m_ncomp; }
42 
59  virtual void PrepareFunctorData ( int i_buffer, bool z_slice_in_domain,
60  amrex::Real current_z_boost,
61  amrex::Box buffer_box, const int k_index_zlab,
62  const int snapshot_full) {
63  amrex::ignore_unused(i_buffer,
64  z_slice_in_domain,
65  current_z_boost, buffer_box,
66  k_index_zlab, snapshot_full);
67  }
68  virtual void InitData() {}
69 
71  amrex::MultiFab& mf_dst, const amrex::MultiFab& mf_src, int dcomp,
72  const amrex::DistributionMapping& dm, bool convertRZmodes2cartesian ) const
73  {
74 #ifdef WARPX_DIM_RZ
75  if (convertRZmodes2cartesian) {
76  // In cylindrical geometry, sum real part of all modes of mf_src in
77  // temporary MultiFab mf_dst_stag, and cell-center it to mf_dst
79  nComp()==1,
80  "The RZ averaging over modes must write into one single component");
81  amrex::MultiFab mf_dst_stag( mf_src.boxArray(), dm, 1, mf_src.nGrowVect() );
82  // Mode 0
83  amrex::MultiFab::Copy( mf_dst_stag, mf_src, 0, 0, 1, mf_src.nGrowVect() );
84  for (int ic=1 ; ic < mf_src.nComp() ; ic += 2) {
85  // Real part of all modes > 0
86  amrex::MultiFab::Add( mf_dst_stag, mf_src, ic, 0, 1, mf_src.nGrowVect() );
87  }
88  ablastr::coarsen::sample::Coarsen( mf_dst, mf_dst_stag, dcomp, 0, nComp(), 0, m_crse_ratio );
89  } else {
90  ablastr::coarsen::sample::Coarsen( mf_dst, mf_src, dcomp, 0, nComp(), 0, m_crse_ratio );
91  }
92 #else
93  // In Cartesian geometry, coarsen and interpolate from temporary MultiFab mf_src
94  // to output diagnostic MultiFab mf_dst
95  ablastr::coarsen::sample::Coarsen(mf_dst, mf_src, dcomp, 0, nComp(), mf_dst.nGrowVect(), m_crse_ratio );
96  amrex::ignore_unused(convertRZmodes2cartesian, dm);
97 #endif
98  }
99 
100 private:
102  int m_ncomp;
103 
104 protected:
107 };
108 
109 #endif // WARPX_COMPUTEDIAGFUNCTOR_H_
#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 &...)