preCICE v3.1.2
Loading...
Searching...
No Matches
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
16/* ****************************************************************************
17 *
18 * A few comments concerning the present design choice.
19 *
20 * All the functions from the base class BAseQNAcceleration are specialized in
21 * the sub classes as needed. This is done vi overwriting the base functions in
22 * the specialized sub classes and calling the respective base function after
23 * performing the specialized stuff (in order to perform the common, generalized
24 * computations i.e. handling of V,W matrices etc.)
25 * However, for the performAcceleration Method we decided (for better readability)
26 * to have this method only in the base class, while introducing a function
27 * performPPSecondaryData that handles all the specialized stuff concerning acceleration
28 * processing for the secondary data in the sub classes.
29 *
30 * Another possibility would have been to introduce a bunch of functions like
31 * initializeSpecialized(), removeMatrixColumnSpecialized(),
32 * iterationsConvergedSpecialized(), etc
33 * and call those function from the base class top down to the sub classes.
34 *
35 * The third possibility was to separate the approximation of the Jacobian from
36 * the common stuff like handling V,W matrices in the acceleration.
37 * Here, we have a class QNAcceleration that handles the V,W stuff an d the basic
38 * scheme of the QN update. Furthermore we have a base class (or rather interface)
39 * JacobianApproximation with sub classes IQNIMVJAPX and IQNAPX that handle all the
40 * specialized stuff like Jacobian approximation, handling of secondary data etc.
41 * However, this approach is not feasible, as we have to call the function
42 * removeMatrixColumn() down in the specialized sub classes IQNIMVJApx and IQNApx.
43 * This is not possible as the function works on the V, W matrices that are
44 * completely treated by QNAcceleration.
45 *
46 * ****************************************************************************
47 */
48
49// ----------------------------------------------------------- CLASS DEFINITION
50
51namespace precice {
52namespace io {
53class TXTReader;
54class TXTWriter;
55} // namespace io
56
57namespace acceleration {
58
64public:
66 double initialRelaxation,
67 bool forceInitialRelaxation,
68 int maxIterationsUsed,
69 int timeWindowsReused,
70 int filter,
71 double singularityLimit,
72 std::vector<int> dataIDs,
73 impl::PtrPreconditioner preconditioner);
74
79 {
80 // not necessary for user, only for developer, if needed, this should be configurable
81 // if (utils::IntraComm::isPrimary() || !utils::IntraComm::isParallel()) {
82 // _infostream.open("precice-accelerationInfo.log", std::ios_base::out);
83 // _infostream << std::setprecision(16);
84 // _infostream << _infostringstream.str();
85 // }
86 }
87
92 {
93 return _dataIDs;
94 }
95
99 virtual void initialize(const DataMap &cplData);
100
106 virtual void performAcceleration(DataMap &cplData);
107
114 virtual void iterationsConverged(const DataMap &cplData);
115
119 virtual void exportState(io::TXTWriter &writer);
120
126 virtual void importState(io::TXTReader &reader);
127
129 virtual int getDeletedColumns() const;
130
132 virtual int getDroppedColumns() const;
133
141 virtual int getLSSystemCols() const;
142
143protected:
144 logging::Logger _log{"acceleration::BaseQNAcceleration"};
145
148
151
154
157
160
163
165 bool _firstIteration = true;
166
167 /* @brief Indicates the first time window, where constant relaxation is used
168 * later, we replace the constant relaxation by a qN-update from last time window.
169 */
170 bool _firstTimeWindow = true;
171
172 /*
173 * @brief True if this process has nodes at the coupling interface
174 */
176
177 /* @brief If true, the QN-scheme always performs a underrelaxation in the first iteration of
178 * a new time window. Otherwise, the LS system from the previous time window is used in the
179 * first iteration.
180 */
182
186 bool _resetLS = false;
187
189 Eigen::VectorXd _oldXTilde;
190
192 Eigen::VectorXd _residuals;
193
196
198 Eigen::MatrixXd _matrixV;
199
201 Eigen::MatrixXd _matrixW;
202
205
210
218
227
232
236
237 int getLSSystemRows();
238
245 virtual void specializedIterationsConverged(const DataMap &cplData) = 0;
246
248 virtual void updateDifferenceMatrices(const DataMap &cplData);
249
251 virtual void splitCouplingData(const DataMap &cplData);
252
254 virtual void applyFilter();
255
257 virtual void computeUnderrelaxationSecondaryData(const DataMap &cplData) = 0;
258
260 virtual void computeQNUpdate(const DataMap &cplData, Eigen::VectorXd &xUpdate) = 0;
261
263 virtual void removeMatrixColumn(int columnIndex);
264
266 void writeInfo(const std::string &s, bool allProcs = false);
267
268 int its = 0, tWindows = 0;
269
270private:
272 Eigen::VectorXd _values;
273
275 Eigen::VectorXd _oldValues;
276
278 Eigen::VectorXd _oldResiduals;
279
284 Eigen::MatrixXd _matrixVBackup;
285 Eigen::MatrixXd _matrixWBackup;
287
289 int _nbDelCols = 0;
290
292 int _nbDropCols = 0;
293};
294} // namespace acceleration
295} // 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)
Eigen::VectorXd _oldResiduals
Difference between solver input and output from last time window.
int _filter
filter method that is used to maintain good conditioning of the least-squares system Either of two ty...
double _initialRelaxation
Constant relaxation factor used for first iteration.
int _timeWindowsReused
Maximum number of old time windows (with data values) kept.
virtual 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.
virtual ~BaseQNAcceleration()
Destructor, empty.
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.
virtual void iterationsConverged(const DataMap &cplData)
Marks a iteration sequence as converged.
Eigen::VectorXd _oldValues
Concatenation of all (old) coupling data involved in the QN system.
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.
std::ostringstream _infostringstream
write some debug/acceleration info to file
bool _resetLS
If true, the LS system has been modified (reset or recomputed) in such a way, that mere updating of m...
BaseQNAcceleration(double initialRelaxation, bool forceInitialRelaxation, int maxIterationsUsed, int timeWindowsReused, int filter, double singularityLimit, std::vector< int > dataIDs, impl::PtrPreconditioner preconditioner)
int _nbDropCols
Number of dropped columns in this time window (old time window out of scope)
virtual int getDroppedColumns() const
how many QN columns were dropped (went out of scope) in this time window
virtual void computeQNUpdate(const DataMap &cplData, Eigen::VectorXd &xUpdate)=0
Computes the quasi-Newton update using the specified pp scheme (IQNIMVJ, IQNILS)
virtual std::vector< int > getDataIDs() const
Returns all IQN involved data IDs.
Eigen::MatrixXd _matrixV
Stores residual deltas.
std::vector< int > _dimOffsets
Stores the local dimensions, i.e., the offsets in _invJacobian for all processors.
virtual void exportState(io::TXTWriter &writer)
Exports the current state of the acceleration to a file.
Eigen::VectorXd _residuals
Current iteration residuals of IQN data. Temporary.
virtual void importState(io::TXTReader &reader)
Imports the last exported state of the acceleration from file.
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 coupling data involved in the QN system.
virtual void splitCouplingData(const DataMap &cplData)
Splits up QN system vector back into the coupling data.
virtual void performAcceleration(DataMap &cplData)
Performs one acceleration step.
double _singularityLimit
Determines sensitivity when two matrix columns are considered equal.
void writeInfo(const std::string &s, bool allProcs=false)
Wwrites 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.
virtual int getDeletedColumns() const
how many QN columns were deleted in this time window
std::map< int, Eigen::VectorXd > _secondaryResiduals
Current iteration residuals of secondary data.
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.
int _maxIterationsUsed
Maximum number of old data iterations kept.
virtual void computeUnderrelaxationSecondaryData(const DataMap &cplData)=0
Computes underrelaxation for the secondary data.
std::vector< int > _secondaryDataIDs
Data IDs of data not involved in IQN coefficient computation.
std::vector< int > _dataIDs
Data IDs of data to be involved in the IQN algorithm.
virtual void initialize(const DataMap &cplData)
Initializes the acceleration.
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:16
File writer for matrix in Matlab V7 ASCII format.
Definition TXTWriter.hpp:14
This class provides a lightweight logger.
Definition Logger.hpp:16
Main namespace of the precice library.