3#include <boost/range/adaptor/map.hpp>
28namespace acceleration {
35 double initialRelaxation,
36 bool forceInitialRelaxation,
37 int maxIterationsUsed,
38 int timeWindowsReused,
40 double singularityLimit,
44 : _preconditioner(
std::move(preconditioner)),
45 _initialRelaxation(initialRelaxation),
46 _maxIterationsUsed(maxIterationsUsed),
47 _timeWindowsReused(timeWindowsReused),
48 _primaryDataIDs(
std::move(dataIDs)),
49 _forceInitialRelaxation(forceInitialRelaxation),
50 _reducedTimeGrid(reducedTimeGrid),
53 _singularityLimit(singularityLimit),
54 _infostringstream(
std::ostringstream::ate)
57 "Initial relaxation factor for QN acceleration has to "
58 "be larger than zero and smaller or equal than one. "
59 "Current initial relaxation is {}",
62 "Maximum number of iterations used in the quasi-Newton acceleration "
63 "scheme has to be larger than zero. "
64 "Current maximum reused iterations is {}",
67 "Number of previous time windows to be reused for "
68 "quasi-Newton acceleration has to be larger than or equal to zero. "
69 "Current number of time windows reused is {}",
84 for (
const DataMap::value_type &pair : cplData) {
85 PRECICE_ASSERT(pair.second->getSize() == pair.second->getPreviousIterationSize(),
"current and previousIteration have to be initialized and of identical size.",
86 pair.second->getSize(), pair.second->getPreviousIterationSize());
91 "Gradient data, which is required by at least one of the configured data mappings, is not yet compatible with quasi-Newton acceleration. This combination might lead to numerical issues. "
92 "Consider switching to a different acceleration scheme or a different data mapping scheme.");
97 for (
const DataMap::value_type &pair : cplData) {
126 "The coupling residual equals almost zero. There is maybe something wrong in your adapter. "
127 "Maybe you always write the same data or you call advance without "
128 "providing new data first or you do not use available read data. "
129 "Or you just converge much further than actually necessary.");
144 "The number of columns in the least squares system exceeded half the number of unknowns at the interface. "
145 "The system will probably become bad or ill-conditioned and the quasi-Newton acceleration may not "
146 "converge. Maybe the number of allowed columns (\"max-used-iterations\") should be limited.");
151 Eigen::VectorXd deltaXTilde =
_values;
161 "Adding a vector with a two-norm of {} to the quasi-Newton V matrix, which will lead to "
162 "ill-conditioning. A filter might delete the column again. Still, this could mean that you are "
163 "converging too tightly, that you reached steady-state, or that you are giving by mistake identical "
164 "data to preCICE in two consecutive iterations.",
169 if (not columnLimitReached && overdetermined) {
233 _timeGrids.value().moveTimeGridToNewWindow(cplData);
266 PRECICE_DEBUG(
" Last time window converged after one iteration. Need to restore the matrices from backup.");
307 applyingFilter.
stop();
317 Eigen::VectorXd xUpdate = Eigen::VectorXd::Zero(
_values.size());
350 "The quasi-Newton update contains NaN values. This means that the quasi-Newton acceleration failed to converge. "
351 "When writing your own adapter this could indicate that you give wrong information to preCICE, such as identical "
352 "data in succeeding iterations. Or you do not properly save and reload checkpoints. "
353 "If you give the correct data this could also mean that the coupled problem is too hard to solve. Try to use a QR "
354 "filter or increase its threshold (larger epsilon).");
376 for (
int i = delIndices.
size() - 1; i >= 0; i--) {
380 PRECICE_DEBUG(
" Filter: removing column with index {} in iteration {} of time window: {}", delIndices[i],
its,
tWindows);
387 const DataMap &cplData,
double windowStart)
392 Eigen::Index offset = 0;
396 auto &couplingData = *cplData.
at(
id);
397 size_t dataSize = couplingData.getSize();
399 Eigen::VectorXd timeGrid =
_timeGrids->getTimeGridAfter(
id, windowStart);
400 couplingData.timeStepsStorage().trimAfter(windowStart);
401 for (
int i = 0; i < timeGrid.size(); i++) {
403 Eigen::VectorXd temp = Eigen::VectorXd::Zero(dataSize);
404 for (
size_t j = 0; j < dataSize; j++) {
409 couplingData.setSampleAtTime(timeGrid(i),
time::Sample(couplingData.getDimensions(), temp));
423 const DataMap &cplData,
double windowStart)
460 stream <<
"Matrix column counters: ";
462 stream << cols <<
", ";
507 for (
int i = 0; i < toRemove; i++) {
542 if (cols > columnIndex) {
635 Eigen::Index offset = 0;
636 for (
int id : dataIDs) {
637 Eigen::Index dataSize = cplData.
at(
id)->getSize();
640 for (
int i = 0; i < timeGrid.size(); i++) {
642 auto current = cplData.
at(
id)->timeStepsStorage().sample(timeGrid(i));
643 auto old = cplData.
at(
id)->getPreviousValuesAtTime(timeGrid(i));
645 PRECICE_ASSERT(data.size() >= offset + dataSize,
"the values were not initialized correctly");
646 PRECICE_ASSERT(oldData.size() >= offset + dataSize,
"the values were not initialized correctly");
648 for (Eigen::Index i = 0; i < dataSize; i++) {
649 data(i + offset) = current(i);
650 oldData(i + offset) = old(i);
664 auto addTimeSliceSize = [&](
size_t sum,
int id,
precice::time::TimeGrids timeGrids) {
return sum + timeGrids.getTimeGridAfter(
id, windowStart).size() * cplData.
at(
id)->getSize(); };
672 _values = Eigen::VectorXd::Zero(dataSize);
695 if (primaryDataSize <= 0) {
709 int accumulatedNumberOfUnknowns = 0;
710 int accumulatedNumberOfPrimaryUnknowns = 0;
713 const auto &offsets = cplData.
at(elem)->getVertexOffsets();
715 accumulatedNumberOfUnknowns += offsets[i] * cplData.
at(elem)->getDimensions() *
_timeGrids.value().getTimeGridAfter(elem, windowStart).size();
718 accumulatedNumberOfPrimaryUnknowns += offsets[i] * cplData.
at(elem)->getDimensions() *
_primaryTimeGrids.value().getTimeGridAfter(elem, windowStart).size();
734 PRECICE_ASSERT(primaryDataSize == primaryUnknowns, primaryDataSize, primaryUnknowns);
#define PRECICE_WARN_IF(condition,...)
#define PRECICE_DEBUG(...)
#define PRECICE_TRACE(...)
#define PRECICE_CHECK(check,...)
#define PRECICE_ASSERT(...)
T back_inserter(T... args)
static void applyRelaxation(double omega, DataMap &cplData, double windowStart)
performs a relaxation given a relaxation factor omega
void checkDataIDs(const DataMap &cplData) const
Checks if all dataIDs are contained in cplData.
static const int QR3FILTER
static const int QR2FILTER
static const int NOFILTER
virtual void updateDifferenceMatrices(const DataMap &cplData)
Updates the V, W matrices (as well as the matrices for the secondary data)
std::deque< int > _matrixColsBackup
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)
Eigen::MatrixXd _matrixWBackup
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.
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.
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...
bool _hasNodesOnInterface
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.
int getPrimaryLSSystemRows() const
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.
int getLSSystemRows() const
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.
bool _forceInitialRelaxation
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.
void setGlobalRows(int gr)
void applyFilter(double singularityLimit, std::vector< int > &delIndices, Eigen::MatrixXd &V)
filters the least squares system, i.e., the decomposition Q*R = V according to the defined filter tec...
void reset()
resets the QR factorization zo zero Q(0:0, 0:0)R(0:0, 0:0)
void popBack()
deletes the column at position _cols-1, i.e., deletes the last column and updates the QR factorizatio...
void pushFront(const Eigen::VectorXd &v)
inserts a new column at position 0, i.e., shifts right and inserts at first position and updates the ...
File reader for matrix/vector in Matlab V7 ASCII format.
File writer for matrix in Matlab V7 ASCII format.
void stop()
Stops a running event.
Interface for storing the time grids in the Quasi-Newton and Aitken methods. A time grid is a ordered...
Eigen::VectorXd getTimeGridAfter(int dataID, double time) const
static int getSize()
Number of ranks. This includes ranks from both participants, e.g. minimal size is 2.
static double l2norm(const Eigen::VectorXd &vec)
The l2 norm of a vector is calculated on distributed data.
static Rank getRank()
Current rank.
static bool isPrimary()
True if this process is running the primary rank.
static bool isParallel()
True if this process is running in parallel.
static com::PtrCommunication & getCommunication()
Intra-participant communication.
constexpr bool equals(const Eigen::MatrixBase< DerivedA > &A, const Eigen::MatrixBase< DerivedB > &B, double tolerance=NUMERICAL_ZERO_DIFFERENCE)
Compares two Eigen::MatrixBase for equality up to tolerance.
static constexpr SynchronizeTag Synchronize
Convenience instance of the SynchronizeTag.
void removeColumnFromMatrix(Eigen::MatrixXd &A, int col)
bool contained(const ELEMENT_T &element, const std::vector< ELEMENT_T > &vec)
Returns true, if given element is in vector, otherwise false.
void appendFront(Eigen::MatrixXd &A, Eigen::VectorXd &v)
void shiftSetFirst(Eigen::MatrixXd &A, const Eigen::VectorXd &v)
Main namespace of the precice library.