IBAMR  IBAMR version 0.19.
Public Member Functions | List of all members
IBAMR::KrylovFreeBodyMobilitySolver Class Reference

#include <ibamr/KrylovFreeBodyMobilitySolver.h>

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

Public Member Functions

 KrylovFreeBodyMobilitySolver (std::string object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, std::string default_options_prefix, SAMRAI::tbox::Pointer< IBAMR::CIBStrategy > cib_strategy, MPI_Comm petsc_comm=PETSC_COMM_WORLD)
 Constructor for \( [T M^{-1} T^{*}]^{-1} \) solver that employs the PETSc KSP solver framework. More...
 
virtual ~KrylovFreeBodyMobilitySolver ()
 Destructor. More...
 
void setMobilitySolver (SAMRAI::tbox::Pointer< IBAMR::CIBMobilitySolver > mobility_solver)
 Set the mobility solver for this class. More...
 
void setInterpScale (const double interp_scale)
 Set scale for interp operator. More...
 
void setSpreadScale (const double spread_scale)
 Set scale for spread operator. More...
 
void setStokesSpecifications (const IBAMR::StokesSpecifications &stokes_spec)
 Set stokes specifications. 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...
 

Functions to access the underlying PETSc objects.

const KSP & getPETScKSP () const
 Get the PETSc KSP object. More...
 
bool solveSystem (Vec x, Vec b)
 Solve the linear system of equations \( Nx = b \) for \( x \). More...
 
void initializeSolverState (Vec x, Vec b)
 Compute hierarchy dependent data required for solving \( Nx = 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...
 
 KrylovFreeBodyMobilitySolver (const KrylovFreeBodyMobilitySolver &from)=delete
 Copy constructor. More...
 
KrylovFreeBodyMobilitySolveroperator= (const KrylovFreeBodyMobilitySolver &that)=delete
 Assignment operator. 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...
 

Static functions for use by PETSc KSP and MatShell objects.

std::string d_object_name
 
std::string d_ksp_type = KSPGMRES
 
std::string d_pc_type = "shell"
 
bool d_is_initialized = false
 
bool d_reinitializing_solver = false
 
Vec d_petsc_b = nullptr
 
Vec d_petsc_temp_v = nullptr
 
Vec d_petsc_temp_f = nullptr
 
std::string d_options_prefix
 
MPI_Comm d_petsc_comm
 
KSP d_petsc_ksp = nullptr
 
Mat d_petsc_mat = nullptr
 
int d_max_iterations = 10000
 
int d_current_iterations
 
double d_abs_residual_tol = 1.0e-50
 
double d_rel_residual_tol = 1.0e-5
 
double d_current_residual_norm
 
bool d_initial_guess_nonzero = true
 
bool d_enable_logging = false
 
SAMRAI::tbox::Pointer< IBAMR::CIBStrategyd_cib_strategy
 
SAMRAI::tbox::Pointer< IBAMR::CIBMobilitySolverd_mobility_solver
 
double d_current_time = std::numeric_limits<double>::signaling_NaN()
 
double d_new_time = std::numeric_limits<double>::signaling_NaN()
 
double d_solution_time = std::numeric_limits<double>::signaling_NaN()
 
double d_dt = std::numeric_limits<double>::signaling_NaN()
 
double d_interp_scale
 
double d_spread_scale
 
double d_mu = 1.0
 
double d_rho = 1.0
 
static PetscErrorCode MatVecMult_KFBMSolver (Mat A, Vec x, Vec y)
 Compute the matrix vector product \(y=Ax\). More...
 
static PetscErrorCode PCApply_KFBMSolver (PC pc, Vec x, Vec y)
 Apply the preconditioner to x and store the result in y. More...
 
static PetscErrorCode monitorKSP (KSP ksp, int it, PetscReal rnorm, void *mctx)
 Set KSP monitoring routine for the KSP. More...
 

Detailed Description

We are trying to solve the problem

\( Nx = [T M^{-1} T^{*}]x = b \); for \( x \),

in this class. Here, \( N \) is the body mobility matrix, \( M \) is the mobility matrix, \( T \) is the rigid body operator, and \( T^{*} \) is the rigid body velocity operator. We employ Krylov solver preconditioned by direct solver to solve the body-mobility equation.

Constructor & Destructor Documentation

◆ KrylovFreeBodyMobilitySolver() [1/2]

IBAMR::KrylovFreeBodyMobilitySolver::KrylovFreeBodyMobilitySolver ( std::string  object_name,
SAMRAI::tbox::Pointer< SAMRAI::tbox::Database input_db,
std::string  default_options_prefix,
SAMRAI::tbox::Pointer< IBAMR::CIBStrategy cib_strategy,
MPI_Comm  petsc_comm = PETSC_COMM_WORLD 
)

◆ ~KrylovFreeBodyMobilitySolver()

virtual IBAMR::KrylovFreeBodyMobilitySolver::~KrylovFreeBodyMobilitySolver ( )
virtual

◆ KrylovFreeBodyMobilitySolver() [2/2]

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

Member Function Documentation

◆ setMobilitySolver()

void IBAMR::KrylovFreeBodyMobilitySolver::setMobilitySolver ( SAMRAI::tbox::Pointer< IBAMR::CIBMobilitySolver mobility_solver)

◆ setInterpScale()

void IBAMR::KrylovFreeBodyMobilitySolver::setInterpScale ( const double  interp_scale)

◆ setSpreadScale()

void IBAMR::KrylovFreeBodyMobilitySolver::setSpreadScale ( const double  spread_scale)

◆ setStokesSpecifications()

void IBAMR::KrylovFreeBodyMobilitySolver::setStokesSpecifications ( const IBAMR::StokesSpecifications stokes_spec)

◆ setKSPType()

void IBAMR::KrylovFreeBodyMobilitySolver::setKSPType ( const std::string &  ksp_type)

◆ setOptionsPrefix()

void IBAMR::KrylovFreeBodyMobilitySolver::setOptionsPrefix ( const std::string &  options_prefix)

◆ getPETScKSP()

const KSP& IBAMR::KrylovFreeBodyMobilitySolver::getPETScKSP ( ) const

◆ solveSystem()

bool IBAMR::KrylovFreeBodyMobilitySolver::solveSystem ( Vec  x,
Vec  b 
)
Parameters
xsolution vector
bright-hand-side vector
Returns
true if the solver converged to the specified tolerances, false otherwise

◆ initializeSolverState()

void IBAMR::KrylovFreeBodyMobilitySolver::initializeSolverState ( Vec  x,
Vec  b 
)
Parameters
xsolution vector
bright-hand-side vector

◆ deallocateSolverState()

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

◆ setSolutionTime()

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

◆ setTimeInterval()

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

◆ operator=()

KrylovFreeBodyMobilitySolver& IBAMR::KrylovFreeBodyMobilitySolver::operator= ( const KrylovFreeBodyMobilitySolver 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.

◆ initializeKSP()

void IBAMR::KrylovFreeBodyMobilitySolver::initializeKSP ( )
private

◆ resetKSPOptions()

void IBAMR::KrylovFreeBodyMobilitySolver::resetKSPOptions ( )
private

◆ resetKSPOperators()

void IBAMR::KrylovFreeBodyMobilitySolver::resetKSPOperators ( )
private

◆ resetKSPPC()

void IBAMR::KrylovFreeBodyMobilitySolver::resetKSPPC ( )
private

◆ destroyKSP()

void IBAMR::KrylovFreeBodyMobilitySolver::destroyKSP ( )
private

◆ MatVecMult_KFBMSolver()

static PetscErrorCode IBAMR::KrylovFreeBodyMobilitySolver::MatVecMult_KFBMSolver ( Mat  A,
Vec  x,
Vec  y 
)
staticprivate

◆ PCApply_KFBMSolver()

static PetscErrorCode IBAMR::KrylovFreeBodyMobilitySolver::PCApply_KFBMSolver ( PC  pc,
Vec  x,
Vec  y 
)
staticprivate

◆ monitorKSP()

static PetscErrorCode IBAMR::KrylovFreeBodyMobilitySolver::monitorKSP ( KSP  ksp,
int  it,
PetscReal  rnorm,
void *  mctx 
)
staticprivate

Member Data Documentation

◆ d_object_name

std::string IBAMR::KrylovFreeBodyMobilitySolver::d_object_name
private

◆ d_ksp_type

std::string IBAMR::KrylovFreeBodyMobilitySolver::d_ksp_type = KSPGMRES
private

◆ d_pc_type

std::string IBAMR::KrylovFreeBodyMobilitySolver::d_pc_type = "shell"
private

◆ d_is_initialized

bool IBAMR::KrylovFreeBodyMobilitySolver::d_is_initialized = false
private

◆ d_reinitializing_solver

bool IBAMR::KrylovFreeBodyMobilitySolver::d_reinitializing_solver = false
private

◆ d_petsc_b

Vec IBAMR::KrylovFreeBodyMobilitySolver::d_petsc_b = nullptr
private

◆ d_petsc_temp_v

Vec IBAMR::KrylovFreeBodyMobilitySolver::d_petsc_temp_v = nullptr
private

◆ d_petsc_temp_f

Vec IBAMR::KrylovFreeBodyMobilitySolver::d_petsc_temp_f = nullptr
private

◆ d_options_prefix

std::string IBAMR::KrylovFreeBodyMobilitySolver::d_options_prefix
private

◆ d_petsc_comm

MPI_Comm IBAMR::KrylovFreeBodyMobilitySolver::d_petsc_comm
private

◆ d_petsc_ksp

KSP IBAMR::KrylovFreeBodyMobilitySolver::d_petsc_ksp = nullptr
private

◆ d_petsc_mat

Mat IBAMR::KrylovFreeBodyMobilitySolver::d_petsc_mat = nullptr
private

◆ d_max_iterations

int IBAMR::KrylovFreeBodyMobilitySolver::d_max_iterations = 10000
private

◆ d_current_iterations

int IBAMR::KrylovFreeBodyMobilitySolver::d_current_iterations
private

◆ d_abs_residual_tol

double IBAMR::KrylovFreeBodyMobilitySolver::d_abs_residual_tol = 1.0e-50
private

◆ d_rel_residual_tol

double IBAMR::KrylovFreeBodyMobilitySolver::d_rel_residual_tol = 1.0e-5
private

◆ d_current_residual_norm

double IBAMR::KrylovFreeBodyMobilitySolver::d_current_residual_norm
private

◆ d_initial_guess_nonzero

bool IBAMR::KrylovFreeBodyMobilitySolver::d_initial_guess_nonzero = true
private

◆ d_enable_logging

bool IBAMR::KrylovFreeBodyMobilitySolver::d_enable_logging = false
private

◆ d_cib_strategy

SAMRAI::tbox::Pointer<IBAMR::CIBStrategy> IBAMR::KrylovFreeBodyMobilitySolver::d_cib_strategy
private

◆ d_mobility_solver

SAMRAI::tbox::Pointer<IBAMR::CIBMobilitySolver> IBAMR::KrylovFreeBodyMobilitySolver::d_mobility_solver
private

◆ d_current_time

double IBAMR::KrylovFreeBodyMobilitySolver::d_current_time = std::numeric_limits<double>::signaling_NaN()
private

◆ d_new_time

double IBAMR::KrylovFreeBodyMobilitySolver::d_new_time = std::numeric_limits<double>::signaling_NaN()
private

◆ d_solution_time

double IBAMR::KrylovFreeBodyMobilitySolver::d_solution_time = std::numeric_limits<double>::signaling_NaN()
private

◆ d_dt

double IBAMR::KrylovFreeBodyMobilitySolver::d_dt = std::numeric_limits<double>::signaling_NaN()
private

◆ d_interp_scale

double IBAMR::KrylovFreeBodyMobilitySolver::d_interp_scale
private

◆ d_spread_scale

double IBAMR::KrylovFreeBodyMobilitySolver::d_spread_scale
private

◆ d_mu

double IBAMR::KrylovFreeBodyMobilitySolver::d_mu = 1.0
private

◆ d_rho

double IBAMR::KrylovFreeBodyMobilitySolver::d_rho = 1.0
private

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