|
IBAMR
IBAMR version 0.19.
|
#include <ibamr/KrylovMobilitySolver.h>

Public Member Functions | |
| KrylovMobilitySolver (std::string object_name, SAMRAI::tbox::Pointer< IBAMR::INSStaggeredHierarchyIntegrator > navier_stokes_integrator, SAMRAI::tbox::Pointer< IBAMR::CIBStrategy > cib_strategy, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, std::string default_options_prefix, MPI_Comm petsc_comm=PETSC_COMM_WORLD) | |
| Constructor for mobility solver that employs the PETSc KSP solver framework. More... | |
| virtual | ~KrylovMobilitySolver () |
| 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... | |
| const KSP & | getPETScKSP () const |
| Get the PETSc KSP object. More... | |
| SAMRAI::tbox::Pointer< IBAMR::StaggeredStokesSolver > | getStokesSolver () const |
| Return the Stokes solver used in the preconditioner of the solver. 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 | 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... | |
| 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... | |
| bool | solveSystem (Vec x, Vec b) |
| Solve the linear system of equations \( Mx=b \) for \( x \). More... | |
| void | initializeSolverState (Vec x, Vec b) |
| Compute hierarchy dependent data required for solving \( Mx = 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... | |
| void | setInterpScale (const double scale_interp) |
| Set scale factor for interp operator. More... | |
| void | setSpreadScale (const double scale_spread) |
| Set scale factor for spread operator. More... | |
| void | setRegularizeMobilityScale (const double scale_reg_mob) |
| Set scale factor for regularizing mobility matrix. More... | |
| void | setNormalizeSpreadForce (const bool normalize_force) |
| Set if the mean of the Lagrangian force is to be subtracted from the Eulerian force variable. More... | |
Private Member Functions | |
| KrylovMobilitySolver (const KrylovMobilitySolver &from)=delete | |
| Copy constructor. More... | |
| KrylovMobilitySolver & | operator= (const KrylovMobilitySolver &that)=delete |
| Assignment operator. More... | |
| void | getFromInput (SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db) |
| Get solver settings from the input file. 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 mobility matrix. 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... | |
| void | destroyKSP () |
| Routine to destroy KSP object. More... | |
We are trying to solve the problem
\( Mx = [J L^{-1} S]x = b \); for \( x \).
Here, \( M \) is the mobility matrix, \( J \) is the interpolation operator, \( L \) is the Stokes operator, and \( S \) is the spreading operator.
| IBAMR::KrylovMobilitySolver::KrylovMobilitySolver | ( | std::string | object_name, |
| SAMRAI::tbox::Pointer< IBAMR::INSStaggeredHierarchyIntegrator > | navier_stokes_integrator, | ||
| SAMRAI::tbox::Pointer< IBAMR::CIBStrategy > | cib_strategy, | ||
| SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > | input_db, | ||
| std::string | default_options_prefix, | ||
| MPI_Comm | petsc_comm = PETSC_COMM_WORLD |
||
| ) |
|
virtual |
|
privatedelete |
| from | The value to copy to this object. |
| void IBAMR::KrylovMobilitySolver::setKSPType | ( | const std::string & | ksp_type | ) |
| void IBAMR::KrylovMobilitySolver::setOptionsPrefix | ( | const std::string & | options_prefix | ) |
| const KSP& IBAMR::KrylovMobilitySolver::getPETScKSP | ( | ) | const |
| SAMRAI::tbox::Pointer<IBAMR::StaggeredStokesSolver> IBAMR::KrylovMobilitySolver::getStokesSolver | ( | ) | const |
| void IBAMR::KrylovMobilitySolver::setVelocityPoissonSpecifications | ( | const SAMRAI::solv::PoissonSpecifications & | u_problem_coefs | ) |
| void IBAMR::KrylovMobilitySolver::setPhysicalBoundaryHelper | ( | SAMRAI::tbox::Pointer< IBAMR::StaggeredStokesPhysicalBoundaryHelper > | bc_helper | ) |
| void IBAMR::KrylovMobilitySolver::setPhysicalBcCoefs | ( | const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > & | u_bc_coefs, |
| SAMRAI::solv::RobinBcCoefStrategy< NDIM > * | p_bc_coef | ||
| ) |
| u_bc_coefs | 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. |
| bool IBAMR::KrylovMobilitySolver::solveSystem | ( | Vec | x, |
| Vec | b | ||
| ) |
| x | solution vector |
| b | right-hand-side vector |
true if the solver converged to the specified tolerances, false otherwise | void IBAMR::KrylovMobilitySolver::initializeSolverState | ( | Vec | x, |
| Vec | b | ||
| ) |
| x | solution vector |
| b | right-hand-side vector |
| void IBAMR::KrylovMobilitySolver::deallocateSolverState | ( | ) |
| void IBAMR::KrylovMobilitySolver::setSolutionTime | ( | double | solution_time | ) |
| void IBAMR::KrylovMobilitySolver::setInterpScale | ( | const double | scale_interp | ) |
| void IBAMR::KrylovMobilitySolver::setSpreadScale | ( | const double | scale_spread | ) |
| void IBAMR::KrylovMobilitySolver::setRegularizeMobilityScale | ( | const double | scale_reg_mob | ) |
| void IBAMR::KrylovMobilitySolver::setNormalizeSpreadForce | ( | const bool | normalize_force | ) |
|
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 |
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 |
|
private |
|
private |
|
private |
|
private |
1.8.17