preCICE v3.1.2
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Private Member Functions | Private Attributes | List of all members
precice::acceleration::impl::SVDFactorization Class Reference

Class that provides functionality to maintain a SVD decomposition of a matrix via successive rank-1 updates and truncation with respect to the truncation threshold eps. More...

#include <SVDFactorization.hpp>

Collaboration diagram for precice::acceleration::impl::SVDFactorization:
[legend]

Public Types

typedef Eigen::MatrixXd Matrix
 
typedef Eigen::VectorXd Vector
 

Public Member Functions

 SVDFactorization (double eps, PtrPreconditioner preconditioner)
 Constructor.
 
virtual ~SVDFactorization ()
 Destructor, empty.
 
template<typename Derived1 , typename Derived2 >
void update (const Eigen::MatrixBase< Derived1 > &A, const Eigen::MatrixBase< Derived2 > &B)
 : updates the SVD decomposition with the rank-1 update A*B^T, i.e., _psi * _sigma * _phi^T + A*B^T and overrides the internal SVD representation. After the update, the SVD is truncated according to the threshold _truncationEps
 
void initialize (PtrParMatrixOps parMatOps, int globalRows)
 : initializes the updated SVD factorization, i.e., sets the object for parallel matrix-matrix operations and the number of global rows.
 
void reset ()
 : resets the SVD factorization
 
MatrixmatrixPsi ()
 : returns a matrix representation of the orthogonal matrix Psi, A = Psi * Sigma * Phi^bs
 
VectorsingularValues ()
 : returns a matrix representation of the orthogonal matrix Sigma, A = Psi * Sigma * Phi^T
 
MatrixmatrixPhi ()
 : returns a matrix representation of the orthogonal matrix Phi, A = Psi * Sigma * Phi^T
 
int cols ()
 : returns the number of columns in the QR-decomposition
 
int rows ()
 : returns the number of rows in the QR-decomposition
 
Rank rank ()
 : returns the rank of the truncated SVD factorization
 
int getWaste ()
 : returns the total number of truncated modes since last call to this method
 
void setThreshold (double eps)
 : sets the threshold for the truncation of the SVD factorization
 
double getThreshold ()
 : returns the truncation threshold for the SVD
 
void setPrecondApplied (bool b)
 : applies the preconditioner to the factorized and truncated representation of the Jacobian matrix
 
void setApplyFilterQR (bool b, double eps=1e-3)
 : enables or disables an additional QR-2 filter for the QR-decomposition
 
bool isSVDinitialized ()
 
void setfstream (std::fstream *stream)
 Optional file-stream for logging output.
 

Private Member Functions

void computeQRdecomposition (Matrix const &A, Matrix &Q, Matrix &R)
 : computes the QR decomposition of a matrix A of type A = PSI^T*A \in R^(rank x n)
 

Private Attributes

logging::Logger _log {"acceleration::SVDFactorization"}
 
PtrPreconditioner _preconditioner
 : preconditioner for least-squares system if vectorial system is used.
 
PtrParMatrixOps _parMatrixOps = nullptr
 : object for parallel matrix operations, i.e., parallel mat-mat/ mat-vec multiplications
 
Matrix _psi
 : SVD factorization of the matrix J = _psi * _sigma * _phi^T
 
Matrix _phi
 
Vector _sigma
 
int _rows = 0
 Number of rows (on each proc, i.e., local)
 
int _cols = 0
 Number of columns, i.e., rank of the truncated svd.
 
int _globalRows = 0
 Number of global rows, i.e., sum of _rows for all procs.
 
int _waste = 0
 
double _truncationEps
 Truncation parameter for the updated SVD decomposition.
 
double _epsQR2 = 1e-3
 Threshold for the QR2 filter for the QR decomposition.
 
bool _preconditionerApplied = false
 true if the preconditioner has been applied appropriate to the updated SVD decomposition
 
bool _initialized = false
 true, if ParallelMatrixOperations object is set, i.e., initialized
 
bool _initialSVD = false
 true, if at least one update has been made, i.e., the number of rows is known and a initial rank is given.
 
bool _applyFilterQR = false
 

Detailed Description

Class that provides functionality to maintain a SVD decomposition of a matrix via successive rank-1 updates and truncation with respect to the truncation threshold eps.

Definition at line 31 of file SVDFactorization.hpp.

Member Typedef Documentation

◆ Matrix

Definition at line 34 of file SVDFactorization.hpp.

◆ Vector

Definition at line 35 of file SVDFactorization.hpp.

Constructor & Destructor Documentation

◆ SVDFactorization()

precice::acceleration::impl::SVDFactorization::SVDFactorization ( double eps,
PtrPreconditioner preconditioner )

Constructor.

Definition at line 12 of file SVDFactorization.cpp.

◆ ~SVDFactorization()

virtual precice::acceleration::impl::SVDFactorization::~SVDFactorization ( )
inlinevirtual

Destructor, empty.

Definition at line 47 of file SVDFactorization.hpp.

Member Function Documentation

◆ cols()

int precice::acceleration::impl::SVDFactorization::cols ( )

: returns the number of columns in the QR-decomposition

Definition at line 286 of file SVDFactorization.cpp.

◆ computeQRdecomposition()

void precice::acceleration::impl::SVDFactorization::computeQRdecomposition ( Matrix const & A,
Matrix & Q,
Matrix & R )
private

: computes the QR decomposition of a matrix A of type A = PSI^T*A \in R^(rank x n)

This method computes a dedicated QR factorization [Q,R] = PSI^T*A. In case of linear dependence, i.e., if the next column cannot be orthogonalized to the columns in Q, this columns is deleted from Q, however not from R. To account for the missing column in Q, the corresponding line in R is deleted, though, the column in R is still filled with the projection weights, such that the dimension (cols) of R_A is aligned with the dimension (cols) of PSI^T*A, from which the matrix K is composed of.

The threshold parameter eps, indicates whether a column is seen to be in the column space of Q via the criterium ||v_orth|| / ||v|| <= eps (cmp. QR2 Filter)

orthogonalize column "col" to columns in Q

  • test if reorthogonalization is necessary - rho0 = |v_init|, t = |r_(i,cols-1)|, rho_orth = |v_orth| rho_orth is small, if the new information incorporated in v is small, i.e., the part of v orthogonal to _Q is small. if rho_orth is very small it is possible, that we are adding (more or less) only round-off errors to the decomposition. Later normalization will scale this new information so that it is equally weighted as the columns in Q. To keep a good orthogonality, the gram-schmidt process is iterated for the bad orthogonalized column v_orth = col'

re-orthogonalize if: ||v_orth|| / ||v|| <= 1/theta

Definition at line 74 of file SVDFactorization.cpp.

Here is the call graph for this function:

◆ getThreshold()

double precice::acceleration::impl::SVDFactorization::getThreshold ( )

: returns the truncation threshold for the SVD

Definition at line 274 of file SVDFactorization.cpp.

◆ getWaste()

int precice::acceleration::impl::SVDFactorization::getWaste ( )

: returns the total number of truncated modes since last call to this method

Definition at line 279 of file SVDFactorization.cpp.

◆ initialize()

void precice::acceleration::impl::SVDFactorization::initialize ( PtrParMatrixOps parMatOps,
int globalRows )

: initializes the updated SVD factorization, i.e., sets the object for parallel matrix-matrix operations and the number of global rows.

Definition at line 20 of file SVDFactorization.cpp.

◆ isSVDinitialized()

bool precice::acceleration::impl::SVDFactorization::isSVDinitialized ( )

Definition at line 264 of file SVDFactorization.cpp.

◆ matrixPhi()

SVDFactorization::Matrix & precice::acceleration::impl::SVDFactorization::matrixPhi ( )

: returns a matrix representation of the orthogonal matrix Phi, A = Psi * Sigma * Phi^T

Definition at line 231 of file SVDFactorization.cpp.

◆ matrixPsi()

SVDFactorization::Matrix & precice::acceleration::impl::SVDFactorization::matrixPsi ( )

: returns a matrix representation of the orthogonal matrix Psi, A = Psi * Sigma * Phi^bs

Definition at line 236 of file SVDFactorization.cpp.

◆ rank()

Rank precice::acceleration::impl::SVDFactorization::rank ( )

: returns the rank of the truncated SVD factorization

Definition at line 296 of file SVDFactorization.cpp.

◆ reset()

void precice::acceleration::impl::SVDFactorization::reset ( )

: resets the SVD factorization

Definition at line 63 of file SVDFactorization.cpp.

◆ rows()

int precice::acceleration::impl::SVDFactorization::rows ( )

: returns the number of rows in the QR-decomposition

Definition at line 291 of file SVDFactorization.cpp.

◆ setApplyFilterQR()

void precice::acceleration::impl::SVDFactorization::setApplyFilterQR ( bool b,
double eps = 1e-3 )

: enables or disables an additional QR-2 filter for the QR-decomposition

Definition at line 251 of file SVDFactorization.cpp.

◆ setfstream()

void precice::acceleration::impl::SVDFactorization::setfstream ( std::fstream * stream)

Optional file-stream for logging output.

◆ setPrecondApplied()

void precice::acceleration::impl::SVDFactorization::setPrecondApplied ( bool b)

: applies the preconditioner to the factorized and truncated representation of the Jacobian matrix

: appplies the inverse preconditioner to the factorized and truncated representation of the Jacobian matrix

Definition at line 246 of file SVDFactorization.cpp.

◆ setThreshold()

void precice::acceleration::impl::SVDFactorization::setThreshold ( double eps)

: sets the threshold for the truncation of the SVD factorization

Definition at line 269 of file SVDFactorization.cpp.

◆ singularValues()

SVDFactorization::Vector & precice::acceleration::impl::SVDFactorization::singularValues ( )

: returns a matrix representation of the orthogonal matrix Sigma, A = Psi * Sigma * Phi^T

Definition at line 241 of file SVDFactorization.cpp.

◆ update()

template<typename Derived1 , typename Derived2 >
void precice::acceleration::impl::SVDFactorization::update ( const Eigen::MatrixBase< Derived1 > & A,
const Eigen::MatrixBase< Derived2 > & B )
inline

: updates the SVD decomposition with the rank-1 update A*B^T, i.e., _psi * _sigma * _phi^T + A*B^T and overrides the internal SVD representation. After the update, the SVD is truncated according to the threshold _truncationEps

updates the truncated svd factorization of the Jacobian with a rank-1 modification

\psi * \sigma * \phi <– \psi * \sigma * \phi + A * B^T

(1): compute orthogonal basis P of (I-\psi\psi^T)A

(2): compute orthogonal basis Q of (I-\phi\phi^T)B

(3) construct matrix $ K \in (K_bar + m -x) x (K_bar +m -y) $ if x .. deleted columns in P -> (m-x) new modes from A (rows of R_A) y .. deleted columns in Q -> (m-y) new modes from B (rows of R_B)

[ \sigma 0] + [ Atil ] * [ Btil ]^T [ 0 0] [ R_A ] [ R_B ] (stored local on each proc).

(4) rotate left and right subspaces

(5) truncation of SVD

Definition at line 55 of file SVDFactorization.hpp.

Here is the call graph for this function:

Member Data Documentation

◆ _applyFilterQR

bool precice::acceleration::impl::SVDFactorization::_applyFilterQR = false
private

Definition at line 293 of file SVDFactorization.hpp.

◆ _cols

int precice::acceleration::impl::SVDFactorization::_cols = 0
private

Number of columns, i.e., rank of the truncated svd.

Definition at line 270 of file SVDFactorization.hpp.

◆ _epsQR2

double precice::acceleration::impl::SVDFactorization::_epsQR2 = 1e-3
private

Threshold for the QR2 filter for the QR decomposition.

Definition at line 282 of file SVDFactorization.hpp.

◆ _globalRows

int precice::acceleration::impl::SVDFactorization::_globalRows = 0
private

Number of global rows, i.e., sum of _rows for all procs.

Definition at line 273 of file SVDFactorization.hpp.

◆ _initialized

bool precice::acceleration::impl::SVDFactorization::_initialized = false
private

true, if ParallelMatrixOperations object is set, i.e., initialized

Definition at line 288 of file SVDFactorization.hpp.

◆ _initialSVD

bool precice::acceleration::impl::SVDFactorization::_initialSVD = false
private

true, if at least one update has been made, i.e., the number of rows is known and a initial rank is given.

Definition at line 291 of file SVDFactorization.hpp.

◆ _log

logging::Logger precice::acceleration::impl::SVDFactorization::_log {"acceleration::SVDFactorization"}
private

Definition at line 253 of file SVDFactorization.hpp.

◆ _parMatrixOps

PtrParMatrixOps precice::acceleration::impl::SVDFactorization::_parMatrixOps = nullptr
private

: object for parallel matrix operations, i.e., parallel mat-mat/ mat-vec multiplications

Definition at line 259 of file SVDFactorization.hpp.

◆ _phi

Matrix precice::acceleration::impl::SVDFactorization::_phi
private

Definition at line 263 of file SVDFactorization.hpp.

◆ _preconditioner

PtrPreconditioner precice::acceleration::impl::SVDFactorization::_preconditioner
private

: preconditioner for least-squares system if vectorial system is used.

Definition at line 256 of file SVDFactorization.hpp.

◆ _preconditionerApplied

bool precice::acceleration::impl::SVDFactorization::_preconditionerApplied = false
private

true if the preconditioner has been applied appropriate to the updated SVD decomposition

Definition at line 285 of file SVDFactorization.hpp.

◆ _psi

Matrix precice::acceleration::impl::SVDFactorization::_psi
private

: SVD factorization of the matrix J = _psi * _sigma * _phi^T

Definition at line 262 of file SVDFactorization.hpp.

◆ _rows

int precice::acceleration::impl::SVDFactorization::_rows = 0
private

Number of rows (on each proc, i.e., local)

Definition at line 267 of file SVDFactorization.hpp.

◆ _sigma

Vector precice::acceleration::impl::SVDFactorization::_sigma
private

Definition at line 264 of file SVDFactorization.hpp.

◆ _truncationEps

double precice::acceleration::impl::SVDFactorization::_truncationEps
private

Truncation parameter for the updated SVD decomposition.

Definition at line 279 of file SVDFactorization.hpp.

◆ _waste

int precice::acceleration::impl::SVDFactorization::_waste = 0
private

Definition at line 276 of file SVDFactorization.hpp.


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