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

Class CIBSaddlePointSolver solves for the fluid velocity $ \vec{v} $, fluid pressure $ p $, and the Lagrange multiplier $ \vec{\lambda} $ maintaining the rigidity constraint. More...

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

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

Public Member Functions

 CIBSaddlePointSolver (std::string object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, SAMRAI::tbox::Pointer< IBAMR::INSStaggeredHierarchyIntegrator > navier_stokes_integrator, SAMRAI::tbox::Pointer< IBAMR::CIBStrategy > cib_strategy, std::string default_options_prefix, MPI_Comm petsc_comm=PETSC_COMM_WORLD)
 Constructor for a saddle-point solver that employs the PETSc KSP solver framework.
 
virtual ~CIBSaddlePointSolver ()
 Destructor.
 
void setKSPType (const std::string &ksp_type)
 Set the KSP type.
 
void setOptionsPrefix (const std::string &options_prefix)
 Set the options prefix used by this PETSc solver object.
 
void setVelocityPoissonSpecifications (const SAMRAI::solv::PoissonSpecifications &u_problem_coefs)
 Set the PoissonSpecifications object used to specify the coefficients for the momentum equation in the incompressible Stokes operator.
 
void setPhysicalBcCoefs (const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > &u_bc_coefs, SAMRAI::solv::RobinBcCoefStrategy< NDIM > *p_bc_coef)
 Set the SAMRAI::solv::RobinBcCoefStrategy objects used to specify physical boundary conditions. More...
 
void setPhysicalBoundaryHelper (SAMRAI::tbox::Pointer< IBAMR::StaggeredStokesPhysicalBoundaryHelper > bc_helper)
 Set the StokesSpecifications object and timestep size used to specify the coefficients for the time-dependent incompressible Stokes operator.
 
SAMRAI::tbox::Pointer< IBTK::LinearOperatorgetA () const
 Return the linear operator for the saddle-point solver.
 
SAMRAI::tbox::Pointer< IBAMR::StaggeredStokesSolvergetStokesSolver () const
 Return the Stokes solver used in the preconditioner of the solver.
 
double getInterpScale () const
 Return the interpolation scale factor.
 
double getSpreadScale () const
 Return the spreading scale factor.
 
double getRegularizationWeight () const
 Return the regularization scale factor.
 
bool getNormalizeSpreadForce () const
 Return if the spread force is to be normalized.
 
bool solveSystem (Vec x, Vec b)
 Solve the linear system of equations $Ax=b$ for $x$. More...
 
void initializeSolverState (Vec x, Vec b)
 Compute hierarchy dependent data required for solving $Ax=b$.
 
void deallocateSolverState ()
 Remove all hierarchy dependent data allocated by initializeSolverState(). More...
 
void setSolutionTime (double solution_time)
 Set the time at which the solution is to be evaluated.
 
void setTimeInterval (double current_time, double new_time)
 Set the current time interval.
 
SAMRAI::tbox::Pointer< IBAMR::CIBMobilitySolvergetCIBMobilitySolver () const
 Get the mobility solver.
 

Detailed Description

Class CIBSaddlePointSolver solves for the fluid velocity $ \vec{v} $, fluid pressure $ p $, and the Lagrange multiplier $ \vec{\lambda} $ maintaining the rigidity constraint.

For free-moving (self-moving bodies) it also solves for the rigid body translational and rotational velocities $ \vec{U} $.

\begin{eqnarray*} A \vec{v} + \nabla p - \gamma S \vec{\lambda} &=& \vec{g} \\ -\nabla \cdot \vec{v} &=& (h = 0) \\ -\beta J \vec{v} + \beta T^{*} \vec{U} -\beta \delta \vec{\lambda} &=& (\vec{w} = -\beta \vec{U}_{\mbox{\scriptsize def}}) \\ T \vec{\lambda} &=& \vec{F} \end{eqnarray*}

\begin{eqnarray*} A \vec{v} + \nabla p - \gamma S \vec{\lambda} &=& \vec{g} \\ -\nabla \cdot \vec{v} &=& (h = 0) \\ -\beta J \vec{v} - \beta \delta \vec{\lambda} &=& (w = -\beta (T^{*}\vec{U}+ \vec{U}_{\mbox{\scriptsize def}}) = -\beta \vec{U}_{\mbox{\scriptsize imposed}}) \end{eqnarray*}

Here, $ A $ is the Helmholtz operator, $ S $ is the spreading operator, $ G $ is the gradient operator, $ J $ is the interpolation operator, and $ T $ is the rigid body operator (with $ T^{*} $ denoting its adjoint), $ \vec{F} $ is the net external force/torque on the body (excluding the hydrodynamic force/torque), $ \vec{U}_{\mbox{\scriptsize def}} $ is the deformational velocity of the body, and $ \vec{U}_{\mbox{\scriptsize imposed}} $ is the velocity imposed on the body. $ \gamma $ and $ \beta $ are two scaling parameters to make the system well-scaled. $ \delta $ is the regularization parameter that regularizes the problem in presence of many IB/control points in a grid cell.

Here, we employ the Krylov solver to solve the above saddle-point problem. We use Schur complement preconditioner to precondition the iterative solver.

Member Function Documentation

◆ deallocateSolverState()

void IBAMR::CIBSaddlePointSolver::deallocateSolverState ( )

Remove all hierarchy dependent data allocated by initializeSolverState().

Note
It is safe to call deallocateSolverState() when the solver state is already deallocated.
See also
initializeSolverState

◆ setPhysicalBcCoefs()

void IBAMR::CIBSaddlePointSolver::setPhysicalBcCoefs ( const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > &  u_bc_coefs,
SAMRAI::solv::RobinBcCoefStrategy< NDIM > *  p_bc_coef 
)

Set the SAMRAI::solv::RobinBcCoefStrategy objects used to specify physical boundary conditions.

Note
Any of the elements of u_bc_coefs may be NULL. In this case, homogeneous Dirichlet boundary conditions are employed for that data depth. p_bc_coef may also be NULL; in that case, homogeneous Neumann boundary conditions are employed for the pressure.
Parameters
u_bc_coefsIBTK::Vector of pointers to objects that can set the Robin boundary condition coefficients for the velocity.
p_bc_coefPointer to object that can set the Robin boundary condition coefficients for the pressure.

◆ solveSystem()

bool IBAMR::CIBSaddlePointSolver::solveSystem ( Vec  x,
Vec  b 
)

Solve the linear system of equations $Ax=b$ for $x$.

Returns
true if the solver converged to the specified tolerances, false otherwise

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