OGS
|
Extrapolator doing a linear least squares extrapolation locally in each element.
The results of the element-wise least squares are averaged at each node.
The described procedure is not least-squares optimal on the whole mesh! But the residuals are computed from the actual interpolation result that is being returned by this class.
The number of integration points in each element must be greater than or equal to the number of nodes of that element. This restriction is due to the use of the least squares which requires an exact or overdetermined equation system.
Definition at line 38 of file LocalLinearLeastSquaresExtrapolator.h.
#include <LocalLinearLeastSquaresExtrapolator.h>
Classes | |
struct | CachedData |
Stores a matrix and its Moore-Penrose pseudo-inverse. More... | |
Public Member Functions | |
LocalLinearLeastSquaresExtrapolator (NumLib::LocalToGlobalIndexMap const &dof_table) | |
void | extrapolate (const unsigned num_components, ExtrapolatableElementCollection const &extrapolatables, const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table) override |
Extrapolates the given property from the given local assemblers. | |
void | calculateResiduals (const unsigned num_components, ExtrapolatableElementCollection const &extrapolatables, const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table) override |
GlobalVector const & | getNodalValues () const override |
GlobalVector const & | getElementResiduals () const override |
Public Member Functions inherited from NumLib::Extrapolator | |
virtual | ~Extrapolator ()=default |
Private Member Functions | |
void | extrapolateElement (std::size_t const element_index, const unsigned num_components, ExtrapolatableElementCollection const &extrapolatables, const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, GlobalVector &counts) |
Extrapolate one element. | |
void | calculateResidualElement (std::size_t const element_index, const unsigned num_components, ExtrapolatableElementCollection const &extrapolatables, const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table) |
Compute the residuals for one element. | |
Private Attributes | |
std::unique_ptr< GlobalVector > | _nodal_values |
extrapolated nodal values | |
std::unique_ptr< GlobalVector > | _residuals |
extrapolation residuals | |
NumLib::LocalToGlobalIndexMap const & | _dof_table_single_component |
DOF table used for writing to global vectors. | |
std::vector< double > | _integration_point_values_cache |
Avoids frequent reallocations. | |
std::map< std::pair< unsigned, unsigned >, CachedData > | _qr_decomposition_cache |
|
explicit |
Constructs a new instance.
dof_table
must point to a d.o.f. table for one single-component variable. Definition at line 25 of file LocalLinearLeastSquaresExtrapolator.cpp.
References NumLib::LocalToGlobalIndexMap::getNumberOfGlobalComponents().
|
private |
Compute the residuals for one element.
Definition at line 281 of file LocalLinearLeastSquaresExtrapolator.cpp.
References _dof_table_single_component, _integration_point_values_cache, _nodal_values, _qr_decomposition_cache, _residuals, NumLib::ExtrapolatableElementCollection::getIntegrationPointValues(), OGS_FATAL, MathLib::LinAlg::setLocalAccessibleVector(), and MathLib::toMatrix().
Referenced by calculateResiduals().
|
overridevirtual |
Computes residuals from the extrapolation of the given property
.
The residuals are computed as element values.
The computed residuals are root-mean-square of the difference between the integration point values obtained from the local assemblers and the extrapolation results when interpolated back to the integration points again.
Implements NumLib::Extrapolator.
Definition at line 107 of file LocalLinearLeastSquaresExtrapolator.cpp.
References _dof_table_single_component, _residuals, calculateResidualElement(), MathLib::LinAlg::finalizeAssembly(), OGS_FATAL, NumLib::ExtrapolatableElementCollection::size(), and NumLib::LocalToGlobalIndexMap::size().
|
overridevirtual |
Extrapolates the given property
from the given local assemblers.
Implements NumLib::Extrapolator.
Definition at line 44 of file LocalLinearLeastSquaresExtrapolator.cpp.
References _dof_table_single_component, _nodal_values, MathLib::LinAlg::componentwiseDivide(), NumLib::LocalToGlobalIndexMap::dofSizeWithoutGhosts(), extrapolateElement(), MathLib::LinAlg::finalizeAssembly(), NumLib::LocalToGlobalIndexMap::getGhostIndices(), and NumLib::ExtrapolatableElementCollection::size().
|
private |
Extrapolate one element.
Definition at line 141 of file LocalLinearLeastSquaresExtrapolator.cpp.
References _dof_table_single_component, _integration_point_values_cache, _nodal_values, _qr_decomposition_cache, MathLib::EigenVector::add(), DBUG(), NumLib::ExtrapolatableElementCollection::getIntegrationPointValues(), NumLib::ExtrapolatableElementCollection::getShapeMatrix(), OGS_FATAL, MathLib::toMatrix(), and MathLib::toVector().
Referenced by extrapolate().
|
inlineoverridevirtual |
Returns the extrapolation residuals.
Implements NumLib::Extrapolator.
Definition at line 77 of file LocalLinearLeastSquaresExtrapolator.h.
References _residuals.
|
inlineoverridevirtual |
Returns the extrapolated nodal values.
Implements NumLib::Extrapolator.
Definition at line 72 of file LocalLinearLeastSquaresExtrapolator.h.
References _nodal_values.
|
private |
DOF table used for writing to global vectors.
Definition at line 104 of file LocalLinearLeastSquaresExtrapolator.h.
Referenced by calculateResidualElement(), calculateResiduals(), extrapolate(), and extrapolateElement().
|
private |
Avoids frequent reallocations.
Definition at line 107 of file LocalLinearLeastSquaresExtrapolator.h.
Referenced by calculateResidualElement(), and extrapolateElement().
|
private |
extrapolated nodal values
Definition at line 100 of file LocalLinearLeastSquaresExtrapolator.h.
Referenced by calculateResidualElement(), extrapolate(), extrapolateElement(), and getNodalValues().
|
private |
Maps (#nodes, #int_pts) to (N_0, QR decomposition), where N_0 is the shape matrix of the first integration point.
Definition at line 129 of file LocalLinearLeastSquaresExtrapolator.h.
Referenced by calculateResidualElement(), and extrapolateElement().
|
private |
extrapolation residuals
Definition at line 101 of file LocalLinearLeastSquaresExtrapolator.h.
Referenced by calculateResidualElement(), calculateResiduals(), and getElementResiduals().