kelly_error_estimator.C
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 // 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"
28 #include "libmesh/error_vector.h"
29 #include "libmesh/fe_base.h"
31 #include "libmesh/elem.h"
32 #include "libmesh/system.h"
33 #include "libmesh/dense_vector.h"
34 #include "libmesh/tensor_tools.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  fine_context->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  std::vector<std::vector<RealGradient>> dphi_coarse = fe_coarse->get_dphi();
105  std::vector<std::vector<RealGradient>> dphi_fine = fe_fine->get_dphi();
106  std::vector<Point> face_normals = fe_fine->get_normals();
107  std::vector<Real> JxW_face = fe_fine->get_JxW();
108 
109  for (unsigned int qp=0; qp != n_qp; ++qp)
110  {
111  // Calculate solution gradients on fine and coarse elements
112  // at this quadrature point
113  Gradient
114  grad_fine = fine_context->side_gradient(var, qp),
115  grad_coarse = coarse_context->side_gradient(var, qp);
116 
117  // Find the jump in the normal derivative
118  // at this quadrature point
119  const Number jump = (grad_fine - grad_coarse)*face_normals[qp];
120  const Real jump2 = TensorTools::norm_sq(jump);
121 
122  // Accumulate the jump integral
123  error += JxW_face[qp] * jump2;
124  }
125 
126  // Add the h-weighted jump integral to each error term
127  fine_error =
128  error * fine_elem.hmax() * error_norm.weight(var);
129  coarse_error =
130  error * coarse_elem.hmax() * error_norm.weight(var);
131 }
132 
133 
134 bool
136 {
137  const Elem & fine_elem = fine_context->get_elem();
138 
139  FEBase * fe_fine = nullptr;
140  fine_context->get_side_fe( var, fe_fine, fine_elem.dim() );
141 
142  const std::string & var_name =
143  fine_context->get_system().variable_name(var);
144 
145  std::vector<std::vector<RealGradient>> dphi_fine = fe_fine->get_dphi();
146  std::vector<Point> face_normals = fe_fine->get_normals();
147  std::vector<Real> JxW_face = fe_fine->get_JxW();
148  std::vector<Point> qface_point = fe_fine->get_xyz();
149 
150  // The reinitialization also recomputes the locations of
151  // the quadrature points on the side. By checking if the
152  // first quadrature point on the side is on a flux boundary
153  // for a particular variable, we will determine if the whole
154  // element is on a flux boundary (assuming quadrature points
155  // are strictly contained in the side).
156  if (this->_bc_function(fine_context->get_system(),
157  qface_point[0], var_name).first)
158  {
159  const Real h = fine_elem.hmax();
160 
161  // The number of quadrature points
162  const unsigned int n_qp = fe_fine->n_quadrature_points();
163 
164  // The error contribution from this face
165  Real error = 1.e-30;
166 
167  // loop over the integration points on the face.
168  for (unsigned int qp=0; qp<n_qp; qp++)
169  {
170  // Value of the imposed flux BC at this quadrature point.
171  const std::pair<bool,Real> flux_bc =
172  this->_bc_function(fine_context->get_system(),
173  qface_point[qp], var_name);
174 
175  // Be sure the BC function still thinks we're on the
176  // flux boundary.
177  libmesh_assert_equal_to (flux_bc.first, true);
178 
179  // The solution gradient from each element
180  Gradient grad_fine = fine_context->side_gradient(var, qp);
181 
182  // The difference between the desired BC and the approximate solution.
183  const Number jump = flux_bc.second - grad_fine*face_normals[qp];
184 
185  // The flux jump squared. If using complex numbers,
186  // TensorTools::norm_sq(z) returns |z|^2, where |z| is the modulus of z.
187  const Real jump2 = TensorTools::norm_sq(jump);
188 
189  // Integrate the error on the face. The error is
190  // scaled by an additional power of h, where h is
191  // the maximum side length for the element. This
192  // arises in the definition of the indicator.
193  error += JxW_face[qp]*jump2;
194 
195  } // End quadrature point loop
196 
197  fine_error = error*h*error_norm.weight(var);
198 
199  return true;
200  } // end if side on flux boundary
201  return false;
202 }
203 
204 
205 
206 void
207 KellyErrorEstimator::attach_flux_bc_function (std::pair<bool,Real> fptr(const System & system,
208  const Point & p,
209  const std::string & var_name))
210 {
211  _bc_function = fptr;
212 
213  // We may be turning boundary side integration on or off
214  if (fptr)
216  else
217  integrate_boundary_sides = false;
218 }
219 
220 } // namespace libMesh
void attach_flux_bc_function(std::pair< bool, Real > fptr(const System &system, const Point &p, const std::string &var_name))
std::pair< bool, Real >(* _bc_function)(const System &system, const Point &p, const std::string &var_name)
const std::vector< Point > & get_xyz() const
Definition: fe_abstract.h:238
std::unique_ptr< FEMContext > fine_context
The base class for all geometric element types.
Definition: elem.h:100
virtual unsigned int n_quadrature_points() const =0
virtual void init_context(FEMContext &c) override
const unsigned int n_vars
Definition: tecplot_io.C:69
virtual Real hmax() const
Definition: elem.C:373
const std::vector< std::vector< OutputGradient > > & get_dphi() const
Definition: fe_base.h:215
const std::vector< Point > & get_normals() const
Definition: fe_abstract.h:385
virtual void internal_side_integration() override
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:92
const std::vector< Real > & get_JxW() const
Definition: fe_abstract.h:245
std::unique_ptr< FEMContext > coarse_context
Real weight(unsigned int var) const
Definition: system_norm.C:132
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual unsigned short dim() const =0
virtual bool boundary_side_integration() override
unsigned int n_vars() const
Definition: diff_context.h:99
A geometric point in (x,y,z) space.
Definition: point.h:38
virtual ErrorEstimatorType type() const override
const std::set< unsigned char > & elem_dimensions() const
Definition: fem_context.h:913