preCICE v3.1.2
Loading...
Searching...
No Matches
IQNILSAcceleration.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <Eigen/Core>
4#include <map>
5#include <vector>
9
10namespace precice {
11namespace acceleration {
12
27public:
29 double initialRelaxation,
30 bool forceInitialRelaxation,
31 int maxIterationsUsed,
32 int pastTimeWindowsReused,
33 int filter,
34 double singularityLimit,
35 std::vector<int> dataIDs,
36 impl::PtrPreconditioner preconditioner);
37
39
41 virtual void initialize(const DataMap &cplData);
42
49 virtual void specializedIterationsConverged(const DataMap &cplData);
50
51private:
54
55 // @brief Secondary data x-tilde deltas.
56 //
57 // Stores x-tilde deltas for data not involved in least-squares computation.
60
62 virtual void updateDifferenceMatrices(const DataMap &cplData);
63
65 virtual void computeQNUpdate(const DataMap &cplData, Eigen::VectorXd &xUpdate);
66
68 virtual void computeUnderrelaxationSecondaryData(const DataMap &cplData);
69
71 virtual void removeMatrixColumn(int columnIndex);
72};
73} // namespace acceleration
74} // namespace precice
Base Class for quasi-Newton acceleration schemes.
Interface quasi-Newton with interface least-squares approximation.
std::map< int, Eigen::MatrixXd > _secondaryMatricesWBackup
virtual void initialize(const DataMap &cplData)
Initializes the acceleration.
virtual void updateDifferenceMatrices(const DataMap &cplData)
updates the V, W matrices (as well as the matrices for the secondary data)
virtual void removeMatrixColumn(int columnIndex)
Removes one iteration from V,W matrices and adapts _matrixCols.
virtual void computeUnderrelaxationSecondaryData(const DataMap &cplData)
computes underrelaxation for the secondary data
IQNILSAcceleration(double initialRelaxation, bool forceInitialRelaxation, int maxIterationsUsed, int pastTimeWindowsReused, int filter, double singularityLimit, std::vector< int > dataIDs, impl::PtrPreconditioner preconditioner)
std::map< int, Eigen::MatrixXd > _secondaryMatricesW
virtual void computeQNUpdate(const DataMap &cplData, Eigen::VectorXd &xUpdate)
computes the IQN-ILS update using QR decomposition
std::map< int, Eigen::VectorXd > _secondaryOldXTildes
Secondary data solver output from last iteration.
virtual void specializedIterationsConverged(const DataMap &cplData)
Marks a iteration sequence as converged.
Main namespace of the precice library.