libMesh::ExactErrorEstimator Class Reference

#include <exact_error_estimator.h>

Inheritance diagram for libMesh::ExactErrorEstimator:

Public Types

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

Public Member Functions

 ExactErrorEstimator ()
 
 ~ExactErrorEstimator ()
 
void attach_exact_values (std::vector< FunctionBase< Number > * > f)
 
void attach_exact_value (unsigned int sys_num, FunctionBase< Number > *f)
 
void attach_exact_value (Number fptr(const Point &p, const Parameters &Parameters, const std::string &sys_name, const std::string &unknown_name))
 
void attach_exact_derivs (std::vector< FunctionBase< Gradient > * > g)
 
void attach_exact_deriv (unsigned int sys_num, FunctionBase< Gradient > *g)
 
void attach_exact_deriv (Gradient gptr(const Point &p, const Parameters &parameters, const std::string &sys_name, const std::string &unknown_name))
 
void attach_exact_hessians (std::vector< FunctionBase< Tensor > * > h)
 
void attach_exact_hessian (unsigned int sys_num, FunctionBase< Tensor > *h)
 
void attach_exact_hessian (Tensor hptr(const Point &p, const Parameters &parameters, const std::string &sys_name, const std::string &unknown_name))
 
void attach_reference_solution (EquationSystems *es_fine)
 
void extra_quadrature_order (const int extraorder)
 
virtual void estimate_error (const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=libmesh_nullptr, bool estimate_parent_error=false) libmesh_override
 
virtual ErrorEstimatorType type () const libmesh_override
 
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=libmesh_nullptr, bool estimate_parent_error=false)
 
virtual void estimate_errors (const EquationSystems &equation_systems, ErrorMap &errors_per_cell, const std::map< const System *, const NumericVector< Number > * > *solution_vectors=libmesh_nullptr, bool estimate_parent_error=false)
 

Public Attributes

SystemNorm error_norm
 

Protected Member Functions

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

Private Member Functions

Real find_squared_element_error (const System &system, const std::string &var_name, const Elem *elem, const DenseVector< Number > &Uelem, FEBase *fe, MeshFunction *fine_values) const
 
void clear_functors ()
 

Private Attributes

Number(* _exact_value )(const Point &p, const Parameters &parameters, const std::string &sys_name, const std::string &unknown_name)
 
Gradient(* _exact_deriv )(const Point &p, const Parameters &parameters, const std::string &sys_name, const std::string &unknown_name)
 
Tensor(* _exact_hessian )(const Point &p, const Parameters &parameters, const std::string &sys_name, const std::string &unknown_name)
 
std::vector< FunctionBase< Number > * > _exact_values
 
std::vector< FunctionBase< Gradient > * > _exact_derivs
 
std::vector< FunctionBase< Tensor > * > _exact_hessians
 
EquationSystems_equation_systems_fine
 
int _extra_order
 

Detailed Description

This class implements an "error estimator" based on the difference between the approximate and exact solution. In theory the quadrature error in this estimate should be much lower than the approximation error in other estimates, so this estimator can be used to calculate effectivity.

Author
Roy Stogner
Date
2006

Definition at line 67 of file exact_error_estimator.h.

Member Typedef Documentation

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 113 of file error_estimator.h.

Constructor & Destructor Documentation

libMesh::ExactErrorEstimator::ExactErrorEstimator ( )
inline

Constructor. Responsible for initializing the _bc_function function pointer to libmesh_nullptr, and defaulting the norm type to H1.

Definition at line 75 of file exact_error_estimator.h.

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

75  :
81  _extra_order(0)
82  { error_norm = H1; }
Number(* _exact_value)(const Point &p, const Parameters &parameters, const std::string &sys_name, const std::string &unknown_name)
const class libmesh_nullptr_t libmesh_nullptr
Tensor(* _exact_hessian)(const Point &p, const Parameters &parameters, const std::string &sys_name, const std::string &unknown_name)
Gradient(* _exact_deriv)(const Point &p, const Parameters &parameters, const std::string &sys_name, const std::string &unknown_name)
libMesh::ExactErrorEstimator::~ExactErrorEstimator ( )
inline

Member Function Documentation

void libMesh::ExactErrorEstimator::attach_exact_deriv ( unsigned int  sys_num,
FunctionBase< Gradient > *  g 
)

Clone and attach an arbitrary functor which computes the exact gradient of the system sys_num solution at any point.

Definition at line 127 of file exact_error_estimator.C.

References _equation_systems_fine, _exact_derivs, _exact_hessian, attach_exact_hessian(), clear_functors(), libMesh::FunctionBase< Output >::clone(), libMesh::libmesh_assert(), and libmesh_nullptr.

Referenced by attach_exact_value(), and ~ExactErrorEstimator().

129 {
130  if (_exact_derivs.size() <= sys_num)
131  _exact_derivs.resize(sys_num+1, libmesh_nullptr);
132 
133  if (g)
134  _exact_derivs[sys_num] = g->clone().release();
135 }
std::vector< FunctionBase< Gradient > * > _exact_derivs
const class libmesh_nullptr_t libmesh_nullptr
void libMesh::ExactErrorEstimator::attach_exact_deriv ( Gradient   gptrconst Point &p, const Parameters &parameters, const std::string &sys_name, const std::string &unknown_name)

Attach an arbitrary function which computes the exact gradient of the solution at any point.

void libMesh::ExactErrorEstimator::attach_exact_derivs ( std::vector< FunctionBase< Gradient > * >  g)

Clone and attach arbitrary functors which compute the exact gradients of the EquationSystems' solutions at any point.

Definition at line 108 of file exact_error_estimator.C.

References _exact_derivs, and libmesh_nullptr.

Referenced by ~ExactErrorEstimator().

109 {
110  // Clear out any previous _exact_derivs entries, then add a new
111  // entry for each system.
112  for (std::size_t i=0; i != _exact_derivs.size(); ++i)
113  delete (_exact_derivs[i]);
114 
115  _exact_derivs.clear();
116  _exact_derivs.resize(g.size(), libmesh_nullptr);
117 
118  // We use clone() to get non-sliced copies of FunctionBase
119  // subclasses, but we don't currently put the resulting UniquePtrs
120  // into an STL container.
121  for (std::size_t i=0; i != g.size(); ++i)
122  if (g[i])
123  _exact_derivs[i] = g[i]->clone().release();
124 }
std::vector< FunctionBase< Gradient > * > _exact_derivs
const class libmesh_nullptr_t libmesh_nullptr
void libMesh::ExactErrorEstimator::attach_exact_hessian ( unsigned int  sys_num,
FunctionBase< Tensor > *  h 
)

Clone and attach an arbitrary functor which computes the exact second derivatives of the system sys_num solution at any point.

Definition at line 175 of file exact_error_estimator.C.

References _exact_hessians, libMesh::FunctionBase< Output >::clone(), and libmesh_nullptr.

Referenced by attach_exact_deriv(), and ~ExactErrorEstimator().

177 {
178  if (_exact_hessians.size() <= sys_num)
179  _exact_hessians.resize(sys_num+1, libmesh_nullptr);
180 
181  if (h)
182  _exact_hessians[sys_num] = h->clone().release();
183 }
std::vector< FunctionBase< Tensor > * > _exact_hessians
const class libmesh_nullptr_t libmesh_nullptr
void libMesh::ExactErrorEstimator::attach_exact_hessian ( Tensor   hptrconst Point &p, const Parameters &parameters, const std::string &sys_name, const std::string &unknown_name)

Attach an arbitrary function which computes the exact second derivatives of the solution at any point.

void libMesh::ExactErrorEstimator::attach_exact_hessians ( std::vector< FunctionBase< Tensor > * >  h)

Clone and attach arbitrary functors which compute the exact second derivatives of the EquationSystems' solutions at any point.

Definition at line 156 of file exact_error_estimator.C.

References _exact_hessians, and libmesh_nullptr.

Referenced by ~ExactErrorEstimator().

157 {
158  // Clear out any previous _exact_hessians entries, then add a new
159  // entry for each system.
160  for (std::size_t i=0; i != _exact_hessians.size(); ++i)
161  delete (_exact_hessians[i]);
162 
163  _exact_hessians.clear();
164  _exact_hessians.resize(h.size(), libmesh_nullptr);
165 
166  // We use clone() to get non-sliced copies of FunctionBase
167  // subclasses, but we don't currently put the resulting UniquePtrs
168  // into an STL container.
169  for (std::size_t i=0; i != h.size(); ++i)
170  if (h[i])
171  _exact_hessians[i] = h[i]->clone().release();
172 }
std::vector< FunctionBase< Tensor > * > _exact_hessians
const class libmesh_nullptr_t libmesh_nullptr
void libMesh::ExactErrorEstimator::attach_exact_value ( unsigned int  sys_num,
FunctionBase< Number > *  f 
)

Clone and attach an arbitrary functor which computes the exact value of the system sys_num solution at any point.

Definition at line 81 of file exact_error_estimator.C.

References _equation_systems_fine, _exact_deriv, _exact_values, attach_exact_deriv(), clear_functors(), libMesh::FunctionBase< Output >::clone(), libMesh::libmesh_assert(), and libmesh_nullptr.

Referenced by ~ExactErrorEstimator().

83 {
84  if (_exact_values.size() <= sys_num)
85  _exact_values.resize(sys_num+1, libmesh_nullptr);
86 
87  if (f)
88  _exact_values[sys_num] = f->clone().release();
89 }
std::vector< FunctionBase< Number > * > _exact_values
const class libmesh_nullptr_t libmesh_nullptr
virtual UniquePtr< FunctionBase< Output > > clone() const =0
void libMesh::ExactErrorEstimator::attach_exact_value ( Number   fptrconst Point &p, const Parameters &Parameters, const std::string &sys_name, const std::string &unknown_name)

Attach an arbitrary function which computes the exact value of the solution at any point.

void libMesh::ExactErrorEstimator::attach_exact_values ( std::vector< FunctionBase< Number > * >  f)

Clone and attach arbitrary functors which compute the exact values of the EquationSystems' solutions at any point.

Definition at line 62 of file exact_error_estimator.C.

References _exact_values, and libmesh_nullptr.

Referenced by ~ExactErrorEstimator().

63 {
64  // Clear out any previous _exact_values entries, then add a new
65  // entry for each system.
66  for (std::size_t i=0; i != _exact_values.size(); ++i)
67  delete (_exact_values[i]);
68 
69  _exact_values.clear();
70  _exact_values.resize(f.size(), libmesh_nullptr);
71 
72  // We use clone() to get non-sliced copies of FunctionBase
73  // subclasses, but we don't currently put the resulting UniquePtrs
74  // into an STL container.
75  for (std::size_t i=0; i != f.size(); ++i)
76  if (f[i])
77  _exact_values[i] = f[i]->clone().release();
78 }
std::vector< FunctionBase< Number > * > _exact_values
const class libmesh_nullptr_t libmesh_nullptr
void libMesh::ExactErrorEstimator::attach_reference_solution ( EquationSystems es_fine)

Attach function similar to system.h which allows the user to attach a second EquationSystems object with a reference fine grid solution.

Definition at line 186 of file exact_error_estimator.C.

References _equation_systems_fine, _exact_deriv, _exact_hessian, _exact_value, clear_functors(), libMesh::libmesh_assert(), and libmesh_nullptr.

Referenced by ~ExactErrorEstimator().

187 {
188  libmesh_assert(es_fine);
189  _equation_systems_fine = es_fine;
190 
191  // If we're using a fine grid solution, we're not using exact value
192  // function pointers or functors.
196 
197  this->clear_functors();
198 }
Number(* _exact_value)(const Point &p, const Parameters &parameters, const std::string &sys_name, const std::string &unknown_name)
const class libmesh_nullptr_t libmesh_nullptr
libmesh_assert(j)
Tensor(* _exact_hessian)(const Point &p, const Parameters &parameters, const std::string &sys_name, const std::string &unknown_name)
Gradient(* _exact_deriv)(const Point &p, const Parameters &parameters, const std::string &sys_name, const std::string &unknown_name)
void libMesh::ExactErrorEstimator::clear_functors ( )
private

Helper method for cleanup

Definition at line 554 of file exact_error_estimator.C.

References _exact_derivs, _exact_hessians, and _exact_values.

Referenced by attach_exact_deriv(), attach_exact_value(), and attach_reference_solution().

555 {
556  // delete will clean up any cloned functors and no-op on any NULL
557  // pointers
558 
559  for (std::size_t i=0; i != _exact_values.size(); ++i)
560  delete (_exact_values[i]);
561 
562  for (std::size_t i=0; i != _exact_derivs.size(); ++i)
563  delete (_exact_derivs[i]);
564 
565  for (std::size_t i=0; i != _exact_hessians.size(); ++i)
566  delete (_exact_hessians[i]);
567 }
std::vector< FunctionBase< Number > * > _exact_values
std::vector< FunctionBase< Tensor > * > _exact_hessians
std::vector< FunctionBase< Gradient > * > _exact_derivs
void libMesh::ExactErrorEstimator::estimate_error ( const System system,
ErrorVector error_per_cell,
const NumericVector< Number > *  solution_vector = libmesh_nullptr,
bool  estimate_parent_error = false 
)
virtual

This function uses the exact solution function to estimate the error on each cell. The estimated error is output in the vector error_per_cell

Implements libMesh::ErrorEstimator.

Definition at line 200 of file exact_error_estimator.C.

References _equation_systems_fine, _exact_derivs, _exact_hessians, _exact_values, _extra_order, libMesh::Elem::active(), libMesh::MeshBase::active_local_elements_begin(), libMesh::MeshBase::active_local_elements_end(), libMesh::FEGenericBase< OutputType >::build(), libMesh::NumericVector< T >::build(), libMesh::Elem::child_ptr(), libMesh::FEGenericBase< OutputType >::coarsened_dof_values(), libMesh::ParallelObject::comm(), libMesh::System::current_local_solution, libMesh::System::current_solution(), libMesh::FEType::default_quadrature_rule(), libMesh::DofMap::dof_indices(), libMesh::ErrorEstimator::error_norm, libMesh::ErrorVectorReal, find_squared_element_error(), libMesh::System::get_dof_map(), libMesh::System::get_mesh(), libMesh::EquationSystems::get_system(), libMesh::DofObject::id(), libMesh::libmesh_ignore(), libMesh::MeshBase::max_elem_id(), mesh, libMesh::MeshBase::mesh_dimension(), libMesh::Elem::n_children(), n_vars, libMesh::System::n_vars(), libMesh::System::name(), libMesh::Elem::parent(), libMesh::ErrorEstimator::reduce_error(), libMesh::SERIAL, libMesh::System::solution, libMesh::NumericVector< T >::swap(), libMesh::sys, libMesh::System::update_global_solution(), libMesh::System::variable_name(), libMesh::System::variable_number(), libMesh::DofMap::variable_type(), and libMesh::SystemNorm::weight().

Referenced by extra_quadrature_order().

204 {
205  // Ignore the fact that this variable is unused when !LIBMESH_ENABLE_AMR
206  libmesh_ignore(estimate_parent_error);
207 
208  // The current mesh
209  const MeshBase & mesh = system.get_mesh();
210 
211  // The dimensionality of the mesh
212  const unsigned int dim = mesh.mesh_dimension();
213 
214  // The number of variables in the system
215  const unsigned int n_vars = system.n_vars();
216 
217  // The DofMap for this system
218  const DofMap & dof_map = system.get_dof_map();
219 
220  // Resize the error_per_cell vector to be
221  // the number of elements, initialize it to 0.
222  error_per_cell.resize (mesh.max_elem_id());
223  std::fill (error_per_cell.begin(), error_per_cell.end(), 0.);
224 
225  // Prepare current_local_solution to localize a non-standard
226  // solution vector if necessary
227  if (solution_vector && solution_vector != system.solution.get())
228  {
229  NumericVector<Number> * newsol =
230  const_cast<NumericVector<Number> *>(solution_vector);
231  System & sys = const_cast<System &>(system);
232  newsol->swap(*sys.solution);
233  sys.update();
234  }
235 
236  // Loop over all the variables in the system
237  for (unsigned int var=0; var<n_vars; var++)
238  {
239  // Possibly skip this variable
240  if (error_norm.weight(var) == 0.0) continue;
241 
242  // The (string) name of this variable
243  const std::string & var_name = system.variable_name(var);
244 
245  // The type of finite element to use for this variable
246  const FEType & fe_type = dof_map.variable_type (var);
247 
248  UniquePtr<FEBase> fe (FEBase::build (dim, fe_type));
249 
250  // Build an appropriate Gaussian quadrature rule
251  UniquePtr<QBase> qrule =
252  fe_type.default_quadrature_rule (dim,
253  _extra_order);
254 
255  fe->attach_quadrature_rule (qrule.get());
256 
257  // Prepare a global solution and a MeshFunction of the fine system if we need one
258  UniquePtr<MeshFunction> fine_values;
259  UniquePtr<NumericVector<Number> > fine_soln = NumericVector<Number>::build(system.comm());
261  {
262  const System & fine_system = _equation_systems_fine->get_system(system.name());
263 
264  std::vector<Number> global_soln;
265  // FIXME - we're assuming that the fine system solution gets
266  // used even when a different vector is used for the coarse
267  // system
268  fine_system.update_global_solution(global_soln);
269  fine_soln->init
270  (cast_int<numeric_index_type>(global_soln.size()), true,
271  SERIAL);
272  (*fine_soln) = global_soln;
273 
274  fine_values = UniquePtr<MeshFunction>
275  (new MeshFunction(*_equation_systems_fine,
276  *fine_soln,
277  fine_system.get_dof_map(),
278  fine_system.variable_number(var_name)));
279  fine_values->init();
280  } else {
281  // Initialize functors if we're using them
282  for (std::size_t i=0; i != _exact_values.size(); ++i)
283  if (_exact_values[i])
284  _exact_values[i]->init();
285 
286  for (std::size_t i=0; i != _exact_derivs.size(); ++i)
287  if (_exact_derivs[i])
288  _exact_derivs[i]->init();
289 
290  for (std::size_t i=0; i != _exact_hessians.size(); ++i)
291  if (_exact_hessians[i])
292  _exact_hessians[i]->init();
293  }
294 
295  // Request the data we'll need to compute with
296  fe->get_JxW();
297  fe->get_phi();
298  fe->get_dphi();
299 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
300  fe->get_d2phi();
301 #endif
302  fe->get_xyz();
303 
304 #ifdef LIBMESH_ENABLE_AMR
305  // If we compute on parent elements, we'll want to do so only
306  // once on each, so we need to keep track of which we've done.
307  std::vector<bool> computed_var_on_parent;
308 
309  if (estimate_parent_error)
310  computed_var_on_parent.resize(error_per_cell.size(), false);
311 #endif
312 
313  // TODO: this ought to be threaded (and using subordinate
314  // MeshFunction objects in each thread rather than a single
315  // master)
316 
317  // Iterate over all the active elements in the mesh
318  // that live on this processor.
319  MeshBase::const_element_iterator
320  elem_it = mesh.active_local_elements_begin();
321  const MeshBase::const_element_iterator
322  elem_end = mesh.active_local_elements_end();
323 
324  for (;elem_it != elem_end; ++elem_it)
325  {
326  // e is necessarily an active element on the local processor
327  const Elem * elem = *elem_it;
328  const dof_id_type e_id = elem->id();
329 
330 #ifdef LIBMESH_ENABLE_AMR
331  // See if the parent of element e has been examined yet;
332  // if not, we may want to compute the estimator on it
333  const Elem * parent = elem->parent();
334 
335  // We only can compute and only need to compute on
336  // parents with all active children
337  bool compute_on_parent = true;
338  if (!parent || !estimate_parent_error)
339  compute_on_parent = false;
340  else
341  for (unsigned int c=0; c != parent->n_children(); ++c)
342  if (!parent->child_ptr(c)->active())
343  compute_on_parent = false;
344 
345  if (compute_on_parent &&
346  !computed_var_on_parent[parent->id()])
347  {
348  computed_var_on_parent[parent->id()] = true;
349 
350  // Compute a projection onto the parent
351  DenseVector<Number> Uparent;
352  FEBase::coarsened_dof_values(*(system.current_local_solution),
353  dof_map, parent, Uparent,
354  var, false);
355 
356  error_per_cell[parent->id()] +=
357  static_cast<ErrorVectorReal>
358  (find_squared_element_error(system, var_name,
359  parent, Uparent,
360  fe.get(),
361  fine_values.get()));
362  }
363 #endif
364 
365  // Get the local to global degree of freedom maps
366  std::vector<dof_id_type> dof_indices;
367  dof_map.dof_indices (elem, dof_indices, var);
368  const unsigned int n_dofs =
369  cast_int<unsigned int>(dof_indices.size());
370  DenseVector<Number> Uelem(n_dofs);
371  for (unsigned int i=0; i != n_dofs; ++i)
372  Uelem(i) = system.current_solution(dof_indices[i]);
373 
374  error_per_cell[e_id] +=
375  static_cast<ErrorVectorReal>
376  (find_squared_element_error(system, var_name, elem,
377  Uelem, fe.get(),
378  fine_values.get()));
379 
380  } // End loop over active local elements
381  } // End loop over variables
382 
383 
384 
385  // Each processor has now computed the error contributions
386  // for its local elements. We need to sum the vector
387  // and then take the square-root of each component. Note
388  // that we only need to sum if we are running on multiple
389  // processors, and we only need to take the square-root
390  // if the value is nonzero. There will in general be many
391  // zeros for the inactive elements.
392 
393  // First sum the vector of estimated error values
394  this->reduce_error(error_per_cell, system.comm());
395 
396  // Compute the square-root of each component.
397  {
398  LOG_SCOPE("std::sqrt()", "ExactErrorEstimator");
399  for (std::size_t i=0; i<error_per_cell.size(); i++)
400  if (error_per_cell[i] != 0.)
401  {
402  libmesh_assert_greater (error_per_cell[i], 0.);
403  error_per_cell[i] = std::sqrt(error_per_cell[i]);
404  }
405  }
406 
407  // If we used a non-standard solution before, now is the time to fix
408  // the current_local_solution
409  if (solution_vector && solution_vector != system.solution.get())
410  {
411  NumericVector<Number> * newsol =
412  const_cast<NumericVector<Number> *>(solution_vector);
413  System & sys = const_cast<System &>(system);
414  newsol->swap(*sys.solution);
415  sys.update();
416  }
417 }
std::vector< FunctionBase< Number > * > _exact_values
ImplicitSystem & sys
Real weight(unsigned int var) const
Definition: system_norm.h:292
std::vector< FunctionBase< Tensor > * > _exact_hessians
std::vector< FunctionBase< Gradient > * > _exact_derivs
MeshBase & mesh
const unsigned int n_vars
Definition: tecplot_io.C:68
DIE A HORRIBLE DEATH HERE typedef float ErrorVectorReal
static UniquePtr< NumericVector< T > > build(const Parallel::Communicator &comm, const SolverPackage solver_package=libMesh::default_solver_package())
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)
Definition: fe_base.C:804
void reduce_error(std::vector< ErrorVectorReal > &error_per_cell, const Parallel::Communicator &comm LIBMESH_CAN_DEFAULT_TO_COMMWORLD) const
void libmesh_ignore(const T &)
virtual void update()
Definition: system.C:423
virtual void get(const std::vector< numeric_index_type > &index, T *values) const
const T_sys & get_system(const std::string &name) const
Real find_squared_element_error(const System &system, const std::string &var_name, const Elem *elem, const DenseVector< Number > &Uelem, FEBase *fe, MeshFunction *fine_values) const
uint8_t dof_id_type
Definition: id_types.h:64
static UniquePtr< FEGenericBase > build(const unsigned int dim, const FEType &type)
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 = libmesh_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 48 of file error_estimator.C.

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

Referenced by libMesh::ErrorEstimator::~ErrorEstimator().

53 {
54  SystemNorm old_error_norm = this->error_norm;
55 
56  // Sum the error values from each system
57  for (unsigned int s = 0; s != equation_systems.n_systems(); ++s)
58  {
59  ErrorVector system_error_per_cell;
60  const System & sys = equation_systems.get_system(s);
61  if (error_norms.find(&sys) == error_norms.end())
62  this->error_norm = old_error_norm;
63  else
64  this->error_norm = error_norms.find(&sys)->second;
65 
66  const NumericVector<Number> * solution_vector = libmesh_nullptr;
67  if (solution_vectors &&
68  solution_vectors->find(&sys) != solution_vectors->end())
69  solution_vector = solution_vectors->find(&sys)->second;
70 
71  this->estimate_error(sys, system_error_per_cell,
72  solution_vector, estimate_parent_error);
73 
74  if (s)
75  {
76  libmesh_assert_equal_to (error_per_cell.size(), system_error_per_cell.size());
77  for (std::size_t i=0; i != error_per_cell.size(); ++i)
78  error_per_cell[i] += system_error_per_cell[i];
79  }
80  else
81  error_per_cell = system_error_per_cell;
82  }
83 
84  // Restore our old state before returning
85  this->error_norm = old_error_norm;
86 }
ImplicitSystem & sys
const class libmesh_nullptr_t libmesh_nullptr
virtual void estimate_error(const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=libmesh_nullptr, bool estimate_parent_error=false)=0
void libMesh::ErrorEstimator::estimate_errors ( const EquationSystems equation_systems,
ErrorMap errors_per_cell,
const std::map< const System *, const NumericVector< Number > * > *  solution_vectors = libmesh_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.

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

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

Reimplemented in libMesh::UniformRefinementEstimator.

Definition at line 94 of file error_estimator.C.

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

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

Increases or decreases the order of the quadrature rule used for numerical integration.

Definition at line 167 of file exact_error_estimator.h.

References _extra_order, estimate_error(), and libmesh_nullptr.

168  { _extra_order = extraorder; }
Real libMesh::ExactErrorEstimator::find_squared_element_error ( const System system,
const std::string &  var_name,
const Elem elem,
const DenseVector< Number > &  Uelem,
FEBase fe,
MeshFunction fine_values 
) const
private

Helper method for calculating on each element

Definition at line 421 of file exact_error_estimator.C.

References _equation_systems_fine, _exact_deriv, _exact_derivs, _exact_hessian, _exact_hessians, _exact_value, _exact_values, libMesh::ErrorEstimator::error_norm, libMesh::FEGenericBase< OutputType >::get_d2phi(), libMesh::FEGenericBase< OutputType >::get_dphi(), libMesh::System::get_equation_systems(), libMesh::FEAbstract::get_JxW(), libMesh::FEGenericBase< OutputType >::get_phi(), libMesh::FEAbstract::get_xyz(), libMesh::MeshFunction::gradient(), libMesh::H1, libMesh::H1_SEMINORM, libMesh::H2, libMesh::H2_SEMINORM, libMesh::MeshFunction::hessian(), libMesh::L2, libMesh::System::name(), libMesh::TensorTools::norm_sq(), libMesh::TypeVector< T >::norm_sq(), libMesh::TypeTensor< T >::norm_sq(), libMesh::System::number(), libMesh::EquationSystems::parameters, libMesh::Real, libMesh::FEAbstract::reinit(), libMesh::DenseVector< T >::size(), libMesh::System::time, libMesh::SystemNorm::type(), libMesh::System::variable_number(), and libMesh::System::variable_scalar_number().

Referenced by estimate_error().

427 {
428  // The (string) name of this system
429  const std::string & sys_name = system.name();
430  const unsigned int sys_num = system.number();
431 
432  const unsigned int var = system.variable_number(var_name);
433  const unsigned int var_component =
434  system.variable_scalar_number(var, 0);
435 
436  const Parameters & parameters = system.get_equation_systems().parameters;
437 
438  // reinitialize the element-specific data
439  // for the current element
440  fe->reinit (elem);
441 
442  // Get the data we need to compute with
443  const std::vector<Real> & JxW = fe->get_JxW();
444  const std::vector<std::vector<Real> > & phi_values = fe->get_phi();
445  const std::vector<std::vector<RealGradient> > & dphi_values = fe->get_dphi();
446  const std::vector<Point> & q_point = fe->get_xyz();
447 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
448  const std::vector<std::vector<RealTensor> > & d2phi_values = fe->get_d2phi();
449 #endif
450 
451  // The number of shape functions
452  const unsigned int n_sf =
453  cast_int<unsigned int>(Uelem.size());
454 
455  // The number of quadrature points
456  const unsigned int n_qp =
457  cast_int<unsigned int>(JxW.size());
458 
459  Real error_val = 0;
460 
461  // Begin the loop over the Quadrature points.
462  //
463  for (unsigned int qp=0; qp<n_qp; qp++)
464  {
465  // Real u_h = 0.;
466  // RealGradient grad_u_h;
467 
468  Number u_h = 0.;
469 
470  Gradient grad_u_h;
471 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
472  Tensor grad2_u_h;
473 #endif
474 
475  // Compute solution values at the current
476  // quadrature point. This requires a sum
477  // over all the shape functions evaluated
478  // at the quadrature point.
479  for (unsigned int i=0; i<n_sf; i++)
480  {
481  // Values from current solution.
482  u_h += phi_values[i][qp]*Uelem(i);
483  grad_u_h += dphi_values[i][qp]*Uelem(i);
484 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
485  grad2_u_h += d2phi_values[i][qp]*Uelem(i);
486 #endif
487  }
488 
489  // Compute the value of the error at this quadrature point
490  if (error_norm.type(var) == L2 ||
491  error_norm.type(var) == H1 ||
492  error_norm.type(var) == H2)
493  {
494  Number val_error = u_h;
495  if (_exact_value)
496  val_error -= _exact_value(q_point[qp],parameters,sys_name,var_name);
497  else if (_exact_values.size() > sys_num && _exact_values[sys_num])
498  val_error -= _exact_values[sys_num]->
499  component(var_component, q_point[qp], system.time);
500  else if (_equation_systems_fine)
501  val_error -= (*fine_values)(q_point[qp]);
502 
503  // Add the squares of the error to each contribution
504  error_val += JxW[qp]*TensorTools::norm_sq(val_error);
505  }
506 
507  // Compute the value of the error in the gradient at this
508  // quadrature point
509  if (error_norm.type(var) == H1 ||
510  error_norm.type(var) == H1_SEMINORM ||
511  error_norm.type(var) == H2)
512  {
513  Gradient grad_error = grad_u_h;
514  if (_exact_deriv)
515  grad_error -= _exact_deriv(q_point[qp],parameters,sys_name,var_name);
516  else if (_exact_derivs.size() > sys_num && _exact_derivs[sys_num])
517  grad_error -= _exact_derivs[sys_num]->
518  component(var_component, q_point[qp], system.time);
519  else if (_equation_systems_fine)
520  grad_error -= fine_values->gradient(q_point[qp]);
521 
522  error_val += JxW[qp]*grad_error.norm_sq();
523  }
524 
525 
526 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
527  // Compute the value of the error in the hessian at this
528  // quadrature point
529  if ((error_norm.type(var) == H2_SEMINORM ||
530  error_norm.type(var) == H2))
531  {
532  Tensor grad2_error = grad2_u_h;
533  if (_exact_hessian)
534  grad2_error -= _exact_hessian(q_point[qp],parameters,sys_name,var_name);
535  else if (_exact_hessians.size() > sys_num && _exact_hessians[sys_num])
536  grad2_error -= _exact_hessians[sys_num]->
537  component(var_component, q_point[qp], system.time);
538  else if (_equation_systems_fine)
539  grad2_error -= fine_values->hessian(q_point[qp]);
540 
541  error_val += JxW[qp]*grad2_error.norm_sq();
542  }
543 #endif
544 
545  } // end qp loop
546 
547  libmesh_assert_greater_equal (error_val, 0.);
548 
549  return error_val;
550 }
Number(* _exact_value)(const Point &p, const Parameters &parameters, const std::string &sys_name, const std::string &unknown_name)
std::vector< FunctionBase< Number > * > _exact_values
std::vector< FunctionBase< Tensor > * > _exact_hessians
std::vector< FunctionBase< Gradient > * > _exact_derivs
virtual unsigned int size() const libmesh_override
Definition: dense_vector.h:87
NumberVectorValue Gradient
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Tensor(* _exact_hessian)(const Point &p, const Parameters &parameters, const std::string &sys_name, const std::string &unknown_name)
NumberTensorValue Tensor
FEMNormType type(unsigned int var) const
Definition: system_norm.h:268
Gradient(* _exact_deriv)(const Point &p, const Parameters &parameters, const std::string &sys_name, const std::string &unknown_name)
void libMesh::ErrorEstimator::reduce_error ( std::vector< ErrorVectorReal > &  error_per_cell,
const Parallel::Communicator &comm  LIBMESH_CAN_DEFAULT_TO_COMMWORLD 
) 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 33 of file error_estimator.C.

References libMesh::Parallel::Communicator::sum().

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

35 {
36  // This function must be run on all processors at once
37  // parallel_object_only();
38 
39  // Each processor has now computed the error contributions
40  // for its local elements. We may need to sum the vector to
41  // recover the error for each element.
42 
43  comm.sum(error_per_cell);
44 }
virtual ErrorEstimatorType libMesh::ExactErrorEstimator::type ( ) const
inlinevirtual
Returns
The type for the ErrorEstimator subclass.

Implements libMesh::ErrorEstimator.

Definition at line 189 of file exact_error_estimator.h.

References libMesh::EXACT.

Member Data Documentation

EquationSystems* libMesh::ExactErrorEstimator::_equation_systems_fine
private

Constant pointer to the EquationSystems object containing a fine grid solution.

Definition at line 243 of file exact_error_estimator.h.

Referenced by attach_exact_deriv(), attach_exact_value(), attach_reference_solution(), estimate_error(), and find_squared_element_error().

Gradient(* libMesh::ExactErrorEstimator::_exact_deriv) (const Point &p, const Parameters &parameters, const std::string &sys_name, const std::string &unknown_name)
private

Function pointer to user-provided function which computes the exact derivative of the solution.

Definition at line 207 of file exact_error_estimator.h.

Referenced by attach_exact_value(), attach_reference_solution(), and find_squared_element_error().

std::vector<FunctionBase<Gradient> *> libMesh::ExactErrorEstimator::_exact_derivs
private

User-provided functors which compute the exact derivative of the solution for each system.

Definition at line 231 of file exact_error_estimator.h.

Referenced by attach_exact_deriv(), attach_exact_derivs(), clear_functors(), estimate_error(), and find_squared_element_error().

Tensor(* libMesh::ExactErrorEstimator::_exact_hessian) (const Point &p, const Parameters &parameters, const std::string &sys_name, const std::string &unknown_name)
private

Function pointer to user-provided function which computes the exact hessian of the solution.

Definition at line 216 of file exact_error_estimator.h.

Referenced by attach_exact_deriv(), attach_reference_solution(), and find_squared_element_error().

std::vector<FunctionBase<Tensor> *> libMesh::ExactErrorEstimator::_exact_hessians
private

User-provided functors which compute the exact hessians of the solution for each system.

Definition at line 237 of file exact_error_estimator.h.

Referenced by attach_exact_hessian(), attach_exact_hessians(), clear_functors(), estimate_error(), and find_squared_element_error().

Number(* libMesh::ExactErrorEstimator::_exact_value) (const Point &p, const Parameters &parameters, const std::string &sys_name, const std::string &unknown_name)
private

Function pointer to user-provided function which computes the exact value of the solution.

Definition at line 198 of file exact_error_estimator.h.

Referenced by attach_reference_solution(), and find_squared_element_error().

std::vector<FunctionBase<Number> *> libMesh::ExactErrorEstimator::_exact_values
private

User-provided functors which compute the exact value of the solution for each system.

Definition at line 225 of file exact_error_estimator.h.

Referenced by attach_exact_value(), attach_exact_values(), clear_functors(), estimate_error(), and find_squared_element_error().

int libMesh::ExactErrorEstimator::_extra_order
private

Extra order to use for quadrature rule

Definition at line 263 of file exact_error_estimator.h.

Referenced by estimate_error(), and extra_quadrature_order().

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 150 of file error_estimator.h.

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


The documentation for this class was generated from the following files: