jump_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_JUMP_ERROR_ESTIMATOR_H
21 #define LIBMESH_JUMP_ERROR_ESTIMATOR_H
22 
23 // Local Includes
24 #include "libmesh/auto_ptr.h" // deprecated
25 #include "libmesh/dense_vector.h"
27 #include "libmesh/fem_context.h"
28 
29 // C++ includes
30 #include <cstddef>
31 #include <string>
32 #include <vector>
33 #include <memory>
34 
35 namespace libMesh
36 {
37 
38 // Forward Declarations
39 class Point;
40 class Elem;
41 
50 {
51 public:
52 
57  : ErrorEstimator(),
58  scale_by_n_flux_faces(false),
60  fine_context(),
62  fine_error(0),
63  coarse_error(0) {}
64 
70  JumpErrorEstimator (const JumpErrorEstimator &) = delete;
72 
78  virtual ~JumpErrorEstimator() = default;
79 
86  virtual void estimate_error (const System & system,
87  ErrorVector & error_per_cell,
88  const NumericVector<Number> * solution_vector = nullptr,
89  bool estimate_parent_error = false) override;
90 
100 
101 protected:
106  void reinit_sides();
107 
112 
117  virtual void init_context(FEMContext & c);
118 
123  virtual void internal_side_integration() = 0;
124 
132  virtual bool boundary_side_integration() { return false; }
133 
139 
144  std::unique_ptr<FEMContext> fine_context, coarse_context;
145 
150 
154  unsigned int var;
155 };
156 
157 
158 } // namespace libMesh
159 
160 #endif // LIBMESH_JUMP_ERROR_ESTIMATOR_H
JumpErrorEstimator & operator=(const JumpErrorEstimator &)=delete
std::unique_ptr< FEMContext > fine_context
virtual void internal_side_integration()=0
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:92
std::unique_ptr< FEMContext > coarse_context
virtual ~JumpErrorEstimator()=default
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void estimate_error(const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=nullptr, bool estimate_parent_error=false) override
virtual void init_context(FEMContext &c)