preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Private Member Functions | Private Attributes | List of all members
precice::acceleration::IQNILSAcceleration Class Reference

Interface quasi-Newton with interface least-squares approximation. More...

#include <IQNILSAcceleration.hpp>

Inheritance diagram for precice::acceleration::IQNILSAcceleration:
[legend]
Collaboration diagram for precice::acceleration::IQNILSAcceleration:
[legend]

Public Member Functions

 IQNILSAcceleration (double initialRelaxation, bool forceInitialRelaxation, int maxIterationsUsed, int pastTimeWindowsReused, int filter, double singularityLimit, std::vector< int > dataIDs, impl::PtrPreconditioner preconditioner, bool reducedTimeGrid)
 
 ~IQNILSAcceleration () override=default
 
void specializedIterationsConverged (const DataMap &cplData) override
 Marks a iteration sequence as converged.
 
- Public Member Functions inherited from precice::acceleration::BaseQNAcceleration
 BaseQNAcceleration (double initialRelaxation, bool forceInitialRelaxation, int maxIterationsUsed, int timeWindowsReused, int filter, double singularityLimit, std::vector< int > dataIDs, impl::PtrPreconditioner preconditioner, bool reducedTimeGrid)
 
 ~BaseQNAcceleration () override
 Destructor, empty.
 
std::vector< int > getPrimaryDataIDs () const final override
 Returns all IQN involved data IDs.
 
void initialize (const DataMap &cplData) final override
 Initializes the acceleration.
 
void performAcceleration (DataMap &cplData, double windowStart, double windowEnd) final override
 Performs one acceleration step.
 
void iterationsConverged (const DataMap &cplData, double windowStart) final override
 Marks a iteration sequence as converged.
 
void exportState (io::TXTWriter &writer) final override
 Exports the current state of the acceleration to a file.
 
void importState (io::TXTReader &reader) final override
 Imports the last exported state of the acceleration from file.
 
int getDeletedColumns () const
 how many QN columns were deleted in this time window
 
int getDroppedColumns () const
 how many QN columns were dropped (went out of scope) in this time window
 
int getLSSystemCols () const
 : computes number of cols in least squares system, i.e, number of cols in _matrixV, _matrixW, _qrV, etc.. This is only necessary if some procs do not have any nodes on the coupling interface. In this case, the matrices are not constructed and we have no information about the number of cols. This info is needed for intra-participant communication. Number of its =! _cols in general.
 
int getMaxUsedIterations () const
 Get the maximum number of reused iterations.
 
int getMaxUsedTimeWindows () const
 Get the maximum number of reused time windows.
 
- Public Member Functions inherited from precice::acceleration::Acceleration
virtual ~Acceleration ()=default
 

Private Member Functions

void updateDifferenceMatrices (const DataMap &cplData) override
 updates the V, W matrices (as well as the matrices for the secondary data)
 
void computeQNUpdate (Eigen::VectorXd &xUpdate) override
 computes the IQN-ILS update using QR decomposition
 
void removeMatrixColumn (int columnIndex) override
 Removes one iteration from V,W matrices and adapts _matrixCols.
 
void specializedInitializeVectorsAndPreconditioner (const DataMap &cplData) final override
 

Private Attributes

std::map< int, Eigen::VectorXd > _secondaryOldXTildes
 Secondary data solver output from last iteration.
 

Additional Inherited Members

- Public Types inherited from precice::acceleration::Acceleration
using DataMap = std::map<int, cplscheme::PtrCouplingData>
 Map from data ID to data values.
 
- Static Public Attributes inherited from precice::acceleration::Acceleration
static const int NOFILTER = 0
 
static const int QR1FILTER = 1
 
static const int QR1FILTER_ABS = 2
 
static const int QR2FILTER = 3
 
static const int PODFILTER = 4
 
static const int QR3FILTER = 5
 
- Protected Member Functions inherited from precice::acceleration::BaseQNAcceleration
int getLSSystemRows () const
 
int getPrimaryLSSystemRows () const
 
virtual void updateCouplingData (const DataMap &cplData, double windowStart)
 Splits up QN system vector back into the waveforms in coupling data.
 
virtual void applyFilter ()
 Applies the filter method for the least-squares system, defined in the configuration.
 
void writeInfo (const std::string &s, bool allProcs=false)
 Writes info to the _infostream (also in parallel)
 
- Protected Member Functions inherited from precice::acceleration::Acceleration
void checkDataIDs (const DataMap &cplData) const
 Checks if all dataIDs are contained in cplData.
 
- Static Protected Member Functions inherited from precice::acceleration::Acceleration
static void applyRelaxation (double omega, DataMap &cplData, double windowStart)
 performs a relaxation given a relaxation factor omega
 
- Protected Attributes inherited from precice::acceleration::BaseQNAcceleration
logging::Logger _log {"acceleration::BaseQNAcceleration"}
 
impl::PtrPreconditioner _preconditioner
 Preconditioner for least-squares system if vectorial system is used.
 
const double _initialRelaxation
 Constant relaxation factor used for first iteration.
 
const int _maxIterationsUsed
 Maximum number of old data iterations kept.
 
const int _timeWindowsReused
 Maximum number of old time windows (with data values) kept.
 
const std::vector< int > _primaryDataIDs
 Data IDs of primary data to be involved in the IQN coefficient computation.
 
std::vector< int > _dataIDs
 Data IDs of all coupling data.
 
bool _firstIteration = true
 Indicates the first iteration, where constant relaxation is used.
 
bool _firstTimeWindow = true
 
bool _hasNodesOnInterface = true
 
bool _forceInitialRelaxation
 
bool _resetLS = false
 If true, the LS system has been modified (reset or recomputed) in such a way, that mere updating of matrices _Wtil, Q, R etc.. is not feasible any more and need to be recomputed.
 
Eigen::VectorXd _oldPrimaryXTilde
 Solver output regarding primary data from last iteration.
 
Eigen::VectorXd _oldXTilde
 Solver output from last iteration.
 
Eigen::VectorXd _primaryResiduals
 Current iteration residuals of primary IQN data. Temporary.
 
Eigen::VectorXd _residuals
 Current iteration residuals of IQN data. Temporary.
 
Eigen::MatrixXd _matrixV
 Stores residual deltas.
 
Eigen::MatrixXd _matrixW
 Stores x tilde deltas, where x tilde are values computed by solvers.
 
const bool _reducedTimeGrid
 if _reducedTimeGrid = false uses the full QN-WI and if _reducedTimeGrid = true uses rQN-WI form the paper https://onlinelibrary.wiley.com/doi/10.1002/nme.6443
 
impl::QRFactorization _qrV
 Stores the current QR decomposition ov _matrixV, can be updated via deletion/insertion of columns.
 
const int _filter
 filter method that is used to maintain good conditioning of the least-squares system Either of two types: QR1FILTER or QR2Filter
 
const double _singularityLimit
 Determines sensitivity when two matrix columns are considered equal.
 
std::deque< int > _matrixCols
 Indices (of columns in W, V matrices) of 1st iterations of time windows.
 
std::vector< int > _dimOffsets
 Stores the local dimensions, i.e., the offsets in _invJacobian for all processors.
 
std::vector< int > _dimOffsetsPrimary
 Stores the local dimensions regarding primary data,.
 
std::ostringstream _infostringstream
 write some debug/acceleration info to file
 
std::fstream _infostream
 
int its = 0
 
int tWindows = 0
 

Detailed Description

Interface quasi-Newton with interface least-squares approximation.

Performs a quasi-Newton to accelerate the convergence of implicit coupling iterations. A least-squares approximation is used to find the best linear combination of old data values. After every coupling iteration, the data values used are enhanced by the new coupling iterates.

If more coupling data is present than used to compute the IQN acceleration, this data is relaxed using the same linear combination as computed for the IQN-related data. The data is called "secondary" henceforth and additional old value and data matrices are needed for it.

Definition at line 25 of file IQNILSAcceleration.hpp.

Constructor & Destructor Documentation

◆ IQNILSAcceleration()

precice::acceleration::IQNILSAcceleration::IQNILSAcceleration ( double initialRelaxation,
bool forceInitialRelaxation,
int maxIterationsUsed,
int pastTimeWindowsReused,
int filter,
double singularityLimit,
std::vector< int > dataIDs,
impl::PtrPreconditioner preconditioner,
bool reducedTimeGrid )

Definition at line 27 of file IQNILSAcceleration.cpp.

◆ ~IQNILSAcceleration()

precice::acceleration::IQNILSAcceleration::~IQNILSAcceleration ( )
overridedefault

Member Function Documentation

◆ computeQNUpdate()

void precice::acceleration::IQNILSAcceleration::computeQNUpdate ( Eigen::VectorXd & xUpdate)
overrideprivatevirtual

computes the IQN-ILS update using QR decomposition

Implements precice::acceleration::BaseQNAcceleration.

Definition at line 49 of file IQNILSAcceleration.cpp.

Here is the call graph for this function:

◆ removeMatrixColumn()

void precice::acceleration::IQNILSAcceleration::removeMatrixColumn ( int columnIndex)
overrideprivatevirtual

Removes one iteration from V,W matrices and adapts _matrixCols.

Reimplemented from precice::acceleration::BaseQNAcceleration.

Definition at line 130 of file IQNILSAcceleration.cpp.

Here is the call graph for this function:

◆ specializedInitializeVectorsAndPreconditioner()

void precice::acceleration::IQNILSAcceleration::specializedInitializeVectorsAndPreconditioner ( const DataMap & cplData)
inlinefinaloverrideprivatevirtual

Implements precice::acceleration::BaseQNAcceleration.

Definition at line 62 of file IQNILSAcceleration.hpp.

◆ specializedIterationsConverged()

void precice::acceleration::IQNILSAcceleration::specializedIterationsConverged ( const DataMap & cplData)
overridevirtual

Marks a iteration sequence as converged.

called by the iterationsConverged() method in the BaseQNAcceleration class handles the acceleration specific action after the convergence of one iteration

Implements precice::acceleration::BaseQNAcceleration.

Definition at line 117 of file IQNILSAcceleration.cpp.

Here is the call graph for this function:

◆ updateDifferenceMatrices()

void precice::acceleration::IQNILSAcceleration::updateDifferenceMatrices ( const DataMap & cplData)
overrideprivatevirtual

updates the V, W matrices (as well as the matrices for the secondary data)

Reimplemented from precice::acceleration::BaseQNAcceleration.

Definition at line 42 of file IQNILSAcceleration.cpp.

Here is the call graph for this function:

Member Data Documentation

◆ _secondaryOldXTildes

std::map<int, Eigen::VectorXd> precice::acceleration::IQNILSAcceleration::_secondaryOldXTildes
private

Secondary data solver output from last iteration.

Definition at line 50 of file IQNILSAcceleration.hpp.


The documentation for this class was generated from the following files: