Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
amrex::EBFArrayBoxFactory Class Reference

#include <AMReX_EBFabFactory.H>

Inheritance diagram for amrex::EBFArrayBoxFactory:
amrex::FabFactory< FArrayBox >

Public Member Functions

 EBFArrayBoxFactory (const EB2::Level &a_level, const Geometry &a_geom, const BoxArray &a_ba, const DistributionMapping &a_dm, const Vector< int > &a_ngrow, EBSupport a_support)
 Build an EB-aware factory backed by EB2 level data.
 
 ~EBFArrayBoxFactory () override=default
 
 EBFArrayBoxFactory (const EBFArrayBoxFactory &)=default
 
 EBFArrayBoxFactory (EBFArrayBoxFactory &&) noexcept=default
 
 EBFArrayBoxFactory ()=delete
 
EBFArrayBoxFactoryoperator= (const EBFArrayBoxFactory &)=delete
 
EBFArrayBoxFactoryoperator= (EBFArrayBoxFactory &&)=delete
 
FArrayBoxcreate (const Box &box, int ncomps, const FabInfo &info, int box_index) const final
 Allocate an EB-aware FAB for box with ncomps components.
 
FArrayBoxcreate_alias (FArrayBox const &rhs, int scomp, int ncomp) const final
 Create an alias view into an existing EB-aware FAB.
 
void destroy (FArrayBox *fab) const final
 Destroy a FAB produced by this factory.
 
EBFArrayBoxFactoryclone () const final
 Produce a deep copy of this factory (excluding EB level ownership).
 
const FabArray< EBCellFlagFab > & getMultiEBCellFlagFab () const noexcept
 EB cell flags for all boxes.
 
const MultiFabgetLevelSet () const noexcept
 Level set field.
 
const MultiFabgetVolFrac () const noexcept
 Volume fraction.
 
const MultiCutFabgetCentroid () const noexcept
 Cell centroid.
 
const MultiCutFabgetBndryCent () const noexcept
 Boundary centroid.
 
const MultiCutFabgetBndryNormal () const noexcept
 Boundary normals.
 
const MultiCutFabgetBndryArea () const noexcept
 Boundary area.
 
Array< const MultiCutFab *, 3 > getAreaFrac () const noexcept
 Area fractions along each axis.
 
Array< const MultiCutFab *, 3 > getFaceCent () const noexcept
 
Array< const MultiCutFab *, 3 > getEdgeCent () const noexcept
 Edge centroids.
 
bool isAllRegular () const noexcept
 
EB2::Level const * getEBLevel () const noexcept
 
EB2::IndexSpace const * getEBIndexSpace () const noexcept
 
int maxCoarseningLevel () const noexcept
 Maximum coarse level supported by the EB hierarchy.
 
const DistributionMappingDistributionMap () const noexcept
 Distribution mapping describing the FABs this factory allocates.
 
const BoxArrayboxArray () const noexcept
 BoxArray describing the FAB layout.
 
const GeometryGeom () const noexcept
 Geometry inherited from the EB level.
 
bool hasEBInfo () const noexcept
 True if EB geometric data beyond flags is available.
 
iMultiFab const * getCutCellMask () const noexcept
 
EBData getEBData (MFIter const &mfi) const noexcept
 Convenience wrapper returning EBData views for the FAB referenced by mfi.
 
EBDataArrays getEBDataArrays () const noexcept
 Return GPU-accessible array views for all FABs owned by this factory.
 
- Public Member Functions inherited from amrex::FabFactory< FArrayBox >
 FabFactory () noexcept=default
 
 FabFactory (FabFactory const &) noexcept=default
 
 FabFactory (FabFactory &&) noexcept=default
 
FabFactoryoperator= (FabFactory const &) noexcept=default
 
FabFactoryoperator= (FabFactory &&) noexcept=default
 
virtual ~FabFactory () noexcept=default
 
virtual Long nBytes (const Box &box, int ncomps, int) const
 

Detailed Description

FabFactory specialization that exposes EB geometry alongside state FABs.

Constructor & Destructor Documentation

◆ EBFArrayBoxFactory() [1/4]

amrex::EBFArrayBoxFactory::EBFArrayBoxFactory ( const EB2::Level a_level,
const Geometry a_geom,
const BoxArray a_ba,
const DistributionMapping a_dm,
const Vector< int > &  a_ngrow,
EBSupport  a_support 
)

Build an EB-aware factory backed by EB2 level data.

Parameters
a_levelEB level to sample.
a_geomGeometry describing the level.
a_baGrid layout.
a_dmDistribution mapping for FAB ownership.
a_ngrowGrow-cell requirements (basic/volume/full).
a_supportRequested EB support level.

◆ ~EBFArrayBoxFactory()

amrex::EBFArrayBoxFactory::~EBFArrayBoxFactory ( )
overridedefault

◆ EBFArrayBoxFactory() [2/4]

amrex::EBFArrayBoxFactory::EBFArrayBoxFactory ( const EBFArrayBoxFactory )
default

◆ EBFArrayBoxFactory() [3/4]

amrex::EBFArrayBoxFactory::EBFArrayBoxFactory ( EBFArrayBoxFactory &&  )
defaultnoexcept

◆ EBFArrayBoxFactory() [4/4]

amrex::EBFArrayBoxFactory::EBFArrayBoxFactory ( )
delete

Member Function Documentation

◆ boxArray()

const BoxArray & amrex::EBFArrayBoxFactory::boxArray ( ) const
noexcept

BoxArray describing the FAB layout.

◆ clone()

EBFArrayBoxFactory * amrex::EBFArrayBoxFactory::clone ( ) const
finalvirtual

Produce a deep copy of this factory (excluding EB level ownership).

Implements amrex::FabFactory< FArrayBox >.

◆ create()

FArrayBox * amrex::EBFArrayBoxFactory::create ( const Box box,
int  ncomps,
const FabInfo info,
int  box_index 
) const
finalvirtual

Allocate an EB-aware FAB for box with ncomps components.

Implements amrex::FabFactory< FArrayBox >.

◆ create_alias()

FArrayBox * amrex::EBFArrayBoxFactory::create_alias ( FArrayBox const &  rhs,
int  scomp,
int  ncomp 
) const
finalvirtual

Create an alias view into an existing EB-aware FAB.

Reimplemented from amrex::FabFactory< FArrayBox >.

◆ destroy()

void amrex::EBFArrayBoxFactory::destroy ( FArrayBox fab) const
finalvirtual

Destroy a FAB produced by this factory.

Implements amrex::FabFactory< FArrayBox >.

◆ DistributionMap()

const DistributionMapping & amrex::EBFArrayBoxFactory::DistributionMap ( ) const
noexcept

Distribution mapping describing the FABs this factory allocates.

◆ Geom()

const Geometry & amrex::EBFArrayBoxFactory::Geom ( ) const
inlinenoexcept

Geometry inherited from the EB level.

◆ getAreaFrac()

Array< const MultiCutFab *, 3 > amrex::EBFArrayBoxFactory::getAreaFrac ( ) const
inlinenoexcept

Area fractions along each axis.

◆ getBndryArea()

const MultiCutFab & amrex::EBFArrayBoxFactory::getBndryArea ( ) const
inlinenoexcept

Boundary area.

◆ getBndryCent()

const MultiCutFab & amrex::EBFArrayBoxFactory::getBndryCent ( ) const
inlinenoexcept

Boundary centroid.

◆ getBndryNormal()

const MultiCutFab & amrex::EBFArrayBoxFactory::getBndryNormal ( ) const
inlinenoexcept

Boundary normals.

◆ getCentroid()

const MultiCutFab & amrex::EBFArrayBoxFactory::getCentroid ( ) const
inlinenoexcept

Cell centroid.

◆ getCutCellMask()

iMultiFab const * amrex::EBFArrayBoxFactory::getCutCellMask ( ) const
inlinenoexcept

Returns nullptr unless this level is built by EB2::addRegularCoarseLevels. One should use getMultiEBCellFlagFab for normal levels.

◆ getEBData()

EBData amrex::EBFArrayBoxFactory::getEBData ( MFIter const &  mfi) const
noexcept

Convenience wrapper returning EBData views for the FAB referenced by mfi.

◆ getEBDataArrays()

EBDataArrays amrex::EBFArrayBoxFactory::getEBDataArrays ( ) const
noexcept

Return GPU-accessible array views for all FABs owned by this factory.

◆ getEBIndexSpace()

EB2::IndexSpace const * amrex::EBFArrayBoxFactory::getEBIndexSpace ( ) const
noexcept
Returns
EB IndexSpace pointer.

◆ getEBLevel()

EB2::Level const * amrex::EBFArrayBoxFactory::getEBLevel ( ) const
inlinenoexcept
Returns
EB level pointer.

◆ getEdgeCent()

Array< const MultiCutFab *, 3 > amrex::EBFArrayBoxFactory::getEdgeCent ( ) const
inlinenoexcept

Edge centroids.

◆ getFaceCent()

Array< const MultiCutFab *, 3 > amrex::EBFArrayBoxFactory::getFaceCent ( ) const
inlinenoexcept

Face centroids along each axis. For y-face-centroid, the two components are for x and then z-directions, respectively.

◆ getLevelSet()

const MultiFab & amrex::EBFArrayBoxFactory::getLevelSet ( ) const
inlinenoexcept

Level set field.

◆ getMultiEBCellFlagFab()

const FabArray< EBCellFlagFab > & amrex::EBFArrayBoxFactory::getMultiEBCellFlagFab ( ) const
inlinenoexcept

EB cell flags for all boxes.

◆ getVolFrac()

const MultiFab & amrex::EBFArrayBoxFactory::getVolFrac ( ) const
inlinenoexcept

Volume fraction.

◆ hasEBInfo()

bool amrex::EBFArrayBoxFactory::hasEBInfo ( ) const
noexcept

True if EB geometric data beyond flags is available.

◆ isAllRegular()

bool amrex::EBFArrayBoxFactory::isAllRegular ( ) const
noexcept
Returns
True if the EB is regular everywhere.

◆ maxCoarseningLevel()

int amrex::EBFArrayBoxFactory::maxCoarseningLevel ( ) const
noexcept

Maximum coarse level supported by the EB hierarchy.

◆ operator=() [1/2]

EBFArrayBoxFactory & amrex::EBFArrayBoxFactory::operator= ( const EBFArrayBoxFactory )
delete

◆ operator=() [2/2]

EBFArrayBoxFactory & amrex::EBFArrayBoxFactory::operator= ( EBFArrayBoxFactory &&  )
delete

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