libMesh
discontinuity_measure.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/discontinuity_measure.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 {
45  error_norm = L2;
46 }
47 
48 
49 
52 {
53  return DISCONTINUITY_MEASURE;
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 values and mapping Jacobians on both sides for
78  // discontinuity computation
79  side_fe->get_phi();
80 
81  // But we only need mapping Jacobians from one side
82  if (&c != coarse_context.get())
83  side_fe->get_JxW();
84  }
85  }
86 }
87 
88 
89 
90 void
92 {
93  const Elem & coarse_elem = coarse_context->get_elem();
94  const Elem & fine_elem = fine_context->get_elem();
95 
96  FEBase * fe_fine = nullptr;
97  fine_context->get_side_fe( var, fe_fine, fine_elem.dim() );
98 
99  FEBase * fe_coarse = nullptr;
100  coarse_context->get_side_fe( var, fe_coarse, fine_elem.dim() );
101 
102  Real error = 1.e-30;
103  unsigned int n_qp = fe_fine->n_quadrature_points();
104 
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 values on fine and coarse elements
110  // at this quadrature point
111  Number
112  u_fine = fine_context->side_value(var, qp),
113  u_coarse = coarse_context->side_value(var, qp);
114 
115  // Find the jump in the value
116  // at this quadrature point
117  const Number jump = u_fine - u_coarse;
118  const Real jump2 = TensorTools::norm_sq(jump);
119  // Accumulate the jump integral
120  error += JxW_face[qp] * jump2;
121  }
122 
123  // Add the h-weighted jump integral to each error term
124  fine_error =
125  error * fine_elem.hmax() * error_norm.weight(var);
126  coarse_error =
127  error * coarse_elem.hmax() * error_norm.weight(var);
128 }
129 
130 
131 bool
133 {
134  const Elem & fine_elem = fine_context->get_elem();
135 
136  FEBase * fe_fine = nullptr;
137  fine_context->get_side_fe( var, fe_fine, fine_elem.dim() );
138 
139  const std::string & var_name =
140  fine_context->get_system().variable_name(var);
141 
142  const std::vector<Real> & JxW_face = fe_fine->get_JxW();
143  const std::vector<Point> & qface_point = fe_fine->get_xyz();
144 
145  // The reinitialization also recomputes the locations of
146  // the quadrature points on the side. By checking if the
147  // first quadrature point on the side is on an essential boundary
148  // for a particular variable, we will determine if the whole
149  // element is on an essential boundary (assuming quadrature points
150  // are strictly contained in the side).
151  if (this->_bc_function(fine_context->get_system(),
152  qface_point[0], var_name).first)
153  {
154  const Real h = fine_elem.hmax();
155 
156  // The number of quadrature points
157  const unsigned int n_qp = fe_fine->n_quadrature_points();
158 
159  // The error contribution from this face
160  Real error = 1.e-30;
161 
162  // loop over the integration points on the face.
163  for (unsigned int qp=0; qp<n_qp; qp++)
164  {
165  // Value of the imposed essential BC at this quadrature point.
166  const std::pair<bool,Real> essential_bc =
167  this->_bc_function(fine_context->get_system(), qface_point[qp],
168  var_name);
169 
170  // Be sure the BC function still thinks we're on the
171  // essential boundary.
172  libmesh_assert_equal_to (essential_bc.first, true);
173 
174  // The solution value on each point
175  Number u_fine = fine_context->side_value(var, qp);
176 
177  // The difference between the desired BC and the approximate solution.
178  const Number jump = essential_bc.second - u_fine;
179 
180  // The flux jump squared. If using complex numbers,
181  // norm_sq(z) returns |z|^2, where |z| is the modulus of z.
182  const Real jump2 = TensorTools::norm_sq(jump);
183 
184  // Integrate the error on the face. The error is
185  // scaled by an additional power of h, where h is
186  // the maximum side length for the element. This
187  // arises in the definition of the indicator.
188  error += JxW_face[qp]*jump2;
189 
190  } // End quadrature point loop
191 
192  fine_error = error*h*error_norm.weight(var);
193 
194  return true;
195  } // end if side on flux boundary
196  return false;
197 }
198 
199 
200 
201 void
203  const Point & p,
204  const std::string & var_name))
205 {
206  _bc_function = fptr;
207 
208  // We may be turning boundary side integration on or off
209  if (fptr)
211  else
212  integrate_boundary_sides = false;
213 }
214 
215 } // namespace libMesh
virtual bool boundary_side_integration() override
The function which calculates a normal derivative jump based error term on a boundary side...
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
The libMesh namespace provides an interface to certain functionality in the library.
virtual Real hmax() const
Definition: elem.C:722
void attach_essential_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 essential BCs.
Number fptr(const Point &p, const Parameters &, const std::string &libmesh_dbg_var(sys_name), const std::string &unknown_name)
Definition: projection.C:80
virtual void init_context(FEMContext &c) override
An initialization function, for requesting specific data from the FE objects.
std::pair< bool, Real >(* _bc_function)(const System &system, const Point &p, const std::string &var_name)
Pointer to function that provides BC information.
unsigned int n_vars
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
virtual void internal_side_integration() override
The function which calculates a normal derivative jump based error term on an internal side...
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 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 ErrorEstimatorType type() const override
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
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
const std::vector< std::vector< OutputShape > > & get_phi() const
Definition: fe_base.h:207