IBAMR  IBAMR version 0.19.
Public Member Functions | Protected Attributes | List of all members
IBTK::FischerGuess Class Reference

#include <ibtk/FischerGuess.h>

Public Member Functions

 FischerGuess ()
 
 FischerGuess (int n_vectors)
 
void submit (const libMesh::NumericVector< double > &solution, const libMesh::NumericVector< double > &rhs)
 
void guess (libMesh::NumericVector< double > &solution, const libMesh::NumericVector< double > &rhs) const
 

Protected Attributes

int d_n_max_vectors = 5
 
int d_n_stored_vectors = 0
 
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > d_correlation_matrix
 
std::vector< std::unique_ptr< libMesh::NumericVector< double > > > d_solutions
 
std::vector< std::unique_ptr< libMesh::NumericVector< double > > > d_rhs
 

Detailed Description

Class implementing a modified version of Fischer's first algorithm from the 1998 manuscript "Projection techniques for iterative solution of A x = b with successive right-hand sides".

In essence, the caller submits pairs of solution and corresponding RHS vectors from solving the same matrix equation several times (i.e., at different time steps). This class uses those pairs to compute an estimated solution vector for a new RHS. In practice, this estimate is a good enough that it cuts the number of required solver iterations by about 50%, while computing the initial guess is relatively inexpensive.

The algorithm used below varies from the one Fischer proposed in a few major ways:

  1. The RHS vectors are not orthogonalized. Instead, the least-squares equation is solved with the SVD, which enables the projection matrix to be singular (or nearly so).
  2. Since this class does not orthogonalize, it can compute the projection matrix incrementally. This means (unlike the PETSc implementation) we can always use the full set of available vectors to compute a guess.
  3. The projection is done with the standard $R^n$ inner product, which may not be optimal, but is much faster than doing sparse matrix-vector multiplications.

Since the systems we solve typically require low iteration counts, these tricks to get a less optimal guess at a lower cost are beneficial.

Constructor & Destructor Documentation

◆ FischerGuess() [1/2]

IBTK::FischerGuess::FischerGuess ( )

Constructor. Sets the number of stored vectors to five.

◆ FischerGuess() [2/2]

IBTK::FischerGuess::FischerGuess ( int  n_vectors)

Constructor.

Parameters
n_vectorsThe number of stored vectors.

Member Function Documentation

◆ submit()

void IBTK::FischerGuess::submit ( const libMesh::NumericVector< double > &  solution,
const libMesh::NumericVector< double > &  rhs 
)

Add a new solution and RHS pair to the stored collection. If the collection is full then the oldest vectors are removed.

◆ guess()

void IBTK::FischerGuess::guess ( libMesh::NumericVector< double > &  solution,
const libMesh::NumericVector< double > &  rhs 
) const

Given a RHS vector, use the stored collection of vectors to compute an estimate of its corresponding solution vector.

Member Data Documentation

◆ d_n_max_vectors

int IBTK::FischerGuess::d_n_max_vectors = 5
protected

◆ d_n_stored_vectors

int IBTK::FischerGuess::d_n_stored_vectors = 0
protected

◆ d_correlation_matrix

Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic> IBTK::FischerGuess::d_correlation_matrix
protected

◆ d_solutions

std::vector<std::unique_ptr<libMesh::NumericVector<double> > > IBTK::FischerGuess::d_solutions
protected

◆ d_rhs

std::vector<std::unique_ptr<libMesh::NumericVector<double> > > IBTK::FischerGuess::d_rhs
protected

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