preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | 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, 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
 

Protected Member Functions

int getLSSystemRows () const
 
int getPrimaryLSSystemRows () const
 
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 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.
 
virtual void computeQNUpdate (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)
 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.
 

Protected Attributes

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
 

Private Member Functions

void initializeVectorsAndPreconditioner (const DataMap &cplData, double windowStart)
 Initializes the vectors, matrices and preconditioner This has to be done after the first iteration of the first time window, since everything in the QN-algorithm is sampled to the timegrid of the first waveform.
 
virtual void specializedInitializeVectorsAndPreconditioner (const DataMap &cplData)=0
 handles the initialization of matrices and vectors in the sub-classes
 
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.
 

Private Attributes

std::optional< time::TimeGrids_timeGrids
 Stores the time grids to which the primary and secondary data involved in the QN system will be interpolated to.
 
std::optional< time::TimeGrids_primaryTimeGrids
 Stores the time grids to which the primary data involved in the QN system will be interpolated to. If _reducedTimeGrids is true then this will only contain the last time stamp of the time window, see https://doi.org/10.1002/nme.6443.
 
Eigen::VectorXd _primaryValues
 Concatenation of all primary data involved in the QN system.
 
Eigen::VectorXd _oldPrimaryValues
 Concatenation of all old primary data involved in the QN system.
 
Eigen::VectorXd _oldPrimaryResiduals
 Difference between solver input and output from last time window regarding primary data.
 
Eigen::VectorXd _values
 Concatenation of all primary and secondary data involved in the QN system.
 
Eigen::VectorXd _oldValues
 Concatenation of all old primary and secondary 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, 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 const int QR3FILTER = 5
 
- 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
 

Detailed Description

Base Class for quasi-Newton acceleration schemes.

Definition at line 64 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,
bool reducedTimeGrid )

Definition at line 34 of file BaseQNAcceleration.cpp.

◆ ~BaseQNAcceleration()

precice::acceleration::BaseQNAcceleration::~BaseQNAcceleration ( )
inlineoverride

Destructor, empty.

Definition at line 80 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 364 of file BaseQNAcceleration.cpp.

Here is the call graph for this function:

◆ computeQNUpdate()

virtual void precice::acceleration::BaseQNAcceleration::computeQNUpdate ( 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.

◆ concatenateCouplingData()

void precice::acceleration::BaseQNAcceleration::concatenateCouplingData ( Eigen::VectorXd & data,
Eigen::VectorXd & oldData,
const DataMap & cplData,
std::vector< int > dataIDs,
precice::time::TimeGrids timeGrids,
double windowStart ) const
private

Samples and concatenates the data and old data in cplData into a long vector.

Definition at line 633 of file BaseQNAcceleration.cpp.

Here is the call graph for this function:

◆ exportState()

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

Exports the current state of the acceleration to a file.

Reimplemented from precice::acceleration::Acceleration.

Definition at line 554 of file BaseQNAcceleration.cpp.

◆ getDeletedColumns()

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

how many QN columns were deleted in this time window

Definition at line 564 of file BaseQNAcceleration.cpp.

◆ getDroppedColumns()

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

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

Definition at line 569 of file BaseQNAcceleration.cpp.

◆ getLSSystemCols()

int precice::acceleration::BaseQNAcceleration::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.

Definition at line 574 of file BaseQNAcceleration.cpp.

Here is the call graph for this function:

◆ getLSSystemRows()

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

Definition at line 598 of file BaseQNAcceleration.cpp.

Here is the call graph for this function:

◆ getMaxUsedIterations()

int precice::acceleration::BaseQNAcceleration::getMaxUsedIterations ( ) const

Get the maximum number of reused iterations.

Definition at line 588 of file BaseQNAcceleration.cpp.

◆ getMaxUsedTimeWindows()

int precice::acceleration::BaseQNAcceleration::getMaxUsedTimeWindows ( ) const

Get the maximum number of reused time windows.

Definition at line 593 of file BaseQNAcceleration.cpp.

◆ getPrimaryDataIDs()

std::vector< int > precice::acceleration::BaseQNAcceleration::getPrimaryDataIDs ( ) const
inlinefinaloverridevirtual

Returns all IQN involved data IDs.

Implements precice::acceleration::Acceleration.

Definition at line 93 of file BaseQNAcceleration.hpp.

◆ getPrimaryLSSystemRows()

int precice::acceleration::BaseQNAcceleration::getPrimaryLSSystemRows ( ) const
protected

Definition at line 606 of file BaseQNAcceleration.cpp.

Here is the call graph for this function:

◆ importState()

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

Imports the last exported state of the acceleration from file.

Is empty at the moment!!!

Reimplemented from precice::acceleration::Acceleration.

Definition at line 559 of file BaseQNAcceleration.cpp.

◆ initialize()

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

Initializes the acceleration.


initialize()

@brief: Initializes all the needed variables and data

Implements precice::acceleration::Acceleration.

Definition at line 79 of file BaseQNAcceleration.cpp.

Here is the call graph for this function:

◆ initializeVectorsAndPreconditioner()

void precice::acceleration::BaseQNAcceleration::initializeVectorsAndPreconditioner ( const DataMap & cplData,
double windowStart )
private

Initializes the vectors, matrices and preconditioner This has to be done after the first iteration of the first time window, since everything in the QN-algorithm is sampled to the timegrid of the first waveform.

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.

Definition at line 657 of file BaseQNAcceleration.cpp.

Here is the call graph for this function:

◆ iterationsConverged()

void precice::acceleration::BaseQNAcceleration::iterationsConverged ( const DataMap & cplData,
double windowStart )
finaloverridevirtual

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

Sample all the data to the corresponding time grid in _timeGrids and concatenate everything into a long vector timeGrids are stored using std::optional, thus the .value() to get the actual object

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 422 of file BaseQNAcceleration.cpp.

Here is the call graph for this function:

◆ performAcceleration()

void precice::acceleration::BaseQNAcceleration::performAcceleration ( DataMap & cplData,
double windowStart,
double windowEnd )
finaloverridevirtual

Performs one acceleration step.

Has to be called after every implicit coupling iteration.


performAcceleration()

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

Sample all the data to the corresponding time grid in _timeGrids and concatenate everything into a long vector timeGrids are stored using std::optional, thus the .value() to get the actual object

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.

Implements precice::acceleration::Acceleration.

Definition at line 212 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 526 of file BaseQNAcceleration.cpp.

Here is the call graph for this function:

◆ specializedInitializeVectorsAndPreconditioner()

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

handles the initialization of matrices and vectors in the sub-classes

called by the initializeVectorsAndPreconditioner method in the BaseQNAcceleration class

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

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

◆ updateCouplingData()

void precice::acceleration::BaseQNAcceleration::updateCouplingData ( const DataMap & cplData,
double windowStart )
protectedvirtual

Splits up QN system vector back into the waveforms in coupling data.

Definition at line 386 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 114 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

Writes info to the _infostream (also in parallel)

Definition at line 614 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 all coupling data.

Definition at line 174 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 249 of file BaseQNAcceleration.hpp.

◆ _dimOffsetsPrimary

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

Stores the local dimensions regarding primary data,.

Definition at line 253 of file BaseQNAcceleration.hpp.

◆ _filter

const 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 227 of file BaseQNAcceleration.hpp.

◆ _firstIteration

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

Indicates the first iteration, where constant relaxation is used.

Definition at line 177 of file BaseQNAcceleration.hpp.

◆ _firstTimeWindow

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

Definition at line 182 of file BaseQNAcceleration.hpp.

◆ _forceInitialRelaxation

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

Definition at line 193 of file BaseQNAcceleration.hpp.

◆ _hasNodesOnInterface

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

Definition at line 187 of file BaseQNAcceleration.hpp.

◆ _infostream

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

Definition at line 257 of file BaseQNAcceleration.hpp.

◆ _infostringstream

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

write some debug/acceleration info to file

Definition at line 256 of file BaseQNAcceleration.hpp.

◆ _initialRelaxation

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

Constant relaxation factor used for first iteration.

Definition at line 162 of file BaseQNAcceleration.hpp.

◆ _log

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

Definition at line 156 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 244 of file BaseQNAcceleration.hpp.

◆ _matrixColsBackup

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

Definition at line 333 of file BaseQNAcceleration.hpp.

◆ _matrixV

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

Stores residual deltas.

Definition at line 213 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 331 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 216 of file BaseQNAcceleration.hpp.

◆ _matrixWBackup

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

Definition at line 332 of file BaseQNAcceleration.hpp.

◆ _maxIterationsUsed

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

Maximum number of old data iterations kept.

Definition at line 165 of file BaseQNAcceleration.hpp.

◆ _nbDelCols

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

Number of filtered out columns in this time window.

Definition at line 336 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 339 of file BaseQNAcceleration.hpp.

◆ _oldPrimaryResiduals

Eigen::VectorXd precice::acceleration::BaseQNAcceleration::_oldPrimaryResiduals
private

Difference between solver input and output from last time window regarding primary data.

Definition at line 319 of file BaseQNAcceleration.hpp.

◆ _oldPrimaryValues

Eigen::VectorXd precice::acceleration::BaseQNAcceleration::_oldPrimaryValues
private

Concatenation of all old primary data involved in the QN system.

Definition at line 316 of file BaseQNAcceleration.hpp.

◆ _oldPrimaryXTilde

Eigen::VectorXd precice::acceleration::BaseQNAcceleration::_oldPrimaryXTilde
protected

Solver output regarding primary data from last iteration.

Definition at line 201 of file BaseQNAcceleration.hpp.

◆ _oldValues

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

Concatenation of all old primary and secondary data involved in the QN system.

Definition at line 325 of file BaseQNAcceleration.hpp.

◆ _oldXTilde

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

Solver output from last iteration.

Definition at line 204 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 159 of file BaseQNAcceleration.hpp.

◆ _primaryDataIDs

const std::vector<int> precice::acceleration::BaseQNAcceleration::_primaryDataIDs
protected

Data IDs of primary data to be involved in the IQN coefficient computation.

Definition at line 171 of file BaseQNAcceleration.hpp.

◆ _primaryResiduals

Eigen::VectorXd precice::acceleration::BaseQNAcceleration::_primaryResiduals
protected

Current iteration residuals of primary IQN data. Temporary.

Definition at line 207 of file BaseQNAcceleration.hpp.

◆ _primaryTimeGrids

std::optional<time::TimeGrids> precice::acceleration::BaseQNAcceleration::_primaryTimeGrids
private

Stores the time grids to which the primary data involved in the QN system will be interpolated to. If _reducedTimeGrids is true then this will only contain the last time stamp of the time window, see https://doi.org/10.1002/nme.6443.

Definition at line 310 of file BaseQNAcceleration.hpp.

◆ _primaryValues

Eigen::VectorXd precice::acceleration::BaseQNAcceleration::_primaryValues
private

Concatenation of all primary data involved in the QN system.

Definition at line 313 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 222 of file BaseQNAcceleration.hpp.

◆ _reducedTimeGrid

const bool precice::acceleration::BaseQNAcceleration::_reducedTimeGrid
protected

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

Definition at line 219 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 198 of file BaseQNAcceleration.hpp.

◆ _residuals

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

Current iteration residuals of IQN data. Temporary.

Definition at line 210 of file BaseQNAcceleration.hpp.

◆ _singularityLimit

const 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 235 of file BaseQNAcceleration.hpp.

◆ _timeGrids

std::optional<time::TimeGrids> precice::acceleration::BaseQNAcceleration::_timeGrids
private

Stores the time grids to which the primary and secondary data involved in the QN system will be interpolated to.

Definition at line 306 of file BaseQNAcceleration.hpp.

◆ _timeWindowsReused

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

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

Definition at line 168 of file BaseQNAcceleration.hpp.

◆ _values

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

Concatenation of all primary and secondary data involved in the QN system.

Definition at line 322 of file BaseQNAcceleration.hpp.

◆ its

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

Definition at line 288 of file BaseQNAcceleration.hpp.

◆ tWindows

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

Definition at line 288 of file BaseQNAcceleration.hpp.


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