|
IBAMR
IBAMR version 0.19.
|
Class CIBSaddlePointSolver solves for the fluid velocity \( \vec{v} \), fluid pressure \( p \), and the Lagrange multiplier \( \vec{\lambda} \) maintaining the rigidity constraint. More...
#include <ibamr/CIBSaddlePointSolver.h>

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. More... | |
| virtual | ~CIBSaddlePointSolver () |
| Destructor. More... | |
| void | setKSPType (const std::string &ksp_type) |
| Set the KSP type. More... | |
| void | setOptionsPrefix (const std::string &options_prefix) |
| Set the options prefix used by this PETSc solver object. More... | |
| 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. More... | |
| 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. More... | |
| SAMRAI::tbox::Pointer< IBTK::LinearOperator > | getA () const |
| Return the linear operator for the saddle-point solver. More... | |
| SAMRAI::tbox::Pointer< IBAMR::StaggeredStokesSolver > | getStokesSolver () const |
| Return the Stokes solver used in the preconditioner of the solver. More... | |
| double | getInterpScale () const |
| Return the interpolation scale factor. More... | |
| double | getSpreadScale () const |
| Return the spreading scale factor. More... | |
| double | getRegularizationWeight () const |
| Return the regularization scale factor. More... | |
| bool | getNormalizeSpreadForce () const |
| Return if the spread force is to be normalized. More... | |
| 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\). More... | |
| 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. More... | |
| void | setTimeInterval (double current_time, double new_time) |
| Set the current time interval. More... | |
| SAMRAI::tbox::Pointer< IBAMR::CIBMobilitySolver > | getCIBMobilitySolver () const |
| Get the mobility solver. More... | |
Private Member Functions | |
| CIBSaddlePointSolver (const CIBSaddlePointSolver &from)=delete | |
| Copy constructor. More... | |
| CIBSaddlePointSolver & | operator= (const CIBSaddlePointSolver &that)=delete |
| Assignment operator. More... | |
| void | getFromInput (SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > db) |
| Get options from input file. More... | |
| void | destroyKSP () |
| Routine to destroy KSP object. More... | |
| void | initializeStokesSolver (const SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &sol_vec, const SAMRAI::solv::SAMRAIVectorReal< NDIM, double > &rhs_vec) |
| Initialize the Stokes solver needed in the preconditioning step. More... | |
| void | initializeKSP () |
| Routine to setup KSP object. More... | |
| void | resetKSPOptions () |
| Reset the values of the convergence tolerances for the PETSc KSP object. More... | |
| void | resetKSPOperators () |
| Reset the KSP operators to correspond to the supplied LinearOperator. More... | |
| void | resetKSPPC () |
| Reset the KSP PC to correspond to the supplied preconditioner. More... | |
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.
| IBAMR::CIBSaddlePointSolver::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 |
||
| ) |
|
virtual |
|
privatedelete |
| from | The value to copy to this object. |
| void IBAMR::CIBSaddlePointSolver::setKSPType | ( | const std::string & | ksp_type | ) |
| void IBAMR::CIBSaddlePointSolver::setOptionsPrefix | ( | const std::string & | options_prefix | ) |
| void IBAMR::CIBSaddlePointSolver::setVelocityPoissonSpecifications | ( | const SAMRAI::solv::PoissonSpecifications & | u_problem_coefs | ) |
| void IBAMR::CIBSaddlePointSolver::setPhysicalBcCoefs | ( | const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > & | u_bc_coefs, |
| SAMRAI::solv::RobinBcCoefStrategy< NDIM > * | p_bc_coef | ||
| ) |
| u_bc_coefs | IBTK::Vector of pointers to objects that can set the Robin boundary condition coefficients for the velocity. |
| p_bc_coef | Pointer to object that can set the Robin boundary condition coefficients for the pressure. |
| void IBAMR::CIBSaddlePointSolver::setPhysicalBoundaryHelper | ( | SAMRAI::tbox::Pointer< IBAMR::StaggeredStokesPhysicalBoundaryHelper > | bc_helper | ) |
| SAMRAI::tbox::Pointer<IBTK::LinearOperator> IBAMR::CIBSaddlePointSolver::getA | ( | ) | const |
| SAMRAI::tbox::Pointer<IBAMR::StaggeredStokesSolver> IBAMR::CIBSaddlePointSolver::getStokesSolver | ( | ) | const |
|
inline |
|
inline |
|
inline |
|
inline |
| bool IBAMR::CIBSaddlePointSolver::solveSystem | ( | Vec | x, |
| Vec | b | ||
| ) |
true if the solver converged to the specified tolerances, false otherwise | void IBAMR::CIBSaddlePointSolver::initializeSolverState | ( | Vec | x, |
| Vec | b | ||
| ) |
| void IBAMR::CIBSaddlePointSolver::deallocateSolverState | ( | ) |
| void IBAMR::CIBSaddlePointSolver::setSolutionTime | ( | double | solution_time | ) |
| SAMRAI::tbox::Pointer<IBAMR::CIBMobilitySolver> IBAMR::CIBSaddlePointSolver::getCIBMobilitySolver | ( | ) | const |
|
privatedelete |
| that | The value to assign to this object. |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
staticprivate |
|
staticprivate |
|
staticprivate |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
This boolean value determines whether the pressure is normalized to have zero mean (i.e., discrete integral) at the end of each timestep.
|
private |
This boolean value determines whether the velocity is normalized to have zero mean (i.e., discrete integral) at the end of each timestep.
This parameter only affects the case in which rho=0 (i.e. the steady Stokes equations).
|
private |
|
private |
1.8.17