libMesh
Public Types | Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | List of all members
libMesh::DiscontinuityMeasure Class Reference

This class measures discontinuities between elements for debugging purposes. More...

#include <discontinuity_measure.h>

Inheritance diagram for libMesh::DiscontinuityMeasure:
[legend]

Public Types

typedef std::map< std::pair< const System *, unsigned int >, ErrorVector * > ErrorMap
 When calculating many error vectors at once, we need a data structure to hold them all. More...
 

Public Member Functions

 DiscontinuityMeasure ()
 Constructor. More...
 
 DiscontinuityMeasure (const DiscontinuityMeasure &)=delete
 This class cannot be (default) copy constructed/assigned because its base class has unique_ptr members. More...
 
DiscontinuityMeasureoperator= (const DiscontinuityMeasure &)=delete
 
 DiscontinuityMeasure (DiscontinuityMeasure &&)=default
 Defaulted move ctor, move assignment operator, and destructor. More...
 
DiscontinuityMeasureoperator= (DiscontinuityMeasure &&)=default
 
virtual ~DiscontinuityMeasure ()=default
 
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. More...
 
virtual ErrorEstimatorType type () const override
 
virtual void estimate_error (const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=nullptr, bool estimate_parent_error=false) override
 This function uses the derived class's jump error estimate formula to estimate the error on each cell. More...
 
virtual void estimate_errors (const EquationSystems &equation_systems, ErrorVector &error_per_cell, const std::map< const System *, SystemNorm > &error_norms, const std::map< const System *, const NumericVector< Number > * > *solution_vectors=nullptr, bool estimate_parent_error=false)
 This virtual function can be redefined in derived classes, but by default computes the sum of the error_per_cell for each system in the equation_systems. More...
 
virtual void estimate_errors (const EquationSystems &equation_systems, ErrorMap &errors_per_cell, const std::map< const System *, const NumericVector< Number > * > *solution_vectors=nullptr, bool estimate_parent_error=false)
 This virtual function can be redefined in derived classes, but by default it calls estimate_error repeatedly to calculate the requested error vectors. More...
 

Public Attributes

bool scale_by_n_flux_faces
 This boolean flag allows you to scale the error indicator result for each element by the number of "flux faces" the element actually has. More...
 
bool use_unweighted_quadrature_rules
 This boolean flag allows you to use "unweighted" quadrature rules (sized to exactly integrate unweighted shape functions in master element space) rather than "default" quadrature rules (sized to exactly integrate polynomials of one higher degree than mass matrix terms). More...
 
SystemNorm error_norm
 When estimating the error in a single system, the error_norm is used to control the scaling and norm choice for each variable. More...
 

Protected Member Functions

virtual void init_context (FEMContext &c) override
 An initialization function, for requesting specific data from the FE objects. More...
 
virtual void internal_side_integration () override
 The function which calculates a normal derivative jump based error term on an internal side. More...
 
virtual bool boundary_side_integration () override
 The function which calculates a normal derivative jump based error term on a boundary side. More...
 
void reinit_sides ()
 A utility function to reinit the finite element data on elements sharing a side. More...
 
float coarse_n_flux_faces_increment ()
 A utility function to correctly increase n_flux_faces for the coarse element. More...
 
void reduce_error (std::vector< ErrorVectorReal > &error_per_cell, const Parallel::Communicator &comm) const
 This method takes the local error contributions in error_per_cell from each processor and combines them to get the global error vector. More...
 

Protected Attributes

std::pair< bool, Real >(* _bc_function )(const System &system, const Point &p, const std::string &var_name)
 Pointer to function that provides BC information. More...
 
bool integrate_boundary_sides
 A boolean flag, by default false, to be set to true if integrations with boundary_side_integration() should be performed. More...
 
std::unique_ptr< FEMContextfine_context
 Context objects for integrating on the fine and coarse elements sharing a face. More...
 
std::unique_ptr< FEMContextcoarse_context
 
Real fine_error
 The fine and coarse error values to be set by each side_integration();. More...
 
Real coarse_error
 
unsigned int var
 The variable number currently being evaluated. More...
 

Detailed Description

This class measures discontinuities between elements for debugging purposes.

It derives from ErrorEstimator just in case someone finds it useful in a DG framework.

Author
Roy H. Stogner
Date
2006

Definition at line 49 of file discontinuity_measure.h.

Member Typedef Documentation

◆ ErrorMap

typedef std::map<std::pair<const System *, unsigned int>, ErrorVector *> libMesh::ErrorEstimator::ErrorMap
inherited

When calculating many error vectors at once, we need a data structure to hold them all.

Definition at line 127 of file error_estimator.h.

Constructor & Destructor Documentation

◆ DiscontinuityMeasure() [1/3]

libMesh::DiscontinuityMeasure::DiscontinuityMeasure ( )

Constructor.

Responsible for initializing the _bc_function function pointer to nullptr. Defaults to L2 norm; changes to system norm are ignored.

Definition at line 41 of file discontinuity_measure.C.

41  :
43  _bc_function(nullptr)
44 {
45  error_norm = L2;
46 }

References libMesh::ErrorEstimator::error_norm, and libMesh::L2.

◆ DiscontinuityMeasure() [2/3]

libMesh::DiscontinuityMeasure::DiscontinuityMeasure ( const DiscontinuityMeasure )
delete

This class cannot be (default) copy constructed/assigned because its base class has unique_ptr members.

◆ DiscontinuityMeasure() [3/3]

libMesh::DiscontinuityMeasure::DiscontinuityMeasure ( DiscontinuityMeasure &&  )
default

Defaulted move ctor, move assignment operator, and destructor.

◆ ~DiscontinuityMeasure()

virtual libMesh::DiscontinuityMeasure::~DiscontinuityMeasure ( )
virtualdefault

Member Function Documentation

◆ attach_essential_bc_function()

void libMesh::DiscontinuityMeasure::attach_essential_bc_function ( std::pair< bool, Real >   fptrconst System &system, const Point &p, const std::string &var_name)

Register a user function to use in computing the essential BCs.

Definition at line 200 of file discontinuity_measure.C.

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 }

References _bc_function, fptr(), and libMesh::JumpErrorEstimator::integrate_boundary_sides.

◆ boundary_side_integration()

bool libMesh::DiscontinuityMeasure::boundary_side_integration ( )
overrideprotectedvirtual

The function which calculates a normal derivative jump based error term on a boundary side.

Returns
true if the flux bc function is in fact defined on the current side.

Reimplemented from libMesh::JumpErrorEstimator.

Definition at line 129 of file discontinuity_measure.C.

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 }

References _bc_function, libMesh::Elem::dim(), libMesh::ErrorEstimator::error_norm, libMesh::JumpErrorEstimator::fine_context, libMesh::JumpErrorEstimator::fine_error, libMesh::FEAbstract::get_JxW(), libMesh::FEGenericBase< OutputType >::get_phi(), libMesh::FEAbstract::get_xyz(), libMesh::Elem::hmax(), libMesh::FEAbstract::n_quadrature_points(), libMesh::TensorTools::norm_sq(), libMesh::Real, libMesh::JumpErrorEstimator::var, and libMesh::SystemNorm::weight().

◆ coarse_n_flux_faces_increment()

float libMesh::JumpErrorEstimator::coarse_n_flux_faces_increment ( )
protectedinherited

A utility function to correctly increase n_flux_faces for the coarse element.

Definition at line 467 of file jump_error_estimator.C.

468 {
469  // Keep track of the number of internal flux sides found on each
470  // element
471  unsigned short dim = coarse_context->get_elem().dim();
472 
473  const unsigned int divisor =
474  1 << (dim-1)*(fine_context->get_elem().level() -
475  coarse_context->get_elem().level());
476 
477  // With a difference of n levels between fine and coarse elements,
478  // we compute a fractional flux face for the coarse element by adding:
479  // 1/2^n in 2D
480  // 1/4^n in 3D
481  // each time. This code will get hit 2^n times in 2D and 4^n
482  // times in 3D so that the final flux face count for the coarse
483  // element will be an integer value.
484 
485  return 1.0f / static_cast<float>(divisor);
486 }

References libMesh::JumpErrorEstimator::coarse_context, dim, and libMesh::JumpErrorEstimator::fine_context.

Referenced by libMesh::JumpErrorEstimator::estimate_error().

◆ estimate_error()

void libMesh::JumpErrorEstimator::estimate_error ( const System system,
ErrorVector error_per_cell,
const NumericVector< Number > *  solution_vector = nullptr,
bool  estimate_parent_error = false 
)
overridevirtualinherited

This function uses the derived class's jump error estimate formula to estimate the error on each cell.

The estimated error is output in the vector error_per_cell

Conventions for assigning the direction of the normal:

  • e & f are global element ids

Case (1.) Elements are at the same level, e<f Compute the flux jump on the face and add it as a contribution to error_per_cell[e] and error_per_cell[f]


| | |

f
e —> n

Case (2.) The neighbor is at a higher level. Compute the flux jump on e's face and add it as a contribution to error_per_cell[e] and error_per_cell[f]


| | | | | | e |—> n | | | | | |--------—| f |


Implements libMesh::ErrorEstimator.

Definition at line 53 of file jump_error_estimator.C.

57 {
58  LOG_SCOPE("estimate_error()", "JumpErrorEstimator");
59 
96  // This parameter is not used when !LIBMESH_ENABLE_AMR.
97  libmesh_ignore(estimate_parent_error);
98 
99  // The current mesh
100  const MeshBase & mesh = system.get_mesh();
101 
102  // The number of variables in the system
103  const unsigned int n_vars = system.n_vars();
104 
105  // The DofMap for this system
106 #ifdef LIBMESH_ENABLE_AMR
107  const DofMap & dof_map = system.get_dof_map();
108 #endif
109 
110  // Resize the error_per_cell vector to be
111  // the number of elements, initialize it to 0.
112  error_per_cell.resize (mesh.max_elem_id());
113  std::fill (error_per_cell.begin(), error_per_cell.end(), 0.);
114 
115  // Declare a vector of floats which is as long as
116  // error_per_cell above, and fill with zeros. This vector will be
117  // used to keep track of the number of edges (faces) on each active
118  // element which are either:
119  // 1) an internal edge
120  // 2) an edge on a Neumann boundary for which a boundary condition
121  // function has been specified.
122  // The error estimator can be scaled by the number of flux edges (faces)
123  // which the element actually has to obtain a more uniform measure
124  // of the error. Use floats instead of ints since in case 2 (above)
125  // f gets 1/2 of a flux face contribution from each of his
126  // neighbors
127  std::vector<float> n_flux_faces;
129  n_flux_faces.resize(error_per_cell.size(), 0);
130 
131  // Prepare current_local_solution to localize a non-standard
132  // solution vector if necessary
133  if (solution_vector && solution_vector != system.solution.get())
134  {
135  NumericVector<Number> * newsol =
136  const_cast<NumericVector<Number> *>(solution_vector);
137  System & sys = const_cast<System &>(system);
138  newsol->swap(*sys.solution);
139  sys.update();
140  }
141 
142  fine_context = libmesh_make_unique<FEMContext>(system);
143  coarse_context = libmesh_make_unique<FEMContext>(system);
144 
145  // Don't overintegrate - we're evaluating differences of FE values,
146  // not products of them.
148  fine_context->use_unweighted_quadrature_rules(system.extra_quadrature_order);
149 
150  // Loop over all the variables we've been requested to find jumps in, to
151  // pre-request
152  for (var=0; var<n_vars; var++)
153  {
154  // Skip variables which aren't part of our norm,
155  // as well as SCALAR variables, which have no jumps
156  if (error_norm.weight(var) == 0.0 ||
157  system.variable_type(var).family == SCALAR)
158  continue;
159 
160  // FIXME: Need to generalize this to vector-valued elements. [PB]
161  FEBase * side_fe = nullptr;
162 
163  const std::set<unsigned char> & elem_dims =
164  fine_context->elem_dimensions();
165 
166  for (const auto & dim : elem_dims)
167  {
168  fine_context->get_side_fe( var, side_fe, dim );
169 
170  side_fe->get_xyz();
171  }
172  }
173 
174  this->init_context(*fine_context);
176 
177  // Iterate over all the active elements in the mesh
178  // that live on this processor.
179  for (const auto & e : mesh.active_local_element_ptr_range())
180  {
181  const dof_id_type e_id = e->id();
182 
183 #ifdef LIBMESH_ENABLE_AMR
184  // See if the parent of element e has been examined yet;
185  // if not, we may want to compute the estimator on it
186  const Elem * parent = e->parent();
187 
188  // We only can compute and only need to compute on
189  // parents with all active children
190  bool compute_on_parent = true;
191  if (!parent || !estimate_parent_error)
192  compute_on_parent = false;
193  else
194  for (auto & child : parent->child_ref_range())
195  if (!child.active())
196  compute_on_parent = false;
197 
198  if (compute_on_parent &&
199  !error_per_cell[parent->id()])
200  {
201  // Compute a projection onto the parent
202  DenseVector<Number> Uparent;
204  (*(system.solution), dof_map, parent, Uparent, false);
205 
206  // Loop over the neighbors of the parent
207  for (auto n_p : parent->side_index_range())
208  {
209  if (parent->neighbor_ptr(n_p) != nullptr) // parent has a neighbor here
210  {
211  // Find the active neighbors in this direction
212  std::vector<const Elem *> active_neighbors;
213  parent->neighbor_ptr(n_p)->
214  active_family_tree_by_neighbor(active_neighbors,
215  parent);
216  // Compute the flux to each active neighbor
217  for (std::size_t a=0,
218  n_active_neighbors = active_neighbors.size();
219  a != n_active_neighbors; ++a)
220  {
221  const Elem * f = active_neighbors[a];
222  // FIXME - what about when f->level <
223  // parent->level()??
224  if (f->level() >= parent->level())
225  {
226  fine_context->pre_fe_reinit(system, f);
227  coarse_context->pre_fe_reinit(system, parent);
228  libmesh_assert_equal_to
229  (coarse_context->get_elem_solution().size(),
230  Uparent.size());
231  coarse_context->get_elem_solution() = Uparent;
232 
233  this->reinit_sides();
234 
235  // Loop over all significant variables in the system
236  for (var=0; var<n_vars; var++)
237  if (error_norm.weight(var) != 0.0 &&
238  system.variable_type(var).family != SCALAR)
239  {
241 
242  error_per_cell[fine_context->get_elem().id()] +=
243  static_cast<ErrorVectorReal>(fine_error);
244  error_per_cell[coarse_context->get_elem().id()] +=
245  static_cast<ErrorVectorReal>(coarse_error);
246  }
247 
248  // Keep track of the number of internal flux
249  // sides found on each element
251  {
252  n_flux_faces[fine_context->get_elem().id()]++;
253  n_flux_faces[coarse_context->get_elem().id()] +=
255  }
256  }
257  }
258  }
259  else if (integrate_boundary_sides)
260  {
261  fine_context->pre_fe_reinit(system, parent);
262  libmesh_assert_equal_to
263  (fine_context->get_elem_solution().size(),
264  Uparent.size());
265  fine_context->get_elem_solution() = Uparent;
266  fine_context->side = cast_int<unsigned char>(n_p);
267  fine_context->side_fe_reinit();
268 
269  // If we find a boundary flux for any variable,
270  // let's just count it as a flux face for all
271  // variables. Otherwise we'd need to keep track of
272  // a separate n_flux_faces and error_per_cell for
273  // every single var.
274  bool found_boundary_flux = false;
275 
276  for (var=0; var<n_vars; var++)
277  if (error_norm.weight(var) != 0.0 &&
278  system.variable_type(var).family != SCALAR)
279  {
280  if (this->boundary_side_integration())
281  {
282  error_per_cell[fine_context->get_elem().id()] +=
283  static_cast<ErrorVectorReal>(fine_error);
284  found_boundary_flux = true;
285  }
286  }
287 
288  if (scale_by_n_flux_faces && found_boundary_flux)
289  n_flux_faces[fine_context->get_elem().id()]++;
290  }
291  }
292  }
293 #endif // #ifdef LIBMESH_ENABLE_AMR
294 
295  // If we do any more flux integration, e will be the fine element
296  fine_context->pre_fe_reinit(system, e);
297 
298  // Loop over the neighbors of element e
299  for (auto n_e : e->side_index_range())
300  {
301  if ((e->neighbor_ptr(n_e) != nullptr) ||
303  {
304  fine_context->side = cast_int<unsigned char>(n_e);
305  fine_context->side_fe_reinit();
306  }
307 
308  if (e->neighbor_ptr(n_e) != nullptr) // e is not on the boundary
309  {
310  const Elem * f = e->neighbor_ptr(n_e);
311  const dof_id_type f_id = f->id();
312 
313  // Compute flux jumps if we are in case 1 or case 2.
314  if ((f->active() && (f->level() == e->level()) && (e_id < f_id))
315  || (f->level() < e->level()))
316  {
317  // f is now the coarse element
318  coarse_context->pre_fe_reinit(system, f);
319 
320  this->reinit_sides();
321 
322  // Loop over all significant variables in the system
323  for (var=0; var<n_vars; var++)
324  if (error_norm.weight(var) != 0.0 &&
325  system.variable_type(var).family != SCALAR)
326  {
328 
329  error_per_cell[fine_context->get_elem().id()] +=
330  static_cast<ErrorVectorReal>(fine_error);
331  error_per_cell[coarse_context->get_elem().id()] +=
332  static_cast<ErrorVectorReal>(coarse_error);
333  }
334 
335  // Keep track of the number of internal flux
336  // sides found on each element
338  {
339  n_flux_faces[fine_context->get_elem().id()]++;
340  n_flux_faces[coarse_context->get_elem().id()] +=
342  }
343  } // end if (case1 || case2)
344  } // if (e->neighbor(n_e) != nullptr)
345 
346  // Otherwise, e is on the boundary. If it happens to
347  // be on a Dirichlet boundary, we need not do anything.
348  // On the other hand, if e is on a Neumann (flux) boundary
349  // with grad(u).n = g, we need to compute the additional residual
350  // (h * \int |g - grad(u_h).n|^2 dS)^(1/2).
351  // We can only do this with some knowledge of the boundary
352  // conditions, i.e. the user must have attached an appropriate
353  // BC function.
354  else if (integrate_boundary_sides)
355  {
356  bool found_boundary_flux = false;
357 
358  for (var=0; var<n_vars; var++)
359  if (error_norm.weight(var) != 0.0 &&
360  system.variable_type(var).family != SCALAR)
361  if (this->boundary_side_integration())
362  {
363  error_per_cell[fine_context->get_elem().id()] +=
364  static_cast<ErrorVectorReal>(fine_error);
365  found_boundary_flux = true;
366  }
367 
368  if (scale_by_n_flux_faces && found_boundary_flux)
369  n_flux_faces[fine_context->get_elem().id()]++;
370  } // end if (e->neighbor_ptr(n_e) == nullptr)
371  } // end loop over neighbors
372  } // End loop over active local elements
373 
374 
375  // Each processor has now computed the error contributions
376  // for its local elements. We need to sum the vector
377  // and then take the square-root of each component. Note
378  // that we only need to sum if we are running on multiple
379  // processors, and we only need to take the square-root
380  // if the value is nonzero. There will in general be many
381  // zeros for the inactive elements.
382 
383  // First sum the vector of estimated error values
384  this->reduce_error(error_per_cell, system.comm());
385 
386  // Compute the square-root of each component.
387  for (auto i : index_range(error_per_cell))
388  if (error_per_cell[i] != 0.)
389  error_per_cell[i] = std::sqrt(error_per_cell[i]);
390 
391 
392  if (this->scale_by_n_flux_faces)
393  {
394  // Sum the vector of flux face counts
395  this->reduce_error(n_flux_faces, system.comm());
396 
397  // Sanity check: Make sure the number of flux faces is
398  // always an integer value
399 #ifdef DEBUG
400  for (const auto & val : n_flux_faces)
401  libmesh_assert_equal_to (val, static_cast<float>(static_cast<unsigned int>(val)));
402 #endif
403 
404  // Scale the error by the number of flux faces for each element
405  for (auto i : index_range(n_flux_faces))
406  {
407  if (n_flux_faces[i] == 0.0) // inactive or non-local element
408  continue;
409 
410  error_per_cell[i] /= static_cast<ErrorVectorReal>(n_flux_faces[i]);
411  }
412  }
413 
414  // If we used a non-standard solution before, now is the time to fix
415  // the current_local_solution
416  if (solution_vector && solution_vector != system.solution.get())
417  {
418  NumericVector<Number> * newsol =
419  const_cast<NumericVector<Number> *>(solution_vector);
420  System & sys = const_cast<System &>(system);
421  newsol->swap(*sys.solution);
422  sys.update();
423  }
424 }

References libMesh::Elem::active(), libMesh::ElemInternal::active_family_tree_by_neighbor(), libMesh::JumpErrorEstimator::boundary_side_integration(), libMesh::Elem::child_ref_range(), libMesh::JumpErrorEstimator::coarse_context, libMesh::JumpErrorEstimator::coarse_error, libMesh::JumpErrorEstimator::coarse_n_flux_faces_increment(), libMesh::FEGenericBase< OutputType >::coarsened_dof_values(), libMesh::ParallelObject::comm(), dim, libMesh::ErrorEstimator::error_norm, libMesh::System::extra_quadrature_order, libMesh::FEType::family, libMesh::JumpErrorEstimator::fine_context, libMesh::JumpErrorEstimator::fine_error, libMesh::System::get_dof_map(), libMesh::System::get_mesh(), libMesh::FEAbstract::get_xyz(), libMesh::DofObject::id(), libMesh::index_range(), libMesh::JumpErrorEstimator::init_context(), libMesh::JumpErrorEstimator::integrate_boundary_sides, libMesh::JumpErrorEstimator::internal_side_integration(), libMesh::Elem::level(), libMesh::libmesh_ignore(), mesh, n_vars, libMesh::System::n_vars(), libMesh::Elem::neighbor_ptr(), libMesh::Elem::parent(), libMesh::ErrorEstimator::reduce_error(), libMesh::JumpErrorEstimator::reinit_sides(), libMesh::SCALAR, libMesh::JumpErrorEstimator::scale_by_n_flux_faces, libMesh::Elem::side_index_range(), libMesh::DenseVector< T >::size(), libMesh::System::solution, std::sqrt(), libMesh::NumericVector< T >::swap(), libMesh::JumpErrorEstimator::use_unweighted_quadrature_rules, libMesh::JumpErrorEstimator::var, libMesh::System::variable_type(), and libMesh::SystemNorm::weight().

Referenced by assemble_and_solve(), and main().

◆ estimate_errors() [1/2]

void libMesh::ErrorEstimator::estimate_errors ( const EquationSystems equation_systems,
ErrorMap errors_per_cell,
const std::map< const System *, const NumericVector< Number > * > *  solution_vectors = nullptr,
bool  estimate_parent_error = false 
)
virtualinherited

This virtual function can be redefined in derived classes, but by default it calls estimate_error repeatedly to calculate the requested error vectors.

FIXME: This is a default implementation - derived classes should reimplement it for efficiency.

Currently this function ignores the error_norm.weight() values because it calculates each variable's error individually, unscaled.

The user selects which errors get computed by filling a map with error vectors: If errors_per_cell[&system][v] exists, it will be filled with the error values in variable v of system

Reimplemented in libMesh::UniformRefinementEstimator.

Definition at line 93 of file error_estimator.C.

97 {
98  SystemNorm old_error_norm = this->error_norm;
99 
100  // Find the requested error values from each system
101  for (auto s : IntRange<unsigned int>(0, equation_systems.n_systems()))
102  {
103  const System & sys = equation_systems.get_system(s);
104 
105  unsigned int n_vars = sys.n_vars();
106 
107  for (unsigned int v = 0; v != n_vars; ++v)
108  {
109  // Only fill in ErrorVectors the user asks for
110  if (errors_per_cell.find(std::make_pair(&sys, v)) ==
111  errors_per_cell.end())
112  continue;
113 
114  // Calculate error in only one variable
115  std::vector<Real> weights(n_vars, 0.0);
116  weights[v] = 1.0;
117  this->error_norm =
118  SystemNorm(std::vector<FEMNormType>(n_vars, old_error_norm.type(v)),
119  weights);
120 
121  const NumericVector<Number> * solution_vector = nullptr;
122  if (solution_vectors &&
123  solution_vectors->find(&sys) != solution_vectors->end())
124  solution_vector = solution_vectors->find(&sys)->second;
125 
126  this->estimate_error
127  (sys, *errors_per_cell[std::make_pair(&sys, v)],
128  solution_vector, estimate_parent_error);
129  }
130  }
131 
132  // Restore our old state before returning
133  this->error_norm = old_error_norm;
134 }

References libMesh::ErrorEstimator::error_norm, libMesh::ErrorEstimator::estimate_error(), libMesh::EquationSystems::get_system(), libMesh::EquationSystems::n_systems(), n_vars, libMesh::System::n_vars(), and libMesh::SystemNorm::type().

◆ estimate_errors() [2/2]

void libMesh::ErrorEstimator::estimate_errors ( const EquationSystems equation_systems,
ErrorVector error_per_cell,
const std::map< const System *, SystemNorm > &  error_norms,
const std::map< const System *, const NumericVector< Number > * > *  solution_vectors = nullptr,
bool  estimate_parent_error = false 
)
virtualinherited

This virtual function can be redefined in derived classes, but by default computes the sum of the error_per_cell for each system in the equation_systems.

Currently this function ignores the error_norm member variable, and uses the function argument error_norms instead.

This function is named estimate_errors instead of estimate_error because otherwise C++ can get confused.

Reimplemented in libMesh::UniformRefinementEstimator.

Definition at line 47 of file error_estimator.C.

52 {
53  SystemNorm old_error_norm = this->error_norm;
54 
55  // Sum the error values from each system
56  for (auto s : IntRange<unsigned int>(0, equation_systems.n_systems()))
57  {
58  ErrorVector system_error_per_cell;
59  const System & sys = equation_systems.get_system(s);
60  if (error_norms.find(&sys) == error_norms.end())
61  this->error_norm = old_error_norm;
62  else
63  this->error_norm = error_norms.find(&sys)->second;
64 
65  const NumericVector<Number> * solution_vector = nullptr;
66  if (solution_vectors &&
67  solution_vectors->find(&sys) != solution_vectors->end())
68  solution_vector = solution_vectors->find(&sys)->second;
69 
70  this->estimate_error(sys, system_error_per_cell,
71  solution_vector, estimate_parent_error);
72 
73  if (s)
74  {
75  libmesh_assert_equal_to (error_per_cell.size(), system_error_per_cell.size());
76  for (auto i : index_range(error_per_cell))
77  error_per_cell[i] += system_error_per_cell[i];
78  }
79  else
80  error_per_cell = system_error_per_cell;
81  }
82 
83  // Restore our old state before returning
84  this->error_norm = old_error_norm;
85 }

References libMesh::ErrorEstimator::error_norm, libMesh::ErrorEstimator::estimate_error(), libMesh::EquationSystems::get_system(), libMesh::index_range(), and libMesh::EquationSystems::n_systems().

◆ init_context()

void libMesh::DiscontinuityMeasure::init_context ( FEMContext c)
overrideprotectedvirtual

An initialization function, for requesting specific data from the FE objects.

Reimplemented from libMesh::JumpErrorEstimator.

Definition at line 59 of file discontinuity_measure.C.

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 }

References dim, libMesh::FEMContext::elem_dimensions(), libMesh::ErrorEstimator::error_norm, libMesh::JumpErrorEstimator::fine_context, libMesh::FEGenericBase< OutputType >::get_phi(), libMesh::DiffContext::n_vars(), n_vars, and libMesh::SystemNorm::weight().

◆ internal_side_integration()

void libMesh::DiscontinuityMeasure::internal_side_integration ( )
overrideprotectedvirtual

The function which calculates a normal derivative jump based error term on an internal side.

Implements libMesh::JumpErrorEstimator.

Definition at line 86 of file discontinuity_measure.C.

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 }

References libMesh::JumpErrorEstimator::coarse_context, libMesh::JumpErrorEstimator::coarse_error, libMesh::Elem::dim(), libMesh::ErrorEstimator::error_norm, libMesh::JumpErrorEstimator::fine_context, libMesh::JumpErrorEstimator::fine_error, libMesh::FEAbstract::get_JxW(), libMesh::FEGenericBase< OutputType >::get_phi(), libMesh::Elem::hmax(), libMesh::FEAbstract::n_quadrature_points(), libMesh::TensorTools::norm_sq(), libMesh::Real, libMesh::JumpErrorEstimator::var, and libMesh::SystemNorm::weight().

◆ operator=() [1/2]

DiscontinuityMeasure& libMesh::DiscontinuityMeasure::operator= ( const DiscontinuityMeasure )
delete

◆ operator=() [2/2]

DiscontinuityMeasure& libMesh::DiscontinuityMeasure::operator= ( DiscontinuityMeasure &&  )
default

◆ reduce_error()

void libMesh::ErrorEstimator::reduce_error ( std::vector< ErrorVectorReal > &  error_per_cell,
const Parallel::Communicator &  comm 
) const
protectedinherited

This method takes the local error contributions in error_per_cell from each processor and combines them to get the global error vector.

Definition at line 32 of file error_estimator.C.

34 {
35  // This function must be run on all processors at once
36  // parallel_object_only();
37 
38  // Each processor has now computed the error contributions
39  // for its local elements. We may need to sum the vector to
40  // recover the error for each element.
41 
42  comm.sum(error_per_cell);
43 }

Referenced by libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::WeightedPatchRecoveryErrorEstimator::estimate_error(), libMesh::PatchRecoveryErrorEstimator::estimate_error(), libMesh::JumpErrorEstimator::estimate_error(), libMesh::AdjointRefinementEstimator::estimate_error(), and libMesh::ExactErrorEstimator::estimate_error().

◆ reinit_sides()

void libMesh::JumpErrorEstimator::reinit_sides ( )
protectedinherited

A utility function to reinit the finite element data on elements sharing a side.

Definition at line 429 of file jump_error_estimator.C.

430 {
431  fine_context->side_fe_reinit();
432 
433  unsigned short dim = fine_context->get_elem().dim();
434  libmesh_assert_equal_to(dim, coarse_context->get_elem().dim());
435 
436  FEBase * fe_fine = nullptr;
437  fine_context->get_side_fe( 0, fe_fine, dim );
438 
439  // Get the physical locations of the fine element quadrature points
440  std::vector<Point> qface_point = fe_fine->get_xyz();
441 
442  // Find the master quadrature point locations on the coarse element
443  FEBase * fe_coarse = nullptr;
444  coarse_context->get_side_fe( 0, fe_coarse, dim );
445 
446  std::vector<Point> qp_coarse;
447 
448  FEMap::inverse_map (coarse_context->get_elem().dim(),
449  &coarse_context->get_elem(), qface_point,
450  qp_coarse);
451 
452  // The number of variables in the system
453  const unsigned int n_vars = fine_context->n_vars();
454 
455  // Calculate all coarse element shape functions at those locations
456  for (unsigned int v=0; v<n_vars; v++)
457  if (error_norm.weight(v) != 0.0 &&
458  fine_context->get_system().variable_type(v).family != SCALAR)
459  {
460  coarse_context->get_side_fe( v, fe_coarse, dim );
461  fe_coarse->reinit (&coarse_context->get_elem(), &qp_coarse);
462  }
463 }

References libMesh::JumpErrorEstimator::coarse_context, dim, libMesh::ErrorEstimator::error_norm, libMesh::JumpErrorEstimator::fine_context, libMesh::FEMap::inverse_map(), n_vars, libMesh::FEAbstract::reinit(), libMesh::SCALAR, and libMesh::SystemNorm::weight().

Referenced by libMesh::JumpErrorEstimator::estimate_error().

◆ type()

ErrorEstimatorType libMesh::DiscontinuityMeasure::type ( ) const
overridevirtual
Returns
The type for the ErrorEstimator subclass.

Implements libMesh::ErrorEstimator.

Definition at line 51 of file discontinuity_measure.C.

52 {
53  return DISCONTINUITY_MEASURE;
54 }

References libMesh::DISCONTINUITY_MEASURE.

Member Data Documentation

◆ _bc_function

std::pair<bool,Real>(* libMesh::DiscontinuityMeasure::_bc_function) (const System &system, const Point &p, const std::string &var_name)
protected

Pointer to function that provides BC information.

Definition at line 108 of file discontinuity_measure.h.

Referenced by attach_essential_bc_function(), and boundary_side_integration().

◆ coarse_context

std::unique_ptr<FEMContext> libMesh::JumpErrorEstimator::coarse_context
protectedinherited

◆ coarse_error

Real libMesh::JumpErrorEstimator::coarse_error
protectedinherited

◆ error_norm

SystemNorm libMesh::ErrorEstimator::error_norm
inherited

When estimating the error in a single system, the error_norm is used to control the scaling and norm choice for each variable.

Not all estimators will support all norm choices. The default scaling is for all variables to be weighted equally. The default norm choice depends on the error estimator.

Part of this functionality was supported via component_scale and sobolev_order in older libMesh versions, and a small part was supported via component_mask in even older versions. Hopefully the encapsulation here will allow us to avoid changing this API again.

Definition at line 164 of file error_estimator.h.

Referenced by libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::AdjointRefinementEstimator::AdjointRefinementEstimator(), boundary_side_integration(), libMesh::KellyErrorEstimator::boundary_side_integration(), DiscontinuityMeasure(), libMesh::JumpErrorEstimator::estimate_error(), libMesh::AdjointResidualErrorEstimator::estimate_error(), libMesh::ExactErrorEstimator::estimate_error(), libMesh::ErrorEstimator::estimate_errors(), libMesh::ExactErrorEstimator::ExactErrorEstimator(), libMesh::ExactErrorEstimator::find_squared_element_error(), libMesh::LaplacianErrorEstimator::init_context(), init_context(), libMesh::KellyErrorEstimator::init_context(), libMesh::LaplacianErrorEstimator::internal_side_integration(), internal_side_integration(), libMesh::KellyErrorEstimator::internal_side_integration(), libMesh::KellyErrorEstimator::KellyErrorEstimator(), libMesh::LaplacianErrorEstimator::LaplacianErrorEstimator(), main(), libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::PatchRecoveryErrorEstimator::PatchRecoveryErrorEstimator(), libMesh::JumpErrorEstimator::reinit_sides(), and libMesh::UniformRefinementEstimator::UniformRefinementEstimator().

◆ fine_context

std::unique_ptr<FEMContext> libMesh::JumpErrorEstimator::fine_context
protectedinherited

◆ fine_error

Real libMesh::JumpErrorEstimator::fine_error
protectedinherited

◆ integrate_boundary_sides

bool libMesh::JumpErrorEstimator::integrate_boundary_sides
protectedinherited

A boolean flag, by default false, to be set to true if integrations with boundary_side_integration() should be performed.

Definition at line 152 of file jump_error_estimator.h.

Referenced by attach_essential_bc_function(), libMesh::KellyErrorEstimator::attach_flux_bc_function(), and libMesh::JumpErrorEstimator::estimate_error().

◆ scale_by_n_flux_faces

bool libMesh::JumpErrorEstimator::scale_by_n_flux_faces
inherited

This boolean flag allows you to scale the error indicator result for each element by the number of "flux faces" the element actually has.

This tends to weight more evenly cells which are on the boundaries and thus have fewer contributions to their flux. The value is initialized to false, simply set it to true if you want to use the feature.

Definition at line 99 of file jump_error_estimator.h.

Referenced by libMesh::JumpErrorEstimator::estimate_error().

◆ use_unweighted_quadrature_rules

bool libMesh::JumpErrorEstimator::use_unweighted_quadrature_rules
inherited

This boolean flag allows you to use "unweighted" quadrature rules (sized to exactly integrate unweighted shape functions in master element space) rather than "default" quadrature rules (sized to exactly integrate polynomials of one higher degree than mass matrix terms).

The results with the former, lower-order rules will be somewhat less accurate in many cases but will be much cheaper to compute.

The value is initialized to false, simply set it to true if you want to use the feature.

Definition at line 113 of file jump_error_estimator.h.

Referenced by libMesh::JumpErrorEstimator::estimate_error(), and main().

◆ var

unsigned int libMesh::JumpErrorEstimator::var
protectedinherited

The documentation for this class was generated from the following files:
libMesh::JumpErrorEstimator::reinit_sides
void reinit_sides()
A utility function to reinit the finite element data on elements sharing a side.
Definition: jump_error_estimator.C:429
libMesh::JumpErrorEstimator::coarse_context
std::unique_ptr< FEMContext > coarse_context
Definition: jump_error_estimator.h:158
libMesh::dof_id_type
uint8_t dof_id_type
Definition: id_types.h:67
libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::JumpErrorEstimator::integrate_boundary_sides
bool integrate_boundary_sides
A boolean flag, by default false, to be set to true if integrations with boundary_side_integration() ...
Definition: jump_error_estimator.h:152
libMesh::JumpErrorEstimator::init_context
virtual void init_context(FEMContext &c)
An initialization function, to give derived classes a chance to request specific data from the FE obj...
Definition: jump_error_estimator.C:47
libMesh::JumpErrorEstimator::JumpErrorEstimator
JumpErrorEstimator()
Constructor.
Definition: jump_error_estimator.h:55
libMesh::ErrorEstimator::error_norm
SystemNorm error_norm
When estimating the error in a single system, the error_norm is used to control the scaling and norm ...
Definition: error_estimator.h:164
n_vars
unsigned int n_vars
Definition: adaptivity_ex3.C:116
libMesh::index_range
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:106
libMesh::JumpErrorEstimator::var
unsigned int var
The variable number currently being evaluated.
Definition: jump_error_estimator.h:168
libMesh::JumpErrorEstimator::fine_error
Real fine_error
The fine and coarse error values to be set by each side_integration();.
Definition: jump_error_estimator.h:163
std::sqrt
MetaPhysicL::DualNumber< T, D > sqrt(const MetaPhysicL::DualNumber< T, D > &in)
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::SystemNorm::weight
Real weight(unsigned int var) const
Definition: system_norm.C:133
libMesh::JumpErrorEstimator::scale_by_n_flux_faces
bool scale_by_n_flux_faces
This boolean flag allows you to scale the error indicator result for each element by the number of "f...
Definition: jump_error_estimator.h:99
libMesh::JumpErrorEstimator::coarse_n_flux_faces_increment
float coarse_n_flux_faces_increment()
A utility function to correctly increase n_flux_faces for the coarse element.
Definition: jump_error_estimator.C:467
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::JumpErrorEstimator::boundary_side_integration
virtual bool boundary_side_integration()
The function, to be implemented by derived classes, which calculates an error term on a boundary side...
Definition: jump_error_estimator.h:146
libMesh::DISCONTINUITY_MEASURE
Definition: enum_error_estimator_type.h:37
libMesh::JumpErrorEstimator::coarse_error
Real coarse_error
Definition: jump_error_estimator.h:163
libMesh::ErrorEstimator::estimate_error
virtual void estimate_error(const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=nullptr, bool estimate_parent_error=false)=0
This pure virtual function must be redefined in derived classes to compute the error for each active ...
libMesh::FEMap::inverse_map
static Point inverse_map(const unsigned int dim, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
Definition: fe_map.C:1622
libMesh::JumpErrorEstimator::fine_context
std::unique_ptr< FEMContext > fine_context
Context objects for integrating on the fine and coarse elements sharing a face.
Definition: jump_error_estimator.h:158
libMesh::FEBase
FEGenericBase< Real > FEBase
Definition: exact_error_estimator.h:39
libMesh::libmesh_ignore
void libmesh_ignore(const Args &...)
Definition: libmesh_common.h:526
libMesh::FEGenericBase::coarsened_dof_values
static void coarsened_dof_values(const NumericVector< Number > &global_vector, const DofMap &dof_map, const Elem *coarse_elem, DenseVector< Number > &coarse_dofs, const unsigned int var, const bool use_old_dof_indices=false)
Creates a local projection on coarse_elem, based on the DoF values in global_vector for it's children...
Definition: fe_base.C:822
libMesh::ElemInternal::active_family_tree_by_neighbor
void active_family_tree_by_neighbor(T elem, std::vector< T > &family, T neighbor_in, bool reset=true)
Definition: elem_internal.h:320
libMesh::DiscontinuityMeasure::_bc_function
std::pair< bool, Real >(* _bc_function)(const System &system, const Point &p, const std::string &var_name)
Pointer to function that provides BC information.
Definition: discontinuity_measure.h:108
libMesh::JumpErrorEstimator::internal_side_integration
virtual void internal_side_integration()=0
The function, to be implemented by derived classes, which calculates an error term on an internal sid...
libMesh::TensorTools::norm_sq
T norm_sq(std::complex< T > a)
Definition: tensor_tools.h:85
libMesh::ErrorEstimator::reduce_error
void reduce_error(std::vector< ErrorVectorReal > &error_per_cell, const Parallel::Communicator &comm) const
This method takes the local error contributions in error_per_cell from each processor and combines th...
Definition: error_estimator.C:32
libMesh::L2
Definition: enum_norm_type.h:36
libMesh::JumpErrorEstimator::use_unweighted_quadrature_rules
bool use_unweighted_quadrature_rules
This boolean flag allows you to use "unweighted" quadrature rules (sized to exactly integrate unweigh...
Definition: jump_error_estimator.h:113
libMesh::SCALAR
Definition: enum_fe_family.h:58
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
fptr
Number fptr(const Point &p, const Parameters &, const std::string &libmesh_dbg_var(sys_name), const std::string &unknown_name)
Definition: projection.C:80