Block-Structured AMR Software Framework
AMReX_EBFabFactory.H
Go to the documentation of this file.
1 #ifndef AMREX_EBFABFACTORY_H_
2 #define AMREX_EBFABFACTORY_H_
3 #include <AMReX_Config.H>
4 
5 #include <AMReX_FabFactory.H>
6 
8 #include <AMReX_Geometry.H>
9 #include <AMReX_EBSupport.H>
10 #include <AMReX_Array.H>
11 
12 namespace amrex
13 {
14 
15 namespace EB2 {
16  class Level;
17  class IndexSpace;
18 }
19 
21  : public FabFactory<FArrayBox>
22 {
23 public:
24 
25  EBFArrayBoxFactory (const EB2::Level& a_level, const Geometry& a_geom,
26  const BoxArray& a_ba, const DistributionMapping& a_dm,
27  const Vector<int>& a_ngrow, EBSupport a_support);
28  ~EBFArrayBoxFactory () override = default;
29 
31  EBFArrayBoxFactory (EBFArrayBoxFactory&&) noexcept = default;
32 
33  EBFArrayBoxFactory () = delete;
34  EBFArrayBoxFactory& operator= (const EBFArrayBoxFactory&) = delete;
35  EBFArrayBoxFactory& operator= (EBFArrayBoxFactory&&) = delete;
36 
38  FArrayBox* create (const Box& box, int ncomps, const FabInfo& info, int box_index) const final;
39 
41  FArrayBox* create_alias (FArrayBox const& rhs, int scomp, int ncomp) const final;
42 
43  void destroy (FArrayBox* fab) const final;
44 
46  EBFArrayBoxFactory* clone () const final;
47 
48  [[nodiscard]] const FabArray<EBCellFlagFab>& getMultiEBCellFlagFab () const noexcept
49  { return m_ebdc->getMultiEBCellFlagFab(); }
50 
51  [[nodiscard]] const MultiFab& getLevelSet () const noexcept { return m_ebdc->getLevelSet(); }
52 
53  [[nodiscard]] const MultiFab& getVolFrac () const noexcept { return m_ebdc->getVolFrac(); }
54 
55  [[nodiscard]] const MultiCutFab& getCentroid () const noexcept { return m_ebdc->getCentroid(); }
56 
57  [[nodiscard]] const MultiCutFab& getBndryCent () const noexcept { return m_ebdc->getBndryCent(); }
58 
59  [[nodiscard]] const MultiCutFab& getBndryNormal () const noexcept { return m_ebdc->getBndryNormal(); }
60 
61  [[nodiscard]] const MultiCutFab& getBndryArea () const noexcept { return m_ebdc->getBndryArea(); }
62 
63  [[nodiscard]] Array<const MultiCutFab*,AMREX_SPACEDIM> getAreaFrac () const noexcept {
64  return m_ebdc->getAreaFrac();
65  }
66 
67  // For y-face-centroid, the two components are for x and then z-directions, respectively.
68  [[nodiscard]] Array<const MultiCutFab*,AMREX_SPACEDIM> getFaceCent () const noexcept {
69  return m_ebdc->getFaceCent();
70  }
71 
72  [[nodiscard]] Array<const MultiCutFab*,AMREX_SPACEDIM> getEdgeCent () const noexcept {
73  return m_ebdc->getEdgeCent();
74  }
75 
76  [[nodiscard]] bool isAllRegular () const noexcept;
77 
78  [[nodiscard]] EB2::Level const* getEBLevel () const noexcept { return m_parent; }
79  [[nodiscard]] EB2::IndexSpace const* getEBIndexSpace () const noexcept;
80  [[nodiscard]] int maxCoarseningLevel () const noexcept;
81 
82  [[nodiscard]] const DistributionMapping& DistributionMap () const noexcept;
83  [[nodiscard]] const BoxArray& boxArray () const noexcept;
84  [[nodiscard]] const Geometry& Geom () const noexcept { return m_geom; }
85 
86  [[nodiscard]] bool hasEBInfo() const noexcept;
87 
90  [[nodiscard]] iMultiFab const* getCutCellMask () const noexcept { return m_ebdc->getCutCellMask(); }
91 
92 private:
93 
96  std::shared_ptr<EBDataCollection> m_ebdc;
97  EB2::Level const* m_parent = nullptr;
98 };
99 
100 std::unique_ptr<EBFArrayBoxFactory>
101 makeEBFabFactory (const Geometry& a_geom,
102  const BoxArray& a_ba,
103  const DistributionMapping& a_dm,
104  const Vector<int>& a_ngrow, EBSupport a_support);
105 
106 std::unique_ptr<EBFArrayBoxFactory>
108  const BoxArray& a_ba,
109  const DistributionMapping& a_dm,
110  const Vector<int>& a_ngrow, EBSupport a_support);
111 
112 std::unique_ptr<EBFArrayBoxFactory>
113 makeEBFabFactory (const EB2::IndexSpace*, const Geometry& a_geom,
114  const BoxArray& a_ba,
115  const DistributionMapping& a_dm,
116  const Vector<int>& a_ngrow, EBSupport a_support);
117 
118 }
119 
120 #endif
#define AMREX_NODISCARD
Definition: AMReX_Extension.H:251
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
Definition: AMReX_EB2.H:26
Definition: AMReX_EB2_Level.H:33
Definition: AMReX_EBCellFlag.H:287
Definition: AMReX_EBFabFactory.H:22
EBSupport m_support
Definition: AMReX_EBFabFactory.H:94
AMREX_NODISCARD FArrayBox * create(const Box &box, int ncomps, const FabInfo &info, int box_index) const final
Definition: AMReX_EBFabFactory.cpp:27
std::shared_ptr< EBDataCollection > m_ebdc
Definition: AMReX_EBFabFactory.H:96
Array< const MultiCutFab *, AMREX_SPACEDIM > getFaceCent() const noexcept
Definition: AMReX_EBFabFactory.H:68
Geometry m_geom
Definition: AMReX_EBFabFactory.H:95
void destroy(FArrayBox *fab) const final
Definition: AMReX_EBFabFactory.cpp:57
const MultiCutFab & getBndryCent() const noexcept
Definition: AMReX_EBFabFactory.H:57
EB2::IndexSpace const * getEBIndexSpace() const noexcept
Definition: AMReX_EBFabFactory.cpp:84
Array< const MultiCutFab *, AMREX_SPACEDIM > getEdgeCent() const noexcept
Definition: AMReX_EBFabFactory.H:72
iMultiFab const * getCutCellMask() const noexcept
Definition: AMReX_EBFabFactory.H:90
bool isAllRegular() const noexcept
Definition: AMReX_EBFabFactory.cpp:78
const BoxArray & boxArray() const noexcept
Definition: AMReX_EBFabFactory.cpp:107
const Geometry & Geom() const noexcept
Definition: AMReX_EBFabFactory.H:84
EB2::Level const * getEBLevel() const noexcept
Definition: AMReX_EBFabFactory.H:78
Array< const MultiCutFab *, AMREX_SPACEDIM > getAreaFrac() const noexcept
Definition: AMReX_EBFabFactory.H:63
const MultiCutFab & getCentroid() const noexcept
Definition: AMReX_EBFabFactory.H:55
int maxCoarseningLevel() const noexcept
Definition: AMReX_EBFabFactory.cpp:90
const DistributionMapping & DistributionMap() const noexcept
Definition: AMReX_EBFabFactory.cpp:101
EB2::Level const * m_parent
Definition: AMReX_EBFabFactory.H:97
const MultiFab & getLevelSet() const noexcept
Definition: AMReX_EBFabFactory.H:51
const FabArray< EBCellFlagFab > & getMultiEBCellFlagFab() const noexcept
Definition: AMReX_EBFabFactory.H:48
const MultiCutFab & getBndryNormal() const noexcept
Definition: AMReX_EBFabFactory.H:59
EBFArrayBoxFactory(EBFArrayBoxFactory &&) noexcept=default
const MultiFab & getVolFrac() const noexcept
Definition: AMReX_EBFabFactory.H:53
const MultiCutFab & getBndryArea() const noexcept
Definition: AMReX_EBFabFactory.H:61
AMREX_NODISCARD FArrayBox * create_alias(FArrayBox const &rhs, int scomp, int ncomp) const final
Definition: AMReX_EBFabFactory.cpp:43
bool hasEBInfo() const noexcept
Definition: AMReX_EBFabFactory.cpp:113
~EBFArrayBoxFactory() override=default
AMREX_NODISCARD EBFArrayBoxFactory * clone() const final
Definition: AMReX_EBFabFactory.cpp:72
EBFArrayBoxFactory(const EBFArrayBoxFactory &)=default
A Fortran Array of REALs.
Definition: AMReX_FArrayBox.H:229
An Array of FortranArrayBox(FAB)-like Objects.
Definition: AMReX_FabArray.H:344
Definition: AMReX_FabFactory.H:50
Rectangular problem domain geometry.
Definition: AMReX_Geometry.H:73
Definition: AMReX_MultiCutFab.H:81
A collection (stored as an array) of FArrayBox objects.
Definition: AMReX_MultiFab.H:38
Definition: AMReX_iMultiFab.H:32
Definition: AMReX_Amr.cpp:49
EBSupport
Definition: AMReX_EBSupport.H:7
std::unique_ptr< EBFArrayBoxFactory > makeEBFabFactory(const Geometry &a_geom, const BoxArray &a_ba, const DistributionMapping &a_dm, const Vector< int > &a_ngrow, EBSupport a_support)
Definition: AMReX_EBFabFactory.cpp:119
std::array< T, N > Array
Definition: AMReX_Array.H:24
Definition: AMReX_FabFactory.H:27