OPAL (Object Oriented Parallel Accelerator Library)  2.2.0
OPAL
Public Types | Public Member Functions | Static Public Attributes | Private Member Functions | Private Attributes | List of all members
BoxLibLayout< T, Dim > Class Template Reference

#include <BoxLibLayout.h>

Inheritance diagram for BoxLibLayout< T, Dim >:
Inheritance graph
[legend]
Collaboration diagram for BoxLibLayout< T, Dim >:
Collaboration graph
[legend]

Public Types

typedef ParticleAmrLayout< T,
Dim >::pair_t 
pair_t
 
typedef ParticleAmrLayout< T,
Dim >::pair_iterator 
pair_iterator
 
typedef ParticleAmrLayout< T,
Dim >::SingleParticlePos_t 
SingleParticlePos_t
 
typedef ParticleAmrLayout< T,
Dim >::Index_t 
Index_t
 
typedef amr::AmrField_t AmrField_t
 
typedef amr::AmrVectorField_t AmrVectorField_t
 
typedef
amr::AmrScalarFieldContainer_t 
AmrScalarFieldContainer_t
 
typedef
amr::AmrVectorFieldContainer_t 
AmrVectorFieldContainer_t
 
typedef ParticleAmrLayout< T,
Dim >::ParticlePos_t 
ParticlePos_t
 
typedef ParticleAttrib< Index_tParticleIndex_t
 
typedef amr::AmrProcMap_t AmrProcMap_t
 
typedef amr::AmrGrid_t AmrGrid_t
 
typedef amr::AmrGeometry_t AmrGeometry_t
 
typedef amr::AmrGeomContainer_t AmrGeomContainer_t
 
typedef amr::AmrGridContainer_t AmrGridContainer_t
 
typedef amr::AmrProcMapContainer_t AmrProcMapContainer_t
 
typedef amr::AmrIntVect_t AmrIntVect_t
 
typedef amr::AmrIntVectContainer_t AmrIntVectContainer_t
 
typedef amr::AmrIntArray_t AmrIntArray_t
 
typedef amr::AmrDomain_t AmrDomain_t
 
typedef amr::AmrBox_t AmrBox_t
 
typedef amr::AmrReal_t AmrReal_t
 
typedef amrex::BaseFab< int > basefab_t
 
typedef amrex::FabArray
< basefab_t
mask_t
 
- Public Types inherited from ParticleAmrLayout< T, Dim >
typedef int pair_t
 
typedef pair_tpair_iterator
 
typedef ParticleLayout< T, Dim >
::SingleParticlePos_t 
SingleParticlePos_t
 
typedef ParticleLayout< T, Dim >
::Index_t 
Index_t
 
typedef ParticleAttrib
< SingleParticlePos_t
ParticlePos_t
 
typedef ParticleAttrib< Index_tParticleIndex_t
 
- Public Types inherited from ParticleLayout< T, Dim >
enum  { Dimension = Dim }
 
enum  UpdateFlags {
  SWAP, BCONDS, NUMFLAGS, OPTDESTROY,
  ALL
}
 
typedef T Position_t
 
typedef unsigned Index_t
 
typedef Vektor< T, DimSingleParticlePos_t
 

Public Member Functions

 BoxLibLayout ()
 
 BoxLibLayout (const BoxLibLayout *layout_p)
 
 BoxLibLayout (int nGridPoints, int maxGridSize)
 
 BoxLibLayout (const AmrGeometry_t &geom, const AmrProcMap_t &dmap, const AmrGrid_t &ba)
 
 BoxLibLayout (const AmrGeomContainer_t &geom, const AmrProcMapContainer_t &dmap, const AmrGridContainer_t &ba, const AmrIntArray_t &rr)
 
void setBoundingBox (double dh)
 
void setDomainRatio (const std::vector< double > &ratio)
 
void update (IpplParticleBase< BoxLibLayout< T, Dim > > &PData, const ParticleAttrib< char > *canSwap=0)
 
void update (AmrParticleBase< BoxLibLayout< T, Dim > > &PData, int lev_min=0, int lev_max=-1, bool isRegrid=false)
 
AmrIntVect_t Index (AmrParticleBase< BoxLibLayout< T, Dim > > &p, const unsigned int ip, int level) const
 
AmrIntVect_t Index (SingleParticlePos_t &R, int lev) const
 
void buildLevelMask (int lev, const int ncells=1)
 
void clearLevelMask (int lev)
 
const std::unique_ptr< mask_t > & getLevelMask (int lev) const
 
void resize (int maxLevel)
 
void define (const AmrGeomContainer_t &geom)
 
void define (const AmrIntVectContainer_t &refRatio)
 
bool LevelDefined (int level) const
 
int finestLevel () const
 
int maxLevel () const
 
AmrIntVect_t refRatio (int level) const
 
int MaxRefRatio (int level) const
 
- Public Member Functions inherited from ParticleAmrLayout< T, Dim >
 ParticleAmrLayout ()
 
void setFinestLevel (int finestLevel)
 
void setMaxLevel (int maxLevel)
 
- Public Member Functions inherited from ParticleLayout< T, Dim >
 ParticleLayout ()
 
 ~ParticleLayout ()
 
void setUpdateFlag (UpdateFlags f, bool val)
 
bool getUpdateFlag (UpdateFlags f) const
 
ParticleBConds< T, Dim > & getBConds ()
 
void setBConds (const ParticleBConds< T, Dim > &bc)
 

Static Public Attributes

static Vector_t lowerBound = - Vector_t(1.0, 1.0, 1.0)
 
static Vector_t upperBound = Vector_t(1.0, 1.0, 1.0)
 

Private Member Functions

void initBaseBox_m (int nGridPoints, int maxGridSize, double dh=0.04)
 
bool Where (AmrParticleBase< BoxLibLayout< T, Dim > > &p, const unsigned int ip, int lev_min=0, int lev_max=-1, int nGrow=0) const
 
bool EnforcePeriodicWhere (AmrParticleBase< BoxLibLayout< T, Dim > > &prt, const unsigned int ip, int lev_min=0, int lev_max=-1) const
 
bool PeriodicShift (SingleParticlePos_t R) const
 
void locateParticle (AmrParticleBase< BoxLibLayout< T, Dim > > &p, const unsigned int ip, int lev_min, int lev_max, int nGrow) const
 

Private Attributes

AmrIntVectContainer_t refRatio_m
 
std::vector< std::unique_ptr
< mask_t > > 
masks_m
 Refinement ratios [0:finest_level-1]. More...
 

Additional Inherited Members

- Protected Member Functions inherited from ParticleLayout< T, Dim >
template<class PPT , class NDI >
void apply_bconds (unsigned n, PPT &R, const ParticleBConds< T, Dim > &bcs, const NDI &nr)
 
- Protected Attributes inherited from ParticleAmrLayout< T, Dim >
int finestLevel_m
 Current finest level of simluation. More...
 
int maxLevel_m
 Maximum level allowed. More...
 

Detailed Description

template<class T, unsigned Dim>
class BoxLibLayout< T, Dim >

In contrast to AMReX, OPAL is optimized for distribution particles to cores. In AMReX the ParGDB object is responsible for the particle to core distribution. This layout is derived from this object and does all important bunch updates. It is the interface for AMReX and Ippl.

In AMReX the geometry, i.e. physical domain, is fixed during the whole computation. Particles leaving the domain would be deleted. In order to prevent this we map the particles onto the domain \([-1, 1]^3\). Furthermore, it makes sure that we have enougth grid points to represent the bunch when its charges are scattered on the grid for the self-field computation.

The self-field computation and the particle-to-core update are performed in the particle mapped domain.

Definition at line 30 of file BoxLibLayout.h.

Member Typedef Documentation

template<class T, unsigned Dim>
typedef amr::AmrBox_t BoxLibLayout< T, Dim >::AmrBox_t

Definition at line 57 of file BoxLibLayout.h.

template<class T, unsigned Dim>
typedef amr::AmrDomain_t BoxLibLayout< T, Dim >::AmrDomain_t

Definition at line 56 of file BoxLibLayout.h.

template<class T, unsigned Dim>
typedef amr::AmrField_t BoxLibLayout< T, Dim >::AmrField_t

Definition at line 40 of file BoxLibLayout.h.

template<class T, unsigned Dim>
typedef amr::AmrGeomContainer_t BoxLibLayout< T, Dim >::AmrGeomContainer_t

Definition at line 50 of file BoxLibLayout.h.

template<class T, unsigned Dim>
typedef amr::AmrGeometry_t BoxLibLayout< T, Dim >::AmrGeometry_t

Definition at line 49 of file BoxLibLayout.h.

template<class T, unsigned Dim>
typedef amr::AmrGrid_t BoxLibLayout< T, Dim >::AmrGrid_t

Definition at line 48 of file BoxLibLayout.h.

template<class T, unsigned Dim>
typedef amr::AmrGridContainer_t BoxLibLayout< T, Dim >::AmrGridContainer_t

Definition at line 51 of file BoxLibLayout.h.

template<class T, unsigned Dim>
typedef amr::AmrIntArray_t BoxLibLayout< T, Dim >::AmrIntArray_t

Definition at line 55 of file BoxLibLayout.h.

template<class T, unsigned Dim>
typedef amr::AmrIntVect_t BoxLibLayout< T, Dim >::AmrIntVect_t

Definition at line 53 of file BoxLibLayout.h.

template<class T, unsigned Dim>
typedef amr::AmrIntVectContainer_t BoxLibLayout< T, Dim >::AmrIntVectContainer_t

Definition at line 54 of file BoxLibLayout.h.

template<class T, unsigned Dim>
typedef amr::AmrProcMap_t BoxLibLayout< T, Dim >::AmrProcMap_t

Definition at line 47 of file BoxLibLayout.h.

template<class T, unsigned Dim>
typedef amr::AmrProcMapContainer_t BoxLibLayout< T, Dim >::AmrProcMapContainer_t

Definition at line 52 of file BoxLibLayout.h.

template<class T, unsigned Dim>
typedef amr::AmrReal_t BoxLibLayout< T, Dim >::AmrReal_t

Definition at line 58 of file BoxLibLayout.h.

template<class T, unsigned Dim>
typedef amr::AmrScalarFieldContainer_t BoxLibLayout< T, Dim >::AmrScalarFieldContainer_t

Definition at line 42 of file BoxLibLayout.h.

template<class T, unsigned Dim>
typedef amr::AmrVectorField_t BoxLibLayout< T, Dim >::AmrVectorField_t

Definition at line 41 of file BoxLibLayout.h.

template<class T, unsigned Dim>
typedef amr::AmrVectorFieldContainer_t BoxLibLayout< T, Dim >::AmrVectorFieldContainer_t

Definition at line 43 of file BoxLibLayout.h.

template<class T, unsigned Dim>
typedef amrex::BaseFab<int> BoxLibLayout< T, Dim >::basefab_t

Definition at line 60 of file BoxLibLayout.h.

template<class T, unsigned Dim>
typedef ParticleAmrLayout<T, Dim>::Index_t BoxLibLayout< T, Dim >::Index_t

Definition at line 38 of file BoxLibLayout.h.

template<class T, unsigned Dim>
typedef amrex::FabArray<basefab_t> BoxLibLayout< T, Dim >::mask_t

Definition at line 61 of file BoxLibLayout.h.

template<class T, unsigned Dim>
typedef ParticleAmrLayout<T, Dim>::pair_iterator BoxLibLayout< T, Dim >::pair_iterator

Definition at line 36 of file BoxLibLayout.h.

template<class T, unsigned Dim>
typedef ParticleAmrLayout<T, Dim>::pair_t BoxLibLayout< T, Dim >::pair_t

Definition at line 35 of file BoxLibLayout.h.

template<class T, unsigned Dim>
typedef ParticleAttrib<Index_t> BoxLibLayout< T, Dim >::ParticleIndex_t

Definition at line 45 of file BoxLibLayout.h.

template<class T, unsigned Dim>
typedef ParticleAmrLayout<T, Dim>::ParticlePos_t BoxLibLayout< T, Dim >::ParticlePos_t

Definition at line 44 of file BoxLibLayout.h.

template<class T, unsigned Dim>
typedef ParticleAmrLayout<T, Dim>::SingleParticlePos_t BoxLibLayout< T, Dim >::SingleParticlePos_t

Definition at line 37 of file BoxLibLayout.h.

Constructor & Destructor Documentation

template<class T , unsigned Dim>
BoxLibLayout< T, Dim >::BoxLibLayout ( )

Initializes default Geometry, DistributionMapping and BoxArray.

Definition at line 21 of file BoxLibLayout.hpp.

References ceil(), IpplInfo::getNodes(), and BoxLibLayout< T, Dim >::initBaseBox_m().

Here is the call graph for this function:

template<class T , unsigned Dim>
BoxLibLayout< T, Dim >::BoxLibLayout ( const BoxLibLayout< T, Dim > *  layout_p)

Given a layout it copies that.

Definition at line 45 of file BoxLibLayout.hpp.

References ParticleAmrLayout< T, Dim >::maxLevel_m, BoxLibLayout< T, Dim >::refRatio(), and BoxLibLayout< T, Dim >::refRatio_m.

Here is the call graph for this function:

template<class T , unsigned Dim>
BoxLibLayout< T, Dim >::BoxLibLayout ( int  nGridPoints,
int  maxGridSize 
)
Parameters
nGridPointsper dimension (nx, ny, nz / nt)
maxGridSizefor all levels.

Definition at line 60 of file BoxLibLayout.hpp.

References BoxLibLayout< T, Dim >::initBaseBox_m().

Here is the call graph for this function:

template<class T , unsigned Dim>
BoxLibLayout< T, Dim >::BoxLibLayout ( const AmrGeometry_t geom,
const AmrProcMap_t dmap,
const AmrGrid_t ba 
)

Single-level constructor.

Parameters
geomspecifies the box domain
dmapis the distribution map for grids
bais the array of boxes for a level

Definition at line 70 of file BoxLibLayout.hpp.

template<class T , unsigned Dim>
BoxLibLayout< T, Dim >::BoxLibLayout ( const AmrGeomContainer_t geom,
const AmrProcMapContainer_t dmap,
const AmrGridContainer_t ba,
const AmrIntArray_t rr 
)

Multi-level constructor.

Parameters
geomis basically the physical domain storing the mesh spacing per level
dmapare all distribution maps of grids to core
baare all boxes of all levels
rris the refinement ratio among the levels (always the ratio from l to l+1)

Definition at line 80 of file BoxLibLayout.hpp.

Member Function Documentation

template<class T , unsigned Dim>
void BoxLibLayout< T, Dim >::buildLevelMask ( int  lev,
const int  ncells = 1 
)

Build mask for a level used for interpolation from grid to particles to reduce spurious self field forces near coarse-fine interfaces.

Definition at line 383 of file BoxLibLayout.hpp.

Referenced by AmrBoxLib::MakeNewLevel(), and AmrBoxLib::RemakeLevel().

template<class T , unsigned Dim>
void BoxLibLayout< T, Dim >::clearLevelMask ( int  lev)

Definition at line 444 of file BoxLibLayout.hpp.

Referenced by AmrBoxLib::postRegrid_m().

template<class T, unsigned Dim>
void BoxLibLayout< T, Dim >::define ( const AmrGeomContainer_t geom)
inline

Set the geometry of the problem. It is called in AmrBoxLib::initBaseLevel_m().

Parameters
geomgeometry of all levels

Definition at line 242 of file BoxLibLayout.h.

Referenced by AmrBoxLib::initBaseLevel_m().

template<class T, unsigned Dim>
void BoxLibLayout< T, Dim >::define ( const AmrIntVectContainer_t refRatio)
inline

Set the refinement ratios. It is called in AmrBoxLib::initBaseLevel_m().

Parameters
refRatioamong levels

Definition at line 254 of file BoxLibLayout.h.

References BoxLibLayout< T, Dim >::refRatio_m.

template<class T, unsigned Dim>
bool BoxLibLayout< T, Dim >::EnforcePeriodicWhere ( AmrParticleBase< BoxLibLayout< T, Dim > > &  prt,
const unsigned int  ip,
int  lev_min = 0,
int  lev_max = -1 
) const
private

Function from AMReX adjusted to work with Ippl AmrParticleBase class Checks/sets whether the particle has crossed a periodic boundary in such a way that it is on levels lev_min and higher.

Parameters
prtis the bunch information
ipis the local (i.e. to a core) particle index
lev_minto check
lev_maxto check
Returns
true if mapped to the other side.

Definition at line 567 of file BoxLibLayout.hpp.

References Attrib::Distribution::R.

template<class T , unsigned Dim>
int BoxLibLayout< T, Dim >::finestLevel ( ) const
inline
Returns
the current finest level

Definition at line 745 of file BoxLibLayout.hpp.

Referenced by AmrYtWriter::writeBunch().

template<class T , unsigned Dim>
const std::unique_ptr< typename BoxLibLayout< T, Dim >::mask_t > & BoxLibLayout< T, Dim >::getLevelMask ( int  lev) const

Definition at line 452 of file BoxLibLayout.hpp.

template<class T, unsigned Dim>
BoxLibLayout< T, Dim >::AmrIntVect_t BoxLibLayout< T, Dim >::Index ( AmrParticleBase< BoxLibLayout< T, Dim > > &  p,
const unsigned int  ip,
int  level 
) const

Get the cell of a particle

Parameters
pis the particle data
ipis the local index of the particle in the container
levelof the particle

Definition at line 356 of file BoxLibLayout.hpp.

Referenced by AmrBoxLib::tagForMaxNumParticles_m(), AmrBoxLib::tagForMinNumParticles_m(), and AmrBoxLib::tagForMomenta_m().

template<class T, unsigned Dim>
BoxLibLayout< T, Dim >::AmrIntVect_t BoxLibLayout< T, Dim >::Index ( SingleParticlePos_t R,
int  lev 
) const

Get the cell of a particle

Parameters
Ris the position of a particle
levis the level

Definition at line 366 of file BoxLibLayout.hpp.

References floor().

Here is the call graph for this function:

template<class T , unsigned Dim>
void BoxLibLayout< T, Dim >::initBaseBox_m ( int  nGridPoints,
int  maxGridSize,
double  dh = 0.04 
)
private

Set up the box for the whole computation. The AMR object owning the bunch is not yet initialized.

Parameters
nGridPointsper dimension (nx, ny, nz / nt)
maxGridSizefor all levels
dhis the mesh enlargement factor

Definition at line 774 of file BoxLibLayout.hpp.

References IpplInfo::getNodes().

Referenced by BoxLibLayout< T, Dim >::BoxLibLayout().

Here is the call graph for this function:

template<class T , unsigned Dim>
bool BoxLibLayout< T, Dim >::LevelDefined ( int  level) const
inline

Check if an AMR level is well defined

Parameters
levelto check

Definition at line 739 of file BoxLibLayout.hpp.

template<class T, unsigned Dim>
void BoxLibLayout< T, Dim >::locateParticle ( AmrParticleBase< BoxLibLayout< T, Dim > > &  p,
const unsigned int  ip,
int  lev_min,
int  lev_max,
int  nGrow 
) const
private

Function from AMReX adjusted to work with Ippl AmrParticleBase class

Parameters
pis basically the bunch
ipis the local particle index
lev_minto check
lev_maxto check
nGrowis the number of ghost cells

Definition at line 686 of file BoxLibLayout.hpp.

template<class T , unsigned Dim>
int BoxLibLayout< T, Dim >::maxLevel ( ) const
inline
Returns
the maximum level of simulation

Definition at line 751 of file BoxLibLayout.hpp.

Referenced by BoxLibLayout< T, Dim >::resize().

template<class T , unsigned Dim>
int BoxLibLayout< T, Dim >::MaxRefRatio ( int  level) const
inline
Parameters
level
Returns
the maximum refinement ratio among all directions for the given level.

Definition at line 765 of file BoxLibLayout.hpp.

References max(), and Hypervolume::n.

Here is the call graph for this function:

template<class T , unsigned Dim>
bool BoxLibLayout< T, Dim >::PeriodicShift ( SingleParticlePos_t  R) const
private

Function from AMReX adjusted to work with Ippl AmrParticleBase class

Move the particle to the opposite side of the domain

Parameters
Ris the particle position
Returns
true if the particle was shifted.

Definition at line 615 of file BoxLibLayout.hpp.

template<class T , unsigned Dim>
BoxLibLayout< T, Dim >::AmrIntVect_t BoxLibLayout< T, Dim >::refRatio ( int  level) const
inline
Parameters
level
Returns
the refinement ratio of this level to the next higher one

Definition at line 758 of file BoxLibLayout.hpp.

Referenced by BoxLibLayout< T, Dim >::BoxLibLayout().

template<class T, unsigned Dim>
void BoxLibLayout< T, Dim >::resize ( int  maxLevel)
inline

The particles live initially on the coarsest level. Furthermore, the order the OPAL input file is parsed does not allow us to know the maximum level of the computation. This is known after that the FieldSolver is initialized. Therefore, we need to update the size of the ParGDB containers.

Parameters
maxLevelis set when the FieldSolver is initialized

Definition at line 224 of file BoxLibLayout.h.

References BoxLibLayout< T, Dim >::maxLevel(), ParticleAmrLayout< T, Dim >::maxLevel_m, and BoxLibLayout< T, Dim >::refRatio_m.

Referenced by AmrBoxLib::AmrBoxLib().

Here is the call graph for this function:

template<class T , unsigned Dim>
void BoxLibLayout< T, Dim >::setBoundingBox ( double  dh)
virtual

This method is used when creating the AMR object. OPAL takes the input argument BBOXINCR that is specified in the field solver command. Up to this point the AMR object is not yet initialized. After that this method shouldn't be called anymore.

Parameters
dhis the mesh enlargement factor

Implements ParticleAmrLayout< T, Dim >.

Definition at line 91 of file BoxLibLayout.hpp.

Referenced by AmrPartBunch::set_meshEnlargement().

template<class T , unsigned Dim>
void BoxLibLayout< T, Dim >::setDomainRatio ( const std::vector< double > &  ratio)

The Poisson computation domain is per default [-1,1]^3. With this method this can be changed in order to account for the different extent per direction.

Parameters
ratio,ife.g. ratio = [1, 2, 0.75], then the box is [-1, 1] x [-2, 2] x [-0.75, 0.75]

Definition at line 102 of file BoxLibLayout.hpp.

References Dim.

Referenced by AmrPartBunch::setAmrDomainRatio().

template<class T, unsigned Dim>
void BoxLibLayout< T, Dim >::update ( IpplParticleBase< BoxLibLayout< T, Dim > > &  PData,
const ParticleAttrib< char > *  canSwap = 0 
)

This method shouldn't be called. Otherwise it throws an exception.

Definition at line 131 of file BoxLibLayout.hpp.

template<class T, unsigned Dim>
void BoxLibLayout< T, Dim >::update ( AmrParticleBase< BoxLibLayout< T, Dim > > &  PData,
int  lev_min = 0,
int  lev_max = -1,
bool  isRegrid = false 
)

The proper update method for AMR.

Parameters
PDatais basically the bunch
lev_minbase level to update
lev_maxfinest level to update (if -1 update all levels starting from lev_min)
isRegridin regrid mode a level might be deleted, thus the particle level counter might be invalidated –> we need to iterate through all particles.

Definition at line 145 of file BoxLibLayout.hpp.

References MsgBuffer::add(), IpplInfo::Comm, Communicate::COMM_ANY_NODE, MsgBuffer::get(), MsgBuffer::getBuffer(), IpplInfo::getComm(), IpplInfo::getNodes(), MsgBuffer::getSize(), IpplInfo::myNode(), TagMaker::next_tag(), P_LAYOUT_CYCLE, P_SPATIAL_TRANSFER_TAG, Communicate::raw_isend(), and Communicate::raw_probe_receive().

Here is the call graph for this function:

template<class T, unsigned Dim>
bool BoxLibLayout< T, Dim >::Where ( AmrParticleBase< BoxLibLayout< T, Dim > > &  p,
const unsigned int  ip,
int  lev_min = 0,
int  lev_max = -1,
int  nGrow = 0 
) const
private

Function from AMReX adjusted to work with Ippl AmrParticleBase class Checks/sets a particles location on levels lev_min and higher.

Parameters
pis the bunch information
ipis the local (i.e. to a core) particle index
lev_minto check
lev_maxto check
nGrowis the number of ghost cells
Returns
false if the particle does not exist on that level.

Definition at line 510 of file BoxLibLayout.hpp.

Member Data Documentation

template<class T, unsigned Dim>
Vector_t BoxLibLayout< T, Dim >::lowerBound = - Vector_t(1.0, 1.0, 1.0)
static

Lower physical domain boundary (each dimension). It has to be smaller than -1 since all particles are within \([-1, 1]^3\). The real computational domain is multiplied with the mesh enlargement factor (in [%]) in BoxLibLayout::initBaseBox_m().

Definition at line 69 of file BoxLibLayout.h.

template<class T, unsigned Dim>
std::vector<std::unique_ptr<mask_t> > BoxLibLayout< T, Dim >::masks_m
private

Refinement ratios [0:finest_level-1].

Definition at line 373 of file BoxLibLayout.h.

template<class T, unsigned Dim>
AmrIntVectContainer_t BoxLibLayout< T, Dim >::refRatio_m
private
template<class T, unsigned Dim>
Vector_t BoxLibLayout< T, Dim >::upperBound = Vector_t(1.0, 1.0, 1.0)
static

Upper physical domain boundary (each dimension). It has to be greater than 1 since all particles are within \([-1, 1]^3\). The real computational domain is multiplied with the mesh enlargement factor (in [%]) in BoxLibLayout::initBaseBox_m().

Definition at line 76 of file BoxLibLayout.h.


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