ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_EB.H
Go to the documentation of this file.
1 #ifndef ERF_EB_H_
2 #define ERF_EB_H_
3 
4 #include <AMReX_Geometry.H>
5 #include <AMReX_DistributionMapping.H>
6 #include <AMReX_BoxArray.H>
7 
8 #include <AMReX_EB2.H>
9 #include <AMReX_EBFabFactory.H>
10 
11 #include <ERF_EBAux.H>
12 
13 class eb_ {
14 
15  public:
16 
17  ~eb_ ();
18 
19  eb_ (amrex::Geometry const& a_geom,
20  amrex::FArrayBox const& terrain_fab,
21  amrex::Gpu::DeviceVector<amrex::Real>& a_dz_stretched,
22  bool is_anelastic);
23  eb_ ();
24 
25  void define (int level,
26  amrex::Geometry const& a_geom,
27  amrex::EB2::Level const* a_eb_level,
28  bool is_anelastic);
29 
30  void make_all_factories ( int level,
31  amrex::Geometry const& a_geom,
32  amrex::BoxArray const& ba,
33  amrex::DistributionMapping const& dm,
34  amrex::EB2::Level const& a_eb_level);
35 
36  void make_cc_factory ( int level,
37  amrex::Geometry const& a_geom,
38  amrex::BoxArray const& ba,
39  amrex::DistributionMapping const& dm,
40  amrex::EB2::Level const& a_eb_level);
41 
42  int nghost_basic () const { return 5; } // nghost_eb_basic ()
43  int nghost_volume () const { return 5; } // nghost_eb_volume ()
44  int nghost_full () const { return 4; } // nghost_eb_full ()
45 
46  const std::unique_ptr<amrex::EBFArrayBoxFactory>& get_const_factory () const noexcept { return m_factory; }
47 
48  void set_connection_flags ();
49 
50  eb_aux_ const* get_u_const_factory() const noexcept { return &m_u_factory; }
51  eb_aux_ const* get_v_const_factory() const noexcept { return &m_v_factory; }
52  eb_aux_ const* get_w_const_factory() const noexcept { return &m_w_factory; }
53 
54  class EBToPVD;
55 
56  private:
57 
58  int m_has_eb;
59  std::string m_type;
60 
61  amrex::EBSupport m_support_level;
62 
64 
65  amrex::FabArray<amrex::EBCellFlagFab>* m_cellflags = nullptr;
66 
67  //! EB level constructed from building GeometryShop
68  amrex::EB2::Level const* m_eb_level;
69 
70  std::unique_ptr<amrex::EBFArrayBoxFactory> m_factory = nullptr;
71 
75 
76  void make_terrain (amrex::Geometry const& a_geom);
77 
78  //! Construct EB levels from Geometry shop.
79  template<class F>
80  void build_level (amrex::Geometry const& a_geom,
81  amrex::EB2::GeometryShop<F> a_gshop)
82  {
83  int const req_lev(0);
84  int const max_lev(2);
85 
86  amrex::EB2::Build(a_gshop, a_geom, req_lev, max_lev);
87  const amrex::EB2::IndexSpace& ebis = amrex::EB2::IndexSpace::top();
88  m_eb_level = &(ebis.getLevel(a_geom));
89  }
90 
91  // Wrapper to get non-const version of EBCellFlagFab using const_cast
92  inline amrex::FabArray<amrex::EBCellFlagFab>&
93  getNonConstEBCellFlags(const amrex::EBFArrayBoxFactory& ebfact)
94  {
95  const amrex::FabArray<amrex::EBCellFlagFab>& flags_const = ebfact.getMultiEBCellFlagFab();
96  return const_cast<amrex::FabArray<amrex::EBCellFlagFab>&>(flags_const);
97  }
98 
99  inline amrex::MultiFab&
100  getNonConstVolFrac(const amrex::EBFArrayBoxFactory& ebfact)
101  {
102  const amrex::MultiFab& vfrac_const = ebfact.getVolFrac();
103  return const_cast<amrex::MultiFab&>(vfrac_const);
104  }
105 
106  inline amrex::MultiCutFab&
107  getNonConstCentroid(const amrex::EBFArrayBoxFactory& ebfact)
108  {
109  const amrex::MultiCutFab& vcent_const = ebfact.getCentroid();
110  return const_cast<amrex::MultiCutFab&>(vcent_const);
111  }
112 
113  inline amrex::Array<amrex::MultiCutFab*, AMREX_SPACEDIM>
114  getNonConstAreaFrac(const amrex::EBFArrayBoxFactory& ebfact)
115  {
116  auto afrac_const = ebfact.getAreaFrac();
117  amrex::Array<amrex::MultiCutFab*, AMREX_SPACEDIM> afrac;
118  for (int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
119  afrac[dir] = const_cast<amrex::MultiCutFab*>(afrac_const[dir]);
120  }
121  return afrac;
122  }
123 
124  inline amrex::Array<amrex::MultiCutFab*, AMREX_SPACEDIM>
125  getNonConstFaceCent(const amrex::EBFArrayBoxFactory& ebfact)
126  {
127  auto fcent_const = ebfact.getFaceCent();
128  amrex::Array<amrex::MultiCutFab*, AMREX_SPACEDIM> fcent;
129  for (int dir = 0; dir < AMREX_SPACEDIM; ++dir) {
130  fcent[dir] = const_cast<amrex::MultiCutFab*>(fcent_const[dir]);
131  }
132  return fcent;
133  }
134 
135 };
136 #endif
Definition: ERF_EB.H:13
amrex::FabArray< amrex::EBCellFlagFab > & getNonConstEBCellFlags(const amrex::EBFArrayBoxFactory &ebfact)
Definition: ERF_EB.H:93
eb_()
Definition: ERF_EB.cpp:26
void define(int level, amrex::Geometry const &a_geom, amrex::EB2::Level const *a_eb_level, bool is_anelastic)
void make_terrain(amrex::Geometry const &a_geom)
eb_aux_ const * get_w_const_factory() const noexcept
Definition: ERF_EB.H:52
void build_level(amrex::Geometry const &a_geom, amrex::EB2::GeometryShop< F > a_gshop)
Construct EB levels from Geometry shop.
Definition: ERF_EB.H:80
amrex::Array< amrex::MultiCutFab *, AMREX_SPACEDIM > getNonConstAreaFrac(const amrex::EBFArrayBoxFactory &ebfact)
Definition: ERF_EB.H:114
amrex::EB2::Level const * m_eb_level
EB level constructed from building GeometryShop.
Definition: ERF_EB.H:68
void make_all_factories(int level, amrex::Geometry const &a_geom, amrex::BoxArray const &ba, amrex::DistributionMapping const &dm, amrex::EB2::Level const &a_eb_level)
Definition: ERF_EB.cpp:33
amrex::FabArray< amrex::EBCellFlagFab > * m_cellflags
Definition: ERF_EB.H:65
const std::unique_ptr< amrex::EBFArrayBoxFactory > & get_const_factory() const noexcept
Definition: ERF_EB.H:46
amrex::MultiFab & getNonConstVolFrac(const amrex::EBFArrayBoxFactory &ebfact)
Definition: ERF_EB.H:100
eb_aux_ m_w_factory
Definition: ERF_EB.H:74
std::unique_ptr< amrex::EBFArrayBoxFactory > m_factory
Definition: ERF_EB.H:70
amrex::EBSupport m_support_level
Definition: ERF_EB.H:61
int nghost_volume() const
Definition: ERF_EB.H:43
amrex::Array< amrex::MultiCutFab *, AMREX_SPACEDIM > getNonConstFaceCent(const amrex::EBFArrayBoxFactory &ebfact)
Definition: ERF_EB.H:125
void make_cc_factory(int level, amrex::Geometry const &a_geom, amrex::BoxArray const &ba, amrex::DistributionMapping const &dm, amrex::EB2::Level const &a_eb_level)
Definition: ERF_EB.cpp:73
amrex::MultiCutFab & getNonConstCentroid(const amrex::EBFArrayBoxFactory &ebfact)
Definition: ERF_EB.H:107
eb_aux_ m_u_factory
Definition: ERF_EB.H:72
~eb_()
Definition: ERF_EB.cpp:21
eb_(amrex::Geometry const &a_geom, amrex::FArrayBox const &terrain_fab, amrex::Gpu::DeviceVector< amrex::Real > &a_dz_stretched, bool is_anelastic)
eb_aux_ const * get_v_const_factory() const noexcept
Definition: ERF_EB.H:51
eb_aux_ m_v_factory
Definition: ERF_EB.H:73
int m_has_eb
Definition: ERF_EB.H:54
int nghost_full() const
Definition: ERF_EB.H:44
int nghost_basic() const
Definition: ERF_EB.H:42
int m_write_eb_surface
Definition: ERF_EB.H:63
eb_aux_ const * get_u_const_factory() const noexcept
Definition: ERF_EB.H:50
std::string m_type
Definition: ERF_EB.H:59
void set_connection_flags()
Definition: ERF_EB.cpp:91
Definition: ERF_EBAux.H:12