18namespace acceleration {
45 double initialRelaxation,
46 bool forceInitialRelaxation,
47 int maxIterationsUsed,
48 int pastTimeWindowsReused,
50 double singularityLimit,
53 bool alwaysBuildJacobian,
56 int RSLSreusedTimeWindows,
57 double RSSVDtruncationEps);
Base Class for quasi-Newton acceleration schemes.
Multi vector quasi-Newton update scheme.
static const int RS_SLIDE
int _RSLSreusedTimeWindows
: Number of reused time windows at restart if restart-mode = RS-LS
void computeNewtonUpdateEfficient(const DataMap &cplData, Eigen::VectorXd &update)
: computes the quasi-Newton update vector based on the same numerics as above. However,...
bool _imvjRestart
: If true, the imvj method is used with the restart chunk based approach that avoids to explicitly bu...
virtual void updateDifferenceMatrices(const DataMap &cplData)
: updates the V, W matrices (as well as the matrices for the secondary data)
static const int NO_RESTART
void pseudoInverse(Eigen::MatrixXd &pseudoInverse)
: computes the pseudo inverse of V multiplied with V^T, i.e., Z = (V^TV)^-1V^T via QR-dec
virtual void initialize(const DataMap &cplData)
Initializes the acceleration.
virtual void computeQNUpdate(const DataMap &cplData, Eigen::VectorXd &xUpdate)
: computes the IQNIMVJ update using QR decomposition of V, furthermore it updates the inverse of the ...
Eigen::MatrixXd _matrixW_RSLS
stores columns from previous _RSLSreusedTimeWindows time windows if RS-LS restart-mode is active
IQNIMVJAcceleration(double initialRelaxation, bool forceInitialRelaxation, int maxIterationsUsed, int pastTimeWindowsReused, int filter, double singularityLimit, std::vector< int > dataIDs, const impl::PtrPreconditioner &preconditioner, bool alwaysBuildJacobian, int imvjRestartType, int chunkSize, int RSLSreusedTimeWindows, double RSSVDtruncationEps)
Constructor.
int _nbRestarts
tracks the number of restarts of IMVJ
int _imvjRestartType
: Indicates the type of the imvj restart-mode:
Eigen::MatrixXd _invJacobian
stores the approximation of the inverse Jacobian of the system at current time window.
virtual void computeUnderrelaxationSecondaryData(const DataMap &cplData)
: computes underrelaxation for the secondary data
virtual void removeMatrixColumn(int columnIndex)
: Removes one iteration from V,W matrices and adapts _matrixCols.
impl::PtrParMatrixOps _parMatrixOps
encapsulates matrix-matrix and matrix-vector multiplications for serial and parallel execution
int _chunkSize
: Number of time windows between restarts for the imvj method in restart mode
Eigen::MatrixXd _Wtil
stores the sub result (W-J_prev*V) for the current iteration
void removeMatrixColumnRSLS(int columnINdex)
: Removes one column form the V_RSLS and W_RSLS matrices and adapts _matrixCols_RSLS
void buildWtil()
: re-computes the matrix _Wtil = ( W - J_prev * V) instead of updating it according to V
bool _alwaysBuildJacobian
If true, the less efficient method to compute the quasi-Newton update is used, that explicitly builds...
virtual void specializedIterationsConverged(const DataMap &cplData)
Marks a iteration sequence as converged.
void restartIMVJ()
: restarts the imvj method, i.e., drops all stored matrices Wtil and Z and computes a initial guess o...
void buildJacobian()
: computes a explicit representation of the Jacobian, i.e., n x n matrix
std::vector< Eigen::MatrixXd > _pseudoInverseChunk
stores all pseudo inverses within the current chunk of the imvj restart mode, disabled if _imvjRestar...
Eigen::MatrixXd _oldInvJacobian
stores the approximation of the inverse Jacobian from the previous time window.
Eigen::MatrixXd _matrixV_RSLS
stores columns from previous _RSLSreusedTimeWindows time windows if RS-LS restart-mode is active
virtual ~IQNIMVJAcceleration()
Destructor, empty.
std::deque< int > _matrixCols_RSLS
number of cols per time window
std::vector< Eigen::MatrixXd > _WtilChunk
stores all Wtil matrices within the current chunk of the imvj restart mode, disabled if _imvjRestart ...
impl::SVDFactorization _svdJ
holds and maintains a truncated SVD decomposition of the Jacobian matrix
void computeNewtonUpdate(const DataMap &cplData, Eigen::VectorXd &update)
: computes the quasi-Newton update vector based on the matrices V and W using a QR decomposition of V...
Class that provides functionality to maintain a SVD decomposition of a matrix via successive rank-1 u...
Main namespace of the precice library.