ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
eb_ Class Reference

#include <ERF_EB.H>

Collaboration diagram for eb_:

Public Member Functions

 ~eb_ ()
 
 eb_ (amrex::Geometry const &a_geom, amrex::FArrayBox const &terrain_fab, amrex::Gpu::DeviceVector< amrex::Real > &a_dz_stretched, bool is_anelastic)
 
 eb_ ()
 
void define (int level, amrex::Geometry const &a_geom, amrex::EB2::Level const *a_eb_level, bool is_anelastic)
 
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)
 
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)
 
int nghost_basic () const
 
int nghost_volume () const
 
int nghost_full () const
 
const std::unique_ptr< amrex::EBFArrayBoxFactory > & get_const_factory () const noexcept
 
void set_connection_flags ()
 
eb_aux_ const * get_u_const_factory () const noexcept
 
eb_aux_ const * get_v_const_factory () const noexcept
 
eb_aux_ const * get_w_const_factory () const noexcept
 

Private Member Functions

void make_terrain (amrex::Geometry const &a_geom)
 
template<class F >
void build_level (amrex::Geometry const &a_geom, amrex::EB2::GeometryShop< F > a_gshop)
 Construct EB levels from Geometry shop. More...
 
amrex::FabArray< amrex::EBCellFlagFab > & getNonConstEBCellFlags (const amrex::EBFArrayBoxFactory &ebfact)
 
amrex::MultiFab & getNonConstVolFrac (const amrex::EBFArrayBoxFactory &ebfact)
 
amrex::MultiCutFab & getNonConstCentroid (const amrex::EBFArrayBoxFactory &ebfact)
 
amrex::Array< amrex::MultiCutFab *, AMREX_SPACEDIM > getNonConstAreaFrac (const amrex::EBFArrayBoxFactory &ebfact)
 
amrex::Array< amrex::MultiCutFab *, AMREX_SPACEDIM > getNonConstFaceCent (const amrex::EBFArrayBoxFactory &ebfact)
 

Private Attributes

int m_has_eb
 
std::string m_type
 
amrex::EBSupport m_support_level
 
int m_write_eb_surface
 
amrex::FabArray< amrex::EBCellFlagFab > * m_cellflags = nullptr
 
amrex::EB2::Level const * m_eb_level
 EB level constructed from building GeometryShop. More...
 
std::unique_ptr< amrex::EBFArrayBoxFactory > m_factory = nullptr
 
eb_aux_ m_u_factory
 
eb_aux_ m_v_factory
 
eb_aux_ m_w_factory
 

Constructor & Destructor Documentation

◆ ~eb_()

eb_::~eb_ ( )
22 {
23  // if (m_factory) { m_factory.reset(nullptr); }
24 }

◆ eb_() [1/2]

eb_::eb_ ( amrex::Geometry const &  a_geom,
amrex::FArrayBox const &  terrain_fab,
amrex::Gpu::DeviceVector< amrex::Real > &  a_dz_stretched,
bool  is_anelastic 
)

◆ eb_() [2/2]

eb_::eb_ ( )
27  : m_has_eb(0),
28  m_support_level(EBSupport::full),
30 { }
amrex::EBSupport m_support_level
Definition: ERF_EB.H:61
int m_has_eb
Definition: ERF_EB.H:54
int m_write_eb_surface
Definition: ERF_EB.H:63

Member Function Documentation

◆ build_level()

template<class F >
void eb_::build_level ( amrex::Geometry const &  a_geom,
amrex::EB2::GeometryShop< F >  a_gshop 
)
inlineprivate

Construct EB levels from Geometry shop.

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  }
amrex::EB2::Level const * m_eb_level
EB level constructed from building GeometryShop.
Definition: ERF_EB.H:68

◆ define()

void eb_::define ( int  level,
amrex::Geometry const &  a_geom,
amrex::EB2::Level const *  a_eb_level,
bool  is_anelastic 
)

◆ get_const_factory()

const std::unique_ptr<amrex::EBFArrayBoxFactory>& eb_::get_const_factory ( ) const
inlinenoexcept
46 { return m_factory; }
std::unique_ptr< amrex::EBFArrayBoxFactory > m_factory
Definition: ERF_EB.H:70

Referenced by compute_gradp(), erf_slow_rhs_pre(), ERF::estTimeStep(), and make_buoyancy().

Here is the caller graph for this function:

◆ get_u_const_factory()

eb_aux_ const* eb_::get_u_const_factory ( ) const
inlinenoexcept
50 { return &m_u_factory; }
eb_aux_ m_u_factory
Definition: ERF_EB.H:72

Referenced by AdvectionSrcForMom_EB(), compute_gradp(), and DiffusionSrcForMom_EB().

Here is the caller graph for this function:

◆ get_v_const_factory()

eb_aux_ const* eb_::get_v_const_factory ( ) const
inlinenoexcept
51 { return &m_v_factory; }
eb_aux_ m_v_factory
Definition: ERF_EB.H:73

Referenced by AdvectionSrcForMom_EB(), compute_gradp(), and DiffusionSrcForMom_EB().

Here is the caller graph for this function:

◆ get_w_const_factory()

eb_aux_ const* eb_::get_w_const_factory ( ) const
inlinenoexcept
52 { return &m_w_factory; }
eb_aux_ m_w_factory
Definition: ERF_EB.H:74

Referenced by AdvectionSrcForMom_EB(), compute_gradp(), and DiffusionSrcForMom_EB().

Here is the caller graph for this function:

◆ getNonConstAreaFrac()

amrex::Array<amrex::MultiCutFab*, AMREX_SPACEDIM> eb_::getNonConstAreaFrac ( const amrex::EBFArrayBoxFactory &  ebfact)
inlineprivate
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  }

◆ getNonConstCentroid()

amrex::MultiCutFab& eb_::getNonConstCentroid ( const amrex::EBFArrayBoxFactory &  ebfact)
inlineprivate
108  {
109  const amrex::MultiCutFab& vcent_const = ebfact.getCentroid();
110  return const_cast<amrex::MultiCutFab&>(vcent_const);
111  }

◆ getNonConstEBCellFlags()

amrex::FabArray<amrex::EBCellFlagFab>& eb_::getNonConstEBCellFlags ( const amrex::EBFArrayBoxFactory &  ebfact)
inlineprivate
94  {
95  const amrex::FabArray<amrex::EBCellFlagFab>& flags_const = ebfact.getMultiEBCellFlagFab();
96  return const_cast<amrex::FabArray<amrex::EBCellFlagFab>&>(flags_const);
97  }

Referenced by set_connection_flags().

Here is the caller graph for this function:

◆ getNonConstFaceCent()

amrex::Array<amrex::MultiCutFab*, AMREX_SPACEDIM> eb_::getNonConstFaceCent ( const amrex::EBFArrayBoxFactory &  ebfact)
inlineprivate
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  }

◆ getNonConstVolFrac()

amrex::MultiFab& eb_::getNonConstVolFrac ( const amrex::EBFArrayBoxFactory &  ebfact)
inlineprivate
101  {
102  const amrex::MultiFab& vfrac_const = ebfact.getVolFrac();
103  return const_cast<amrex::MultiFab&>(vfrac_const);
104  }

◆ make_all_factories()

void eb_::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 
)
38 {
39  Print() << "making EB factory\n";
40  m_factory = std::make_unique<EBFArrayBoxFactory>(a_eb_level, a_geom, ba, dm,
42 
43  // Correct cell connectivity
45 
46  { int const idim(0);
47  Print() << "making EB staggered u-factory\n";
48  //m_u_factory.set_verbose();
49  m_u_factory.define(level, idim, a_geom, ba, dm,
50  Vector<int>{nghost_basic(), nghost_volume(), nghost_full()},
51  m_factory.get());
52  }
53 
54  { int const idim(1);
55  Print() << "making EB staggered v-factory\n";
56  //m_v_factory.set_verbose();
57  m_v_factory.define(level, idim, a_geom, ba, dm,
58  Vector<int>{nghost_basic(), nghost_volume(), nghost_full()},
59  m_factory.get());
60  }
61 
62  { int const idim(2);
63  Print() << "making EB staggered w-factory\n";
64  //m_w_factory.set_verbose();
65  m_w_factory.define(level, idim, a_geom, ba, dm,
66  Vector<int>{nghost_basic(), nghost_volume(), nghost_full()},
67  m_factory.get());
68  }
69  Print() << "\nDone making EB factory at level = " << level << ".\n\n";
70 }
int nghost_volume() const
Definition: ERF_EB.H:43
int nghost_full() const
Definition: ERF_EB.H:44
int nghost_basic() const
Definition: ERF_EB.H:42
void set_connection_flags()
Definition: ERF_EB.cpp:91
void define(int const &a_level, int const &a_idim, amrex::Geometry const &a_geom, amrex::BoxArray const &a_grids, amrex::DistributionMapping const &a_dmap, amrex::Vector< int > const &a_ngrow, amrex::EBFArrayBoxFactory const *a_factory)
Definition: ERF_EBAux.cpp:25
Here is the call graph for this function:

◆ make_cc_factory()

void eb_::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 
)
78 {
79  Print() << "making EB factory\n";
80  m_factory = std::make_unique<EBFArrayBoxFactory>(a_eb_level, a_geom, ba, dm,
82 
83  Print() << "\nDone making EB factory at level " << level << ".\n\n";
84 }
Here is the call graph for this function:

◆ make_terrain()

void eb_::make_terrain ( amrex::Geometry const &  a_geom)
private

◆ nghost_basic()

int eb_::nghost_basic ( ) const
inline
42 { return 5; } // nghost_eb_basic ()

Referenced by make_all_factories(), and make_cc_factory().

Here is the caller graph for this function:

◆ nghost_full()

int eb_::nghost_full ( ) const
inline
44 { return 4; } // nghost_eb_full ()

Referenced by make_all_factories(), and make_cc_factory().

Here is the caller graph for this function:

◆ nghost_volume()

int eb_::nghost_volume ( ) const
inline
43 { return 5; } // nghost_eb_volume ()

Referenced by make_all_factories(), and make_cc_factory().

Here is the caller graph for this function:

◆ set_connection_flags()

void eb_::set_connection_flags ( )
92 {
93  // Get non-const reference to EBCellFlagFab FabArray
94  FabArray<EBCellFlagFab>& cellflag = getNonConstEBCellFlags(*m_factory);
95 
96  const MultiFab& volfrac = m_factory->getVolFrac();
97 
98  for (MFIter mfi(cellflag, false); mfi.isValid(); ++mfi) {
99  const Box& bx = mfi.validbox();
100  const Box gbx = amrex::grow(bx, cellflag.nGrow()-1); // Leave one cell layer
101 
102  Array4<EBCellFlag> const& flag = cellflag.array(mfi);
103  Array4<Real const> const& vfrac = volfrac.const_array(mfi);
104 
105  ParallelFor(gbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
106  {
107  for(int kk(-1); kk<=1; kk++) {
108  for(int jj(-1); jj<=1; jj++) {
109  for(int ii(-1); ii<=1; ii++)
110  {
111  if (vfrac(i+ii,j+jj,k+kk) == 0.0) {
112  flag(i,j,k).setDisconnected(ii,jj,kk);
113  }
114  }}}
115  });
116 
117  ParallelFor(gbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
118  {
119  if (vfrac(i,j,k)==0.0) {
120  flag(i,j,k).setCovered();
121  }
122  });
123 
124  }
125 }
amrex::FabArray< amrex::EBCellFlagFab > & getNonConstEBCellFlags(const amrex::EBFArrayBoxFactory &ebfact)
Definition: ERF_EB.H:93

Referenced by make_all_factories().

Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ m_cellflags

amrex::FabArray<amrex::EBCellFlagFab>* eb_::m_cellflags = nullptr
private

◆ m_eb_level

amrex::EB2::Level const* eb_::m_eb_level
private

EB level constructed from building GeometryShop.

Referenced by build_level().

◆ m_factory

std::unique_ptr<amrex::EBFArrayBoxFactory> eb_::m_factory = nullptr
private

◆ m_has_eb

int eb_::m_has_eb
private

◆ m_support_level

amrex::EBSupport eb_::m_support_level
private

◆ m_type

std::string eb_::m_type
private

◆ m_u_factory

eb_aux_ eb_::m_u_factory
private

◆ m_v_factory

eb_aux_ eb_::m_v_factory
private

◆ m_w_factory

eb_aux_ eb_::m_w_factory
private

◆ m_write_eb_surface

int eb_::m_write_eb_surface
private

The documentation for this class was generated from the following files: