Proto  3.2
Public Member Functions | Private Attributes | List of all members
Proto::MBLevelMap< MAP, MEM > Class Template Reference

Single Level Mapped Multiblock Map. More...

#include <Proto_MBLevelBoxData.H>

Public Member Functions

 MBLevelMap ()
 Trivial Constructor. More...
 
 MBLevelMap (const MBDisjointBoxLayout &a_layout, const Array< Point, DIM+1 > &a_ghost, unsigned int a_focalBlock=0)
 Non-Trivial Constructor. More...
 
 MBLevelMap (const MBDisjointBoxLayout &a_layout, Point a_ghost, unsigned int a_focalBlock=0)
 Non-Trivial Constructor. More...
 
void define (const MBDisjointBoxLayout &a_layout, const Array< Point, DIM+1 > &a_ghost, unsigned int a_focalBlock=0)
 Define. More...
 
void define (const MBDisjointBoxLayout &a_layout, Point a_ghost, unsigned int a_focalBlock=0)
 
void initialize () const
 
MBLevelBoxData< double, DIM, MEM, PR_NODE > & map ()
 Access Cached Coordinate Values. More...
 
std::shared_ptr< MBLevelBoxData< double, DIM, MEM, PR_NODE > > mapData ()
 
const MBLevelBoxData< double, DIM, MEM, PR_NODE > & map () const
 Access Cached Coordinate Values (Const Overload) More...
 
MBLevelBoxData< double, 1, MEM, PR_CELL > & jacobian ()
 Jacobian Access. More...
 
std::shared_ptr< MBLevelBoxData< double, 1, MEM, PR_CELL > > jacobianData ()
 
const MBLevelBoxData< double, 1, MEM, PR_CELL > & jacobian () const
 Jacobian Access (Const Overload) More...
 
void apply (BoxData< double, DIM, MEM > &a_X, BoxData< double, 1, MEM > &a_J, unsigned int a_block) const
 Compute Map. More...
 
void apply (BoxData< double, DIM, MEM > &a_X, BoxData< double, 1, MEM > &a_J, FluxBoxData< double, DIM, MEM > &a_NT, unsigned int a_block) const
 Compute Map (With Metrics) More...
 
void doApply (BoxData< double, DIM, MEM > &a_X, BoxData< double, 1, MEM > &a_J, unsigned int a_computeBlock, unsigned int a_outBlock) const
 Compute Map With Rotation. More...
 
BoxData< double, DIM, MEM > cellCentered (const Box &a_box, unsigned int a_computeBlock, unsigned int a_outBlock) const
 Compute Cell Centered Coordinates. More...
 
BoxData< double, DIM, MEM > cellAveraged (const Box &a_box, unsigned int a_computeBlock, unsigned int a_outBlock) const
 Compute Cell Centered Coordinates. More...
 
Array< double, DIM > mappedCoords (const Point &point, const BlockIndex &block) const
 
Array< double, DIM > mappedCoords (const Point &point, const BlockIndex &block, const Array< double, DIM > &offset) const
 
Array< double, DIM > apply (const MBDataPoint &a_point) const
 Compute Mapped Coordinates. More...
 
Array< double, DIM > apply (const MBDataPoint &a_point, const Array< double, DIM > &a_offset) const
 Analytic Map (DataPoint With Offset) More...
 
Array< double, DIM > cellCentered (const MBDataPoint &a_point) const
 Compute Cell Centered Coordinates. More...
 
Array< double, DIM > analyticCellCentered (const MBDataPoint &a_point) const
 Compute Cell Centered Coordinates, analytic case. More...
 
Array< double, DIM > cellAveraged (const MBDataPoint &a_point) const
 Compute Cell Averaged Coordinates. More...
 
const MBDisjointBoxLayoutlayout () const
 Get Layout. More...
 
const Array< double, DIM > dx (BlockIndex a_block) const
 Get Mapped Grid Spacing. More...
 
MAP & op (BlockIndex a_block) const
 Get Operator. More...
 
Array< Point, DIM+1 > ghost () const
 Get Ghost Sizes. More...
 
BoxData< double, DIM, MEM > X (const Box &a_box, const Array< double, DIM > &a_dx) const
 Get Mapped Coordinate Values. More...
 
MAP & operator[] (BlockIndex block)
 
Point convertPoint (Point srcPoint, BlockIndex srcBlock, BlockIndex dstBlock) const
 Convert a Point from one block's index space to another by using forward / inverse mapping. More...
 
std::optional< std::pair< MBIndex, Point > > whichNeighborContains (const MBIndex &localIndex, const Point &localPoint, bool checkLocalInteriorGhostCells=true) const
 

Private Attributes

Array< Point, DIM+1 > m_ghost
 
unsigned int m_focalBlock
 
Stencil< double > m_c2c
 
std::vector< std::shared_ptr< MAP > > m_ops
 Map operators. More...
 
std::vector< Array< double, DIM > > m_dx
 Mapped space grid spacing. More...
 
MBDisjointBoxLayout m_layout
 
bool m_initialized
 
std::shared_ptr< MBLevelBoxData< double, DIM, MEM, PR_NODE > > m_X
 Cached coodinate values. More...
 
std::shared_ptr< MBLevelBoxData< double, 1, MEM, PR_CELL > > m_J
 Cached Jacobian values. More...
 

Detailed Description

template<typename MAP, MemType MEM = MEMTYPE_DEFAULT>
class Proto::MBLevelMap< MAP, MEM >

Single Level Mapped Multiblock Map.

This class provides an interface for constructing maps for use in single level mapped multiblock applications.

When a map is constructed, coordinate and Jacobian data are cached for all valid regions of the supplied MBDisjointBoxLayout as well as proscribed ghost regions. However, the map's apply function can be used to compute these quantities on any domain where the associated functions are defined.

USAGE: To use this interface, create a class which inherits from MBLevelMap. The map itself is implemented by defining the one of the virtual "apply" functions (the user may supply either one or both depending on whether or not they need and/or would like to define how the metric tensor is computed). The user may also implement the virtual init function which will be called during construction of the Map. The init function is useful for caching constant data such as operators, stencils, or metrics in order to prevent computing these quantities each time apply is called.

Constructor & Destructor Documentation

◆ MBLevelMap() [1/3]

template<typename MAP, MemType MEM = MEMTYPE_DEFAULT>
Proto::MBLevelMap< MAP, MEM >::MBLevelMap ( )
inline

Trivial Constructor.

◆ MBLevelMap() [2/3]

template<typename MAP , MemType MEM>
Proto::MBLevelMap< MAP, MEM >::MBLevelMap ( const MBDisjointBoxLayout a_layout,
const Array< Point, DIM+1 > &  a_ghost,
unsigned int  a_focalBlock = 0 
)
inline

Non-Trivial Constructor.

◆ MBLevelMap() [3/3]

template<typename MAP , MemType MEM>
Proto::MBLevelMap< MAP, MEM >::MBLevelMap ( const MBDisjointBoxLayout a_layout,
Point  a_ghost,
unsigned int  a_focalBlock = 0 
)
inline

Non-Trivial Constructor.

Member Function Documentation

◆ define() [1/2]

template<typename MAP , MemType MEM>
void Proto::MBLevelMap< MAP, MEM >::define ( const MBDisjointBoxLayout a_layout,
const Array< Point, DIM+1 > &  a_ghost,
unsigned int  a_focalBlock = 0 
)
inline

Define.

Delayed construction.

Parameters
a_layoutA mapped multiblock layout defining the domain
a_ghostAn array of ghost data sizes. This ghost array is used To define how much ghost data should be allocated for cached arrays of the coordinates and Jacobian.
a_focalBlockThe "focus" block id associated with this map. Often unused.

◆ define() [2/2]

template<typename MAP , MemType MEM>
void Proto::MBLevelMap< MAP, MEM >::define ( const MBDisjointBoxLayout a_layout,
Point  a_ghost,
unsigned int  a_focalBlock = 0 
)
inline

◆ initialize()

template<typename MAP , MemType MEM>
void Proto::MBLevelMap< MAP, MEM >::initialize ( ) const
inline

◆ map() [1/2]

template<typename MAP , MemType MEM>
MBLevelBoxData< double, DIM, MEM, PR_NODE > & Proto::MBLevelMap< MAP, MEM >::map ( )
inline

Access Cached Coordinate Values.

Returns an MBLevelBoxData containing the cached coordinate data. The returned data has NODE centering.

◆ mapData()

template<typename MAP, MemType MEM = MEMTYPE_DEFAULT>
std::shared_ptr<MBLevelBoxData<double, DIM, MEM, PR_NODE> > Proto::MBLevelMap< MAP, MEM >::mapData ( )
inline

◆ map() [2/2]

template<typename MAP , MemType MEM>
const MBLevelBoxData< double, DIM, MEM, PR_NODE > & Proto::MBLevelMap< MAP, MEM >::map ( ) const
inline

Access Cached Coordinate Values (Const Overload)

◆ jacobian() [1/2]

template<typename MAP , MemType MEM>
MBLevelBoxData< double, 1, MEM, PR_CELL > & Proto::MBLevelMap< MAP, MEM >::jacobian ( )
inline

Jacobian Access.

Returns an MBLevelBoxData containing the cached Jacobian data. The returned data has CELL centering and is assumed to be cell averaged.

Referenced by Proto::MBLevelMap< MAP, MEM >::jacobianData().

◆ jacobianData()

template<typename MAP, MemType MEM = MEMTYPE_DEFAULT>
std::shared_ptr<MBLevelBoxData<double, 1, MEM, PR_CELL> > Proto::MBLevelMap< MAP, MEM >::jacobianData ( )
inline

◆ jacobian() [2/2]

template<typename MAP , MemType MEM>
const MBLevelBoxData< double, 1, MEM, PR_CELL > & Proto::MBLevelMap< MAP, MEM >::jacobian ( ) const
inline

Jacobian Access (Const Overload)

◆ apply() [1/4]

template<typename MAP , MemType MEM>
void Proto::MBLevelMap< MAP, MEM >::apply ( BoxData< double, DIM, MEM > &  a_X,
BoxData< double, 1, MEM > &  a_J,
unsigned int  a_block 
) const
inline

Compute Map.

Users should override this function when defining their map if they do not want to explicitly define the metric tensor (NT). The domains of the coordinate and Jacobian quantities are assumed to be defined before they are input to this function. Note that the ghost values used to define the map during construction have absolutely no bearing on the valid domains for this function, so long as the relevant analytic functions are defined.

Parameters
a_XCoordinates at nodes [input/output]
a_JCell averaged jacobian [input/output]
a_indexMBIndex corresponding to the patch where the operator is being applied

◆ apply() [2/4]

template<typename MAP , MemType MEM>
void Proto::MBLevelMap< MAP, MEM >::apply ( BoxData< double, DIM, MEM > &  a_X,
BoxData< double, 1, MEM > &  a_J,
FluxBoxData< double, DIM, MEM > &  a_NT,
unsigned int  a_block 
) const
inline

Compute Map (With Metrics)

Users should override this function when defining their map if they want to explicitly define the metric tensor (NT). The domains of the coordinate and Jacobian quantities are assumed to be defined before they are input to this function. Note that the ghost values used to define the map during construction have absolutely no bearing on the valid domains for this function, so long as the relevant analytic functions are defined.

Parameters
a_XCoordinates at nodes [input/output]
a_JCell averaged jacobian [input/output]
a_NTFace averaged metric terms in each coordinate direction [output]
a_indexMBIndex corresponding to the patch where the operator is being applied

◆ doApply()

template<typename MAP , MemType MEM>
void Proto::MBLevelMap< MAP, MEM >::doApply ( BoxData< double, DIM, MEM > &  a_X,
BoxData< double, 1, MEM > &  a_J,
unsigned int  a_computeBlock,
unsigned int  a_outBlock 
) const
inline

Compute Map With Rotation.

Used internally for the case where X and J are defined on domains in a different coordinate system from the block of the map computation (defined by a_index). This function creates a temporary version of X and J with domains in the block associated with a_index, computes X and J, and copies the result with rotation into the outputs.

Not recommended for public use.

◆ cellCentered() [1/2]

template<typename MAP , MemType MEM>
BoxData< double, DIM, MEM > Proto::MBLevelMap< MAP, MEM >::cellCentered ( const Box a_box,
unsigned int  a_computeBlock,
unsigned int  a_outBlock 
) const
inline

Compute Cell Centered Coordinates.

Returns the cell-centered coordinates on a specified range.

Parameters
a_boxOutput range
a_computeBlockBlock index associated with the coordinates being computed
a_outBlockBlock index associated with the input range

◆ cellAveraged() [1/2]

template<typename MAP , MemType MEM>
BoxData< double, DIM, MEM > Proto::MBLevelMap< MAP, MEM >::cellAveraged ( const Box a_box,
unsigned int  a_computeBlock,
unsigned int  a_outBlock 
) const
inline

Compute Cell Centered Coordinates.

Returns the cell-averaged coordinates on a specified range.

Parameters
a_boxOutput range
a_computeBlockBlock index associated with the coordinates being computed
a_outBlockBlock index associated with the input range

◆ mappedCoords() [1/2]

template<typename MAP , MemType MEM>
Array< double, DIM > Proto::MBLevelMap< MAP, MEM >::mappedCoords ( const Point point,
const BlockIndex block 
) const
inline

◆ mappedCoords() [2/2]

template<typename MAP , MemType MEM>
Array< double, DIM > Proto::MBLevelMap< MAP, MEM >::mappedCoords ( const Point point,
const BlockIndex block,
const Array< double, DIM > &  offset 
) const
inline

◆ apply() [3/4]

template<typename MAP , MemType MEM>
Array< double, DIM > Proto::MBLevelMap< MAP, MEM >::apply ( const MBDataPoint a_point) const
inline

Compute Mapped Coordinates.

Returns coordinates in mapped space of a given data point. The Coordinates are of the cell center of the point plus the offset which is in units of the grid spacing. Analytic Map (DataPoint Cell Center)

◆ apply() [4/4]

template<typename MAP , MemType MEM>
Array< double, DIM > Proto::MBLevelMap< MAP, MEM >::apply ( const MBDataPoint a_point,
const Array< double, DIM > &  a_offset 
) const
inline

Analytic Map (DataPoint With Offset)

◆ cellCentered() [2/2]

template<typename MAP , MemType MEM>
Array< double, DIM > Proto::MBLevelMap< MAP, MEM >::cellCentered ( const MBDataPoint a_point) const
inline

Compute Cell Centered Coordinates.

Computes the cell-centered coodinate values at a specified MBDataPoint

◆ analyticCellCentered()

template<typename MAP , MemType MEM>
Array< double, DIM > Proto::MBLevelMap< MAP, MEM >::analyticCellCentered ( const MBDataPoint a_point) const
inline

Compute Cell Centered Coordinates, analytic case.

Computes the cell-centered coodinate valuesas a specified analytic function,

auto computeBlock = a_point.srcBlock(); Array<double,DIM> X; Array<double,DIM> thisDx = this->dx(a_point.srcBlock()); Point pt = a_point.point(); f_cubeSphereMapPoint<MEM>(pt,X,thisDx,CUBED_SPHERE_SHELL_R0, CUBED_SPHERE_SHELL_R1,a_point.srcBlock(), CUBED_SPHERE_SHELL_RADIAL_COORD); return X; }

◆ cellAveraged() [2/2]

template<typename MAP , MemType MEM>
Array< double, DIM > Proto::MBLevelMap< MAP, MEM >::cellAveraged ( const MBDataPoint a_point) const
inline

Compute Cell Averaged Coordinates.

Computes the cell-averaged coodinate values at a specified MBDataPoint

◆ layout()

template<typename MAP , MemType MEM>
const MBDisjointBoxLayout & Proto::MBLevelMap< MAP, MEM >::layout ( ) const
inline

Get Layout.

◆ dx()

template<typename MAP , MemType MEM>
const Array< double, DIM > Proto::MBLevelMap< MAP, MEM >::dx ( BlockIndex  a_block) const
inline

Get Mapped Grid Spacing.

◆ op()

template<typename MAP , MemType MEM>
MAP & Proto::MBLevelMap< MAP, MEM >::op ( BlockIndex  a_block) const
inline

Get Operator.

◆ ghost()

template<typename MAP, MemType MEM = MEMTYPE_DEFAULT>
Array<Point, DIM+1> Proto::MBLevelMap< MAP, MEM >::ghost ( ) const
inline

Get Ghost Sizes.

References Proto::MBLevelMap< MAP, MEM >::m_ghost.

◆ X()

template<typename MAP , MemType MEM>
BoxData< double, DIM, MEM > Proto::MBLevelMap< MAP, MEM >::X ( const Box a_box,
const Array< double, DIM > &  a_dx 
) const
inline

Get Mapped Coordinate Values.

Given a node-centered box and grid spacing, generate the node centered mapped coordinates. This assumes that the origin corresponds to the point (0,0,...,0) and that the input box already accounts for node centering (e.g. Box::grow(PR_NODE) is assumed to already have been called.)

/param a_box A box that already accounts for node centering /param a_dx Grid spacing

◆ operator[]()

template<typename MAP , MemType MEM>
MAP & Proto::MBLevelMap< MAP, MEM >::operator[] ( BlockIndex  block)
inline

◆ convertPoint()

template<typename MAP , MemType MEM>
Point Proto::MBLevelMap< MAP, MEM >::convertPoint ( Point  srcPoint,
BlockIndex  srcBlock,
BlockIndex  dstBlock 
) const
inline

Convert a Point from one block's index space to another by using forward / inverse mapping.

Parameters
srcPointPoint in the source block's index space
srcBlockSource block
dstBlockDestination block
Returns

◆ whichNeighborContains()

template<typename MAP , MemType MEM>
std::optional< std::pair< MBIndex, Point > > Proto::MBLevelMap< MAP, MEM >::whichNeighborContains ( const MBIndex localIndex,
const Point localPoint,
bool  checkLocalInteriorGhostCells = true 
) const
inline

Member Data Documentation

◆ m_ghost

template<typename MAP, MemType MEM = MEMTYPE_DEFAULT>
Array<Point, DIM+1> Proto::MBLevelMap< MAP, MEM >::m_ghost
private

◆ m_focalBlock

template<typename MAP, MemType MEM = MEMTYPE_DEFAULT>
unsigned int Proto::MBLevelMap< MAP, MEM >::m_focalBlock
private

◆ m_c2c

template<typename MAP, MemType MEM = MEMTYPE_DEFAULT>
Stencil<double> Proto::MBLevelMap< MAP, MEM >::m_c2c
private

◆ m_ops

template<typename MAP, MemType MEM = MEMTYPE_DEFAULT>
std::vector<std::shared_ptr<MAP> > Proto::MBLevelMap< MAP, MEM >::m_ops
mutableprivate

Map operators.

◆ m_dx

template<typename MAP, MemType MEM = MEMTYPE_DEFAULT>
std::vector<Array<double, DIM> > Proto::MBLevelMap< MAP, MEM >::m_dx
private

Mapped space grid spacing.

◆ m_layout

template<typename MAP, MemType MEM = MEMTYPE_DEFAULT>
MBDisjointBoxLayout Proto::MBLevelMap< MAP, MEM >::m_layout
private

◆ m_initialized

template<typename MAP, MemType MEM = MEMTYPE_DEFAULT>
bool Proto::MBLevelMap< MAP, MEM >::m_initialized
mutableprivate

◆ m_X

template<typename MAP, MemType MEM = MEMTYPE_DEFAULT>
std::shared_ptr<MBLevelBoxData<double, DIM, MEM, PR_NODE> > Proto::MBLevelMap< MAP, MEM >::m_X
mutableprivate

Cached coodinate values.

Referenced by Proto::MBLevelMap< MAP, MEM >::mapData().

◆ m_J

template<typename MAP, MemType MEM = MEMTYPE_DEFAULT>
std::shared_ptr<MBLevelBoxData<double, 1, MEM, PR_CELL> > Proto::MBLevelMap< MAP, MEM >::m_J
mutableprivate

Cached Jacobian values.

Referenced by Proto::MBLevelMap< MAP, MEM >::jacobianData().


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