IBAMR  IBAMR version 0.19.
Public Member Functions | Private Member Functions | Private Attributes | List of all members
SAMRAI::geom::CartesianGridGeometry< DIM > Class Template Reference

#include <ibamr/Wall.h>

Inheritance diagram for SAMRAI::geom::CartesianGridGeometry< DIM >:
Inheritance graph
[legend]

Public Member Functions

 CartesianGridGeometry (const std::string &object_name, tbox::Pointer< tbox::Database > input_db, bool register_for_restart=true)
 
 CartesianGridGeometry (const std::string &object_name, const double *x_lo, const double *x_up, const hier::BoxArray< DIM > &domain, bool register_for_restart=true)
 
virtual ~CartesianGridGeometry ()
 
tbox::Pointer< hier::GridGeometry< DIM > > makeRefinedGridGeometry (const std::string &fine_geom_name, const hier::IntVector< DIM > &refine_ratio, bool register_for_restart) const
 
tbox::Pointer< hier::GridGeometry< DIM > > makeCoarsenedGridGeometry (const std::string &coarse_geom_name, const hier::IntVector< DIM > &coarsen_ratio, bool register_for_restart) const
 
void setGeometryDataOnPatch (hier::Patch< DIM > &patch, const hier::IntVector< DIM > &ratio_to_level_zero, const tbox::Array< tbox::Array< bool > > &touches_regular_bdry, const tbox::Array< tbox::Array< bool > > &touches_periodic_bdry) const
 
void setGeometryData (const double *x_lo, const double *x_up, const hier::BoxArray< DIM > &domain)
 
const doublegetDx () const
 
const doublegetXLower () const
 
const doublegetXUpper () const
 
virtual void printClassData (std::ostream &os) const
 
virtual void putToDatabase (tbox::Pointer< tbox::Database > db)
 
virtual void addSpatialCoarsenOperator (tbox::Pointer< CoarsenOperator< DIM > > coarsen_op)
 
virtual void addSpatialRefineOperator (tbox::Pointer< RefineOperator< DIM > > refine_op)
 
virtual void addTimeInterpolateOperator (tbox::Pointer< TimeInterpolateOperator< DIM > > time_op)
 
virtual tbox::Pointer< CoarsenOperator< DIM > > lookupCoarsenOperator (const tbox::Pointer< hier::Variable< DIM > > &var, const std::string &op_name) const
 
virtual tbox::Pointer< RefineOperator< DIM > > lookupRefineOperator (const tbox::Pointer< hier::Variable< DIM > > &var, const std::string &op_name) const
 
virtual tbox::Pointer< TimeInterpolateOperator< DIM > > lookupTimeInterpolateOperator (const tbox::Pointer< hier::Variable< DIM > > &var, const std::string &op_name="STD_LINEAR_TIME_INTERPOLATE") const
 

Private Member Functions

void getFromInput (tbox::Pointer< tbox::Database > db, bool is_from_restart)
 
void getFromRestart ()
 
void makeStandardOperators ()
 

Private Attributes

std::string d_object_name
 
bool d_registered_for_restart
 
double d_dx [DIM]
 
double d_x_lo [DIM]
 
double d_x_up [DIM]
 
hier::Box< DIM > d_domain_box
 
bool d_using_original_locations
 
tbox::List< tbox::Pointer< CoarsenOperator< DIM > > > d_coarsen_operators
 
tbox::List< tbox::Pointer< RefineOperator< DIM > > > d_refine_operators
 
tbox::List< tbox::Pointer< TimeInterpolateOperator< DIM > > > d_time_operators
 

Functions for computing boundary boxes

BoxArray< DIM > d_physical_domain
 
bool d_domain_is_single_box
 
IntVector< DIM > d_periodic_shift
 
IntVector< DIM > d_max_data_ghost_width
 
void findPatchesTouchingBoundaries (tbox::Array< tbox::Array< tbox::Array< bool > > > &touches_regular_bdry, tbox::Array< tbox::Array< tbox::Array< bool > > > &touches_periodic_bdry, const PatchLevel< DIM > &level, const IntVector< DIM > &periodic_shift, const BoxArray< DIM > &domain) const
 Determine for every patch on a level if it touches a regular physical boundary or a periodic boundary. More...
 
virtual void setGeometryOnPatches (hier::PatchLevel< DIM > &level, const hier::IntVector< DIM > &ratio_to_level_zero, tbox::Array< tbox::Array< tbox::Array< bool > > > &touches_regular_bdry, tbox::Array< tbox::Array< tbox::Array< bool > > > &touches_periodic_bdry, bool defer_boundary_box_creation)
 Pass the arrays holding the boundary information to be stored in the concrete geometry classes, and construct boundary boxes if required. More...
 
void setBoundaryBoxes (hier::PatchLevel< DIM > &level)
 Construct the boundary boxes for each patch and set them on the patch geometries. More...
 
void computeShiftsForLevel (tbox::Array< tbox::List< IntVector< DIM > > > &shifts, const PatchLevel< DIM > &level, const BoxArray< DIM > &physical_domain) const
 
void computePhysicalDomain (BoxArray< DIM > &domain, const IntVector< DIM > &ratio_to_level_zero) const
 
void setPhysicalDomain (const BoxArray< DIM > &domain)
 
const BoxArray< DIM > & getPhysicalDomain () const
 
bool getDomainIsSingleBox () const
 
void initializePeriodicShift (const IntVector< DIM > &directions)
 
IntVector< DIM > getPeriodicShift (const IntVector< DIM > &ratio_to_level_zero=IntVector< DIM >(1)) const
 
IntVector< DIM > computeMaxGhostWidth (tbox::Pointer< PatchDescriptor< DIM > > descriptor)
 Compute the maximum ghost width of all of the components associated with the patch descriptor. More...
 
void computeBoundaryBoxesOnLevel (tbox::Array< BoundaryBox< DIM > > boundaries[], const PatchLevel< DIM > &level, const IntVector< DIM > &periodic_shift, const IntVector< DIM > &ghost_width, const BoxArray< DIM > &domain, bool do_all_patches=false) const
 Compute boundary boxes for each patch in patch level and assign them to the array of boundary box arrays, assumed to be of length DIM * (num patches). More...
 
void getBoundaryBoxes (tbox::Array< BoundaryBox< DIM > > boundaries[DIM], const Box< DIM > &box, const BoxArray< DIM > &domain_boxes, const IntVector< DIM > &ghosts, const IntVector< DIM > &periodic_shift) const
 Compute boundary boxes for patch. More...
 
bool checkPeriodicValidity (const BoxArray< DIM > &domain)
 Check that the domain is valid for periodic boundary conditions. More...
 
bool checkBoundaryBox (const BoundaryBox< DIM > &boundary_box, const Patch< DIM > &patch, const BoxArray< DIM > &domain, const int num_per_dirs, const IntVector< DIM > &max_data_ghost_width) const
 Check on each BoundaryBox when it is created. More...
 
void computeShiftsForPatch (tbox::List< IntVector< DIM > > &shifts, const Box< DIM > &box, const BoxArray< DIM > &domain, const IntVector< DIM > &periodic_shift) const
 Find if a box is on a periodic boundary and compute its shifts. More...
 

Detailed Description

template<int DIM>
class SAMRAI::geom::CartesianGridGeometry< DIM >

Class CartesianGridGeometry<DIM> provides simple Cartesian mesh geometry management on an AMR hierarchy. The mesh coordinates on each hierarchy level are limited to mesh increments specified as DIM-tuple (dx[0],...,dx[DIM-1]) and spatial coordinates of the lower and upper corners of the smallest parallelepiped bounding the entire computational domain. The mesh increments on each level are defined with respect to the coarsest hierarchy level and multiplying those values by the proper refinement ratio. This class sets geometry information on each patch in an AMR hierarchy. This class is derived from the xfer::Geometry<DIM> base class which is further derived from the hier::GridGeometry<DIM> base class.

An object of this class requires numerous parameters to be read from input. Also, data must be written to and read from files for restart. The input and restart data are summarized as follows:

Required input keys and data types:

 - \b    domain_boxes   
    tbox::Array of boxes representing the index space for the entire
    domain (on the coarsest refinement level).

 - \b    x_lo   
    tbox::Array of double values representing the spatial coordinates of
    the lower corner of the physical domain.

 - \b    x_up   
    tbox::Array of double values representing the spatial coordinates of
    the upper corner of the physical domain.

Optional input keys, data types, and defaults:

 - \b    periodic_dimension   
    tbox::Array of integer values representing the directions in which
    the physical domain is periodic.  A non-zero value indicates
    that the direction is periodic.  A zero value indicates that  
    the direction is not periodic.  If no values are specified, then
    the array is initialized to all zeros (no periodic directions).

 - \b    use_original_location_indices 
    Boolean argument to handle backward compatibility with new
    location index scheme for codimension 2 in 3 dimensions.
    Set to true to use the location index scheme that was used in
    SAMRAI v. 1.4 and earlier versions, and false to use the new
    scheme.  In 3 dimensions, this key defaults to true if not
    present in the input.  In all other dimensions, this key
    is irrelevant.

No input values can overwrite restart values.

A sample input file for a two-dimensional problem might look like:

*
*    domain_boxes = [(0,0) , (49,39)]
*    x_lo = 0.0 , 0.0
*    x_up = 50.0 , 40.0
*    periodic_dimension = 0, 1  // periodic in y only
*
* 

This generates a two-dimensional rectangular domain periodic in the y-direction, and having 50 cells in the x-direction and 40 cells in the y-direction, with the cell size 1 unit in each direction.

See also
xfer::Geometry
hier::GridGeometry

Constructor & Destructor Documentation

◆ CartesianGridGeometry() [1/2]

template<int DIM>
SAMRAI::geom::CartesianGridGeometry< DIM >::CartesianGridGeometry ( const std::string &  object_name,
tbox::Pointer< tbox::Database input_db,
bool  register_for_restart = true 
)

Constructor for CartesianGridGeometry<DIM> initializes data members based on parameters read from the specified input database or from the restart database corresponding to the specified object name. The constructor also registers this object for restart using the specified object name when the boolean argument is true. Whether object will write its state to restart files during program execution is determined by this argument.
Note that it has a default state of true.

Errors: passing in a null database pointer or an empty string will result in an unrecoverable assertion.

◆ CartesianGridGeometry() [2/2]

template<int DIM>
SAMRAI::geom::CartesianGridGeometry< DIM >::CartesianGridGeometry ( const std::string &  object_name,
const double x_lo,
const double x_up,
const hier::BoxArray< DIM > &  domain,
bool  register_for_restart = true 
)

Constructor for CartesianGridGeometry<DIM> sets data members based on arguments. The constructor also registers this object for restart using the specified object name when the boolean argument is true. Whether object will write its state to restart files during program execution is determined by this argument.
Note that it has a default state of true.

Errors: passing in an empty string, or null data pointers will result in an unrecoverable assertion.

◆ ~CartesianGridGeometry()

template<int DIM>
virtual SAMRAI::geom::CartesianGridGeometry< DIM >::~CartesianGridGeometry ( )
virtual

Destructor for CartesianGridGeometry<DIM> deallocates data describing grid geometry and unregisters the object with the restart manager if previously registered.

Member Function Documentation

◆ makeRefinedGridGeometry()

template<int DIM>
tbox::Pointer<hier::GridGeometry<DIM> > SAMRAI::geom::CartesianGridGeometry< DIM >::makeRefinedGridGeometry ( const std::string &  fine_geom_name,
const hier::IntVector< DIM > &  refine_ratio,
bool  register_for_restart 
) const
virtual

Create and return a pointer to a refined version of this Cartesian grid geometry object. This function is pure virtual in the hier::GridGeometry<DIM> base class.

Implements SAMRAI::hier::GridGeometry< DIM >.

◆ makeCoarsenedGridGeometry()

template<int DIM>
tbox::Pointer<hier::GridGeometry<DIM> > SAMRAI::geom::CartesianGridGeometry< DIM >::makeCoarsenedGridGeometry ( const std::string &  coarse_geom_name,
const hier::IntVector< DIM > &  coarsen_ratio,
bool  register_for_restart 
) const
virtual

Create and return a pointer to a coarsened version of this Cartesian grid geometry object. This function is pure virtual in the hier::GridGeometry<DIM> base class.

Implements SAMRAI::hier::GridGeometry< DIM >.

◆ setGeometryDataOnPatch()

template<int DIM>
void SAMRAI::geom::CartesianGridGeometry< DIM >::setGeometryDataOnPatch ( hier::Patch< DIM > &  patch,
const hier::IntVector< DIM > &  ratio_to_level_zero,
const tbox::Array< tbox::Array< bool > > &  touches_regular_bdry,
const tbox::Array< tbox::Array< bool > > &  touches_periodic_bdry 
) const
virtual

◆ setGeometryData()

template<int DIM>
void SAMRAI::geom::CartesianGridGeometry< DIM >::setGeometryData ( const double x_lo,
const double x_up,
const hier::BoxArray< DIM > &  domain 
)

Set data members for this CartesianGridGeometry<DIM> object.

◆ getDx()

template<int DIM>
const double* SAMRAI::geom::CartesianGridGeometry< DIM >::getDx ( ) const

Return const pointer to dx array for reference level in hierarchy.

◆ getXLower()

template<int DIM>
const double* SAMRAI::geom::CartesianGridGeometry< DIM >::getXLower ( ) const

Return const pointer to lower spatial coordinate for reference level in hierarchy.

◆ getXUpper()

template<int DIM>
const double* SAMRAI::geom::CartesianGridGeometry< DIM >::getXUpper ( ) const

Return const pointer to upper spatial coordinate for reference level in hierarchy.

◆ printClassData()

template<int DIM>
virtual void SAMRAI::geom::CartesianGridGeometry< DIM >::printClassData ( std::ostream &  os) const
virtual

Print class data representation.

Reimplemented from SAMRAI::xfer::Geometry< DIM >.

◆ putToDatabase()

template<int DIM>
virtual void SAMRAI::geom::CartesianGridGeometry< DIM >::putToDatabase ( tbox::Pointer< tbox::Database db)
virtual

Writes the state of the CartesianGridGeometry object to the database.

When assertion checking is active, db cannot be a null database pointer.

Implements SAMRAI::tbox::Serializable.

◆ getFromInput()

template<int DIM>
void SAMRAI::geom::CartesianGridGeometry< DIM >::getFromInput ( tbox::Pointer< tbox::Database db,
bool  is_from_restart 
)
private

◆ getFromRestart()

template<int DIM>
void SAMRAI::geom::CartesianGridGeometry< DIM >::getFromRestart ( )
private

◆ makeStandardOperators()

template<int DIM>
void SAMRAI::geom::CartesianGridGeometry< DIM >::makeStandardOperators ( )
private

◆ addSpatialCoarsenOperator()

template<int DIM>
virtual void SAMRAI::xfer::Geometry< DIM >::addSpatialCoarsenOperator ( tbox::Pointer< CoarsenOperator< DIM > >  coarsen_op)
virtualinherited

Add concrete spatial coarsening operator instance to appropriate lookup list. Note that each concrete operator must implement a lookup function through which it can be identified.

◆ addSpatialRefineOperator()

template<int DIM>
virtual void SAMRAI::xfer::Geometry< DIM >::addSpatialRefineOperator ( tbox::Pointer< RefineOperator< DIM > >  refine_op)
virtualinherited

Add concrete spatial refinement operator instance to appropriate lookup list. Note that each concrete operator must implement a lookup function through which it can be identified.

◆ addTimeInterpolateOperator()

template<int DIM>
virtual void SAMRAI::xfer::Geometry< DIM >::addTimeInterpolateOperator ( tbox::Pointer< TimeInterpolateOperator< DIM > >  time_op)
virtualinherited

Add concrete time interpolation operator instance to appropriate lookup list. Note that each concrete operator must implement a lookup function through which it can be identified.

◆ lookupCoarsenOperator()

template<int DIM>
virtual tbox::Pointer< CoarsenOperator<DIM> > SAMRAI::xfer::Geometry< DIM >::lookupCoarsenOperator ( const tbox::Pointer< hier::Variable< DIM > > &  var,
const std::string &  op_name 
) const
virtualinherited

Search list for the spatial coarsening operator matching the request for the given variable. If the operator is found, a pointer to it will be returned. Otherwise, an unrecoverable error will result and the program will abort.

◆ lookupRefineOperator()

template<int DIM>
virtual tbox::Pointer< RefineOperator<DIM> > SAMRAI::xfer::Geometry< DIM >::lookupRefineOperator ( const tbox::Pointer< hier::Variable< DIM > > &  var,
const std::string &  op_name 
) const
virtualinherited

Search list for the spatial refinement operator matching the request for the given variable. If the operator is found, a pointer to it will be returned. Otherwise, an unrecoverable error will result and the program will abort.

◆ lookupTimeInterpolateOperator()

template<int DIM>
virtual tbox::Pointer< TimeInterpolateOperator<DIM> > SAMRAI::xfer::Geometry< DIM >::lookupTimeInterpolateOperator ( const tbox::Pointer< hier::Variable< DIM > > &  var,
const std::string &  op_name = "STD_LINEAR_TIME_INTERPOLATE" 
) const
virtualinherited

Search list for the time interpolation operator matching the request for the given variable. If the operator is found, a pointer to it will be returned. Otherwise, an unrecoverable error will result and the program will abort.

◆ findPatchesTouchingBoundaries()

template<int DIM>
void SAMRAI::hier::GridGeometry< DIM >::findPatchesTouchingBoundaries ( tbox::Array< tbox::Array< tbox::Array< bool > > > &  touches_regular_bdry,
tbox::Array< tbox::Array< tbox::Array< bool > > > &  touches_periodic_bdry,
const PatchLevel< DIM > &  level,
const IntVector< DIM > &  periodic_shift,
const BoxArray< DIM > &  domain 
) const
inherited

This routine loops through all of the patches on the given level and determines which kinds of boundaries each patch touches. The 3-dimensional boolean arrays are set to store for each path whether it touches a regular boundary, a periodic boundary, both, or neither.

The array arguments should be uninitialized when they are passed into this function.

Parameters
touches_regular_bdryArray to store which patches touch non-periodic boundaries.
touches_periodic_bdryArray to store which patches touch periodic boundaries.
levelcontaining the patches to be checked
periodic_shiftperiodic shift for the level (see getPeriodicShift)
domainPhysical domain (at the same level of refinement as level)

◆ setGeometryOnPatches()

template<int DIM>
virtual void SAMRAI::hier::GridGeometry< DIM >::setGeometryOnPatches ( hier::PatchLevel< DIM > &  level,
const hier::IntVector< DIM > &  ratio_to_level_zero,
tbox::Array< tbox::Array< tbox::Array< bool > > > &  touches_regular_bdry,
tbox::Array< tbox::Array< tbox::Array< bool > > > &  touches_periodic_bdry,
bool  defer_boundary_box_creation 
)
virtualinherited

This routine will pass the arrays containing the information about which patches touch which boundaries to the concrete grid geometry class. Also, if defer_boundary_box_creation is false, this routine will call a routine to construct all of the boundary boxes for the patches on level and will set them in the patch geometry for each patch.

Parameters
levelcontaining the patches to be checked.
ratio_to_level_zeroratio to the coarsest level.
touches_regular_bdryArray storing which patches touch non-periodic boundaries.
touches_periodic_bdryArray storing which patches touch periodic boundaries.
defer_boundary_box_creationBoundary boxes will be created here if false, not if true.

◆ setBoundaryBoxes()

template<int DIM>
void SAMRAI::hier::GridGeometry< DIM >::setBoundaryBoxes ( hier::PatchLevel< DIM > &  level)
inherited

This routine constructs the boundary boxes for every patch in the level. Once constructed, the boundary boxes are set on each patch's PatchGeometry object.

Parameters
levelThe level for which boundary boxes are constructed.

◆ computeShiftsForLevel()

template<int DIM>
void SAMRAI::hier::GridGeometry< DIM >::computeShiftsForLevel ( tbox::Array< tbox::List< IntVector< DIM > > > &  shifts,
const PatchLevel< DIM > &  level,
const BoxArray< DIM > &  physical_domain 
) const
inherited

Compute the valid periodic shifts for each patch on a level. The shifts array will store a list of IntVectors for each patch. Each list will contain the valid possible periodic shifts for each particular patch. If there are no periodic boundary conditions or a patch does not touch a periodic boundary, the list for a patch will be empty.

The patch geometry object for each patch on the level must be properly initialized before calling this routine. This is typically done in the patch level constructor.

When assertion checking is active, the array of shifts in the argument list must have the same length as the number of patches on the level.

◆ computePhysicalDomain()

template<int DIM>
void SAMRAI::hier::GridGeometry< DIM >::computePhysicalDomain ( BoxArray< DIM > &  domain,
const IntVector< DIM > &  ratio_to_level_zero 
) const
inherited

Compute physical domain box array describing the index space of the physical domain managed by this geometry object. If any entry of ratio vector is negative, the index space is coarsened with respect to the physical domain description. Otherwise, the index space is refined.

◆ setPhysicalDomain()

template<int DIM>
void SAMRAI::hier::GridGeometry< DIM >::setPhysicalDomain ( const BoxArray< DIM > &  domain)
inherited

Set physical domain to input box array and determine whether domain is a single box.

◆ getPhysicalDomain()

template<int DIM>
const BoxArray<DIM>& SAMRAI::hier::GridGeometry< DIM >::getPhysicalDomain ( ) const
inherited

Return const reference to physical domain description for level 0.

◆ getDomainIsSingleBox()

template<int DIM>
bool SAMRAI::hier::GridGeometry< DIM >::getDomainIsSingleBox ( ) const
inherited

Return boolean value indicating whether the physical domain can be represented as a single box.

◆ initializePeriodicShift()

template<int DIM>
void SAMRAI::hier::GridGeometry< DIM >::initializePeriodicShift ( const IntVector< DIM > &  directions)
inherited

Initialize the periodic shift on the coarsest level. The IntVector argument should be set to 1 for periodic directions and 0 for all other directions. The shift will be calculated to be the number of cells in the periodic direction and zero in all other directions.

◆ getPeriodicShift()

template<int DIM>
IntVector<DIM> SAMRAI::hier::GridGeometry< DIM >::getPeriodicShift ( const IntVector< DIM > &  ratio_to_level_zero = IntVector< DIM >(1)) const
inherited

Return IntVector<DIM> containing the periodic shift in each direction for a domain represented by a refinement of the reference physical domain (i.e. level zero) by the given ratio vector. tbox::Array entries will be zero for non-periodic directions. By default (i.e., when no argument is passed, the function returns the periodic shift for level zero in the hierarchy.

◆ computeMaxGhostWidth()

template<int DIM>
IntVector<DIM> SAMRAI::hier::GridGeometry< DIM >::computeMaxGhostWidth ( tbox::Pointer< PatchDescriptor< DIM > >  descriptor)
inherited

Calculates the maximum ghost width for all the variables associated with the patch descriptor. This must only be called after all of the variables have been registered with the VariableDatabase. If a variable is added that changes the maximum ghost width, then an assertion failure will result.

◆ computeBoundaryBoxesOnLevel()

template<int DIM>
void SAMRAI::hier::GridGeometry< DIM >::computeBoundaryBoxesOnLevel ( tbox::Array< BoundaryBox< DIM > >  boundaries[],
const PatchLevel< DIM > &  level,
const IntVector< DIM > &  periodic_shift,
const IntVector< DIM > &  ghost_width,
const BoxArray< DIM > &  domain,
bool  do_all_patches = false 
) const
inherited

The DIM arrays of boundary boxes for each patch will be stored in groups of DIM. For example, in 3d with n patches on the level, the array For example, in 3d with n patches on the level, the array of boundary box arrays will be ordered as follows:

* (patch 0 face array, patch 0 edge array, patch 0 node array,
*  patch 1 face array, patch 1 edge array, patch 1 node array, . . . ,
*  patch n-1 face array, patch n-1 edge array, patch n-1 node array)
* 

The optional argument do_all_patches defaults to false, in which case the boundary box computation is executed only on patches that touch a non-periodic boundary. When this routine is called during patch level construction to describe a physical boundary, it is known that only patches that touch a non-periodic boundary will have non-empty sets of boundary boxes, so for efficiency's sake the boundary box box computation is supressed for all other patches. When this routine is called to create boundary boxes that describe a coarse-fine boundary, the computation must occur for every patch, so do_all_patches mush be set to true.

Parameters
boundariesoutput boundary description
levellevel on which to generate boundaries
periodic_shiftperiodic shift for the level (see getPeriodicShift)
ghost_widthghost width to compute geometry for
domainPhysical domain (in index space of level) for computing boundary boxes.
do_all_patchesExecute boundary box computation on all patches, even those known to not touch a boundary

◆ getBoundaryBoxes()

template<int DIM>
void SAMRAI::hier::GridGeometry< DIM >::getBoundaryBoxes ( tbox::Array< BoundaryBox< DIM > >  boundaries[DIM],
const Box< DIM > &  box,
const BoxArray< DIM > &  domain_boxes,
const IntVector< DIM > &  ghosts,
const IntVector< DIM > &  periodic_shift 
) const
inherited

Decompose patch boundary region into pieces depending on spatial dim. Boxes are extended along the boundary to the edge of the ghost layer if necessary.

◆ checkPeriodicValidity()

template<int DIM>
bool SAMRAI::hier::GridGeometry< DIM >::checkPeriodicValidity ( const BoxArray< DIM > &  domain)
privateinherited

◆ checkBoundaryBox()

template<int DIM>
bool SAMRAI::hier::GridGeometry< DIM >::checkBoundaryBox ( const BoundaryBox< DIM > &  boundary_box,
const Patch< DIM > &  patch,
const BoxArray< DIM > &  domain,
const int  num_per_dirs,
const IntVector< DIM > &  max_data_ghost_width 
) const
privateinherited

This is a check performed on each BoundaryBox when it is created. It returns true when a BoundaryBox has a width of 1 in at least one direction, is adjacent to the patch boundary (possible extended into the patch's ghost region) and is outside the physical domain.

◆ computeShiftsForPatch()

template<int DIM>
void SAMRAI::hier::GridGeometry< DIM >::computeShiftsForPatch ( tbox::List< IntVector< DIM > > &  shifts,
const Box< DIM > &  box,
const BoxArray< DIM > &  domain,
const IntVector< DIM > &  periodic_shift 
) const
privateinherited

If box is located on a periodic boundary, all of its possible shifts will be computed and stored in shifts. If box is not on a periodic boundary, shifts will be an empty list.

Member Data Documentation

◆ d_object_name

template<int DIM>
std::string SAMRAI::geom::CartesianGridGeometry< DIM >::d_object_name
private

◆ d_registered_for_restart

template<int DIM>
bool SAMRAI::geom::CartesianGridGeometry< DIM >::d_registered_for_restart
private

◆ d_dx

template<int DIM>
double SAMRAI::geom::CartesianGridGeometry< DIM >::d_dx[DIM]
private

◆ d_x_lo

template<int DIM>
double SAMRAI::geom::CartesianGridGeometry< DIM >::d_x_lo[DIM]
private

◆ d_x_up

template<int DIM>
double SAMRAI::geom::CartesianGridGeometry< DIM >::d_x_up[DIM]
private

◆ d_domain_box

template<int DIM>
hier::Box<DIM> SAMRAI::geom::CartesianGridGeometry< DIM >::d_domain_box
private

◆ d_using_original_locations

template<int DIM>
bool SAMRAI::geom::CartesianGridGeometry< DIM >::d_using_original_locations
private

◆ d_coarsen_operators

template<int DIM>
tbox::List< tbox::Pointer< CoarsenOperator<DIM> > > SAMRAI::xfer::Geometry< DIM >::d_coarsen_operators
privateinherited

◆ d_refine_operators

template<int DIM>
tbox::List< tbox::Pointer< RefineOperator<DIM> > > SAMRAI::xfer::Geometry< DIM >::d_refine_operators
privateinherited

◆ d_time_operators

template<int DIM>
tbox::List< tbox::Pointer< TimeInterpolateOperator<DIM> > > SAMRAI::xfer::Geometry< DIM >::d_time_operators
privateinherited

◆ d_physical_domain

template<int DIM>
BoxArray<DIM> SAMRAI::hier::GridGeometry< DIM >::d_physical_domain
privateinherited

Box array defining computational domain on coarsest level and boolean flag that is true when domain is a single box.

◆ d_domain_is_single_box

template<int DIM>
bool SAMRAI::hier::GridGeometry< DIM >::d_domain_is_single_box
privateinherited

◆ d_periodic_shift

template<int DIM>
IntVector<DIM> SAMRAI::hier::GridGeometry< DIM >::d_periodic_shift
privateinherited

Integer array vector describing periodic shift coarsest level. An entry of zero means direction is not periodic.

◆ d_max_data_ghost_width

template<int DIM>
IntVector<DIM> SAMRAI::hier::GridGeometry< DIM >::d_max_data_ghost_width
privateinherited

Current maximum ghost cell width over all patch data objects known to the patch descriptor. This is used to compute boundary boxes.


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