libMesh
kelly_error_estimator.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 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 // C++ includes
20 #include <algorithm> // for std::fill
21 #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC
22 #include <cmath> // for sqrt
23 
24 
25 // Local Includes
26 #include "libmesh/libmesh_common.h"
27 #include "libmesh/kelly_error_estimator.h"
28 #include "libmesh/error_vector.h"
29 #include "libmesh/fe_base.h"
30 #include "libmesh/libmesh_logging.h"
31 #include "libmesh/elem.h"
32 #include "libmesh/system.h"
33 #include "libmesh/dense_vector.h"
34 #include "libmesh/tensor_tools.h"
35 #include "libmesh/enum_error_estimator_type.h"
36 #include "libmesh/enum_norm_type.h"
37 
38 namespace libMesh
39 {
40 
43  _bc_function(nullptr)
44 {
46 }
47 
48 
49 
52 {
53  return KELLY;
54 }
55 
56 
57 
58 void
60 {
61  const unsigned int n_vars = c.n_vars();
62  for (unsigned int v=0; v<n_vars; v++)
63  {
64  // Possibly skip this variable
65  if (error_norm.weight(v) == 0.0) continue;
66 
67  // FIXME: Need to generalize this to vector-valued elements. [PB]
68  FEBase * side_fe = nullptr;
69 
70  const std::set<unsigned char> & elem_dims =
71  c.elem_dimensions();
72 
73  for (const auto & dim : elem_dims)
74  {
75  c.get_side_fe( v, side_fe, dim );
76 
77  // We'll need gradients on both sides for flux jump computation
78  side_fe->get_dphi();
79 
80  // But we only need normal vectors from one side
81  if (&c != coarse_context.get())
82  side_fe->get_normals();
83  }
84  }
85 }
86 
87 
88 
89 void
91 {
92  const Elem & coarse_elem = coarse_context->get_elem();
93  const Elem & fine_elem = fine_context->get_elem();
94 
95  FEBase * fe_fine = nullptr;
96  fine_context->get_side_fe( var, fe_fine, fine_elem.dim() );
97 
98  FEBase * fe_coarse = nullptr;
99  coarse_context->get_side_fe( var, fe_coarse, fine_elem.dim() );
100 
101  Real error = 1.e-30;
102  unsigned int n_qp = fe_fine->n_quadrature_points();
103 
104  const std::vector<Point> & face_normals = fe_fine->get_normals();
105  const std::vector<Real> & JxW_face = fe_fine->get_JxW();
106 
107  for (unsigned int qp=0; qp != n_qp; ++qp)
108  {
109  // Calculate solution gradients on fine and coarse elements
110  // at this quadrature point
111  Gradient
112  grad_fine = fine_context->side_gradient(var, qp),
113  grad_coarse = coarse_context->side_gradient(var, qp);
114 
115  // Find the jump in the normal derivative
116  // at this quadrature point
117  const Number jump = (grad_fine - grad_coarse)*face_normals[qp];
118  const Real jump2 = TensorTools::norm_sq(jump);
119 
120  // Accumulate the jump integral
121  error += JxW_face[qp] * jump2;
122  }
123 
124  // Add the h-weighted jump integral to each error term
125  fine_error =
126  error * fine_elem.hmax() * error_norm.weight(var);
127  coarse_error =
128  error * coarse_elem.hmax() * error_norm.weight(var);
129 }
130 
131 
132 bool
134 {
135  const Elem & fine_elem = fine_context->get_elem();
136 
137  FEBase * fe_fine = nullptr;
138  fine_context->get_side_fe( var, fe_fine, fine_elem.dim() );
139 
140  const std::string & var_name =
141  fine_context->get_system().variable_name(var);
142 
143  const std::vector<Point> & face_normals = fe_fine->get_normals();
144  const std::vector<Real> & JxW_face = fe_fine->get_JxW();
145  const std::vector<Point> & qface_point = fe_fine->get_xyz();
146 
147  // The reinitialization also recomputes the locations of
148  // the quadrature points on the side. By checking if the
149  // first quadrature point on the side is on a flux boundary
150  // for a particular variable, we will determine if the whole
151  // element is on a flux boundary (assuming quadrature points
152  // are strictly contained in the side).
153  if (this->_bc_function(fine_context->get_system(),
154  qface_point[0], var_name).first)
155  {
156  const Real h = fine_elem.hmax();
157 
158  // The number of quadrature points
159  const unsigned int n_qp = fe_fine->n_quadrature_points();
160 
161  // The error contribution from this face
162  Real error = 1.e-30;
163 
164  // loop over the integration points on the face.
165  for (unsigned int qp=0; qp<n_qp; qp++)
166  {
167  // Value of the imposed flux BC at this quadrature point.
168  const std::pair<bool,Real> flux_bc =
169  this->_bc_function(fine_context->get_system(),
170  qface_point[qp], var_name);
171 
172  // Be sure the BC function still thinks we're on the
173  // flux boundary.
174  libmesh_assert_equal_to (flux_bc.first, true);
175 
176  // The solution gradient from each element
177  Gradient grad_fine = fine_context->side_gradient(var, qp);
178 
179  // The difference between the desired BC and the approximate solution.
180  const Number jump = flux_bc.second - grad_fine*face_normals[qp];
181 
182  // The flux jump squared. If using complex numbers,
183  // TensorTools::norm_sq(z) returns |z|^2, where |z| is the modulus of z.
184  const Real jump2 = TensorTools::norm_sq(jump);
185 
186  // Integrate the error on the face. The error is
187  // scaled by an additional power of h, where h is
188  // the maximum side length for the element. This
189  // arises in the definition of the indicator.
190  error += JxW_face[qp]*jump2;
191 
192  } // End quadrature point loop
193 
194  fine_error = error*h*error_norm.weight(var);
195 
196  return true;
197  } // end if side on flux boundary
198  return false;
199 }
200 
201 
202 
203 void
204 KellyErrorEstimator::attach_flux_bc_function (std::pair<bool,Real> fptr(const System & system,
205  const Point & p,
206  const std::string & var_name))
207 {
208  _bc_function = fptr;
209 
210  // We may be turning boundary side integration on or off
211  if (fptr)
213  else
214  integrate_boundary_sides = false;
215 }
216 
217 } // namespace libMesh
void attach_flux_bc_function(std::pair< bool, Real > fptr(const System &system, const Point &p, const std::string &var_name))
Register a user function to use in computing the flux BCs.
std::pair< bool, Real >(* _bc_function)(const System &system, const Point &p, const std::string &var_name)
Pointer to function that provides BC information.
void get_side_fe(unsigned int var, FEGenericBase< OutputShape > *&fe) const
Accessor for edge/face (2D/3D) finite element object for variable var for the largest dimension in th...
Definition: fem_context.h:317
SystemNorm error_norm
When estimating the error in a single system, the error_norm is used to control the scaling and norm ...
std::unique_ptr< FEMContext > fine_context
Context objects for integrating on the fine and coarse elements sharing a face.
unsigned int dim
ErrorEstimatorType
Defines an enum for the different types of error estimators which are available.
This abstract base class implements utility functions for error estimators which are based on integra...
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
virtual void init_context(FEMContext &c) override
An initialization function, for requesting specific data from the FE objects.
This class defines a vector in LIBMESH_DIM dimensional Real or Complex space.
The libMesh namespace provides an interface to certain functionality in the library.
virtual Real hmax() const
Definition: elem.C:722
Number fptr(const Point &p, const Parameters &, const std::string &libmesh_dbg_var(sys_name), const std::string &unknown_name)
Definition: projection.C:80
const std::vector< std::vector< OutputGradient > > & get_dphi() const
Definition: fe_base.h:230
unsigned int n_vars
virtual void internal_side_integration() override
The function which calculates a normal derivative jump based error term on an internal side...
virtual unsigned int n_quadrature_points() const
Definition: fe_abstract.C:1279
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:96
virtual_for_inffe const std::vector< Real > & get_JxW() const
Definition: fe_abstract.h:303
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:62
virtual_for_inffe const std::vector< Point > & get_xyz() const
Definition: fe_abstract.h:280
std::unique_ptr< FEMContext > coarse_context
Real fine_error
The fine and coarse error values to be set by each side_integration();.
Real weight(unsigned int var) const
Definition: system_norm.C:134
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual_for_inffe const std::vector< Point > & get_normals() const
Definition: fe_abstract.h:459
virtual unsigned short dim() const =0
auto norm_sq(const T &a) -> decltype(std::norm(a))
Definition: tensor_tools.h:104
unsigned int var
The variable number currently being evaluated.
virtual bool boundary_side_integration() override
The function which calculates a normal derivative jump based error term on a boundary side...
unsigned int n_vars() const
Number of variables in solution.
Definition: diff_context.h:100
bool integrate_boundary_sides
A boolean flag, by default false, to be set to true if integrations with boundary_side_integration() ...
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
virtual ErrorEstimatorType type() const override
This class forms the foundation from which generic finite elements may be derived.
const std::set< unsigned char > & elem_dimensions() const
Definition: fem_context.h:951