IBAMR  IBAMR version 0.19.
Public Member Functions | Protected Attributes | List of all members
IBAMR::INSStaggeredVelocityBcCoef Class Referenceabstract

Class INSStaggeredVelocityBcCoef is a concrete StokesBcCoefStrategy that is used to specify velocity boundary conditions for the staggered grid incompressible Navier-Stokes solver. More...

#include <ibamr/INSStaggeredVelocityBcCoef.h>

Inheritance diagram for IBAMR::INSStaggeredVelocityBcCoef:
Inheritance graph
[legend]

Public Member Functions

 INSStaggeredVelocityBcCoef (unsigned int comp_idx, const INSStaggeredHierarchyIntegrator *fluid_solver, const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > &bc_coefs, TractionBcType traction_bc_type, bool homogeneous_bc=false)
 Constructor. More...
 
 ~INSStaggeredVelocityBcCoef ()=default
 Destructor. More...
 
void setPhysicalBcCoefs (const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > &bc_coefs)
 Set the SAMRAI::solv::RobinBcCoefStrategy objects used to specify physical boundary conditions. More...
 
void setSolutionTime (double solution_time)
 Set the time at which the solution is to be evaluated. More...
 
void setTimeInterval (double current_time, double new_time)
 Set the current time interval. More...
 
virtual void setTractionBcType (TractionBcType bc_type)
 Set the type of traction boundary conditions. Supported options are: TRACTION_BOUNDARY_CONDITIONS and PSEUDO_TRACTION_BOUNDARY_CONDITIONS. More...
 
virtual TractionBcType getTractionBcType () const
 Get the type of traction boundary conditions. More...
 
StokesBcCoefStrategy interface.
void setStokesSpecifications (const StokesSpecifications *problem_coefs) override
 Set the StokesSpecifications object used by this boundary condition specification object. More...
 
void setTargetVelocityPatchDataIndex (int u_target_data_idx) override
 Set the target velocity data index to use when setting physical boundary conditions and the time at which it is defined. More...
 
void clearTargetVelocityPatchDataIndex () override
 Clear the target velocity data index used when setting physical boundary conditions. More...
 
void setTargetPressurePatchDataIndex (int p_target_data_idx) override
 Set the target pressure data index to use when setting physical boundary conditions and the time at which it is defined. More...
 
void clearTargetPressurePatchDataIndex () override
 Clear the target pressure data index used when setting physical boundary conditions. More...
 
Extended SAMRAI::solv::RobinBcCoefStrategy interface.
void setTargetPatchDataIndex (int target_idx) override
 Set the target data index. More...
 
void clearTargetPatchDataIndex () override
 Clear the target data index. More...
 
void setHomogeneousBc (bool homogeneous_bc) override
 Set whether the class is filling homogeneous or inhomogeneous boundary conditions. More...
 

Protected Attributes

const StokesSpecificationsd_problem_coefs
 
int d_u_target_data_idx = IBTK::invalid_index
 
int d_p_target_data_idx = IBTK::invalid_index
 
TractionBcType d_traction_bc_type = TRACTION
 

Implementation of SAMRAI::solv::RobinBcCoefStrategy interface.

const unsigned int d_comp_idx
 
const INSStaggeredHierarchyIntegratord_fluid_solver
 
std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > d_bc_coefs
 
void setBcCoefs (SAMRAI::tbox::Pointer< SAMRAI::pdat::ArrayData< NDIM, double > > &acoef_data, SAMRAI::tbox::Pointer< SAMRAI::pdat::ArrayData< NDIM, double > > &bcoef_data, SAMRAI::tbox::Pointer< SAMRAI::pdat::ArrayData< NDIM, double > > &gcoef_data, const SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > &variable, const SAMRAI::hier::Patch< NDIM > &patch, const SAMRAI::hier::BoundaryBox< NDIM > &bdry_box, double fill_time=0.0) const override
 Function to fill arrays of Robin boundary condition coefficients at a patch boundary. More...
 
SAMRAI::hier::IntVector< NDIM > numberOfExtensionsFillable () const override
 
 INSStaggeredVelocityBcCoef ()=delete
 Default constructor. More...
 
 INSStaggeredVelocityBcCoef (const INSStaggeredVelocityBcCoef &from)=delete
 Copy constructor. More...
 
INSStaggeredVelocityBcCoefoperator= (const INSStaggeredVelocityBcCoef &that)=delete
 Assignment operator. More...
 

Detailed Description

This class interprets pure Dirichlet boundary conditions on the velocity as prescribed velocity boundary conditions, whereas pure Neumann boundary conditions are interpreted as prescribed traction (stress) boundary conditions. These are translated into Dirichlet and generalized Neumann boundary conditions, respectively, for the velocity.

Dirichlet, true traction, and pseudo-traction boundary conditions are all supported.

Constructor & Destructor Documentation

◆ INSStaggeredVelocityBcCoef() [1/3]

IBAMR::INSStaggeredVelocityBcCoef::INSStaggeredVelocityBcCoef ( unsigned int  comp_idx,
const INSStaggeredHierarchyIntegrator fluid_solver,
const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > &  bc_coefs,
TractionBcType  traction_bc_type,
bool  homogeneous_bc = false 
)

◆ ~INSStaggeredVelocityBcCoef()

IBAMR::INSStaggeredVelocityBcCoef::~INSStaggeredVelocityBcCoef ( )
default

◆ INSStaggeredVelocityBcCoef() [2/3]

IBAMR::INSStaggeredVelocityBcCoef::INSStaggeredVelocityBcCoef ( )
privatedelete
Note
This constructor is not implemented and should not be used.

◆ INSStaggeredVelocityBcCoef() [3/3]

IBAMR::INSStaggeredVelocityBcCoef::INSStaggeredVelocityBcCoef ( const INSStaggeredVelocityBcCoef from)
privatedelete
Note
This constructor is not implemented and should not be used.
Parameters
fromThe value to copy to this object.

Member Function Documentation

◆ setPhysicalBcCoefs()

void IBAMR::INSStaggeredVelocityBcCoef::setPhysicalBcCoefs ( const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > &  bc_coefs)
Parameters
bc_coefsIBTK::Vector of boundary condition specification objects

◆ setSolutionTime()

void IBAMR::INSStaggeredVelocityBcCoef::setSolutionTime ( double  solution_time)

◆ setTimeInterval()

void IBAMR::INSStaggeredVelocityBcCoef::setTimeInterval ( double  current_time,
double  new_time 
)

◆ setStokesSpecifications()

void IBAMR::INSStaggeredVelocityBcCoef::setStokesSpecifications ( const StokesSpecifications problem_coefs)
overridevirtual
Parameters
problem_coefsProblem coefficients

Reimplemented from IBAMR::StokesBcCoefStrategy.

◆ setTargetVelocityPatchDataIndex()

void IBAMR::INSStaggeredVelocityBcCoef::setTargetVelocityPatchDataIndex ( int  u_target_data_idx)
overridevirtual

Reimplemented from IBAMR::StokesBcCoefStrategy.

◆ clearTargetVelocityPatchDataIndex()

void IBAMR::INSStaggeredVelocityBcCoef::clearTargetVelocityPatchDataIndex ( )
overridevirtual

Reimplemented from IBAMR::StokesBcCoefStrategy.

◆ setTargetPressurePatchDataIndex()

void IBAMR::INSStaggeredVelocityBcCoef::setTargetPressurePatchDataIndex ( int  p_target_data_idx)
overridevirtual

Reimplemented from IBAMR::StokesBcCoefStrategy.

◆ clearTargetPressurePatchDataIndex()

void IBAMR::INSStaggeredVelocityBcCoef::clearTargetPressurePatchDataIndex ( )
overridevirtual

Reimplemented from IBAMR::StokesBcCoefStrategy.

◆ setTargetPatchDataIndex()

void IBAMR::INSStaggeredVelocityBcCoef::setTargetPatchDataIndex ( int  target_idx)
overridevirtual

Reimplemented from IBTK::ExtendedRobinBcCoefStrategy.

◆ clearTargetPatchDataIndex()

void IBAMR::INSStaggeredVelocityBcCoef::clearTargetPatchDataIndex ( )
overridevirtual

Reimplemented from IBTK::ExtendedRobinBcCoefStrategy.

◆ setHomogeneousBc()

void IBAMR::INSStaggeredVelocityBcCoef::setHomogeneousBc ( bool  homogeneous_bc)
overridevirtual

Reimplemented from IBTK::ExtendedRobinBcCoefStrategy.

◆ setBcCoefs() [1/2]

void IBAMR::INSStaggeredVelocityBcCoef::setBcCoefs ( SAMRAI::tbox::Pointer< SAMRAI::pdat::ArrayData< NDIM, double > > &  acoef_data,
SAMRAI::tbox::Pointer< SAMRAI::pdat::ArrayData< NDIM, double > > &  bcoef_data,
SAMRAI::tbox::Pointer< SAMRAI::pdat::ArrayData< NDIM, double > > &  gcoef_data,
const SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > &  variable,
const SAMRAI::hier::Patch< NDIM > &  patch,
const SAMRAI::hier::BoundaryBox< NDIM > &  bdry_box,
double  fill_time = 0.0 
) const
override
Note
In the original SAMRAI::solv::RobinBcCoefStrategy interface, it was assumed that \( b = (1-a) \). In the new interface, \(a\) and \(b\) are independent.
See also
SAMRAI::solv::RobinBcCoefStrategy::setBcCoefs()
Parameters
acoef_dataBoundary coefficient data. The array will have been defined to include index range for corresponding to the boundary box bdry_box and appropriate for the alignment of the given variable. If this is a null pointer, then the calling function is not interested in a, and you can disregard it.
bcoef_dataBoundary coefficient data. This array is exactly like acoef_data, except that it is to be filled with the b coefficient.
gcoef_dataBoundary coefficient data. This array is exactly like acoef_data, except that it is to be filled with the g coefficient.
variableVariable to set the coefficients for. If implemented for multiple variables, this parameter can be used to determine which variable's coefficients are being sought.
patchPatch requiring bc coefficients.
bdry_boxBoundary box showing where on the boundary the coefficient data is needed.
fill_timeSolution time corresponding to filling, for use when coefficients are time-dependent.

◆ numberOfExtensionsFillable() [1/2]

SAMRAI::hier::IntVector<NDIM> IBAMR::INSStaggeredVelocityBcCoef::numberOfExtensionsFillable ( ) const
override

◆ operator=()

INSStaggeredVelocityBcCoef& IBAMR::INSStaggeredVelocityBcCoef::operator= ( const INSStaggeredVelocityBcCoef that)
privatedelete
Note
This operator is not implemented and should not be used.
Parameters
thatThe value to assign to this object.
Returns
A reference to this object.

◆ setTractionBcType()

virtual void IBAMR::StokesBcCoefStrategy::setTractionBcType ( TractionBcType  bc_type)
virtualinherited

The default is TRACTION_BOUNDARY_CONDITIONS.

◆ getTractionBcType()

virtual TractionBcType IBAMR::StokesBcCoefStrategy::getTractionBcType ( ) const
virtualinherited

◆ setBcCoefs() [2/2]

virtual void SAMRAI::solv::RobinBcCoefStrategy< DIM >::setBcCoefs ( tbox::Pointer< pdat::ArrayData< DIM, double > > &  acoef_data,
tbox::Pointer< pdat::ArrayData< DIM, double > > &  bcoef_data,
tbox::Pointer< pdat::ArrayData< DIM, double > > &  gcoef_data,
const tbox::Pointer< hier::Variable< DIM > > &  variable,
const hier::Patch< DIM > &  patch,
const hier::BoundaryBox< DIM > &  bdry_box,
const double  fill_time = 0.0 
) const
pure virtualinherited

This class specifies the Robin boundary condition coefficients at discrete locations on the patch boundary. Though these locations are defined by boundary box object, they do not necessarily coincide with the centers of the cells referred to by those boxes. These locations typically coincide with the nodes or face centers which do lie on the patch boundary. Accordingly, you use this function to provide the boundary coefficients at those locations by filling an array at indices corresponding to those locations.

When setting the values of the boundary condition coefficients it is useful to note that for any cell (i,j,k), the indices of the sides, edges and nodes are easily determined. The index on the lower side of the cell is the same as the index of the cell, whereas the index on the upper side of the cell has the next higher value. In 2D, the cell and its surrounding nodes and faces has the following indices:

*
*       (i,j+1)----(i,j+1)---(i+1,j+1)
*          |                     |
*          |                     |
*          |                     |
*          |                     |
*        (i,j)      (i,j)     (i+1,j)
*          |                     |
*          |                     |
*          |                     |
*          |                     |
*        (i,j)------(i,j)-----(i+1,j)
*
* 

Once this is understood, translation between the index in the boundary box index space to the index of things on the boundary is simple.

The boundary condition coefficients should be placed in the pdat::ArrayData<DIM> objects, acoef_data and gcoef_data (see argument list), which are dimensioned to contain the indices of the points alligned with variable and lying on the the boundary defined by bdry_box.

This function is only used with type-1 boundary boxes, such as faces in 3D. Other types of boundaries do not have a well-defined surface normal.

The parameter variable is passed through to tell the implementation of this function what variable to set the coefficients for. You may wish to ignore it if your implementation is intended for a specific variable.

Parameters
acoef_databoundary coefficient data. The array will have been defined to include index range for corresponding to the boundary box bdry_box and appropriate for the alignment of the given variable. If this is a null pointer, then the calling function is not interested in a, and you can disregard it.
bcoef_databoundary coefficient data. This array is exactly like acoef_data, except that it is to be filled with the b coefficient.
gcoef_databoundary coefficient data. This array is exactly like acoef_data, except that it is to be filled with the g coefficient.
variablevariable to set the coefficients for. If implemented for multiple variables, this parameter can be used to determine which variable's coefficients are being sought.
patchpatch requiring bc coefficients
bdry_boxboundary box showing where on the boundary the coefficient data is needed.
fill_timesolution time corresponding to filling, for use when coefficients are time-dependent.

◆ numberOfExtensionsFillable() [2/2]

virtual hier::IntVector<DIM> SAMRAI::solv::RobinBcCoefStrategy< DIM >::numberOfExtensionsFillable ( ) const
pure virtualinherited

Member Data Documentation

◆ d_comp_idx

const unsigned int IBAMR::INSStaggeredVelocityBcCoef::d_comp_idx
private

◆ d_fluid_solver

const INSStaggeredHierarchyIntegrator* IBAMR::INSStaggeredVelocityBcCoef::d_fluid_solver
private

◆ d_bc_coefs

std::vector<SAMRAI::solv::RobinBcCoefStrategy<NDIM>*> IBAMR::INSStaggeredVelocityBcCoef::d_bc_coefs
private

◆ d_problem_coefs

const StokesSpecifications* IBAMR::StokesBcCoefStrategy::d_problem_coefs
protectedinherited

◆ d_u_target_data_idx

int IBAMR::StokesBcCoefStrategy::d_u_target_data_idx = IBTK::invalid_index
protectedinherited

Patch data indices.

◆ d_p_target_data_idx

int IBAMR::StokesBcCoefStrategy::d_p_target_data_idx = IBTK::invalid_index
protectedinherited

◆ d_traction_bc_type

TractionBcType IBAMR::StokesBcCoefStrategy::d_traction_bc_type = TRACTION
protectedinherited

◆ d_target_data_idx

int IBTK::ExtendedRobinBcCoefStrategy::d_target_data_idx = IBTK::invalid_index
protectedinherited

◆ d_homogeneous_bc

bool IBTK::ExtendedRobinBcCoefStrategy::d_homogeneous_bc = false
protectedinherited

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