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