exact_error_estimator.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #ifndef LIBMESH_EXACT_ERROR_ESTIMATOR_H
21 #define LIBMESH_EXACT_ERROR_ESTIMATOR_H
22 
23 // Local Includes
25 #include "libmesh/function_base.h"
26 #include "libmesh/vector_value.h"
27 #include "libmesh/tensor_value.h"
28 
29 // C++ includes
30 #include <cstddef>
31 #include <string>
32 #include <vector>
33 
34 namespace libMesh
35 {
36 
37 // Forward Declarations
38 class Elem;
39 template <typename T> class FEGenericBase;
41 class MeshFunction;
42 class Point;
43 class Parameters;
44 template <typename T> class DenseVector;
45 
58 {
59 public:
60 
66 
72  ExactErrorEstimator (const ExactErrorEstimator &) = delete;
74 
80  virtual ~ExactErrorEstimator() = default;
81 
86  void attach_exact_values (std::vector<FunctionBase<Number> *> f);
87 
92  void attach_exact_value (unsigned int sys_num,
94 
99  typedef Number (*ValueFunctionPointer)(const Point & p,
100  const Parameters & Parameters,
101  const std::string & sys_name,
102  const std::string & unknown_name);
104 
109  void attach_exact_derivs (std::vector<FunctionBase<Gradient> *> g);
110 
115  void attach_exact_deriv (unsigned int sys_num,
117 
122  typedef Gradient (*GradientFunctionPointer)(const Point & p,
123  const Parameters & parameters,
124  const std::string & sys_name,
125  const std::string & unknown_name);
127 
132  void attach_exact_hessians (std::vector<FunctionBase<Tensor> *> h);
133 
138  void attach_exact_hessian (unsigned int sys_num,
140 
145  typedef Tensor (*HessianFunctionPointer)(const Point & p,
146  const Parameters & parameters,
147  const std::string & sys_name,
148  const std::string & unknown_name);
150 
157 
158 
163  void extra_quadrature_order (const int extraorder)
164  { _extra_order = extraorder; }
165 
166 
167  // Bring the base class functionality into the name lookup
168  // procedure. This allows for alternative calling formats
169  // defined in the base class. Thanks Wolfgang.
170  // GCC 2.95.3 cannot compile such code. Since it was not really
171  // essential to the functioning of this class, it's been removed.
172  // using ErrorEstimator::estimate_error;
173 
180  virtual void estimate_error (const System & system,
181  ErrorVector & error_per_cell,
182  const NumericVector<Number> * solution_vector = nullptr,
183  bool estimate_parent_error = false) override;
184 
185  virtual ErrorEstimatorType type() const override;
186 
187 private:
188 
194 
200 
206 
211  std::vector<std::unique_ptr<FunctionBase<Number>>> _exact_values;
212 
217  std::vector<std::unique_ptr<FunctionBase<Gradient>>> _exact_derivs;
218 
223  std::vector<std::unique_ptr<FunctionBase<Tensor>>> _exact_hessians;
224 
230 
234  Real find_squared_element_error (const System & system,
235  const std::string & var_name,
236  const Elem * elem,
237  const DenseVector<Number> & Uelem,
238  FEBase * fe,
239  MeshFunction * fine_values) const;
240 
244  void clear_functors ();
245 
250 };
251 
252 
253 } // namespace libMesh
254 
255 #endif // LIBMESH_EXACT_ERROR_ESTIMATOR_H
HessianFunctionPointer _exact_hessian
Manages multiples systems of equations.
void attach_exact_hessian(unsigned int sys_num, FunctionBase< Tensor > *h)
std::vector< std::unique_ptr< FunctionBase< Tensor > > > _exact_hessians
void attach_exact_hessians(std::vector< FunctionBase< Tensor > *> h)
The base class for all geometric element types.
Definition: elem.h:100
void extra_quadrature_order(const int extraorder)
void attach_exact_value(unsigned int sys_num, FunctionBase< Number > *f)
std::vector< std::unique_ptr< FunctionBase< Gradient > > > _exact_derivs
virtual void estimate_error(const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=nullptr, bool estimate_parent_error=false) override
virtual ErrorEstimatorType type() const override
Tensor(* HessianFunctionPointer)(const Point &p, const Parameters &parameters, const std::string &sys_name, const std::string &unknown_name)
virtual ~ExactErrorEstimator()=default
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:92
NumberVectorValue Gradient
FEGenericBase< Real > FEBase
void attach_exact_derivs(std::vector< FunctionBase< Gradient > *> g)
Gradient(* GradientFunctionPointer)(const Point &p, const Parameters &parameters, const std::string &sys_name, const std::string &unknown_name)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
GradientFunctionPointer _exact_deriv
NumberTensorValue Tensor
Number(* ValueFunctionPointer)(const Point &p, const Parameters &Parameters, const std::string &sys_name, const std::string &unknown_name)
Real find_squared_element_error(const System &system, const std::string &var_name, const Elem *elem, const DenseVector< Number > &Uelem, FEBase *fe, MeshFunction *fine_values) const
void attach_exact_values(std::vector< FunctionBase< Number > *> f)
void attach_reference_solution(EquationSystems *es_fine)
void attach_exact_deriv(unsigned int sys_num, FunctionBase< Gradient > *g)
A geometric point in (x,y,z) space.
Definition: point.h:38
ExactErrorEstimator & operator=(const ExactErrorEstimator &)=delete
std::vector< std::unique_ptr< FunctionBase< Number > > > _exact_values