|
preCICE v3.3.0
|
Interface for preconditioner variants that can be applied to quasi-Newton acceleration schemes. More...
#include <Preconditioner.hpp>
Public Member Functions | |
| Preconditioner (int maxNonConstTimeWindows) | |
| virtual | ~Preconditioner ()=default |
| Destructor, empty. | |
| virtual void | initialize (std::vector< size_t > &svs) |
| initialize the preconditioner | |
| void | apply (Eigen::MatrixXd &M, bool transpose) |
| Apply preconditioner to matrix. | |
| void | revert (Eigen::MatrixXd &M, bool transpose) |
| Apply inverse preconditioner to matrix. | |
| void | apply (Eigen::MatrixXd &M) |
| To transform physical values to balanced values. Matrix version. | |
| void | apply (Eigen::VectorXd &v) |
| To transform physical values to balanced values. Vector version. | |
| void | revert (Eigen::MatrixXd &M) |
| To transform balanced values back to physical values. Matrix version. | |
| void | revert (Eigen::VectorXd &v) |
| To transform balanced values back to physical values. Vector version. | |
| void | update (bool timeWindowComplete, const Eigen::VectorXd &oldValues, const Eigen::VectorXd &res) |
| Update the scaling after every FSI iteration and require a new QR decomposition (if necessary) | |
| bool | requireNewQR () |
| returns true if a QR decomposition from scratch is necessary | |
| void | newQRfulfilled () |
| to tell the preconditioner that QR-decomposition has been recomputed | |
| std::vector< double > & | getWeights () |
| bool | isConst () |
Protected Member Functions | |
| virtual void | _update_ (bool timeWindowComplete, const Eigen::VectorXd &oldValues, const Eigen::VectorXd &res)=0 |
| Update the scaling after every FSI iteration and require a new QR decomposition (if necessary) | |
Protected Attributes | |
| std::vector< double > | _weights |
| Weights used to scale the matrix V and the residual. | |
| std::vector< double > | _invWeights |
| Inverse weights (for efficiency reasons) | |
| std::vector< size_t > | _subVectorSizes |
| Sizes of each sub-vector, i.e. each coupling data. | |
| std::vector< size_t > | _subVectorOffsets |
| Offsets of each sub-vector in concatenated data, i.e. each coupling data. | |
| int | _maxNonConstTimeWindows |
| maximum number of non-const time windows, i.e., after this number of time windows, the preconditioner is frozen with the current weights and becomes a constant preconditioner | |
| int | _nbNonConstTimeWindows = 0 |
| Counts the number of completed time windows with a non-const weighting. | |
| bool | _requireNewQR = false |
| True if a QR decomposition from scratch is necessary. | |
| bool | _frozen = false |
| True if _nbNonConstTimeWindows >= _maxNonConstTimeWindows, i.e., preconditioner is not updated any more. | |
Private Attributes | |
| logging::Logger | _log {"acceleration::Preconditioner"} |
Interface for preconditioner variants that can be applied to quasi-Newton acceleration schemes.
Preconditioning (formerly also known as scaling) improves the balance between several sub-vectors for parallel or multi coupling.
apply() applies the weighting, i.e. transforms from physical values to balanced values. revert() reverts the weighting, i.e. transforms from balanced values back to physical values. update() updates the preconditioner, after every FSI iteration (though some variants might only be updated after a complete time window)
Definition at line 23 of file Preconditioner.hpp.
|
inline |
Definition at line 25 of file Preconditioner.hpp.
|
virtualdefault |
Destructor, empty.
|
protectedpure virtual |
Update the scaling after every FSI iteration and require a new QR decomposition (if necessary)
| [in] | timeWindowComplete | True if this FSI iteration also completed a time windows |
Implemented in precice::acceleration::impl::ConstantPreconditioner, precice::acceleration::impl::ResidualPreconditioner, precice::acceleration::impl::ResidualSumPreconditioner, and precice::acceleration::impl::ValuePreconditioner.
|
inline |
To transform physical values to balanced values. Matrix version.
Definition at line 112 of file Preconditioner.hpp.
|
inline |
Apply preconditioner to matrix.
| transpose | false = from left, true = from right |
Definition at line 58 of file Preconditioner.hpp.
|
inline |
To transform physical values to balanced values. Vector version.
Definition at line 127 of file Preconditioner.hpp.
|
inline |
Definition at line 204 of file Preconditioner.hpp.
|
inlinevirtual |
initialize the preconditioner
| size | of the pp system (e.g. rows of V) |
Reimplemented in precice::acceleration::impl::ConstantPreconditioner, and precice::acceleration::impl::ResidualSumPreconditioner.
Definition at line 37 of file Preconditioner.hpp.
|
inline |
Definition at line 209 of file Preconditioner.hpp.
|
inline |
to tell the preconditioner that QR-decomposition has been recomputed
Definition at line 199 of file Preconditioner.hpp.
|
inline |
returns true if a QR decomposition from scratch is necessary
Definition at line 192 of file Preconditioner.hpp.
|
inline |
To transform balanced values back to physical values. Matrix version.
Definition at line 140 of file Preconditioner.hpp.
|
inline |
Apply inverse preconditioner to matrix.
| transpose | false = from left, true = from right |
Definition at line 86 of file Preconditioner.hpp.
|
inline |
To transform balanced values back to physical values. Vector version.
Definition at line 155 of file Preconditioner.hpp.
|
inline |
Update the scaling after every FSI iteration and require a new QR decomposition (if necessary)
| [in] | timeWindowComplete | True if this FSI iteration also completed a time window |
Definition at line 172 of file Preconditioner.hpp.
|
protected |
True if _nbNonConstTimeWindows >= _maxNonConstTimeWindows, i.e., preconditioner is not updated any more.
Definition at line 239 of file Preconditioner.hpp.
|
protected |
Inverse weights (for efficiency reasons)
Definition at line 219 of file Preconditioner.hpp.
|
private |
Definition at line 249 of file Preconditioner.hpp.
|
protected |
maximum number of non-const time windows, i.e., after this number of time windows, the preconditioner is frozen with the current weights and becomes a constant preconditioner
Definition at line 230 of file Preconditioner.hpp.
|
protected |
Counts the number of completed time windows with a non-const weighting.
Definition at line 233 of file Preconditioner.hpp.
|
protected |
True if a QR decomposition from scratch is necessary.
Definition at line 236 of file Preconditioner.hpp.
|
protected |
Offsets of each sub-vector in concatenated data, i.e. each coupling data.
Definition at line 225 of file Preconditioner.hpp.
|
protected |
Sizes of each sub-vector, i.e. each coupling data.
Definition at line 222 of file Preconditioner.hpp.
|
protected |
Weights used to scale the matrix V and the residual.
Definition at line 216 of file Preconditioner.hpp.