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)
 
void make_factory (amrex::Geometry const &a_geom, amrex::DistributionMapping const &a_dmap, amrex::BoxArray const &a_grids)
 
int nghost_basic () const
 
int nghost_volume () const
 
int nghost_full () const
 
amrex::EBFArrayBoxFactory const * get_const_factory () noexcept
 
amrex::EB2::Level const * get_level () const noexcept
 

Private Member Functions

void make_box (amrex::Geometry const &a_geom)
 
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...
 

Private Attributes

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

Constructor & Destructor Documentation

◆ ~eb_()

eb_::~eb_ ( )
13 {
14  if (m_factory) { m_factory.reset(nullptr); }
15 }
std::unique_ptr< amrex::EBFArrayBoxFactory > m_factory
Definition: ERF_EB.H:50

◆ eb_()

eb_::eb_ ( amrex::Geometry const &  a_geom,
amrex::FArrayBox const &  terrain_fab,
amrex::Gpu::DeviceVector< amrex::Real > &  a_dz_stretched,
bool  is_anelastic 
)
20  : m_has_eb(0),
21  m_support_level(EBSupport::full),
23 {
24  m_type = "terrain";
25 
26  int max_coarsening_level;
27  if (is_anelastic) {
28  max_coarsening_level = 100;
29  } else {
30  max_coarsening_level = 0;
31  }
32 
33  int max_level_here = 0;
34 
35  if (m_type == "terrain")
36  {
37  Print() << "\nBuilding EB geometry based on terrain.\n";
38 
39  TerrainIF ebterrain(terrain_fab, a_geom, a_dz_stretched);
40 
41  auto gshop = EB2::makeShop(ebterrain);
42 
43  EB2::Build(gshop, a_geom, max_level_here, max_level_here+max_coarsening_level);
44 
45  m_has_eb = 1;
46 
47  } else if (m_type == "box") {
48 
49  Print() << "\nBuilding box geometry.\n";
50  make_box(a_geom);
51  m_has_eb = 1;
52 
53  } else {
54 
55  EB2::AllRegularIF regular;
56  auto gshop = EB2::makeShop(regular);
57  build_level(a_geom, gshop);
58  }
59 }
Definition: ERF_EBIFTerrain.H:14
void build_level(amrex::Geometry const &a_geom, amrex::EB2::GeometryShop< F > a_gshop)
Construct EB levels from Geometry shop.
Definition: ERF_EB.H:61
amrex::EBSupport m_support_level
Definition: ERF_EB.H:43
int m_has_eb
Definition: ERF_EB.H:40
int m_write_eb_surface
Definition: ERF_EB.H:45
void make_box(amrex::Geometry const &a_geom)
Definition: ERF_EBBox.cpp:22
std::string m_type
Definition: ERF_EB.H:41
Here is the call graph for this function:

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.

63  {
64  int const req_lev(0);
65  int const max_lev(2);
66 
67  amrex::EB2::Build(a_gshop, a_geom, req_lev, max_lev);
68  const amrex::EB2::IndexSpace& ebis = amrex::EB2::IndexSpace::top();
69  m_eb_level = &(ebis.getLevel(a_geom));
70  }
amrex::EB2::Level const * m_eb_level
EB level constructed from building GeometryShop.
Definition: ERF_EB.H:48

Referenced by eb_().

Here is the caller graph for this function:

◆ get_const_factory()

amrex::EBFArrayBoxFactory const* eb_::get_const_factory ( )
inlinenoexcept
33  { return m_factory.get(); }

◆ get_level()

amrex::EB2::Level const* eb_::get_level ( ) const
inlinenoexcept
36  { return m_eb_level; }

◆ make_box()

void eb_::make_box ( amrex::Geometry const &  a_geom)
private
23 {
24  ParmParse pp_box("eb.box");
25 
26  bool inside = true;
27 
28  Vector<Real> boxLo(AMREX_SPACEDIM), boxHi(AMREX_SPACEDIM);
29 
30  Real const* dx = a_geom.CellSize();
31 
32  if ( !almostEqual(dx[0],dx[1]) || !almostEqual(dx[1],dx[2])) {
33  amrex::Error(" EB Error: Mesh spacing must be uniform!.\n");
34  }
35 
36  Real offset = 0.01*dx[0];
37 
38  for (int i = 0; i < AMREX_SPACEDIM; i++) {
39  boxLo[i] = a_geom.ProbLo(i);
40  boxHi[i] = a_geom.ProbHi(i);
41  }
42 
43  pp_box.queryarr("lo", boxLo, 0, AMREX_SPACEDIM);
44  pp_box.queryarr("hi", boxHi, 0, AMREX_SPACEDIM);
45 
46  pp_box.query("internal_flow", inside);
47  pp_box.query("offset", offset);
48 
49  Real xlo = boxLo[0] + offset;
50  Real xhi = boxHi[0] - offset;
51 
52  // This ensures that the walls won't even touch the ghost cells. By
53  // putting them one domain width away
54  if (a_geom.isPeriodic(0))
55  {
56  xlo = 2.0*a_geom.ProbLo(0) - a_geom.ProbHi(0);
57  xhi = 2.0*a_geom.ProbHi(0) - a_geom.ProbLo(0);
58  }
59 
60 
61  Real ylo = boxLo[1] + offset;
62  Real yhi = boxHi[1] - offset;
63 
64  // This ensures that the walls won't even touch the ghost cells. By
65  // putting them one domain width away
66  if (a_geom.isPeriodic(1))
67  {
68  ylo = 2.0*a_geom.ProbLo(1) - a_geom.ProbHi(1);
69  yhi = 2.0*a_geom.ProbHi(1) - a_geom.ProbLo(1);
70  }
71 
72  Real zlo = boxLo[2] + offset;
73  Real zhi = boxHi[2] - offset;
74 
75  // This ensures that the walls won't even touch the ghost cells. By
76  // putting them one domain width away
77  if (a_geom.isPeriodic(2))
78  {
79  zlo = 2.0*a_geom.ProbLo(2) - a_geom.ProbHi(2);
80  zhi = 2.0*a_geom.ProbHi(2) - a_geom.ProbLo(2);
81  }
82 
83  Array<Real,3> point_lox{ xlo, 0.0, 0.0};
84  Array<Real,3> normal_lox{-1.0, 0.0, 0.0};
85  Array<Real,3> point_hix{ xhi, 0.0, 0.0};
86  Array<Real,3> normal_hix{ 1.0, 0.0, 0.0};
87 
88  Array<Real,3> point_loy{0.0, ylo, 0.0};
89  Array<Real,3> normal_loy{0.0,-1.0, 0.0};
90  Array<Real,3> point_hiy{0.0, yhi, 0.0};
91  Array<Real,3> normal_hiy{0.0, 1.0, 0.0};
92 
93  Array<Real,3> point_loz{0.0, 0.0, zlo};
94  Array<Real,3> normal_loz{0.0, 0.0,-1.0};
95  Array<Real,3> point_hiz{0.0, 0.0, zhi};
96  Array<Real,3> normal_hiz{0.0, 0.0, 1.0};
97 
98  EB2::PlaneIF plane_lox(point_lox,normal_lox);
99  EB2::PlaneIF plane_hix(point_hix,normal_hix);
100 
101  EB2::PlaneIF plane_loy(point_loy,normal_loy);
102  EB2::PlaneIF plane_hiy(point_hiy,normal_hiy);
103 
104  EB2::PlaneIF plane_loz(point_loz,normal_loz);
105  EB2::PlaneIF plane_hiz(point_hiz,normal_hiz);
106 
107  auto box = EB2::makeUnion(plane_lox, plane_hix, plane_loy,
108  plane_hiy, plane_loz, plane_hiz );
109 
110 
111  if( inside ) {
112 
113  auto gshop = EB2::makeShop(box);
114  build_level(a_geom, gshop);
115 
116  } else {
117 
118  auto xob = EB2::makeComplement(box);
119  auto gshop = EB2::makeShop(xob);
120 
121  build_level(a_geom, gshop);
122  }
123 
124 }
AMREX_FORCE_INLINE IntVect offset(const int face_dir, const int normal)
Definition: ERF_ReadBndryPlanes.cpp:28

Referenced by eb_().

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

◆ make_factory()

void eb_::make_factory ( amrex::Geometry const &  a_geom,
amrex::DistributionMapping const &  a_dmap,
amrex::BoxArray const &  a_grids 
)
66 {
67  Print() << "making EB factory\n";
68  m_factory = std::make_unique<EBFArrayBoxFactory>(*m_eb_level, a_geom, a_grids, a_dmap,
70 
71  if (m_write_eb_surface) {
72  WriteEBSurface(a_grids, a_dmap, a_geom, m_factory.get());
73  }
74 
75  { int const idim(0);
76 
77  Print() << "making EB staggered u-factory\n";
78  //m_u_factory.set_verbose();
79  m_u_factory.define(idim, a_geom, a_grids, a_dmap,
80  Vector<int>{nghost_basic(), nghost_volume(), nghost_full()},
81  m_factory.get());
82  }
83 
84  { int const idim(1);
85  Print() << "making EB staggered v-factory\n";
86  //m_v_factory.set_verbose();
87  m_v_factory.define(idim, a_geom, a_grids, a_dmap,
88  Vector<int>{nghost_basic(), nghost_volume(), nghost_full()},
89  m_factory.get());
90  }
91 
92  { int const idim(2);
93  Print() << "making EB staggered w-factory\n";
94  //m_w_factory.set_verbose();
95  m_w_factory.define(idim, a_geom, a_grids, a_dmap,
96  Vector<int>{nghost_basic(), nghost_volume(), nghost_full()},
97  m_factory.get());
98  }
99 
100  Print() << "\nDone making EB factory.\n\n";
101 }
eb_aux_ m_w_factory
Definition: ERF_EB.H:54
int nghost_volume() const
Definition: ERF_EB.H:29
eb_aux_ m_u_factory
Definition: ERF_EB.H:52
eb_aux_ m_v_factory
Definition: ERF_EB.H:53
int nghost_full() const
Definition: ERF_EB.H:30
int nghost_basic() const
Definition: ERF_EB.H:28
void define(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:19
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
28 { return 2; }

Referenced by make_factory().

Here is the caller graph for this function:

◆ nghost_full()

int eb_::nghost_full ( ) const
inline
30 { return 2; }

Referenced by make_factory().

Here is the caller graph for this function:

◆ nghost_volume()

int eb_::nghost_volume ( ) const
inline
29 { return 2; }

Referenced by make_factory().

Here is the caller graph for this function:

Member Data Documentation

◆ m_eb_level

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

EB level constructed from building GeometryShop.

Referenced by build_level(), get_level(), and make_factory().

◆ m_factory

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

Referenced by get_const_factory(), and make_factory().

◆ m_has_eb

int eb_::m_has_eb
private

Referenced by eb_().

◆ m_support_level

amrex::EBSupport eb_::m_support_level
private

Referenced by make_factory().

◆ m_type

std::string eb_::m_type
private

Referenced by eb_().

◆ m_u_factory

eb_aux_ eb_::m_u_factory
private

Referenced by make_factory().

◆ m_v_factory

eb_aux_ eb_::m_v_factory
private

Referenced by make_factory().

◆ m_w_factory

eb_aux_ eb_::m_w_factory
private

Referenced by make_factory().

◆ m_write_eb_surface

int eb_::m_write_eb_surface
private

Referenced by make_factory().


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