54 void apply(Eigen::MatrixXd &M,
bool transpose)
61 for (
int i = 0; i < validCols; i++) {
62 for (
int j = 0; j < M.rows(); j++) {
70 for (
int i = 0; i < M.cols(); i++) {
71 for (
int j = 0; j < validRows; j++) {
82 void revert(Eigen::MatrixXd &M,
bool transpose)
90 for (
int i = 0; i < validCols; i++) {
91 for (
int j = 0; j < M.rows(); j++) {
99 for (
int i = 0; i < M.cols(); i++) {
100 for (
int j = 0; j < validRows; j++) {
114 int validRows =
std::min(
static_cast<int>(M.rows()), (
int)
_weights.size());
115 for (
int i = 0; i < M.cols(); i++) {
116 for (
int j = 0; j < validRows; j++) {
129 int validSize =
std::min(
static_cast<int>(v.size()), (
int)
_weights.size());
130 for (
int j = 0; j < validSize; j++) {
142 int validRows =
std::min(
static_cast<int>(M.rows()), (
int)
_weights.size());
143 for (
int i = 0; i < M.cols(); i++) {
144 for (
int j = 0; j < validRows; j++) {
157 int validSize =
std::min(
static_cast<int>(v.size()), (
int)
_weights.size());
158 for (
int j = 0; j < validSize; j++) {
168 void update(
bool timeWindowComplete,
const Eigen::VectorXd &oldValues,
const Eigen::VectorXd &res)
177 if (timeWindowComplete) {
184 _update_(timeWindowComplete, oldValues, res);
239 virtual void _update_(
bool timeWindowComplete,
const Eigen::VectorXd &oldValues,
const Eigen::VectorXd &res) = 0;
#define PRECICE_DEBUG_IF(condition,...)
#define PRECICE_TRACE(...)
void apply(Eigen::MatrixXd &M, bool transpose)
Apply preconditioner to matrix.
void revert(Eigen::MatrixXd &M, bool transpose)
Apply inverse preconditioner to matrix.
int _nbNonConstTimeWindows
Counts the number of completed time windows with a non-const weighting.
std::vector< size_t > _subVectorSizes
Sizes of each sub-vector, i.e. each coupling data.
void newQRfulfilled()
to tell the preconditioner that QR-decomposition has been recomputed
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)
Preconditioner(int maxNonConstTimeWindows)
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)
virtual void initialize(std::vector< size_t > &svs)
initialize the preconditioner
std::vector< double > & getWeights()
virtual ~Preconditioner()=default
Destructor, empty.
void revert(Eigen::MatrixXd &M)
To transform balanced values back to physical values. Matrix version.
std::vector< double > _invWeights
Inverse weights (for efficiency reasons)
bool requireNewQR()
returns true if a QR decomposition from scratch is necessary
int _maxNonConstTimeWindows
maximum number of non-const time windows, i.e., after this number of time windows,...
std::vector< double > _weights
Weights used to scale the matrix V and the residual.
void revert(Eigen::VectorXd &v)
To transform balanced values back to physical values. Vector version.
bool _requireNewQR
True if a QR decomposition from scratch is necessary.
void apply(Eigen::VectorXd &v)
To transform physical values to balanced values. Vector version.
void apply(Eigen::MatrixXd &M)
To transform physical values to balanced values. Matrix version.
bool _frozen
True if _nbNonConstTimeWindows >= _maxNonConstTimeWindows, i.e., preconditioner is not updated any mo...
This class provides a lightweight logger.