Block-Structured AMR Software Framework
AMReX_MultiFab.H
Go to the documentation of this file.
1 
2 #ifndef BL_MULTIFAB_H
3 #define BL_MULTIFAB_H
4 #include <AMReX_Config.H>
5 
6 #include <AMReX_BLassert.H>
7 #include <AMReX_BaseFab.H>
8 #include <AMReX_FArrayBox.H>
9 #include <AMReX_FabArray.H>
10 #include <AMReX_FabArrayUtility.H>
11 #include <AMReX_Periodicity.H>
12 #include <AMReX_NonLocalBC.H>
13 
14 #include <cstdint>
15 
16 namespace amrex
17 {
18 
20 
21 class iMultiFab;
22 
35 class MultiFab
36  :
37  public FabArray<FArrayBox>
38 {
39 public:
40 
47  MultiFab () noexcept;
48 
57  explicit MultiFab (Arena* a) noexcept;
58 
73  MultiFab (const BoxArray& bxs,
74  const DistributionMapping& dm,
75  int ncomp,
76  int ngrow,
77 #ifdef AMREX_STRICT_MODE
78  const MFInfo& info,
79  const FabFactory<FArrayBox>& factory);
80 #else
81  const MFInfo& info = MFInfo(),
82  const FabFactory<FArrayBox>& factory = FArrayBoxFactory());
83 #endif
84 
85  MultiFab (const BoxArray& bxs,
86  const DistributionMapping& dm,
87  int ncomp,
88  const IntVect& ngrow,
89 #ifdef AMREX_STRICT_MODE
90  const MFInfo& info,
91  const FabFactory<FArrayBox>& factory);
92 #else
93  const MFInfo& info = MFInfo(),
94  const FabFactory<FArrayBox>& factory = FArrayBoxFactory());
95 #endif
96 
104  MultiFab (const MultiFab& rhs, MakeType maketype, int scomp, int ncomp);
105 
106  ~MultiFab ();
107 
108  MultiFab (MultiFab&& rhs) noexcept;
109  MultiFab& operator= (MultiFab&& rhs) noexcept = default;
110 
111  MultiFab (const MultiFab& rhs) = delete;
112  MultiFab& operator= (const MultiFab& rhs) = delete;
113 
114  void define (const BoxArray& bxs,
115  const DistributionMapping& dm,
116  int nvar,
117  int ngrow,
118 #ifdef AMREX_STRICT_MODE
119  const MFInfo& info,
120  const FabFactory<FArrayBox>& factory);
121 #else
122  const MFInfo& info = MFInfo(),
123  const FabFactory<FArrayBox>& factory = FArrayBoxFactory());
124 #endif
125 
126  void define (const BoxArray& bxs,
127  const DistributionMapping& dm,
128  int nvar,
129  const IntVect& ngrow,
130 #ifdef AMREX_STRICT_MODE
131  const MFInfo& info,
132  const FabFactory<FArrayBox>& factory);
133 #else
134  const MFInfo& info = MFInfo(),
135  const FabFactory<FArrayBox>& factory = FArrayBoxFactory());
136 #endif
137 
138  MultiFab& operator= (Real r);
139 
148  [[nodiscard]] Real min (int comp,
149  int nghost = 0,
150  bool local = false) const;
155  [[nodiscard]] Real min (const Box& region,
156  int comp,
157  int nghost = 0,
158  bool local = false) const;
167  [[nodiscard]] Real max (int comp,
168  int nghost = 0,
169  bool local = false) const;
174  [[nodiscard]] Real max (const Box& region,
175  int comp,
176  int nghost = 0,
177  bool local = false) const;
182  [[nodiscard]] Real norm0 (int comp = 0, int nghost = 0, bool local = false, bool ignore_covered = false ) const;
183  [[nodiscard]] Real norminf (int comp = 0, int nghost = 0, bool local = false, bool ignore_covered = false ) const {
184  return norm0(comp,nghost,local,ignore_covered);
185  }
186 
187  [[nodiscard]] Real norm0 (const iMultiFab& mask, int comp = 0, int nghost = 0, bool local = false) const;
188  [[nodiscard]] Real norminf (const iMultiFab& mask, int comp = 0, int nghost = 0, bool local = false) const {
189  return norm0(mask,comp,nghost,local);
190  }
191 
192  [[nodiscard]] Real norm0 (int comp, int ncomp, IntVect const& nghost, bool local = false,
193  bool ignore_covered = false) const;
194 
196 
201  [[nodiscard]] Vector<Real> norm0 (const Vector<int>& comps, int nghost = 0, bool local = false, bool ignore_covered = false ) const;
202  [[nodiscard]] Vector<Real> norminf (const Vector<int>& comps, int nghost = 0, bool local = false, bool ignore_covered = false) const {
203  return norm0(comps,nghost,local,ignore_covered);
204  }
205 
211  [[nodiscard]] Real norm1 (int comp, const Periodicity& period, bool ignore_covered = false) const;
216  [[nodiscard]] Real norm1 (int comp = 0, int ngrow = 0, bool local = false) const;
221  [[nodiscard]] Vector<Real> norm1 (const Vector<int>& comps, int ngrow = 0, bool local = false) const;
226  [[nodiscard]] Real norm2 (int comp = 0) const;
231  [[nodiscard]] Real norm2 (int comp, const Periodicity& period) const;
236  [[nodiscard]] Vector<Real> norm2 (const Vector<int>& comps) const;
240  [[nodiscard]] Real sum (int comp = 0, bool local = false) const;
244  [[nodiscard]] Real sum (Box const& region, int comp = 0, bool local = false) const;
245 
247 
252  [[nodiscard]] Real sum_unique (int comp = 0,
253  bool local = false,
254  const Periodicity& period = Periodicity::NonPeriodic()) const;
261  [[nodiscard]] Real sum_unique (Box const& region, int comp = 0, bool local = false) const;
271  void plus (Real val,
272  int comp,
273  int num_comp,
274  int nghost = 0);
280  void plus (Real val,
281  const Box& region,
282  int comp,
283  int num_comp,
284  int nghost = 0);
291  void plus (Real val,
292  int nghost);
300  void plus (Real val,
301  const Box& region,
302  int nghost);
311  void mult (Real val,
312  int comp,
313  int num_comp,
314  int nghost = 0);
322  void mult (Real val,
323  const Box& region,
324  int comp,
325  int num_comp,
326  int nghost = 0);
333  void mult (Real val,
334  int nghost = 0);
342  void mult (Real val,
343  const Box& region,
344  int nghost = 0);
353  void invert (Real numerator,
354  int comp,
355  int num_comp,
356  int nghost = 0);
364  void invert (Real numerator,
365  const Box& region,
366  int comp,
367  int num_comp,
368  int nghost = 0);
375  void invert (Real numerator,
376  int nghost);
384  void invert (Real numerator,
385  const Box& region,
386  int nghost);
394  void negate (int comp,
395  int num_comp,
396  int nghost = 0);
402  void negate (const Box& region,
403  int comp,
404  int num_comp,
405  int nghost = 0);
411  void negate (int nghost = 0);
418  void negate (const Box& region,
419  int nghost = 0);
420 
421  [[nodiscard]] IntVect minIndex (int comp,
422  int nghost = 0) const;
423 
424  [[nodiscard]] IntVect maxIndex (int comp,
425  int nghost = 0) const;
435  void plus (const MultiFab& mf,
436  int strt_comp,
437  int num_comp,
438  int nghost);
448  void minus (const MultiFab& mf,
449  int strt_comp,
450  int num_comp,
451  int nghost);
462  void divide (const MultiFab& mf,
463  int strt_comp,
464  int num_comp,
465  int nghost);
469  static Real Dot (const MultiFab& x, int xcomp,
470  const MultiFab& y, int ycomp,
471  int numcomp, int nghost, bool local = false);
472 
476  static Real Dot (const MultiFab& x, int xcomp,
477  int numcomp, int nghost, bool local = false);
478 
479  static Real Dot (const iMultiFab& mask,
480  const MultiFab& x, int xcomp,
481  const MultiFab& y, int ycomp,
482  int numcomp, int nghost, bool local = false);
487  static void Add (MultiFab& dst,
488  const MultiFab& src,
489  int srccomp,
490  int dstcomp,
491  int numcomp,
492  int nghost);
493 
494  static void Add (MultiFab& dst,
495  const MultiFab& src,
496  int srccomp,
497  int dstcomp,
498  int numcomp,
499  const IntVect& nghost);
500 
506  [[nodiscard]] MultiFab deepCopy () const;
507 
513  static void Copy (MultiFab& dst,
514  const MultiFab& src,
515  int srccomp,
516  int dstcomp,
517  int numcomp,
518  int nghost);
519 
520  static void Copy (MultiFab& dst,
521  const MultiFab& src,
522  int srccomp,
523  int dstcomp,
524  int numcomp,
525  const IntVect& nghost);
526 
532  static void Swap (MultiFab& dst,
533  MultiFab& src,
534  int srccomp,
535  int dstcomp,
536  int numcomp,
537  int nghost);
538 
539  static void Swap (MultiFab& dst,
540  MultiFab& src,
541  int srccomp,
542  int dstcomp,
543  int numcomp,
544  const IntVect& nghost);
545 
550  static void Subtract (MultiFab& dst,
551  const MultiFab& src,
552  int srccomp,
553  int dstcomp,
554  int numcomp,
555  int nghost);
556 
557  static void Subtract (MultiFab& dst,
558  const MultiFab& src,
559  int srccomp,
560  int dstcomp,
561  int numcomp,
562  const IntVect& nghost);
567  static void Multiply (MultiFab& dst,
568  const MultiFab& src,
569  int srccomp,
570  int dstcomp,
571  int numcomp,
572  int nghost);
573 
574  static void Multiply (MultiFab& dst,
575  const MultiFab& src,
576  int srccomp,
577  int dstcomp,
578  int numcomp,
579  const IntVect& nghost);
584  static void Divide (MultiFab& dst,
585  const MultiFab& src,
586  int srccomp,
587  int dstcomp,
588  int numcomp,
589  int nghost);
590 
591  static void Divide (MultiFab& dst,
592  const MultiFab& src,
593  int srccomp,
594  int dstcomp,
595  int numcomp,
596  const IntVect& nghost);
600  static void Saxpy (MultiFab& dst,
601  Real a,
602  const MultiFab& src,
603  int srccomp,
604  int dstcomp,
605  int numcomp,
606  int nghost);
607 
609 
613  static void Xpay (MultiFab& dst,
614  Real a,
615  const MultiFab& src,
616  int srccomp,
617  int dstcomp,
618  int numcomp,
619  int nghost);
620 
622 
626  static void LinComb (MultiFab& dst,
627  Real a,
628  const MultiFab& x,
629  int xcomp,
630  Real b,
631  const MultiFab& y,
632  int ycomp,
633  int dstcomp,
634  int numcomp,
635  int nghost);
636 
638 
642  static void AddProduct (MultiFab& dst,
643  const MultiFab& src1,
644  int comp1,
645  const MultiFab& src2,
646  int comp2,
647  int dstcomp,
648  int numcomp,
649  int nghost);
650 
651  static void AddProduct (MultiFab& dst,
652  const MultiFab& src1,
653  int comp1,
654  const MultiFab& src2,
655  int comp2,
656  int dstcomp,
657  int numcomp,
658  const IntVect& nghost);
659 
667  [[nodiscard]] bool is_finite (bool local=false) const;
668  [[nodiscard]] bool is_finite (int scomp, int ncomp, int ngrow = 0, bool local=false) const;
669  [[nodiscard]] bool is_finite (int scomp, int ncomp, const IntVect& ngrow, bool local=false) const;
670 
678  [[nodiscard]] bool contains_nan (bool local=false) const;
679  [[nodiscard]] bool contains_nan (int scomp, int ncomp, int ngrow = 0, bool local=false) const;
680  [[nodiscard]] bool contains_nan (int scomp, int ncomp, const IntVect& ngrow, bool local=false) const;
687  [[nodiscard]] bool contains_inf (bool local=false) const;
688  [[nodiscard]] bool contains_inf (int scomp, int ncomp, int ngrow = 0, bool local=false) const;
689  [[nodiscard]] bool contains_inf (int scomp, int ncomp, const IntVect& ngrow, bool local=false) const;
690 
694  [[nodiscard]] std::unique_ptr<MultiFab> OverlapMask (const Periodicity& period = Periodicity::NonPeriodic()) const;
696  [[nodiscard]] std::unique_ptr<iMultiFab> OwnerMask (const Periodicity& period = Periodicity::NonPeriodic()) const;
697 
699  void AverageSync (const Periodicity& period = Periodicity::NonPeriodic());
701  void WeightedSync (const MultiFab& wgt, const Periodicity& period = Periodicity::NonPeriodic());
703  void OverrideSync (const iMultiFab& msk, const Periodicity& period = Periodicity::NonPeriodic());
704 
706 
707  static void Initialize ();
708  static void Finalize ();
709 
710 private:
711  //
715 
716  void initVal ();
717 };
718 
719 #ifndef _MSC_VER
720 inline void GccPlacaterMF ()
721 {
722  std::allocator<MultiFab*> a1;
723  std::allocator<MultiFab const*> a2;
724  std::allocator<FabArray<FArrayBox>*> a3;
725  std::allocator<FabArray<FArrayBox> const*> a4;
726 
731 }
732 #endif
733 
734 }
735 
736 #endif /*BL_MULTIFAB_H*/
Array4< int const > mask
Definition: AMReX_InterpFaceRegister.cpp:93
A virtual base class for objects that manage their own dynamic memory allocation.
Definition: AMReX_Arena.H:100
A collection of Boxes stored in an Array.
Definition: AMReX_BoxArray.H:550
Calculates the distribution of FABs to MPI processes.
Definition: AMReX_DistributionMapping.H:41
CopyComTag::CopyComTagsContainer CopyComTagsContainer
Definition: AMReX_FabArrayBase.H:219
CopyComTag::MapOfCopyComTagContainers MapOfCopyComTagContainers
Definition: AMReX_FabArrayBase.H:220
An Array of FortranArrayBox(FAB)-like Objects.
Definition: AMReX_FabArray.H:344
A collection (stored as an array) of FArrayBox objects.
Definition: AMReX_MultiFab.H:38
static void LinComb(MultiFab &dst, Real a, const MultiFab &x, int xcomp, Real b, const MultiFab &y, int ycomp, int dstcomp, int numcomp, int nghost)
dst = a*x + b*y
Definition: AMReX_MultiFab.cpp:359
void define(const BoxArray &bxs, const DistributionMapping &dm, int nvar, int ngrow, const MFInfo &info=MFInfo(), const FabFactory< FArrayBox > &factory=FArrayBoxFactory())
Definition: AMReX_MultiFab.cpp:572
void invert(Real numerator, int comp, int num_comp, int nghost=0)
Replaces the value of each cell in the specified subregion of the MultiFab with its reciprocal multip...
Definition: AMReX_MultiFab.cpp:1476
void minus(const MultiFab &mf, int strt_comp, int num_comp, int nghost)
This function subtracts the values of the cells in mf from the corresponding cells of this MultiFab....
Definition: AMReX_MultiFab.cpp:1417
Real norm0(int comp=0, int nghost=0, bool local=false, bool ignore_covered=false) const
Returns the maximum absolute value contained in component comp of the MultiFab.
Definition: AMReX_MultiFab.cpp:1084
static void Swap(MultiFab &dst, MultiFab &src, int srccomp, int dstcomp, int numcomp, int nghost)
Swap from src to dst including nghost ghost cells. The two MultiFabs MUST have the same underlying Bo...
Definition: AMReX_MultiFab.cpp:213
std::unique_ptr< iMultiFab > OwnerMask(const Periodicity &period=Periodicity::NonPeriodic()) const
Owner is the grid with the lowest grid number containing the data.
Definition: AMReX_MultiFab.cpp:1580
static void Add(MultiFab &dst, const MultiFab &src, int srccomp, int dstcomp, int numcomp, int nghost)
Add src to dst including nghost ghost cells. The two MultiFabs MUST have the same underlying BoxArray...
Definition: AMReX_MultiFab.cpp:151
bool contains_nan(bool local=false) const
Are there any NaNs in the MF? This may return false, even if the MF contains NaNs,...
Definition: AMReX_MultiFab.cpp:713
static void Saxpy(MultiFab &dst, Real a, const MultiFab &src, int srccomp, int dstcomp, int numcomp, int nghost)
dst += a*src
Definition: AMReX_MultiFab.cpp:345
static Real Dot(const MultiFab &x, int xcomp, const MultiFab &y, int ycomp, int numcomp, int nghost, bool local=false)
Returns the dot product of two MultiFabs.
Definition: AMReX_MultiFab.cpp:37
Real min(int comp, int nghost=0, bool local=false) const
Returns the minimum value contained in component comp of the MultiFab.
Definition: AMReX_MultiFab.cpp:774
Real norm1(int comp, const Periodicity &period, bool ignore_covered=false) const
Returns the L1 norm of component comp over the MultiFab.
Definition: AMReX_MultiFab.cpp:1182
static void Subtract(MultiFab &dst, const MultiFab &src, int srccomp, int dstcomp, int numcomp, int nghost)
Subtract src from dst including nghost ghost cells. The two MultiFabs MUST have the same underlying B...
Definition: AMReX_MultiFab.cpp:285
Real sum(int comp=0, bool local=false) const
Returns the sum of component "comp" over the MultiFab – no ghost cells are included.
Definition: AMReX_MultiFab.cpp:1262
static void Multiply(MultiFab &dst, const MultiFab &src, int srccomp, int dstcomp, int numcomp, int nghost)
Multiply dst by src including nghost ghost cells. The two MultiFabs MUST have the same underlying Box...
Definition: AMReX_MultiFab.cpp:305
IntVect maxIndex(int comp, int nghost=0) const
Definition: AMReX_MultiFab.cpp:1070
static void AddProduct(MultiFab &dst, const MultiFab &src1, int comp1, const MultiFab &src2, int comp2, int dstcomp, int numcomp, int nghost)
dst += src1*src2
Definition: AMReX_MultiFab.cpp:368
bool is_finite(bool local=false) const
Are the numbers in the MF finite (i.e., neither nan nor inf)? This may return true,...
Definition: AMReX_MultiFab.cpp:658
static void Initialize()
Definition: AMReX_MultiFab.cpp:475
void AverageSync(const Periodicity &period=Periodicity::NonPeriodic())
Sync up nodal data via averaging.
Definition: AMReX_MultiFab.cpp:1586
Vector< Real > norminf(const Vector< int > &comps, int nghost=0, bool local=false, bool ignore_covered=false) const
Definition: AMReX_MultiFab.H:202
bool contains_inf(bool local=false) const
Are there any Infs in the MF? This may return false, even if the MF contains Infs,...
Definition: AMReX_MultiFab.cpp:768
Real norminf(int comp=0, int nghost=0, bool local=false, bool ignore_covered=false) const
Definition: AMReX_MultiFab.H:183
static void Xpay(MultiFab &dst, Real a, const MultiFab &src, int srccomp, int dstcomp, int numcomp, int nghost)
dst = src + a*dst
Definition: AMReX_MultiFab.cpp:352
Real norm2(int comp=0) const
Returns the L2 norm of component comp over the MultiFab. No ghost cells are used.
Definition: AMReX_MultiFab.cpp:1114
MultiFab & operator=(MultiFab &&rhs) noexcept=default
void initVal()
Definition: AMReX_MultiFab.cpp:596
static void Divide(MultiFab &dst, const MultiFab &src, int srccomp, int dstcomp, int numcomp, int nghost)
Divide dst by src including nghost ghost cells. The two MultiFabs MUST have the same underlying BoxAr...
Definition: AMReX_MultiFab.cpp:325
void mult(Real val, int comp, int num_comp, int nghost=0)
Scales the value of each cell in the specified subregion of the MultiFab by the scalar val (a[i] <- a...
Definition: AMReX_MultiFab.cpp:1456
void plus(Real val, int comp, int num_comp, int nghost=0)
Adds the scalar value val to the value of each cell in the specified subregion of the MultiFab.
Definition: AMReX_MultiFab.cpp:1429
IntVect minIndex(int comp, int nghost=0) const
Definition: AMReX_MultiFab.cpp:1062
static void Copy(MultiFab &dst, const MultiFab &src, int srccomp, int dstcomp, int numcomp, int nghost)
Copy from src to dst including nghost ghost cells. The two MultiFabs MUST have the same underlying Bo...
Definition: AMReX_MultiFab.cpp:193
Real max(int comp, int nghost=0, bool local=false) const
Returns the maximum value contained in component comp of the MultiFab.
Definition: AMReX_MultiFab.cpp:903
Real norminf(const iMultiFab &mask, int comp=0, int nghost=0, bool local=false) const
Definition: AMReX_MultiFab.H:188
MultiFab() noexcept
Constructs an empty MultiFab.
Definition: AMReX_MultiFab.cpp:496
void divide(const MultiFab &mf, int strt_comp, int num_comp, int nghost)
This function divides the values of the cells in mf from the corresponding cells of this MultiFab....
Definition: AMReX_MultiFab.cpp:1423
void negate(int comp, int num_comp, int nghost=0)
Negates the value of each cell in the specified subregion of the MultiFab. The subregion consists of ...
Definition: AMReX_MultiFab.cpp:1496
MultiFab deepCopy() const
Definition: AMReX_MultiFab.cpp:171
void WeightedSync(const MultiFab &wgt, const Periodicity &period=Periodicity::NonPeriodic())
Sync up nodal data with weights.
Definition: AMReX_MultiFab.cpp:1597
MultiFab(const MultiFab &rhs)=delete
~MultiFab()
Definition: AMReX_MultiFab.cpp:557
static void Finalize()
Definition: AMReX_MultiFab.cpp:491
std::unique_ptr< MultiFab > OverlapMask(const Periodicity &period=Periodicity::NonPeriodic()) const
Return a mask indicating how many duplicates are in each point.
Definition: AMReX_MultiFab.cpp:1514
void OverrideSync(const iMultiFab &msk, const Periodicity &period=Periodicity::NonPeriodic())
Sync up nodal data with owners overriding non-owners.
Definition: AMReX_MultiFab.cpp:1617
Real sum_unique(int comp=0, bool local=false, const Periodicity &period=Periodicity::NonPeriodic()) const
Same as sum with local =false, but for non-cell-centered data, this only adds non-unique points that ...
Definition: AMReX_MultiFab.cpp:1310
This provides length of period for periodic domains. 0 means it is not periodic in that direction....
Definition: AMReX_Periodicity.H:17
static const Periodicity & NonPeriodic() noexcept
Definition: AMReX_Periodicity.cpp:52
Definition: AMReX_iMultiFab.H:32
Definition: AMReX_Amr.cpp:49
MakeType
Definition: AMReX_MakeType.H:7
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition: AMReX.H:111
DefaultFabFactory< FArrayBox > FArrayBoxFactory
Definition: AMReX_FArrayBox.H:494
void GccPlacaterMF()
Definition: AMReX_MultiFab.H:720
FabArray memory allocation information.
Definition: AMReX_FabArray.H:66