preCICE
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::acceleration {
18
33public:
34 static const int NO_RESTART = 0;
35 static const int RS_ZERO = 1;
36 static const int RS_LS = 2;
37 static const int RS_SVD = 3;
38 static const int RS_SLIDE = 4;
39
44 double initialRelaxation,
45 bool forceInitialRelaxation,
46 int maxIterationsUsed,
47 int pastTimeWindowsReused,
48 int filter,
49 double singularityLimit,
50 std::vector<int> dataIDs,
51 OnBoundViolation onBoundViolation,
52 const impl::PtrPreconditioner &preconditioner,
53 bool alwaysBuildJacobian,
54 int imvjRestartType,
55 int chunkSize,
56 int RSLSreusedTimeWindows,
57 double RSSVDtruncationEps,
58 bool reducedTimeGrid);
59
64
71 void specializedIterationsConverged(const DataMap &cplData) override;
72
73private:
75 Eigen::MatrixXd _invJacobian;
76
78 Eigen::MatrixXd _oldInvJacobian;
79
81 Eigen::MatrixXd _Wtil;
82
84 std::vector<Eigen::MatrixXd> _WtilChunk;
85
87 std::vector<Eigen::MatrixXd> _pseudoInverseChunk;
88
90 Eigen::MatrixXd _matrixV_RSLS;
91
93 Eigen::MatrixXd _matrixW_RSLS;
94
96 std::deque<int> _matrixCols_RSLS;
97
100
103
109
117
122 const bool _imvjRestart;
123
126
129
132
133 // DEBUG
134 // std::fstream _info2;
135 double _avgRank;
136
140 void computeQNUpdate(Eigen::VectorXd &xUpdate) override;
141
143 void updateDifferenceMatrices(const DataMap &cplData) override;
144
152 void computeNewtonUpdate(Eigen::VectorXd &update);
153
162 void computeNewtonUpdateEfficient(Eigen::VectorXd &update);
163
166 void pseudoInverse(Eigen::MatrixXd &pseudoInverse);
167
170 void buildJacobian();
171
174 void buildWtil();
175
182 void restartIMVJ();
183
185 void removeMatrixColumn(int columnIndex) override;
186
188 void specializedInitializeVectorsAndPreconditioner(const DataMap &cplData) final override;
189
191 void removeMatrixColumnRSLS(int columnINdex);
192};
193} // namespace precice::acceleration
194
195#endif /* PRECICE_NO_MPI */
OnBoundViolation
Options for handling bound violations.
std::map< int, cplscheme::PtrCouplingData > DataMap
Map from data ID to data values.
BaseQNAcceleration(double initialRelaxation, bool forceInitialRelaxation, int maxIterationsUsed, int timeWindowsReused, int filter, double singularityLimit, std::vector< int > dataIDs, OnBoundViolation onBoundViolation, impl::PtrPreconditioner preconditioner, bool reducedTimeGrid)
const bool _imvjRestart
: If true, the imvj method is used with the restart chunk based approach that avoids to explicitly bu...
int _RSLSreusedTimeWindows
: Number of reused time windows at restart if restart-mode = RS-LS
const int _imvjRestartType
: Indicates the type of the imvj restart-mode:
void updateDifferenceMatrices(const DataMap &cplData) override
: 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
~IQNIMVJAcceleration() override
Destructor, empty.
Eigen::MatrixXd _matrixW_RSLS
stores columns from previous _RSLSreusedTimeWindows time windows if RS-LS restart-mode is active
int _nbRestarts
tracks the number of restarts of IMVJ
void computeQNUpdate(Eigen::VectorXd &xUpdate) override
: computes the IQNIMVJ update using QR decomposition of V, furthermore it updates the inverse of the ...
Eigen::MatrixXd _invJacobian
stores the approximation of the inverse Jacobian of the system at current time window.
void removeMatrixColumn(int columnIndex) override
: Removes one iteration from V,W matrices and adapts _matrixCols.
void specializedInitializeVectorsAndPreconditioner(const DataMap &cplData) final override
IQNIMVJAcceleration(double initialRelaxation, bool forceInitialRelaxation, int maxIterationsUsed, int pastTimeWindowsReused, int filter, double singularityLimit, std::vector< int > dataIDs, OnBoundViolation onBoundViolation, const impl::PtrPreconditioner &preconditioner, bool alwaysBuildJacobian, int imvjRestartType, int chunkSize, int RSLSreusedTimeWindows, double RSSVDtruncationEps, bool reducedTimeGrid)
Constructor.
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
void computeNewtonUpdate(Eigen::VectorXd &update)
: computes the quasi-Newton update vector based on the matrices V and W using a QR decomposition of V...
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
void specializedIterationsConverged(const DataMap &cplData) override
Marks a iteration sequence as converged.
bool _alwaysBuildJacobian
If true, the less efficient method to compute the quasi-Newton update is used, that explicitly builds...
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
void computeNewtonUpdateEfficient(Eigen::VectorXd &update)
: computes the quasi-Newton update vector based on the same numerics as above. However,...
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
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
Class that provides functionality to maintain a SVD decomposition of a matrix via successive rank-1 u...
std::shared_ptr< Preconditioner > PtrPreconditioner
std::shared_ptr< ParallelMatrixOperations > PtrParMatrixOps
contains implementations of acceleration schemes.