Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
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
7#include <AMReX_EBData.H>
9#include <AMReX_Geometry.H>
10#include <AMReX_EBSupport.H>
11#include <AMReX_Array.H>
12#include <AMReX_MFIter.H>
13
14namespace amrex
15{
16
17namespace EB2 {
18 class Level;
19 class IndexSpace;
20}
21
24 : public FabFactory<FArrayBox>
25{
26public:
27
28 EBFArrayBoxFactory (const EB2::Level& a_level, const Geometry& a_geom,
29 const BoxArray& a_ba, const DistributionMapping& a_dm,
30 const Vector<int>& a_ngrow, EBSupport a_support);
31 ~EBFArrayBoxFactory () override = default;
32
34 EBFArrayBoxFactory (EBFArrayBoxFactory&&) noexcept = default;
35
36 EBFArrayBoxFactory () = delete;
37 EBFArrayBoxFactory& operator= (const EBFArrayBoxFactory&) = delete;
38 EBFArrayBoxFactory& operator= (EBFArrayBoxFactory&&) = delete;
39
41 FArrayBox* create (const Box& box, int ncomps, const FabInfo& info, int box_index) const final;
42
44 FArrayBox* create_alias (FArrayBox const& rhs, int scomp, int ncomp) const final;
45
46 void destroy (FArrayBox* fab) const final;
47
49 EBFArrayBoxFactory* clone () const final;
50
51 [[nodiscard]] const FabArray<EBCellFlagFab>& getMultiEBCellFlagFab () const noexcept
52 { return m_ebdc->getMultiEBCellFlagFab(); }
53
54 [[nodiscard]] const MultiFab& getLevelSet () const noexcept { return m_ebdc->getLevelSet(); }
55
56 [[nodiscard]] const MultiFab& getVolFrac () const noexcept { return m_ebdc->getVolFrac(); }
57
58 [[nodiscard]] const MultiCutFab& getCentroid () const noexcept { return m_ebdc->getCentroid(); }
59
60 [[nodiscard]] const MultiCutFab& getBndryCent () const noexcept { return m_ebdc->getBndryCent(); }
61
62 [[nodiscard]] const MultiCutFab& getBndryNormal () const noexcept { return m_ebdc->getBndryNormal(); }
63
64 [[nodiscard]] const MultiCutFab& getBndryArea () const noexcept { return m_ebdc->getBndryArea(); }
65
66 [[nodiscard]] Array<const MultiCutFab*,AMREX_SPACEDIM> getAreaFrac () const noexcept {
67 return m_ebdc->getAreaFrac();
68 }
69
70 // For y-face-centroid, the two components are for x and then z-directions, respectively.
71 [[nodiscard]] Array<const MultiCutFab*,AMREX_SPACEDIM> getFaceCent () const noexcept {
72 return m_ebdc->getFaceCent();
73 }
74
75 [[nodiscard]] Array<const MultiCutFab*,AMREX_SPACEDIM> getEdgeCent () const noexcept {
76 return m_ebdc->getEdgeCent();
77 }
78
79 [[nodiscard]] bool isAllRegular () const noexcept;
80
81 [[nodiscard]] EB2::Level const* getEBLevel () const noexcept { return m_parent; }
82 [[nodiscard]] EB2::IndexSpace const* getEBIndexSpace () const noexcept;
83 [[nodiscard]] int maxCoarseningLevel () const noexcept;
84
85 [[nodiscard]] const DistributionMapping& DistributionMap () const noexcept;
86 [[nodiscard]] const BoxArray& boxArray () const noexcept;
87 [[nodiscard]] const Geometry& Geom () const noexcept { return m_geom; }
88
89 [[nodiscard]] bool hasEBInfo() const noexcept;
90
93 [[nodiscard]] iMultiFab const* getCutCellMask () const noexcept { return m_ebdc->getCutCellMask(); }
94
95 [[nodiscard]] EBData getEBData (MFIter const& mfi) const noexcept;
96
97private:
98
99 EBSupport m_support;
100 Geometry m_geom;
101 std::shared_ptr<EBDataCollection> m_ebdc;
102 EB2::Level const* m_parent = nullptr;
104};
105
107std::unique_ptr<EBFArrayBoxFactory>
108makeEBFabFactory (const Geometry& a_geom,
109 const BoxArray& a_ba,
110 const DistributionMapping& a_dm,
111 const Vector<int>& a_ngrow, EBSupport a_support);
112
114std::unique_ptr<EBFArrayBoxFactory>
115makeEBFabFactory (const EB2::Level*,
116 const BoxArray& a_ba,
117 const DistributionMapping& a_dm,
118 const Vector<int>& a_ngrow, EBSupport a_support);
119
121std::unique_ptr<EBFArrayBoxFactory>
122makeEBFabFactory (const EB2::IndexSpace*, const Geometry& a_geom,
123 const BoxArray& a_ba,
124 const DistributionMapping& a_dm,
125 const Vector<int>& a_ngrow, EBSupport a_support);
126
127}
128
129#endif
#define AMREX_NODISCARD
Definition AMReX_Extension.H:251
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:567
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
Definition AMReX_EB2.H:28
Definition AMReX_EB2_Level.H:40
Definition AMReX_EBCellFlag.H:287
Definition AMReX_EBFabFactory.H:25
const MultiCutFab & getBndryArea() const noexcept
Definition AMReX_EBFabFactory.H:64
void destroy(FArrayBox *fab) const final
Definition AMReX_EBFabFactory.cpp:127
EBFArrayBoxFactory * clone() const final
Definition AMReX_EBFabFactory.cpp:142
EB2::IndexSpace const * getEBIndexSpace() const noexcept
Definition AMReX_EBFabFactory.cpp:154
EB2::Level const * getEBLevel() const noexcept
Definition AMReX_EBFabFactory.H:81
const MultiCutFab & getBndryNormal() const noexcept
Definition AMReX_EBFabFactory.H:62
const FabArray< EBCellFlagFab > & getMultiEBCellFlagFab() const noexcept
Definition AMReX_EBFabFactory.H:51
bool isAllRegular() const noexcept
Definition AMReX_EBFabFactory.cpp:148
const BoxArray & boxArray() const noexcept
Definition AMReX_EBFabFactory.cpp:177
EBData getEBData(MFIter const &mfi) const noexcept
Definition AMReX_EBFabFactory.cpp:189
iMultiFab const * getCutCellMask() const noexcept
Definition AMReX_EBFabFactory.H:93
int maxCoarseningLevel() const noexcept
Definition AMReX_EBFabFactory.cpp:160
const DistributionMapping & DistributionMap() const noexcept
Definition AMReX_EBFabFactory.cpp:171
FArrayBox * create_alias(FArrayBox const &rhs, int scomp, int ncomp) const final
Definition AMReX_EBFabFactory.cpp:113
const MultiCutFab & getCentroid() const noexcept
Definition AMReX_EBFabFactory.H:58
Array< const MultiCutFab *, 3 > getEdgeCent() const noexcept
Definition AMReX_EBFabFactory.H:75
const MultiCutFab & getBndryCent() const noexcept
Definition AMReX_EBFabFactory.H:60
EBFArrayBoxFactory(EBFArrayBoxFactory &&) noexcept=default
Array< const MultiCutFab *, 3 > getAreaFrac() const noexcept
Definition AMReX_EBFabFactory.H:66
Array< const MultiCutFab *, 3 > getFaceCent() const noexcept
Definition AMReX_EBFabFactory.H:71
bool hasEBInfo() const noexcept
Definition AMReX_EBFabFactory.cpp:183
~EBFArrayBoxFactory() override=default
const MultiFab & getLevelSet() const noexcept
Definition AMReX_EBFabFactory.H:54
const MultiFab & getVolFrac() const noexcept
Definition AMReX_EBFabFactory.H:56
const Geometry & Geom() const noexcept
Definition AMReX_EBFabFactory.H:87
FArrayBox * create(const Box &box, int ncomps, const FabInfo &info, int box_index) const final
Definition AMReX_EBFabFactory.cpp:97
EBFArrayBoxFactory(const EBFArrayBoxFactory &)=default
A Fortran Array of REALs.
Definition AMReX_FArrayBox.H:231
An Array of FortranArrayBox(FAB)-like Objects.
Definition AMReX_FabArray.H:347
Definition AMReX_FabFactory.H:50
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:74
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:85
Definition AMReX_MultiCutFab.H:81
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:40
Dynamically allocated vector for trivially copyable data.
Definition AMReX_PODVector.H:308
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
A Collection of IArrayBoxes.
Definition AMReX_iMultiFab.H:34
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:202
std::array< T, N > Array
Definition AMReX_Array.H:26
Definition AMReX_Amr.cpp:49
EBSupport
Definition AMReX_EBSupport.H:7
Definition AMReX_EBData.H:26
Definition AMReX_FabFactory.H:27