ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
eb_ Class Reference

#include <ERF_EB.H>

Collaboration diagram for eb_:

Classes

class  EBToPVD
 

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_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 WriteEBSurface (const amrex::BoxArray &ba, const amrex::DistributionMapping &dmap, const amrex::Geometry &geom, const amrex::EBFArrayBoxFactory *ebf, const int level)
 
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_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::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
 
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:59
int m_has_eb
Definition: ERF_EB.H:52
int m_write_eb_surface
Definition: ERF_EB.H:61

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.

81  {
82  int const req_lev(0);
83  int const max_lev(2);
84 
85  amrex::EB2::Build(a_gshop, a_geom, req_lev, max_lev);
86  const amrex::EB2::IndexSpace& ebis = amrex::EB2::IndexSpace::top();
87  m_eb_level = &(ebis.getLevel(a_geom));
88  }
amrex::EB2::Level const * m_eb_level
EB level constructed from building GeometryShop.
Definition: ERF_EB.H:66

◆ 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
40 { return m_factory; }
std::unique_ptr< amrex::EBFArrayBoxFactory > m_factory
Definition: ERF_EB.H:68

Referenced by ERF::estTimeStep().

Here is the caller graph for this function:

◆ get_u_const_factory()

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

◆ get_v_const_factory()

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

◆ get_w_const_factory()

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

◆ 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
void build_level(amrex::Geometry const &a_geom, amrex::EB2::GeometryShop< F > a_gshop)
Construct EB levels from Geometry shop.
Definition: ERF_EB.H:79
Here is the call graph for this function:

◆ make_factory()

void eb_::make_factory ( int  level,
amrex::Geometry const &  a_geom,
amrex::BoxArray const &  ba,
amrex::DistributionMapping const &  dm,
amrex::EB2::Level const &  a_eb_level 
)
38 {
39 
40  Print() << "making EB factory\n";
41  m_factory = std::make_unique<EBFArrayBoxFactory>(a_eb_level, a_geom, ba, dm,
43 
44  eb_::WriteEBSurface(ba, dm, a_geom, m_factory.get(), level);
45 
46  { int const idim(0);
47 
48  Print() << "making EB staggered u-factory\n";
49  //m_u_factory.set_verbose();
50  m_u_factory.define(idim, a_geom, ba, dm,
51  Vector<int>{nghost_basic(), nghost_volume(), nghost_full()},
52  m_factory.get());
53  }
54 
55  { int const idim(1);
56  Print() << "making EB staggered v-factory\n";
57  //m_v_factory.set_verbose();
58  m_v_factory.define(idim, a_geom, ba, dm,
59  Vector<int>{nghost_basic(), nghost_volume(), nghost_full()},
60  m_factory.get());
61  }
62 
63  { int const idim(2);
64  Print() << "making EB staggered w-factory\n";
65  //m_w_factory.set_verbose();
66  m_w_factory.define(idim, a_geom, ba, dm,
67  Vector<int>{nghost_basic(), nghost_volume(), nghost_full()},
68  m_factory.get());
69  }
70 
71  Print() << "\nDone making EB factory.\n\n";
72 }
void WriteEBSurface(const amrex::BoxArray &ba, const amrex::DistributionMapping &dmap, const amrex::Geometry &geom, const amrex::EBFArrayBoxFactory *ebf, const int level)
Definition: ERF_EB.cpp:77
int nghost_volume() const
Definition: ERF_EB.H:37
int nghost_full() const
Definition: ERF_EB.H:38
int nghost_basic() const
Definition: ERF_EB.H:36
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
36 { return 5; } // nghost_eb_basic ()

Referenced by make_factory().

Here is the caller graph for this function:

◆ nghost_full()

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

Referenced by make_factory().

Here is the caller graph for this function:

◆ nghost_volume()

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

Referenced by make_factory().

Here is the caller graph for this function:

◆ WriteEBSurface()

void eb_::WriteEBSurface ( const amrex::BoxArray &  ba,
const amrex::DistributionMapping &  dmap,
const amrex::Geometry &  geom,
const amrex::EBFArrayBoxFactory *  ebf,
const int  level 
)
82 {
83 
84  EBToPVD eb_to_pvd;
85 
86  const Real* dx = geom.CellSize();
87  const Real* problo = geom.ProbLo();
88 
89  MultiFab mf_ba(ba, dmap, 1, 0, MFInfo(), *ebf);
90 
91  for (MFIter mfi(mf_ba); mfi.isValid(); ++mfi) {
92 
93  const auto & sfab = static_cast<EBFArrayBox const &>(mf_ba[mfi]);
94  const auto & my_flag = sfab.getEBCellFlagFab();
95  const auto * my_flag_ptr = &my_flag;
96 
97  const Box & bx = mfi.validbox();
98 
99  if (my_flag.getType(bx) == FabType::covered ||
100  my_flag.getType(bx) == FabType::regular) { continue; }
101 
102  std::array<const CutFab *, AMREX_SPACEDIM> areafrac;
103  const CutFab * bndrycent;
104 
105  for (int d = 0; d < AMREX_SPACEDIM; ++d) {
106  areafrac[d] = &(*ebf->getAreaFrac()[d])[mfi];
107  }
108  bndrycent = &(ebf->getBndryCent()[mfi]);
109 
110 #ifdef AMREX_USE_GPU
111  std::unique_ptr<EBCellFlagFab> host_flag;
112  if (my_flag.arena()->isManaged() || my_flag.arena()->isDevice()) {
113  host_flag = std::make_unique<EBCellFlagFab>(my_flag.box(), my_flag.nComp(),
114  The_Pinned_Arena());
115  Gpu::dtoh_memcpy_async(host_flag->dataPtr(), my_flag.dataPtr(),
116  host_flag->nBytes());
117  Gpu::streamSynchronize();
118  my_flag_ptr = host_flag.get();
119  }
120 
121  std::array<std::unique_ptr<CutFab>, AMREX_SPACEDIM> areafrac_h;
122  for (int d = 0; d < AMREX_SPACEDIM; ++d) {
123  if (areafrac[d]->arena()->isManaged() || areafrac[d]->arena()->isDevice()) {
124  areafrac_h[d] = std::make_unique<CutFab>(areafrac[d]->box(), areafrac[d]->nComp(),
125  The_Pinned_Arena());
126  Gpu::dtoh_memcpy_async(areafrac_h[d]->dataPtr(), areafrac[d]->dataPtr(),
127  areafrac[d]->size()*sizeof(Real));
128  Gpu::streamSynchronize();
129  areafrac[d] = areafrac_h[d].get();
130  }
131  }
132 
133  std::unique_ptr<CutFab> bndrycent_h;
134  if (bndrycent->arena()->isManaged() || bndrycent->arena()->isDevice()) {
135  bndrycent_h = std::make_unique<CutFab>(bndrycent->box(), bndrycent->nComp(),
136  The_Pinned_Arena());
137  Gpu::dtoh_memcpy_async(bndrycent_h->dataPtr(), bndrycent->dataPtr(),
138  bndrycent->size()*sizeof(Real));
139  Gpu::streamSynchronize();
140  bndrycent = bndrycent_h.get();
141  }
142 #endif
143 
144  eb_to_pvd.EBToPolygon(
145  problo, dx,
146  bx, my_flag_ptr->const_array(),
147  bndrycent->const_array(),
148  areafrac[0]->const_array(),
149  areafrac[1]->const_array(),
150  areafrac[2]->const_array());
151  }
152 
153  int cpu = ParallelDescriptor::MyProc();
154  int nProcs = ParallelDescriptor::NProcs();
155 
156  eb_to_pvd.WriteEBVTP(cpu, level);
157 
158  if(ParallelDescriptor::IOProcessor()) {
159  EBToPVD::WritePVTP(nProcs);
160  }
161 
162  for (MFIter mfi(mf_ba); mfi.isValid(); ++mfi) {
163 
164  const auto & sfab = static_cast<EBFArrayBox const &>(mf_ba[mfi]);
165  const auto & my_flag = sfab.getEBCellFlagFab();
166  const auto * my_flag_ptr = &my_flag;
167 
168  const Box & bx = mfi.validbox();
169 
170  if (my_flag.getType(bx) == FabType::covered ||
171  my_flag.getType(bx) == FabType::regular) { continue; }
172 
173 #ifdef AMREX_USE_GPU
174  std::unique_ptr<EBCellFlagFab> host_flag;
175  if (my_flag.arena()->isManaged() || my_flag.arena()->isDevice()) {
176  host_flag = std::make_unique<EBCellFlagFab>(my_flag.box(), my_flag.nComp(),
177  The_Pinned_Arena());
178  Gpu::dtoh_memcpy_async(host_flag->dataPtr(), my_flag.dataPtr(),
179  host_flag->nBytes());
180  Gpu::streamSynchronize();
181  my_flag_ptr = host_flag.get();
182  }
183 #endif
184 
185  eb_to_pvd.EBGridCoverage(cpu, problo, dx, bx, my_flag_ptr->const_array());
186  }
187 }
static void WritePVTP(int nProcs)
Definition: ERF_EBToPVD.cpp:225

Referenced by make_factory().

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
private

Referenced by get_const_factory(), and make_factory().

◆ m_has_eb

int eb_::m_has_eb
private

◆ m_support_level

amrex::EBSupport eb_::m_support_level
private

Referenced by make_factory().

◆ 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: