IBAMR  IBAMR version 0.19.
Public Member Functions | List of all members
SAMRAI::solv::CartesianRobinBcHelper< DIM > Class Template Reference

Helper utility for setting Robin boundary conditions. More...

#include <CartesianRobinBcHelper.h>

Inheritance diagram for SAMRAI::solv::CartesianRobinBcHelper< DIM >:
Inheritance graph
[legend]

Public Member Functions

 CartesianRobinBcHelper (std::string object_name=std::string(), RobinBcCoefStrategy< DIM > *coef_strategy=NULL)
 Constructor using. More...
 
virtual ~CartesianRobinBcHelper ()
 Destructor. More...
 
virtual void setPhysicalBoundaryConditions (hier::Patch< DIM > &patch, const double fill_time, const hier::IntVector< DIM > &ghost_width_to_fill)
 
hier::IntVector< DIM > getRefineOpStencilWidth () const
 
virtual void preprocessRefineBoxes (hier::Patch< DIM > &fine, const hier::Patch< DIM > &coarse, const hier::BoxList< DIM > &fine_boxes, const hier::IntVector< DIM > &ratio)
 
virtual void preprocessRefine (hier::Patch< DIM > &fine, const hier::Patch< DIM > &coarse, const hier::Box< DIM > &fine_box, const hier::IntVector< DIM > &ratio)
 
virtual void postprocessRefineBoxes (hier::Patch< DIM > &fine, const hier::Patch< DIM > &coarse, const hier::BoxList< DIM > &fine_boxes, const hier::IntVector< DIM > &ratio)
 
virtual void postprocessRefine (hier::Patch< DIM > &fine, const hier::Patch< DIM > &coarse, const hier::Box< DIM > &fine_box, const hier::IntVector< DIM > &ratio)
 
Functions to set boundary condition values
void setBoundaryValuesInCells (hier::Patch< DIM > &patch, const double fill_time, const hier::IntVector< DIM > &ghost_width_to_fill, int target_data_id, bool homogeneous_bc=false) const
 Set the physical boundary condition by setting the value of the first ghost cells. More...
 
void setBoundaryValuesInCells (hier::PatchLevel< DIM > &level, const double fill_time, const hier::IntVector< DIM > &ghost_width_to_fill, int target_data_id, bool homogeneous_bc=false) const
 Set ghost cells for an entire level. More...
 
void setBoundaryValuesAtNodes (hier::Patch< DIM > &patch, const double fill_time, int target_data_id, bool homogeneous_bc=false) const
 Set the physical boundary condition by setting the value of the boundary nodes. More...
 

Ways to provide the Robin bc coefficients

std::string d_object_name
 
const RobinBcCoefStrategy< DIM > * d_coef_strategy
 Coefficient strategy giving a way to get to Robin bc coefficients. More...
 
int d_target_data_id
 hier::Index of target patch data when filling ghosts. More...
 
bool d_homogeneous_bc
 Whether to assumg g=0 when filling ghosts. More...
 
tbox::Pointer< tbox::Timert_set_boundary_values_in_cells
 Timers for performance measurement. More...
 
tbox::Pointer< tbox::Timert_use_set_bc_coefs
 
void setCoefImplementation (const RobinBcCoefStrategy< DIM > *coef_strategy)
 Provide an implementation of the RobinBcCoefStrategy<DIM> for determining the boundary coefficients. More...
 
void setTargetDataId (int target_data_id)
 Set the data id that should be filled when setting physical boundary conditions. More...
 
void setHomogeneousBc (bool homogeneous_bc)
 Set whether boundary filling should assume homogeneous conditions. More...
 
hier::BoundaryBox< DIM > trimBoundaryBox (const hier::BoundaryBox< DIM > &boundary_box, const hier::Box< DIM > &limit_box) const
 Trim a boundary box so that it does not stick out past a given box. More...
 
hier::Box< DIM > makeNodeBoundaryBox (const hier::BoundaryBox< DIM > &boundary_box) const
 Return box describing the index space of boundary nodes defined by a boundary box. More...
 
hier::Box< DIM > makeFaceBoundaryBox (const hier::BoundaryBox< DIM > &boundary_box) const
 Return box describing the index space of faces defined by a boundary box. More...
 

Detailed Description

template<int DIM>
class SAMRAI::solv::CartesianRobinBcHelper< DIM >

This class is intended as a helper for performing the tedious task of setting boundary values for scalar quantities for the general case of boundary conditions known as the Robin boundary condition.

It uses the Robin boundary condition coefficients specified by a RobinBcCoefStrategy<DIM> object to determine the boundary value to set. The exact value set depends on the allignment of the data and is derived from various discrete approximations of the Robin formula. This class currently supports cell-centered alignment and will support node-centered alignment in the future.

See RobinBcCoefStrategy<DIM> for the description of the Robin boundary condition.

This class inherits and implements virtual functions from xfer::RefinePatchStrategy<DIM> so it may be used to help create communication schedules if desired.

Constructor & Destructor Documentation

◆ CartesianRobinBcHelper()

template<int DIM>
SAMRAI::solv::CartesianRobinBcHelper< DIM >::CartesianRobinBcHelper ( std::string  object_name = std::string(),
RobinBcCoefStrategy< DIM > *  coef_strategy = NULL 
)

Requires a concrete implementation of RobinBcCoefStrategy<DIM>.

Parameters
object_nameName of the object, for general referencing.
coef_strategyCoefficients strategy being helped.

◆ ~CartesianRobinBcHelper()

template<int DIM>
virtual SAMRAI::solv::CartesianRobinBcHelper< DIM >::~CartesianRobinBcHelper ( )
virtual

Member Function Documentation

◆ setPhysicalBoundaryConditions()

template<int DIM>
virtual void SAMRAI::solv::CartesianRobinBcHelper< DIM >::setPhysicalBoundaryConditions ( hier::Patch< DIM > &  patch,
const double  fill_time,
const hier::IntVector< DIM > &  ghost_width_to_fill 
)
virtual

Pure virtual function to set data associated with the given list of patch data indices at patch boundaries that intersect the physical domain boundary. The specific boundary conditions are determined by the user. The patch data components set in this routine correspond to the "scratch" components specified in calls to the registerRefine() function in the RefineAlgorithm<DIM> class.

Parameters
patchhier::Patch on which to fill boundary data.
fill_timeDouble simulation time for boundary filling.
ghost_width_to_fillInteger vector describing maximum ghost width to fill over all registered scratch components.

Implements SAMRAI::xfer::RefinePatchStrategy< DIM >.

◆ getRefineOpStencilWidth()

template<int DIM>
hier::IntVector<DIM> SAMRAI::solv::CartesianRobinBcHelper< DIM >::getRefineOpStencilWidth ( ) const
virtual

Pure virtual function to return maximum stencil width needed over user-defined data interpolation operations. This is needed to determine the correct interpolation data dependencies.

Implements SAMRAI::xfer::RefinePatchStrategy< DIM >.

◆ preprocessRefineBoxes()

template<int DIM>
virtual void SAMRAI::solv::CartesianRobinBcHelper< DIM >::preprocessRefineBoxes ( hier::Patch< DIM > &  fine,
const hier::Patch< DIM > &  coarse,
const hier::BoxList< DIM > &  fine_boxes,
const hier::IntVector< DIM > &  ratio 
)
virtual

Virtual function to perform user-defined refine operations. This member function is called before standard refining operations (expressed using concrete subclasses of the RefineOperator<DIM> base class). The preprocess function must refine data from the scratch components of the coarse patch into the scratch components of the fine patch on the specified fine box regions.

Typically, only the pure virtual members of this class are implemented in user-defined subclasses of this base class. This version of the preprocess function operates on an entire box list. By default, this version simply loops over the box list and calls the single-box version, which is a pure virtual function.

Parameters
fineFine patch containing destination data.
coarseCoarse patch containing source data.
fine_boxestbox::List of box regions on fine patch into which data is refined.
ratioInteger vector containing ratio relating index space between coarse and fine patches.

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

◆ preprocessRefine()

template<int DIM>
virtual void SAMRAI::solv::CartesianRobinBcHelper< DIM >::preprocessRefine ( hier::Patch< DIM > &  fine,
const hier::Patch< DIM > &  coarse,
const hier::Box< DIM > &  fine_box,
const hier::IntVector< DIM > &  ratio 
)
virtual

Pure virtual function to perform user-defined preprocess data refine operations. This member function is called before standard refine operations (expressed using concrete subclasses of the RefineOperator<DIM> base class). The preprocess function must refine data from the scratch components of the coarse patch into the scratch components of the fine patch on the specified fine box region. Recall that the scratch components are specified in calls to the registerRefine() function in the RefineAlgorithm<DIM> class.

Parameters
fineFine patch containing destination data.
coarseCoarse patch containing source data.
fine_boxhier::Box region on fine patch into which data is refined.
ratioInteger vector containing ratio relating index space between coarse and fine patches.

Implements SAMRAI::xfer::RefinePatchStrategy< DIM >.

◆ postprocessRefineBoxes()

template<int DIM>
virtual void SAMRAI::solv::CartesianRobinBcHelper< DIM >::postprocessRefineBoxes ( hier::Patch< DIM > &  fine,
const hier::Patch< DIM > &  coarse,
const hier::BoxList< DIM > &  fine_boxes,
const hier::IntVector< DIM > &  ratio 
)
virtual

Virtual function to perform user-defined refine operations. This member function is called after standard refining operations (expressed using concrete subclasses of the RefineOperator<DIM> base class). The postprocess function must refine data from the scratch components of the coarse patch into the scratch components of the fine patch on the specified fine box regions.

Typically, only the pure virtual members of this class are implemented in user-defined subclasses of this base class. This version of the postprocess function operates on an entire box list. By default, this version simply loops over the box list and calls the single-box version, which is a pure virtual function.

Parameters
fineFine patch containing destination data.
coarseCoarse patch containing source data.
fine_boxestbox::List of box regions on fine patch into which data is refined.
ratioInteger vector containing ratio relating index space between coarse and fine patches.

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

◆ postprocessRefine()

template<int DIM>
virtual void SAMRAI::solv::CartesianRobinBcHelper< DIM >::postprocessRefine ( hier::Patch< DIM > &  fine,
const hier::Patch< DIM > &  coarse,
const hier::Box< DIM > &  fine_box,
const hier::IntVector< DIM > &  ratio 
)
virtual

Pure virtual function to perform user-defined postprocess data refine operations. This member function is called after standard refine operations (expressed using concrete subclasses of the RefineOperator<DIM> base class). The postprocess function must refine data from the scratch components of the coarse patch into the scratch components of the fine patch on the specified fine box region. Recall that the scratch components are specified in calls to the registerRefine() function in the RefineAlgorithm<DIM> class.

Parameters
fineFine patch containing destination data.
coarseCoarse patch containing source data.
fine_boxhier::Box region on fine patch into which data is refined.
ratioInteger vector containing ratio relating index space between coarse and fine patches.

Implements SAMRAI::xfer::RefinePatchStrategy< DIM >.

◆ setBoundaryValuesInCells() [1/2]

template<int DIM>
void SAMRAI::solv::CartesianRobinBcHelper< DIM >::setBoundaryValuesInCells ( hier::Patch< DIM > &  patch,
const double  fill_time,
const hier::IntVector< DIM > &  ghost_width_to_fill,
int  target_data_id,
bool  homogeneous_bc = false 
) const

This function has an interface similar to the virtual function xfer::RefinePatchStrategy<DIM>::setPhysicalBoundaryConditions(), and it may be used to help implement that function, but it does not serve the same purpose. The primary differences are:

  1. It specializes to cell-centered variables.
  2. Only one ghost cell width is filled. Setting a Robin boundary condition for cell-centered quantities requires only one ghost cell to be set. (More ghost cells can be filled by continuing the linear distribution of data beyond the first cell, but that is not implemented at this time.)
  3. User must specify the index of the data whose ghost cells need to be filled. This index is used to determine the variable for which to set the boundary coefficients and to get the data to be set.

This function calls RobinBcStrategy::setBcCoefs() to get the coefficients, then it sets the values in the first ghost cell on the boundary.

To determine the value for the ghost cell, a linear approximation in the direction normal to the boundary is assumed. We write the following discrete approximations:

\[ u_b = \frac{ u_i + u_o }{2} \]

\[ [u_n]_b = \frac{ u_o - u_i }{h} \]

where the subscript b stands for the the boundary, i stands for the first cell inside the boundary and o stands for the first cell outside the boundary and h is the grid spacing normal to the boundary. Applying this to the Robin formula gives

\[ u_o = \frac{ h\gamma + u_i( \beta - \frac{h}{2} \alpha ) } { \beta + \frac{h}{2} \alpha } \]

or equivalently

\[ u_o = \frac{ hg + u_i (1-a(1+\frac{h}{2})) }{ 1-a(1-\frac{h}{2}) } \]

After setting the edge (face in 3D) boundary conditions, linear approximations are used to set the boundary conditions of higher boundary types (nodes in 2D, edges and nodes in 3D).

In some cases, the calling function wants to set the boundary condition homogeneously, with g=0. This is useful in problems where the the solution of the homogeneous problem is required in solving the inhomogeneous problem. This function respects such requests specified through the argument homogeneous_bc.

◆ setBoundaryValuesInCells() [2/2]

template<int DIM>
void SAMRAI::solv::CartesianRobinBcHelper< DIM >::setBoundaryValuesInCells ( hier::PatchLevel< DIM > &  level,
const double  fill_time,
const hier::IntVector< DIM > &  ghost_width_to_fill,
int  target_data_id,
bool  homogeneous_bc = false 
) const

Loop through all patches on the given level and call setBoundaryValuesInCells(hier::Patch<DIM> &patch, const double fill_time , const hier::IntVector<DIM> &ghost_width_to_fill , int target_data_id , bool homogeneous_bc=false ) for each.

Parameters
levelPatchLevel on which to set boundary condition
fill_timeSolution time corresponding to filling
ghost_width_to_fillMax ghost width requiring fill
target_data_idhier::Patch data index of data to be set. This data must be a cell-centered double.
homogeneous_bcSet a homogeneous boundary condition. This means g=0 for the boundary.

◆ setBoundaryValuesAtNodes()

template<int DIM>
void SAMRAI::solv::CartesianRobinBcHelper< DIM >::setBoundaryValuesAtNodes ( hier::Patch< DIM > &  patch,
const double  fill_time,
int  target_data_id,
bool  homogeneous_bc = false 
) const

This function is not yet implemented!

There are some decisions that must be made before the implementation can be written.

  1. Do we set the values on the boundary or one cell away from the boundary?
  2. What is the discrete formulation we should use to compute the value to be set?

This function has an interface similar to the virtual function xfer::RefinePatchStrategy<DIM>::setPhysicalBoundaryConditions(), and it may be used to help implement that function, but it does not serve the same purpose. The primary differences are:

  1. It specializes to node-centered variables.
  2. User must specify the index of the data whose ghost cells need to be filled. This index is used to determine the variable for which to set the boundary coefficients and to get the data to be set.

This function calls RobinBcStrategy::setBcCoefs() to get the coefficients, then it sets the values at the boundary nodes.

In some cases, the calling function wants to set the boundary condition homogeneously, with g=0. This is useful in problems where the the solution of the homogeneous problem is required to solving the inhomogeneous problem. This function respects such requests specified through the argument homogeneous_bc.

Parameters
patchhier::Patch on which to set boundary condition
fill_timeSolution time corresponding to filling
target_data_idhier::Patch data index of data to be set.
homogeneous_bcSet a homogeneous boundary condition. This means g=0 for the boundary.

◆ setCoefImplementation()

template<int DIM>
void SAMRAI::solv::CartesianRobinBcHelper< DIM >::setCoefImplementation ( const RobinBcCoefStrategy< DIM > *  coef_strategy)

Provide the implementation that can be used to set the Robin bc coefficients.

Parameters
coef_strategytbox::Pointer to a concrete inmplementation of the coefficient strategy.

◆ setTargetDataId()

template<int DIM>
void SAMRAI::solv::CartesianRobinBcHelper< DIM >::setTargetDataId ( int  target_data_id)

When setPhysicalBoundaryConditions is called, the data specified will be set. This information is required because the it is not passed in through the argument list of setPhysicalBounaryConditions.

◆ setHomogeneousBc()

template<int DIM>
void SAMRAI::solv::CartesianRobinBcHelper< DIM >::setHomogeneousBc ( bool  homogeneous_bc)

In certain circumstances, only the value of a is needed, while the value of g is temporarily not required and taken to be zero. (An example is in setting the boundary condition for error value in an iterative method.) In such cases, use this function to set a flag that will cause a null pointer to be given to setBcCoefs() to indicate that fact.

◆ trimBoundaryBox()

template<int DIM>
hier::BoundaryBox<DIM> SAMRAI::solv::CartesianRobinBcHelper< DIM >::trimBoundaryBox ( const hier::BoundaryBox< DIM > &  boundary_box,
const hier::Box< DIM > &  limit_box 
) const
private

Certain boundary-related operations occur on patch data that do not or cannot extend past the edgr or corner of a patch. This function is used to trim down the parts of the boundary box that extend past those points so that a suitable index range is achieved.

The boundary box trimmed must be of type 1 or 2.

Parameters
boundary_boxBoundary box to be trimmed.
limit_boxhier::Box to not stick past
Returns
New trimmed boundary box.

◆ makeNodeBoundaryBox()

template<int DIM>
hier::Box<DIM> SAMRAI::solv::CartesianRobinBcHelper< DIM >::makeNodeBoundaryBox ( const hier::BoundaryBox< DIM > &  boundary_box) const
private

Define a box describing the indices of the nodes corresponding to the input boundary box. These nodes lie on the boundary itself.

The input boundary_box must be of type 1 (see hier::BoundaryBox::getBoundaryType()).

Parameters
boundary_boxinput boundary box
Returns
a box to define the node indices corresponding to boundary_box

◆ makeFaceBoundaryBox()

template<int DIM>
hier::Box<DIM> SAMRAI::solv::CartesianRobinBcHelper< DIM >::makeFaceBoundaryBox ( const hier::BoundaryBox< DIM > &  boundary_box) const
private

Define a box describing the indices of the codimension 1 surface corresponding to the input boundary box.

The input boundary_box must be of type 1 (see hier::BoundaryBox::getBoundaryType()).

This is a utility function for working with the indices coresponding to a boundary box but coincide with the patch boundary.

Parameters
boundary_boxinput boundary box
Returns
a box to define the face indices corresponding to boundary_box

Member Data Documentation

◆ d_object_name

template<int DIM>
std::string SAMRAI::solv::CartesianRobinBcHelper< DIM >::d_object_name
private

◆ d_coef_strategy

template<int DIM>
const RobinBcCoefStrategy<DIM>* SAMRAI::solv::CartesianRobinBcHelper< DIM >::d_coef_strategy
private

◆ d_target_data_id

template<int DIM>
int SAMRAI::solv::CartesianRobinBcHelper< DIM >::d_target_data_id
private

◆ d_homogeneous_bc

template<int DIM>
bool SAMRAI::solv::CartesianRobinBcHelper< DIM >::d_homogeneous_bc
private

◆ t_set_boundary_values_in_cells

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::solv::CartesianRobinBcHelper< DIM >::t_set_boundary_values_in_cells
private

◆ t_use_set_bc_coefs

template<int DIM>
tbox::Pointer<tbox::Timer> SAMRAI::solv::CartesianRobinBcHelper< DIM >::t_use_set_bc_coefs
private

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