preCICE v3.1.2
Loading...
Searching...
No Matches
IQNIMVJAcceleration.hpp
Go to the documentation of this file.
1#ifndef PRECICE_NO_MPI
2
3#pragma once
4
5#include <Eigen/Core>
6#include <deque>
7#include <vector>
13#include "com/SharedPointer.hpp"
14
15// ----------------------------------------------------------- CLASS DEFINITION
16
17namespace precice {
18namespace acceleration {
19
34public:
35 static const int NO_RESTART = 0;
36 static const int RS_ZERO = 1;
37 static const int RS_LS = 2;
38 static const int RS_SVD = 3;
39 static const int RS_SLIDE = 4;
40
45 double initialRelaxation,
46 bool forceInitialRelaxation,
47 int maxIterationsUsed,
48 int pastTimeWindowsReused,
49 int filter,
50 double singularityLimit,
51 std::vector<int> dataIDs,
52 const impl::PtrPreconditioner &preconditioner,
53 bool alwaysBuildJacobian,
54 int imvjRestartType,
55 int chunkSize,
56 int RSLSreusedTimeWindows,
57 double RSSVDtruncationEps);
58
63
67 virtual void initialize(const DataMap &cplData);
68
75 virtual void specializedIterationsConverged(const DataMap &cplData);
76
77private:
79 Eigen::MatrixXd _invJacobian;
80
82 Eigen::MatrixXd _oldInvJacobian;
83
85 Eigen::MatrixXd _Wtil;
86
89
92
94 Eigen::MatrixXd _matrixV_RSLS;
95
97 Eigen::MatrixXd _matrixW_RSLS;
98
101
104
107
113
121
127
130
133
136
137 // DEBUG
138 // std::fstream _info2;
139 double _avgRank;
140
144 virtual void computeQNUpdate(const DataMap &cplData, Eigen::VectorXd &xUpdate);
145
147 virtual void updateDifferenceMatrices(const DataMap &cplData);
148
150 virtual void computeUnderrelaxationSecondaryData(const DataMap &cplData);
151
159 void computeNewtonUpdate(const DataMap &cplData, Eigen::VectorXd &update);
160
169 void computeNewtonUpdateEfficient(const DataMap &cplData, Eigen::VectorXd &update);
170
173 void pseudoInverse(Eigen::MatrixXd &pseudoInverse);
174
177 void buildJacobian();
178
181 void buildWtil();
182
189 void restartIMVJ();
190
192 virtual void removeMatrixColumn(int columnIndex);
193
195 void removeMatrixColumnRSLS(int columnINdex);
196};
197} // namespace acceleration
198} // namespace precice
199
200#endif /* PRECICE_NO_MPI */
Base Class for quasi-Newton acceleration schemes.
Multi vector quasi-Newton update scheme.
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)
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.