OGS
LocalLinearLeastSquaresExtrapolator.h
Go to the documentation of this file.
1 
11 #pragma once
12 
13 #include <map>
14 
15 #include "Extrapolator.h"
18 
19 namespace NumLib
20 {
39 {
40 public:
48  NumLib::LocalToGlobalIndexMap const& dof_table);
49 
50  void extrapolate(const unsigned num_components,
51  ExtrapolatableElementCollection const& extrapolatables,
52  const double t,
53  std::vector<GlobalVector*> const& x,
54  std::vector<NumLib::LocalToGlobalIndexMap const*> const&
55  dof_table) override;
56 
64  void calculateResiduals(
65  const unsigned num_components,
66  ExtrapolatableElementCollection const& extrapolatables,
67  const double t,
68  std::vector<GlobalVector*> const& x,
69  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table)
70  override;
71 
72  GlobalVector const& getNodalValues() const override
73  {
74  return *_nodal_values;
75  }
76 
77  GlobalVector const& getElementResiduals() const override
78  {
79  return *_residuals;
80  }
81 
82 private:
84  void extrapolateElement(
85  std::size_t const element_index, const unsigned num_components,
86  ExtrapolatableElementCollection const& extrapolatables, const double t,
87  std::vector<GlobalVector*> const& x,
88  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table,
89  GlobalVector& counts);
90 
93  std::size_t const element_index,
94  const unsigned num_components,
95  ExtrapolatableElementCollection const& extrapolatables,
96  const double t,
97  std::vector<GlobalVector*> const& x,
98  std::vector<NumLib::LocalToGlobalIndexMap const*> const& dof_table);
99 
100  std::unique_ptr<GlobalVector> _nodal_values;
101  std::unique_ptr<GlobalVector> _residuals;
102 
105 
107  std::vector<double> _integration_point_values_cache;
108 
110  struct CachedData
111  {
113  Eigen::MatrixXd A;
114 
116  Eigen::MatrixXd A_pinv;
117  };
118 
129  std::map<std::pair<unsigned, unsigned>, CachedData> _qr_decomposition_cache;
130 };
131 
132 } // namespace NumLib
Global vector based on Eigen vector.
Definition: EigenVector.h:26
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.