preCICE v3.1.2
Loading...
Searching...
No Matches
Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes | List of all members
precice::acceleration::BaseQNAcceleration Class Referenceabstract

Base Class for quasi-Newton acceleration schemes. More...

#include <BaseQNAcceleration.hpp>

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

Public Member Functions

 BaseQNAcceleration (double initialRelaxation, bool forceInitialRelaxation, int maxIterationsUsed, int timeWindowsReused, int filter, double singularityLimit, std::vector< int > dataIDs, impl::PtrPreconditioner preconditioner)
 
virtual ~BaseQNAcceleration ()
 Destructor, empty.
 
virtual std::vector< int > getDataIDs () const
 Returns all IQN involved data IDs.
 
virtual void initialize (const DataMap &cplData)
 Initializes the acceleration.
 
virtual void performAcceleration (DataMap &cplData)
 Performs one acceleration step.
 
virtual void iterationsConverged (const DataMap &cplData)
 Marks a iteration sequence as converged.
 
virtual void exportState (io::TXTWriter &writer)
 Exports the current state of the acceleration to a file.
 
virtual void importState (io::TXTReader &reader)
 Imports the last exported state of the acceleration from file.
 
virtual int getDeletedColumns () const
 how many QN columns were deleted in this time window
 
virtual int getDroppedColumns () const
 how many QN columns were dropped (went out of scope) in this time window
 
virtual 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.
 
- Public Member Functions inherited from precice::acceleration::Acceleration
virtual ~Acceleration ()=default
 

Protected Member Functions

int getLSSystemRows ()
 
virtual void specializedIterationsConverged (const DataMap &cplData)=0
 Marks a iteration sequence as converged.
 
virtual void updateDifferenceMatrices (const DataMap &cplData)
 Updates the V, W matrices (as well as the matrices for the secondary data)
 
virtual void splitCouplingData (const DataMap &cplData)
 Splits up QN system vector back into the coupling data.
 
virtual void applyFilter ()
 Applies the filter method for the least-squares system, defined in the configuration.
 
virtual void computeUnderrelaxationSecondaryData (const DataMap &cplData)=0
 Computes underrelaxation for the secondary data.
 
virtual void computeQNUpdate (const DataMap &cplData, Eigen::VectorXd &xUpdate)=0
 Computes the quasi-Newton update using the specified pp scheme (IQNIMVJ, IQNILS)
 
virtual void removeMatrixColumn (int columnIndex)
 Removes one iteration from V,W matrices and adapts _matrixCols.
 
void writeInfo (const std::string &s, bool allProcs=false)
 Wwrites 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.
 
void concatenateCouplingData (const DataMap &cplData, const std::vector< DataID > &dataIDs, Eigen::VectorXd &targetValues, Eigen::VectorXd &targetOldValues) const
 Concatenates all coupling data involved into a single vector.
 

Protected Attributes

logging::Logger _log {"acceleration::BaseQNAcceleration"}
 
impl::PtrPreconditioner _preconditioner
 Preconditioner for least-squares system if vectorial system is used.
 
double _initialRelaxation
 Constant relaxation factor used for first iteration.
 
int _maxIterationsUsed
 Maximum number of old data iterations kept.
 
int _timeWindowsReused
 Maximum number of old time windows (with data values) kept.
 
std::vector< int > _dataIDs
 Data IDs of data to be involved in the IQN algorithm.
 
std::vector< int > _secondaryDataIDs
 Data IDs of data not involved in IQN coefficient computation.
 
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 _oldXTilde
 Solver output from last iteration.
 
Eigen::VectorXd _residuals
 Current iteration residuals of IQN data. Temporary.
 
std::map< int, Eigen::VectorXd > _secondaryResiduals
 Current iteration residuals of secondary data.
 
Eigen::MatrixXd _matrixV
 Stores residual deltas.
 
Eigen::MatrixXd _matrixW
 Stores x tilde deltas, where x tilde are values computed by solvers.
 
impl::QRFactorization _qrV
 Stores the current QR decomposition ov _matrixV, can be updated via deletion/insertion of columns.
 
int _filter
 filter method that is used to maintain good conditioning of the least-squares system Either of two types: QR1FILTER or QR2Filter
 
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::ostringstream _infostringstream
 write some debug/acceleration info to file
 
std::fstream _infostream
 
int its = 0
 
int tWindows = 0
 

Private Attributes

Eigen::VectorXd _values
 Concatenation of all coupling data involved in the QN system.
 
Eigen::VectorXd _oldValues
 Concatenation of all (old) coupling data involved in the QN system.
 
Eigen::VectorXd _oldResiduals
 Difference between solver input and output from last time window.
 
Eigen::MatrixXd _matrixVBackup
 backup of the V,W and matrixCols data structures. Needed for the skipping of initial relaxation, if previous time window converged within one iteration i.e., V and W are empty – in this case restore V and W with time window t-2.
 
Eigen::MatrixXd _matrixWBackup
 
std::deque< int > _matrixColsBackup
 
int _nbDelCols = 0
 Number of filtered out columns in this time window.
 
int _nbDropCols = 0
 Number of dropped columns in this time window (old time window out of scope)
 

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 Protected Member Functions inherited from precice::acceleration::Acceleration
static void applyRelaxation (double omega, DataMap &cplData)
 performs a relaxation given a relaxation factor omega
 

Detailed Description

Base Class for quasi-Newton acceleration schemes.

Definition at line 63 of file BaseQNAcceleration.hpp.

Constructor & Destructor Documentation

◆ BaseQNAcceleration()

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

Definition at line 33 of file BaseQNAcceleration.cpp.

◆ ~BaseQNAcceleration()

virtual precice::acceleration::BaseQNAcceleration::~BaseQNAcceleration ( )
inlinevirtual

Destructor, empty.

Definition at line 78 of file BaseQNAcceleration.hpp.

Member Function Documentation

◆ applyFilter()

void precice::acceleration::BaseQNAcceleration::applyFilter ( )
protectedvirtual

Applies the filter method for the least-squares system, defined in the configuration.

Definition at line 416 of file BaseQNAcceleration.cpp.

Here is the call graph for this function:

◆ computeQNUpdate()

virtual void precice::acceleration::BaseQNAcceleration::computeQNUpdate ( const DataMap & cplData,
Eigen::VectorXd & xUpdate )
protectedpure virtual

Computes the quasi-Newton update using the specified pp scheme (IQNIMVJ, IQNILS)

Implemented in precice::acceleration::IQNILSAcceleration, and precice::acceleration::IQNIMVJAcceleration.

◆ computeUnderrelaxationSecondaryData()

virtual void precice::acceleration::BaseQNAcceleration::computeUnderrelaxationSecondaryData ( const DataMap & cplData)
protectedpure virtual

Computes underrelaxation for the secondary data.

Implemented in precice::acceleration::IQNILSAcceleration, and precice::acceleration::IQNIMVJAcceleration.

◆ exportState()

void precice::acceleration::BaseQNAcceleration::exportState ( io::TXTWriter & writer)
virtual

Exports the current state of the acceleration to a file.

Reimplemented from precice::acceleration::Acceleration.

Definition at line 580 of file BaseQNAcceleration.cpp.

◆ getDataIDs()

virtual std::vector< int > precice::acceleration::BaseQNAcceleration::getDataIDs ( ) const
inlinevirtual

Returns all IQN involved data IDs.

Implements precice::acceleration::Acceleration.

Definition at line 91 of file BaseQNAcceleration.hpp.

◆ getDeletedColumns()

int precice::acceleration::BaseQNAcceleration::getDeletedColumns ( ) const
virtual

how many QN columns were deleted in this time window

Reimplemented from precice::acceleration::Acceleration.

Definition at line 590 of file BaseQNAcceleration.cpp.

◆ getDroppedColumns()

int precice::acceleration::BaseQNAcceleration::getDroppedColumns ( ) const
virtual

how many QN columns were dropped (went out of scope) in this time window

Reimplemented from precice::acceleration::Acceleration.

Definition at line 595 of file BaseQNAcceleration.cpp.

◆ getLSSystemCols()

int precice::acceleration::BaseQNAcceleration::getLSSystemCols ( ) const
virtual

: 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.

Reimplemented from precice::acceleration::Acceleration.

Definition at line 600 of file BaseQNAcceleration.cpp.

Here is the call graph for this function:

◆ getLSSystemRows()

int precice::acceleration::BaseQNAcceleration::getLSSystemRows ( )
protected

Definition at line 614 of file BaseQNAcceleration.cpp.

Here is the call graph for this function:

◆ importState()

void precice::acceleration::BaseQNAcceleration::importState ( io::TXTReader & reader)
virtual

Imports the last exported state of the acceleration from file.

Is empty at the moment!!!

Reimplemented from precice::acceleration::Acceleration.

Definition at line 585 of file BaseQNAcceleration.cpp.

◆ initialize()

void precice::acceleration::BaseQNAcceleration::initialize ( const DataMap & cplData)
virtual

Initializes the acceleration.


initialize()

@brief: Initializes all the needed variables and data

make dimensions public to all procs, last entry _dimOffsets[IntraComm::getSize()] holds the global dimension, global,n

provide vertex offset information for all processors mesh->getVertexOffsets() provides an array that stores the number of mesh vertices on each processor This information needs to be gathered for all meshes. To get the number of respective unknowns of a specific processor we need to multiply the number of vertices with the dimensionality of the vector-valued data for each coupling data.

Implements precice::acceleration::Acceleration.

Reimplemented in precice::acceleration::IQNILSAcceleration, and precice::acceleration::IQNIMVJAcceleration.

Definition at line 76 of file BaseQNAcceleration.cpp.

Here is the call graph for this function:

◆ iterationsConverged()

void precice::acceleration::BaseQNAcceleration::iterationsConverged ( const DataMap & cplData)
virtual

Marks a iteration sequence as converged.

: Is called when the convergence criterion for the coupling is fulfilled and finalizes the quasi Newton acceleration. Stores new differences in F and C, clears or

Since convergence measurements are done outside the acceleration, this method has to be used to signalize convergence to the acceleration.


iterationsConverged()

updates F and C according to the number of reused time windows

pending deletion (after first iteration of next time window Using the matrices from the old time window for the first iteration is better than doing underrelaxation as first iteration of every time window

Implements precice::acceleration::Acceleration.

Definition at line 462 of file BaseQNAcceleration.cpp.

Here is the call graph for this function:

◆ performAcceleration()

void precice::acceleration::BaseQNAcceleration::performAcceleration ( DataMap & cplData)
virtual

Performs one acceleration step.

Has to be called after every implicit coupling iteration.


performAcceleration()

@brief: performs one iteration of the quasi Newton acceleration.

update the difference matrices V,W includes: scaling of values computation of residuals appending the difference matrices

=== update and apply preconditioner ===

The preconditioner is only applied to the matrix V and the columns that are inserted into the QR-decomposition of V.

compute quasi-Newton update PRECONDITION: All objects are unscaled, except the matrices within the QR-dec of V. Thus, the pseudo inverse needs to be reverted before using it.

apply quasiNewton update

Implements precice::acceleration::Acceleration.

Definition at line 275 of file BaseQNAcceleration.cpp.

Here is the call graph for this function:

◆ removeMatrixColumn()

void precice::acceleration::BaseQNAcceleration::removeMatrixColumn ( int columnIndex)
protectedvirtual

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


removeMatrixColumn()

@brief: removes a column from the least squares system, i. e., from the matrices F and C

Reimplemented in precice::acceleration::IQNILSAcceleration, and precice::acceleration::IQNIMVJAcceleration.

Definition at line 552 of file BaseQNAcceleration.cpp.

Here is the call graph for this function:

◆ specializedIterationsConverged()

virtual void precice::acceleration::BaseQNAcceleration::specializedIterationsConverged ( const DataMap & cplData)
protectedpure virtual

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

Implemented in precice::acceleration::IQNILSAcceleration, and precice::acceleration::IQNIMVJAcceleration.

◆ splitCouplingData()

void precice::acceleration::BaseQNAcceleration::splitCouplingData ( const DataMap & cplData)
protectedvirtual

Splits up QN system vector back into the coupling data.

Definition at line 438 of file BaseQNAcceleration.cpp.

Here is the call graph for this function:

◆ updateDifferenceMatrices()

void precice::acceleration::BaseQNAcceleration::updateDifferenceMatrices ( const DataMap & cplData)
protectedvirtual

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

: computes the current residual and stores it, computes the differences and


updateDifferenceMatrices()

updates the difference matrices F and C.

Reimplemented in precice::acceleration::IQNILSAcceleration, and precice::acceleration::IQNIMVJAcceleration.

Definition at line 180 of file BaseQNAcceleration.cpp.

Here is the call graph for this function:

◆ writeInfo()

void precice::acceleration::BaseQNAcceleration::writeInfo ( const std::string & s,
bool allProcs = false )
protected

Wwrites info to the _infostream (also in parallel)

Definition at line 622 of file BaseQNAcceleration.cpp.

Here is the call graph for this function:

Member Data Documentation

◆ _dataIDs

std::vector<int> precice::acceleration::BaseQNAcceleration::_dataIDs
protected

Data IDs of data to be involved in the IQN algorithm.

Definition at line 159 of file BaseQNAcceleration.hpp.

◆ _dimOffsets

std::vector<int> precice::acceleration::BaseQNAcceleration::_dimOffsets
protected

Stores the local dimensions, i.e., the offsets in _invJacobian for all processors.

Definition at line 231 of file BaseQNAcceleration.hpp.

◆ _filter

int precice::acceleration::BaseQNAcceleration::_filter
protected

filter method that is used to maintain good conditioning of the least-squares system Either of two types: QR1FILTER or QR2Filter

Definition at line 209 of file BaseQNAcceleration.hpp.

◆ _firstIteration

bool precice::acceleration::BaseQNAcceleration::_firstIteration = true
protected

Indicates the first iteration, where constant relaxation is used.

Definition at line 165 of file BaseQNAcceleration.hpp.

◆ _firstTimeWindow

bool precice::acceleration::BaseQNAcceleration::_firstTimeWindow = true
protected

Definition at line 170 of file BaseQNAcceleration.hpp.

◆ _forceInitialRelaxation

bool precice::acceleration::BaseQNAcceleration::_forceInitialRelaxation
protected

Definition at line 181 of file BaseQNAcceleration.hpp.

◆ _hasNodesOnInterface

bool precice::acceleration::BaseQNAcceleration::_hasNodesOnInterface = true
protected

Definition at line 175 of file BaseQNAcceleration.hpp.

◆ _infostream

std::fstream precice::acceleration::BaseQNAcceleration::_infostream
protected

Definition at line 235 of file BaseQNAcceleration.hpp.

◆ _infostringstream

std::ostringstream precice::acceleration::BaseQNAcceleration::_infostringstream
protected

write some debug/acceleration info to file

Definition at line 234 of file BaseQNAcceleration.hpp.

◆ _initialRelaxation

double precice::acceleration::BaseQNAcceleration::_initialRelaxation
protected

Constant relaxation factor used for first iteration.

Definition at line 150 of file BaseQNAcceleration.hpp.

◆ _log

logging::Logger precice::acceleration::BaseQNAcceleration::_log {"acceleration::BaseQNAcceleration"}
protected

Definition at line 144 of file BaseQNAcceleration.hpp.

◆ _matrixCols

std::deque<int> precice::acceleration::BaseQNAcceleration::_matrixCols
protected

Indices (of columns in W, V matrices) of 1st iterations of time windows.

When old time windows are reused (_timeWindowsReused > 0), the indices of the first iteration of each time window needs to be stored, such that, e.g., all iterations of the last time window, or one specific iteration that leads to a singular matrix in the QR decomposition can be removed and tracked.

Definition at line 226 of file BaseQNAcceleration.hpp.

◆ _matrixColsBackup

std::deque<int> precice::acceleration::BaseQNAcceleration::_matrixColsBackup
private

Definition at line 286 of file BaseQNAcceleration.hpp.

◆ _matrixV

Eigen::MatrixXd precice::acceleration::BaseQNAcceleration::_matrixV
protected

Stores residual deltas.

Definition at line 198 of file BaseQNAcceleration.hpp.

◆ _matrixVBackup

Eigen::MatrixXd precice::acceleration::BaseQNAcceleration::_matrixVBackup
private

backup of the V,W and matrixCols data structures. Needed for the skipping of initial relaxation, if previous time window converged within one iteration i.e., V and W are empty – in this case restore V and W with time window t-2.

Definition at line 284 of file BaseQNAcceleration.hpp.

◆ _matrixW

Eigen::MatrixXd precice::acceleration::BaseQNAcceleration::_matrixW
protected

Stores x tilde deltas, where x tilde are values computed by solvers.

Definition at line 201 of file BaseQNAcceleration.hpp.

◆ _matrixWBackup

Eigen::MatrixXd precice::acceleration::BaseQNAcceleration::_matrixWBackup
private

Definition at line 285 of file BaseQNAcceleration.hpp.

◆ _maxIterationsUsed

int precice::acceleration::BaseQNAcceleration::_maxIterationsUsed
protected

Maximum number of old data iterations kept.

Definition at line 153 of file BaseQNAcceleration.hpp.

◆ _nbDelCols

int precice::acceleration::BaseQNAcceleration::_nbDelCols = 0
private

Number of filtered out columns in this time window.

Definition at line 289 of file BaseQNAcceleration.hpp.

◆ _nbDropCols

int precice::acceleration::BaseQNAcceleration::_nbDropCols = 0
private

Number of dropped columns in this time window (old time window out of scope)

Definition at line 292 of file BaseQNAcceleration.hpp.

◆ _oldResiduals

Eigen::VectorXd precice::acceleration::BaseQNAcceleration::_oldResiduals
private

Difference between solver input and output from last time window.

Definition at line 278 of file BaseQNAcceleration.hpp.

◆ _oldValues

Eigen::VectorXd precice::acceleration::BaseQNAcceleration::_oldValues
private

Concatenation of all (old) coupling data involved in the QN system.

Definition at line 275 of file BaseQNAcceleration.hpp.

◆ _oldXTilde

Eigen::VectorXd precice::acceleration::BaseQNAcceleration::_oldXTilde
protected

Solver output from last iteration.

Definition at line 189 of file BaseQNAcceleration.hpp.

◆ _preconditioner

impl::PtrPreconditioner precice::acceleration::BaseQNAcceleration::_preconditioner
protected

Preconditioner for least-squares system if vectorial system is used.

Definition at line 147 of file BaseQNAcceleration.hpp.

◆ _qrV

impl::QRFactorization precice::acceleration::BaseQNAcceleration::_qrV
protected

Stores the current QR decomposition ov _matrixV, can be updated via deletion/insertion of columns.

Definition at line 204 of file BaseQNAcceleration.hpp.

◆ _resetLS

bool precice::acceleration::BaseQNAcceleration::_resetLS = false
protected

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.

Definition at line 186 of file BaseQNAcceleration.hpp.

◆ _residuals

Eigen::VectorXd precice::acceleration::BaseQNAcceleration::_residuals
protected

Current iteration residuals of IQN data. Temporary.

Definition at line 192 of file BaseQNAcceleration.hpp.

◆ _secondaryDataIDs

std::vector<int> precice::acceleration::BaseQNAcceleration::_secondaryDataIDs
protected

Data IDs of data not involved in IQN coefficient computation.

Definition at line 162 of file BaseQNAcceleration.hpp.

◆ _secondaryResiduals

std::map<int, Eigen::VectorXd> precice::acceleration::BaseQNAcceleration::_secondaryResiduals
protected

Current iteration residuals of secondary data.

Definition at line 195 of file BaseQNAcceleration.hpp.

◆ _singularityLimit

double precice::acceleration::BaseQNAcceleration::_singularityLimit
protected

Determines sensitivity when two matrix columns are considered equal.

When during the QR decomposition of the V matrix a pivot element smaller than the singularity limit is found, the matrix is considered to be singular and the corresponding (older) iteration is removed.

Definition at line 217 of file BaseQNAcceleration.hpp.

◆ _timeWindowsReused

int precice::acceleration::BaseQNAcceleration::_timeWindowsReused
protected

Maximum number of old time windows (with data values) kept.

Definition at line 156 of file BaseQNAcceleration.hpp.

◆ _values

Eigen::VectorXd precice::acceleration::BaseQNAcceleration::_values
private

Concatenation of all coupling data involved in the QN system.

Definition at line 272 of file BaseQNAcceleration.hpp.

◆ its

int precice::acceleration::BaseQNAcceleration::its = 0
protected

Definition at line 268 of file BaseQNAcceleration.hpp.

◆ tWindows

int precice::acceleration::BaseQNAcceleration::tWindows = 0
protected

Definition at line 268 of file BaseQNAcceleration.hpp.


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