OGS
LocalLinearLeastSquaresExtrapolator.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
6#include <map>
7
8#include "Extrapolator.h"
11
12namespace NumLib
13{
32{
33public:
41 NumLib::LocalToGlobalIndexMap const& dof_table);
42
43 void extrapolate(const unsigned num_components,
44 ExtrapolatableElementCollection const& extrapolatables,
45 const double t,
46 std::vector<GlobalVector*> const& x,
47 std::vector<NumLib::LocalToGlobalIndexMap const*> const&
48 dof_table) override;
49
58 const unsigned num_components,
59 ExtrapolatableElementCollection const& extrapolatables,
60 const double t,
61 std::vector<GlobalVector*> const& x,
62 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table)
63 override;
64
65 GlobalVector const& getNodalValues() const override
66 {
67 return *_nodal_values;
68 }
69
70 GlobalVector const& getElementResiduals() const override
71 {
72 return *_residuals;
73 }
74
75private:
78 std::size_t const element_index, const unsigned num_components,
79 ExtrapolatableElementCollection const& extrapolatables, const double t,
80 std::vector<GlobalVector*> const& x,
81 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
82 GlobalVector& counts);
83
86 std::size_t const element_index,
87 const unsigned num_components,
88 ExtrapolatableElementCollection const& extrapolatables,
89 const double t,
90 std::vector<GlobalVector*> const& x,
91 std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table);
92
93 std::unique_ptr<GlobalVector> _nodal_values;
94 std::unique_ptr<GlobalVector> _residuals;
95
98
101
104 {
106 Eigen::MatrixXd A;
107
109 Eigen::MatrixXd A_pinv;
110 };
111
122 std::map<std::pair<unsigned, unsigned>, CachedData> _qr_decomposition_cache;
123};
124
125} // namespace NumLib
MathLib::EigenVector GlobalVector
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.
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.
LocalLinearLeastSquaresExtrapolator(NumLib::LocalToGlobalIndexMap const &dof_table)
NumLib::LocalToGlobalIndexMap const & _dof_table_single_component
DOF table used for writing to global vectors.
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.
std::unique_ptr< GlobalVector > _nodal_values
extrapolated nodal values
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
std::unique_ptr< GlobalVector > _residuals
extrapolation residuals
std::vector< double > _integration_point_values_cache
Avoids frequent reallocations.
std::map< std::pair< unsigned, unsigned >, CachedData > _qr_decomposition_cache
Stores a matrix and its Moore-Penrose pseudo-inverse.