IBAMR
An adaptive and distributed-memory parallel implementation of the immersed boundary (IB) method
Classes | Public Types | Public Member Functions | List of all members
IBAMR::BrinkmanAdvDiffBcHelper Class Reference

BrinkmanAdvDiffBcHelper is an abstract class that provides an interface to implement Brinkman penalization body force in the advection-diffusion equation in order to enforce Dirichlet, Neumann and Robin boundary conditions on surfaces of rigid immersed bodies. A single instance of this class is meant to handle all of the Brinkman penalization zones for multiple transported quantities with various boundary conditions. More...

#include </home/runner/work/IBAMR/IBAMR/include/ibamr/BrinkmanAdvDiffBcHelper.h>

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

Public Types

using BrinkmanInhomogeneousBCsFcnPtr = void(*)(int B_idx, SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > ls_var, SAMRAI::tbox::Pointer< IBTK::HierarchyMathOps > hier_math_ops, double time, void *ctx)
 Function specifying the optional forcing function for inhomogeneous boundary conditions $ \zeta q + \kappa n \dot \nabla q = \g $. More...
 

Public Member Functions

 BrinkmanAdvDiffBcHelper (std::string object_name, SAMRAI::tbox::Pointer< IBAMR::AdvDiffHierarchyIntegrator > adv_diff_solver)
 Constructor of the class.
 
 ~BrinkmanAdvDiffBcHelper ()=default
 Destructor of the class.
 
void setTimeInterval (double current_time, double new_time)
 Set the time interval in which Brinkman forcing is computed.
 
void preprocessBrinkmanAdvDiffBcHelper (double current_time, double new_time, int num_cycles)
 Preprocess routine before computing Brinkman penalization terms.
 
void postprocessBrinkmanAdvDiffBcHelper (double current_time, double new_time, int num_cycles)
 Postprocess routine after computing Brinkman penalization terms.
 
void setPenaltyCoefficient (double eta_penalty_coef)
 Set Brinkman penalization penalty factor for all level sets.
 
void setNumInterfaceCells (double num_interface_cells)
 Set the number of interface cells for all level sets.
 
const std::stringgetName () const
 Get the name of the object.
 
std::pair< double, double > getCurrentTimeInterval () const
 Get the current time interval $ [t^{n+1}, t^n] $ in which Brinkman velocity is computed.
 
bool hasBrinkmanBoundaryCondition (SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > Q_var) const
 Function to determine if a transported quantity has Brinkman boundary conditions associated with it.
 
void registerHomogeneousBC (SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > Q_var, SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > ls_solid_var, std::string bc_type, std::string indicator_func_type="SMOOTH", double num_interface_cells=2.0, double eta_penalty_coef=1.0e-8)
 Register a transported quantity with this object, along with the solid level set variable for which to apply a homogeneous boundary condition. More...
 
void registerInhomogeneousBC (SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > Q_var, SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > ls_solid_var, std::string bc_type, BrinkmanInhomogeneousBCsFcnPtr callback=nullptr, void *ctx=nullptr, std::string indicator_func_type="SMOOTH", double num_interface_cells=2.0, double eta_penalty_coef=1.0e-8, double bc_val=std::numeric_limits< double >::signaling_NaN())
 Register a transported quantity with this object, along with the solid level set variable on which to apply inhomogeneous boundary conditions. More...
 
void computeDampingCoefficient (int C_idx, SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > Q_var)
 Function to compute the cell-centered coefficient to the damping linear operator and RHS of the advection-diffusion equation for a specified transported quantity Q_var. More...
 
void computeDiffusionCoefficient (int D_idx, SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > Q_var, int kappa_idx, double kappa)
 Function to compute the side-centered coefficient to the diffusion linear operator and RHS of the advection-diffusion equation for a specified transported quantity Q_var with diffusion coefficient kappa. More...
 
void computeForcing (int F_idx, SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > Q_var)
 Function to compute the Brinkman forcing contribution to the RHS of the advection-diffusion solver for a specified transported quantity Q_var. More...
 
void maskForcingTerm (int N_idx, SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > Q_var, const bool mask_smeared_region=false)
 Function to mask the additional forcing terms on the RHS of the advection-diffusion solver e.g. $ u \dot \grad Q$ and body forces. More...
 

Detailed Description

BrinkmanAdvDiffBcHelper is an abstract class that provides an interface to implement Brinkman penalization body force in the advection-diffusion equation in order to enforce Dirichlet, Neumann and Robin boundary conditions on surfaces of rigid immersed bodies. A single instance of this class is meant to handle all of the Brinkman penalization zones for multiple transported quantities with various boundary conditions.

BrinkmanAdvDiffBcHelper provides an implementation of a volume penalized body force and linear operator modifications required to impose Dirichlet, Neumann and Robin boundary conditions to scalar quantities maintained by BrinkmanSemiImplicitAdvDiffHierarchyIntegrator.

Boundary conditions can be applied to multiple interfaces, which are demarcated using level set variables. This class assumes that the penalized region coincides with negative values of the level set. The sign convention of the level set variable is specified by the user.

Reference Sakurai, T., Yoshimatsu, K., Okamoto N. and Schneider K.,Volume penalization for inhomogeneous Neumann boundary conditions modeling scalar flux in complicated geometry

Member Typedef Documentation

◆ BrinkmanInhomogeneousBCsFcnPtr

Function specifying the optional forcing function for inhomogeneous boundary conditions $ \zeta q + \kappa n \dot \nabla q = \g $.

The user must set the patch data B_idx such that $ n \dot B = g $. The applied forcing term is then computed internally as $ \nabla \dot (\chi B) - \chi \nabla \dot B $. Note that B_idx contains side-centered patch data.

Member Function Documentation

◆ computeDampingCoefficient()

void IBAMR::BrinkmanAdvDiffBcHelper::computeDampingCoefficient ( int  C_idx,
SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  Q_var 
)

Function to compute the cell-centered coefficient to the damping linear operator and RHS of the advection-diffusion equation for a specified transported quantity Q_var.

Note
It is assumed that the physical damping coefficient $\lambda$ is zero.

The functional form of the Brinkman damping coefficient is $ C = \sum_{Dirichlet} \chi_i/\eta + \sum_{Robin} \nabla \dot (\chi_i n_i) - \chi_i \nabla \dot n_i$ where the sum is taken over all level sets with Dirichlet and Robin BCs. Here $\chi_i = 1-H_i$. Note that Neumann BCs do not contribute anything to this term.

◆ computeDiffusionCoefficient()

void IBAMR::BrinkmanAdvDiffBcHelper::computeDiffusionCoefficient ( int  D_idx,
SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  Q_var,
int  kappa_idx,
double  kappa 
)

Function to compute the side-centered coefficient to the diffusion linear operator and RHS of the advection-diffusion equation for a specified transported quantity Q_var with diffusion coefficient kappa.

Note
This function is able to handle both constant and variable kappa.

The functional form of the Brinkman diffusion coefficient is $ D = \kappa + \sum_{Neumann} (-\chi_i + \eta \chi_i) + \sum_{Robin} (-\chi_i + \eta \chi_i)$ where $\chi_i = 1-H_i$ and the sum is taken over all level sets with Neumann and Robin BCs . Note that Dirichlet BCs do not contribute anything to this term.

◆ computeForcing()

void IBAMR::BrinkmanAdvDiffBcHelper::computeForcing ( int  F_idx,
SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  Q_var 
)

Function to compute the Brinkman forcing contribution to the RHS of the advection-diffusion solver for a specified transported quantity Q_var.

For Inhomogeneous Dirichlet BCs, $ F_i = \chi_i/\eta Q_{bc}$ where $\chi_i = 1-H_i$. For inhomogeneous Neumann and Robin BCs, $ F_i = \nabla \dot (\chi_i B) - \chi_i \nabla \dot B_i $, with a user defined $B_i$. The overall functional form of the Brinkman body force is $F = \sum_{i} F_i$.

◆ maskForcingTerm()

void IBAMR::BrinkmanAdvDiffBcHelper::maskForcingTerm ( int  N_idx,
SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  Q_var,
const bool  mask_smeared_region = false 
)

Function to mask the additional forcing terms on the RHS of the advection-diffusion solver e.g. $ u \dot \grad Q$ and body forces.

The functional form of the Brinkman masking term is $ N = (1-\sum_{Neumann and Robin} \chi_i) N$ where $\chi_i = 1-H_i$ and the sum is taken over all level sets with Neumann and Robin BCs. Note that Dirichlet BCs do not mask this term presently.

◆ registerHomogeneousBC()

void IBAMR::BrinkmanAdvDiffBcHelper::registerHomogeneousBC ( SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  Q_var,
SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  ls_solid_var,
std::string  bc_type,
std::string  indicator_func_type = "SMOOTH",
double  num_interface_cells = 2.0,
double  eta_penalty_coef = 1.0e-8 
)

Register a transported quantity with this object, along with the solid level set variable for which to apply a homogeneous boundary condition.

Note
This function can only be used to register homogeneous Dirichlet, Neumann and Robin BCs.

◆ registerInhomogeneousBC()

void IBAMR::BrinkmanAdvDiffBcHelper::registerInhomogeneousBC ( SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  Q_var,
SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > >  ls_solid_var,
std::string  bc_type,
BrinkmanInhomogeneousBCsFcnPtr  callback = nullptr,
void *  ctx = nullptr,
std::string  indicator_func_type = "SMOOTH",
double  num_interface_cells = 2.0,
double  eta_penalty_coef = 1.0e-8,
double  bc_val = std::numeric_limits<double>::signaling_NaN() 
)

Register a transported quantity with this object, along with the solid level set variable on which to apply inhomogeneous boundary conditions.

Note
Inhomogeneous BCs are treated uniquely within this class and require additional user callback inputs.

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