preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
BaseQNAcceleration.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <Eigen/Core>
4#include <algorithm>
5#include <deque>
6#include <fstream>
7#include <map>
8#include <sstream>
9#include <string>
10#include <vector>
14#include "logging/Logger.hpp"
15#include "time/TimeGrids.hpp"
16
17/* ****************************************************************************
18 *
19 * A few comments concerning the present design choice.
20 *
21 * All the functions from the base class BAseQNAcceleration are specialized in
22 * the sub classes as needed. This is done vi overwriting the base functions in
23 * the specialized sub classes and calling the respective base function after
24 * performing the specialized stuff (in order to perform the common, generalized
25 * computations i.e. handling of V,W matrices etc.)
26 * However, for the performAcceleration Method we decided (for better readability)
27 * to have this method only in the base class, while introducing a function
28 * performPPSecondaryData that handles all the specialized stuff concerning acceleration
29 * processing for the secondary data in the sub classes.
30 *
31 * Another possibility would have been to introduce a bunch of functions like
32 * initializeSpecialized(), removeMatrixColumnSpecialized(),
33 * iterationsConvergedSpecialized(), etc
34 * and call those function from the base class top down to the sub classes.
35 *
36 * The third possibility was to separate the approximation of the Jacobian from
37 * the common stuff like handling V,W matrices in the acceleration.
38 * Here, we have a class QNAcceleration that handles the V,W stuff and the basic
39 * scheme of the QN update. Furthermore we have a base class (or rather interface)
40 * JacobianApproximation with sub classes IQNIMVJAPX and IQNAPX that handle all the
41 * specialized stuff like Jacobian approximation, handling of secondary data etc.
42 * However, this approach is not feasible, as we have to call the function
43 * removeMatrixColumn() down in the specialized sub classes IQNIMVJApx and IQNApx.
44 * This is not possible as the function works on the V, W matrices that are
45 * completely treated by QNAcceleration.
46 *
47 * ****************************************************************************
48 */
49
50// ----------------------------------------------------------- CLASS DEFINITION
51
52namespace precice {
53namespace io {
54class TXTReader;
55class TXTWriter;
56} // namespace io
57
58namespace acceleration {
59
65public:
67 double initialRelaxation,
68 bool forceInitialRelaxation,
69 int maxIterationsUsed,
70 int timeWindowsReused,
71 int filter,
72 double singularityLimit,
73 std::vector<int> dataIDs,
74 impl::PtrPreconditioner preconditioner,
75 bool reducedTimeGrid);
76
81 {
82 // not necessary for user, only for developer, if needed, this should be configurable
83 // if (utils::IntraComm::isPrimary() || !utils::IntraComm::isParallel()) {
84 // _infostream.open("precice-accelerationInfo.log", std::ios_base::out);
85 // _infostream << std::setprecision(16);
86 // _infostream << _infostringstream.str();
87 // }
88 }
89
93 std::vector<int> getPrimaryDataIDs() const final override
94 {
95 return _primaryDataIDs;
96 }
97
101 void initialize(const DataMap &cplData) final override;
102
108 void performAcceleration(DataMap &cplData, double windowStart, double windowEnd) final override;
109
116 void iterationsConverged(const DataMap &cplData, double windowStart) final override;
117
121 void exportState(io::TXTWriter &writer) final override;
122
128 void importState(io::TXTReader &reader) final override;
129
131 int getDeletedColumns() const;
132
134 int getDroppedColumns() const;
135
143 int getLSSystemCols() const;
144
148 int getMaxUsedIterations() const;
149
153 int getMaxUsedTimeWindows() const;
154
155protected:
156 logging::Logger _log{"acceleration::BaseQNAcceleration"};
157
160
162 const double _initialRelaxation;
163
166
169
172
175
177 bool _firstIteration = true;
178
179 /* @brief Indicates the first time window, where constant relaxation is used
180 * later, we replace the constant relaxation by a qN-update from last time window.
181 */
182 bool _firstTimeWindow = true;
183
184 /*
185 * @brief True if this process has nodes at the coupling interface
186 */
188
189 /* @brief If true, the QN-scheme always performs a underrelaxation in the first iteration of
190 * a new time window. Otherwise, the LS system from the previous time window is used in the
191 * first iteration.
192 */
194
198 bool _resetLS = false;
199
201 Eigen::VectorXd _oldPrimaryXTilde;
202
204 Eigen::VectorXd _oldXTilde;
205
207 Eigen::VectorXd _primaryResiduals;
208
210 Eigen::VectorXd _residuals; // @todo is this member still needed? Potential refactoring.
211
213 Eigen::MatrixXd _matrixV;
214
216 Eigen::MatrixXd _matrixW;
217
220
223
227 const int _filter;
228
235 const double _singularityLimit;
236
245
250
254
258
259 int getLSSystemRows() const;
260 int getPrimaryLSSystemRows() const;
261
268 virtual void specializedIterationsConverged(const DataMap &cplData) = 0;
269
271 virtual void updateDifferenceMatrices(const DataMap &cplData);
272
274 virtual void updateCouplingData(const DataMap &cplData, double windowStart);
275
277 virtual void applyFilter();
278
280 virtual void computeQNUpdate(Eigen::VectorXd &xUpdate) = 0;
281
283 virtual void removeMatrixColumn(int columnIndex);
284
286 void writeInfo(const std::string &s, bool allProcs = false);
287
288 int its = 0, tWindows = 0;
289
290private:
293 void initializeVectorsAndPreconditioner(const DataMap &cplData, double windowStart);
294
301
303 void concatenateCouplingData(Eigen::VectorXd &data, Eigen::VectorXd &oldData, const DataMap &cplData, std::vector<int> dataIDs, precice::time::TimeGrids timeGrids, double windowStart) const;
304
307
311
313 Eigen::VectorXd _primaryValues;
314
316 Eigen::VectorXd _oldPrimaryValues;
317
319 Eigen::VectorXd _oldPrimaryResiduals;
320
322 Eigen::VectorXd _values;
323
325 Eigen::VectorXd _oldValues;
326
331 Eigen::MatrixXd _matrixVBackup;
332 Eigen::MatrixXd _matrixWBackup;
334
336 int _nbDelCols = 0;
337
339 int _nbDropCols = 0;
340};
341} // namespace acceleration
342} // namespace precice
Base Class for quasi-Newton acceleration schemes.
virtual void updateDifferenceMatrices(const DataMap &cplData)
Updates the V, W matrices (as well as the matrices for the secondary data)
const double _initialRelaxation
Constant relaxation factor used for first iteration.
virtual void computeQNUpdate(Eigen::VectorXd &xUpdate)=0
Computes the quasi-Newton update using the specified pp scheme (IQNIMVJ, IQNILS)
std::vector< int > _dimOffsetsPrimary
Stores the local dimensions regarding primary data,.
Eigen::VectorXd _oldPrimaryValues
Concatenation of all old primary data involved in the QN system.
std::vector< int > getPrimaryDataIDs() const final override
Returns all IQN involved data IDs.
int getLSSystemCols() const
: computes number of cols in least squares system, i.e, number of cols in _matrixV,...
int _nbDelCols
Number of filtered out columns in this time window.
void concatenateCouplingData(Eigen::VectorXd &data, Eigen::VectorXd &oldData, const DataMap &cplData, std::vector< int > dataIDs, precice::time::TimeGrids timeGrids, double windowStart) const
Samples and concatenates the data and old data in cplData into a long vector.
virtual void specializedIterationsConverged(const DataMap &cplData)=0
Marks a iteration sequence as converged.
virtual void applyFilter()
Applies the filter method for the least-squares system, defined in the configuration.
Eigen::VectorXd _oldValues
Concatenation of all old primary and secondary data involved in the QN system.
int getMaxUsedTimeWindows() const
Get the maximum number of reused time windows.
~BaseQNAcceleration() override
Destructor, empty.
int getMaxUsedIterations() const
Get the maximum number of reused iterations.
Eigen::MatrixXd _matrixVBackup
backup of the V,W and matrixCols data structures. Needed for the skipping of initial relaxation,...
Eigen::VectorXd _oldXTilde
Solver output from last iteration.
Eigen::VectorXd _primaryResiduals
Current iteration residuals of primary IQN data. Temporary.
std::ostringstream _infostringstream
write some debug/acceleration info to file
virtual void updateCouplingData(const DataMap &cplData, double windowStart)
Splits up QN system vector back into the waveforms in coupling data.
bool _resetLS
If true, the LS system has been modified (reset or recomputed) in such a way, that mere updating of m...
void importState(io::TXTReader &reader) final override
Imports the last exported state of the acceleration from file.
Eigen::VectorXd _primaryValues
Concatenation of all primary data involved in the QN system.
int _nbDropCols
Number of dropped columns in this time window (old time window out of scope)
int getDroppedColumns() const
how many QN columns were dropped (went out of scope) in this time window
const int _filter
filter method that is used to maintain good conditioning of the least-squares system Either of two ty...
Eigen::MatrixXd _matrixV
Stores residual deltas.
std::vector< int > _dimOffsets
Stores the local dimensions, i.e., the offsets in _invJacobian for all processors.
const double _singularityLimit
Determines sensitivity when two matrix columns are considered equal.
const int _maxIterationsUsed
Maximum number of old data iterations kept.
void exportState(io::TXTWriter &writer) final override
Exports the current state of the acceleration to a file.
Eigen::VectorXd _residuals
Current iteration residuals of IQN data. Temporary.
std::optional< time::TimeGrids > _primaryTimeGrids
Stores the time grids to which the primary data involved in the QN system will be interpolated to....
virtual void specializedInitializeVectorsAndPreconditioner(const DataMap &cplData)=0
handles the initialization of matrices and vectors in the sub-classes
virtual void removeMatrixColumn(int columnIndex)
Removes one iteration from V,W matrices and adapts _matrixCols.
Eigen::MatrixXd _matrixW
Stores x tilde deltas, where x tilde are values computed by solvers.
Eigen::VectorXd _values
Concatenation of all primary and secondary data involved in the QN system.
std::optional< time::TimeGrids > _timeGrids
Stores the time grids to which the primary and secondary data involved in the QN system will be inter...
void initialize(const DataMap &cplData) final override
Initializes the acceleration.
const int _timeWindowsReused
Maximum number of old time windows (with data values) kept.
BaseQNAcceleration(double initialRelaxation, bool forceInitialRelaxation, int maxIterationsUsed, int timeWindowsReused, int filter, double singularityLimit, std::vector< int > dataIDs, impl::PtrPreconditioner preconditioner, bool reducedTimeGrid)
const std::vector< int > _primaryDataIDs
Data IDs of primary data to be involved in the IQN coefficient computation.
void writeInfo(const std::string &s, bool allProcs=false)
Writes info to the _infostream (also in parallel)
std::deque< int > _matrixCols
Indices (of columns in W, V matrices) of 1st iterations of time windows.
bool _firstIteration
Indicates the first iteration, where constant relaxation is used.
int getDeletedColumns() const
how many QN columns were deleted in this time window
const bool _reducedTimeGrid
if _reducedTimeGrid = false uses the full QN-WI and if _reducedTimeGrid = true uses rQN-WI form the p...
impl::PtrPreconditioner _preconditioner
Preconditioner for least-squares system if vectorial system is used.
impl::QRFactorization _qrV
Stores the current QR decomposition ov _matrixV, can be updated via deletion/insertion of columns.
std::vector< int > _dataIDs
Data IDs of all coupling data.
void iterationsConverged(const DataMap &cplData, double windowStart) final override
Marks a iteration sequence as converged.
void performAcceleration(DataMap &cplData, double windowStart, double windowEnd) final override
Performs one acceleration step.
void initializeVectorsAndPreconditioner(const DataMap &cplData, double windowStart)
Initializes the vectors, matrices and preconditioner This has to be done after the first iteration of...
Eigen::VectorXd _oldPrimaryResiduals
Difference between solver input and output from last time window regarding primary data.
Eigen::VectorXd _oldPrimaryXTilde
Solver output regarding primary data from last iteration.
Class that provides functionality for a dynamic QR-decomposition, that can be updated in O(mn) flops ...
File reader for matrix/vector in Matlab V7 ASCII format.
Definition TXTReader.hpp:15
File writer for matrix in Matlab V7 ASCII format.
Definition TXTWriter.hpp:13
This class provides a lightweight logger.
Definition Logger.hpp:17
Interface for storing the time grids in the Quasi-Newton and Aitken methods. A time grid is a ordered...
Definition TimeGrids.hpp:21
Main namespace of the precice library.