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

This class implements an "error estimator" based on the difference between the approximate and exact solution. More...

#include <exact_error_estimator.h>

Inheritance diagram for libMesh::ExactErrorEstimator:
[legend]

Public Types

typedef Number(* ValueFunctionPointer) (const 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. More...
 
typedef Gradient(* GradientFunctionPointer) (const 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. More...
 
typedef Tensor(* HessianFunctionPointer) (const 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. More...
 
typedef std::map< std::pair< const System *, unsigned int >, std::unique_ptr< ErrorVector > > ErrorMap
 When calculating many error vectors at once, we need a data structure to hold them all. More...
 

Public Member Functions

 ExactErrorEstimator ()
 Constructor. More...
 
 ExactErrorEstimator (const ExactErrorEstimator &)=delete
 This class cannot be (default) copy constructed/assigned because it has containers of unique_ptrs. More...
 
ExactErrorEstimatoroperator= (const ExactErrorEstimator &)=delete
 
 ExactErrorEstimator (ExactErrorEstimator &&)=default
 Defaulted move ctor, move assignment operator, and destructor. More...
 
ExactErrorEstimatoroperator= (ExactErrorEstimator &&)=default
 
virtual ~ExactErrorEstimator ()=default
 
void 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. More...
 
void 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. More...
 
void attach_exact_value (ValueFunctionPointer fptr)
 
void 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. More...
 
void 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. More...
 
void attach_exact_deriv (GradientFunctionPointer gptr)
 
void 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. More...
 
void 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. More...
 
void attach_exact_hessian (HessianFunctionPointer hptr)
 
void 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. More...
 
void extra_quadrature_order (const int extraorder)
 Increases or decreases the order of the quadrature rule used for numerical integration. More...
 
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 exact solution function to estimate the error on each cell. More...
 
virtual ErrorEstimatorType type () const 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=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

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

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...
 

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
 Helper method for calculating on each element. More...
 
void clear_functors ()
 Helper method for cleanup. More...
 

Private Attributes

ValueFunctionPointer _exact_value
 Function pointer to user-provided function which computes the exact value of the solution. More...
 
GradientFunctionPointer _exact_deriv
 Function pointer to user-provided function which computes the exact derivative of the solution. More...
 
HessianFunctionPointer _exact_hessian
 Function pointer to user-provided function which computes the exact hessian of the solution. More...
 
std::vector< std::unique_ptr< FunctionBase< Number > > > _exact_values
 User-provided functors which compute the exact value of the solution for each system. More...
 
std::vector< std::unique_ptr< FunctionBase< Gradient > > > _exact_derivs
 User-provided functors which compute the exact derivative of the solution for each system. More...
 
std::vector< std::unique_ptr< FunctionBase< Tensor > > > _exact_hessians
 User-provided functors which compute the exact hessians of the solution for each system. More...
 
EquationSystems_equation_systems_fine
 Constant pointer to the EquationSystems object containing a fine grid solution. More...
 
int _extra_order
 Extra order to use for quadrature rule. More...
 

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 57 of file exact_error_estimator.h.

Member Typedef Documentation

◆ ErrorMap

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

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

Definition at line 121 of file error_estimator.h.

◆ GradientFunctionPointer

typedef Gradient(* libMesh::ExactErrorEstimator::GradientFunctionPointer) (const 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.

Definition at line 122 of file exact_error_estimator.h.

◆ HessianFunctionPointer

typedef Tensor(* libMesh::ExactErrorEstimator::HessianFunctionPointer) (const 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.

Definition at line 145 of file exact_error_estimator.h.

◆ ValueFunctionPointer

typedef Number(* libMesh::ExactErrorEstimator::ValueFunctionPointer) (const 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.

Definition at line 99 of file exact_error_estimator.h.

Constructor & Destructor Documentation

◆ ExactErrorEstimator() [1/3]

libMesh::ExactErrorEstimator::ExactErrorEstimator ( )

Constructor.

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

Definition at line 49 of file exact_error_estimator.C.

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

49  :
51  _exact_value(nullptr),
52  _exact_deriv(nullptr),
53  _exact_hessian(nullptr),
54  _equation_systems_fine(nullptr),
55  _extra_order(1)
56 {
57  error_norm = H1;
58 }
HessianFunctionPointer _exact_hessian
Function pointer to user-provided function which computes the exact hessian of the solution...
SystemNorm error_norm
When estimating the error in a single system, the error_norm is used to control the scaling and norm ...
int _extra_order
Extra order to use for quadrature rule.
ValueFunctionPointer _exact_value
Function pointer to user-provided function which computes the exact value of the solution.
ErrorEstimator()=default
Constructor.
GradientFunctionPointer _exact_deriv
Function pointer to user-provided function which computes the exact derivative of the solution...
EquationSystems * _equation_systems_fine
Constant pointer to the EquationSystems object containing a fine grid solution.

◆ ExactErrorEstimator() [2/3]

libMesh::ExactErrorEstimator::ExactErrorEstimator ( const ExactErrorEstimator )
delete

This class cannot be (default) copy constructed/assigned because it has containers of unique_ptrs.

Explicitly deleting these functions is the best way to document this fact.

◆ ExactErrorEstimator() [3/3]

libMesh::ExactErrorEstimator::ExactErrorEstimator ( ExactErrorEstimator &&  )
default

Defaulted move ctor, move assignment operator, and destructor.

◆ ~ExactErrorEstimator()

virtual libMesh::ExactErrorEstimator::~ExactErrorEstimator ( )
virtualdefault

Member Function Documentation

◆ attach_exact_deriv() [1/2]

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 126 of file exact_error_estimator.C.

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

Referenced by main().

128 {
129  if (_exact_derivs.size() <= sys_num)
130  _exact_derivs.resize(sys_num+1);
131 
132  if (g)
133  _exact_derivs[sys_num] = g->clone();
134 }
std::vector< std::unique_ptr< FunctionBase< Gradient > > > _exact_derivs
User-provided functors which compute the exact derivative of the solution for each system...
virtual std::unique_ptr< FunctionBase< Output > > clone() const =0

◆ attach_exact_deriv() [2/2]

void libMesh::ExactErrorEstimator::attach_exact_deriv ( GradientFunctionPointer  gptr)

Definition at line 102 of file exact_error_estimator.C.

References _equation_systems_fine, _exact_deriv, clear_functors(), gptr(), and libMesh::libmesh_assert().

103 {
105  _exact_deriv = gptr;
106 
107  // We're not using a fine grid solution
108  _equation_systems_fine = nullptr;
109 
110  // We're not using user-provided functors
111  this->clear_functors();
112 }
libmesh_assert(ctx)
void clear_functors()
Helper method for cleanup.
GradientFunctionPointer _exact_deriv
Function pointer to user-provided function which computes the exact derivative of the solution...
EquationSystems * _equation_systems_fine
Constant pointer to the EquationSystems object containing a fine grid solution.
Gradient gptr(const Point &p, const Parameters &, const std::string &libmesh_dbg_var(sys_name), const std::string &unknown_name)
Definition: projection.C:95

◆ attach_exact_derivs()

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 115 of file exact_error_estimator.C.

References _exact_derivs.

116 {
117  // Automatically delete any previous _exact_derivs entries, then add a new
118  // entry for each system.
119  _exact_derivs.clear();
120 
121  for (auto ptr : g)
122  _exact_derivs.emplace_back(ptr ? ptr->clone() : nullptr);
123 }
std::vector< std::unique_ptr< FunctionBase< Gradient > > > _exact_derivs
User-provided functors which compute the exact derivative of the solution for each system...

◆ attach_exact_hessian() [1/2]

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 163 of file exact_error_estimator.C.

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

165 {
166  if (_exact_hessians.size() <= sys_num)
167  _exact_hessians.resize(sys_num+1);
168 
169  if (h)
170  _exact_hessians[sys_num] = h->clone();
171 }
std::vector< std::unique_ptr< FunctionBase< Tensor > > > _exact_hessians
User-provided functors which compute the exact hessians of the solution for each system.

◆ attach_exact_hessian() [2/2]

void libMesh::ExactErrorEstimator::attach_exact_hessian ( HessianFunctionPointer  hptr)

Definition at line 139 of file exact_error_estimator.C.

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

140 {
141  libmesh_assert(hptr);
142  _exact_hessian = hptr;
143 
144  // We're not using a fine grid solution
145  _equation_systems_fine = nullptr;
146 
147  // We're not using user-provided functors
148  this->clear_functors();
149 }
HessianFunctionPointer _exact_hessian
Function pointer to user-provided function which computes the exact hessian of the solution...
libmesh_assert(ctx)
void clear_functors()
Helper method for cleanup.
EquationSystems * _equation_systems_fine
Constant pointer to the EquationSystems object containing a fine grid solution.

◆ attach_exact_hessians()

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 152 of file exact_error_estimator.C.

References _exact_hessians.

153 {
154  // Automatically delete any previous _exact_hessians entries, then add a new
155  // entry for each system.
156  _exact_hessians.clear();
157 
158  for (auto ptr : h)
159  _exact_hessians.emplace_back(ptr ? ptr->clone() : nullptr);
160 }
std::vector< std::unique_ptr< FunctionBase< Tensor > > > _exact_hessians
User-provided functors which compute the exact hessians of the solution for each system.

◆ attach_exact_value() [1/2]

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 91 of file exact_error_estimator.C.

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

Referenced by main().

93 {
94  if (_exact_values.size() <= sys_num)
95  _exact_values.resize(sys_num+1);
96 
97  if (f)
98  _exact_values[sys_num] = f->clone();
99 }
virtual std::unique_ptr< FunctionBase< Output > > clone() const =0
std::vector< std::unique_ptr< FunctionBase< Number > > > _exact_values
User-provided functors which compute the exact value of the solution for each system.

◆ attach_exact_value() [2/2]

void libMesh::ExactErrorEstimator::attach_exact_value ( ValueFunctionPointer  fptr)

Definition at line 67 of file exact_error_estimator.C.

References _equation_systems_fine, _exact_value, clear_functors(), fptr(), and libMesh::libmesh_assert().

68 {
71 
72  // We're not using a fine grid solution
73  _equation_systems_fine = nullptr;
74 
75  // We're not using user-provided functors
76  this->clear_functors();
77 }
ValueFunctionPointer _exact_value
Function pointer to user-provided function which computes the exact value of the solution.
Number fptr(const Point &p, const Parameters &, const std::string &libmesh_dbg_var(sys_name), const std::string &unknown_name)
Definition: projection.C:80
libmesh_assert(ctx)
void clear_functors()
Helper method for cleanup.
EquationSystems * _equation_systems_fine
Constant pointer to the EquationSystems object containing a fine grid solution.

◆ attach_exact_values()

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 80 of file exact_error_estimator.C.

References _exact_values.

81 {
82  // Automatically delete any previous _exact_values entries, then add a new
83  // entry for each system.
84  _exact_values.clear();
85 
86  for (auto ptr : f)
87  _exact_values.emplace_back(ptr ? ptr->clone() : nullptr);
88 }
std::vector< std::unique_ptr< FunctionBase< Number > > > _exact_values
User-provided functors which compute the exact value of the solution for each system.

◆ attach_reference_solution()

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 174 of file exact_error_estimator.C.

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

175 {
176  libmesh_assert(es_fine);
177  _equation_systems_fine = es_fine;
178 
179  // If we're using a fine grid solution, we're not using exact value
180  // function pointers or functors.
181  _exact_value = nullptr;
182  _exact_deriv = nullptr;
183  _exact_hessian = nullptr;
184 
185  this->clear_functors();
186 }
HessianFunctionPointer _exact_hessian
Function pointer to user-provided function which computes the exact hessian of the solution...
ValueFunctionPointer _exact_value
Function pointer to user-provided function which computes the exact value of the solution.
libmesh_assert(ctx)
void clear_functors()
Helper method for cleanup.
GradientFunctionPointer _exact_deriv
Function pointer to user-provided function which computes the exact derivative of the solution...
EquationSystems * _equation_systems_fine
Constant pointer to the EquationSystems object containing a fine grid solution.

◆ clear_functors()

void libMesh::ExactErrorEstimator::clear_functors ( )
private

Helper method for cleanup.

Definition at line 535 of file exact_error_estimator.C.

References _exact_derivs, _exact_hessians, and _exact_values.

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

536 {
537  _exact_values.clear();
538  _exact_derivs.clear();
539  _exact_hessians.clear();
540 }
std::vector< std::unique_ptr< FunctionBase< Tensor > > > _exact_hessians
User-provided functors which compute the exact hessians of the solution for each system.
std::vector< std::unique_ptr< FunctionBase< Gradient > > > _exact_derivs
User-provided functors which compute the exact derivative of the solution for each system...
std::vector< std::unique_ptr< FunctionBase< Number > > > _exact_values
User-provided functors which compute the exact value of the solution for each system.

◆ estimate_error()

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

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 188 of file exact_error_estimator.C.

References _equation_systems_fine, _exact_derivs, _exact_hessians, _exact_values, _extra_order, libMesh::FEGenericBase< OutputType >::build(), libMesh::NumericVector< T >::build(), libMesh::Elem::child_ref_range(), libMesh::FEGenericBase< OutputType >::coarsened_dof_values(), libMesh::ParallelObject::comm(), libMesh::System::current_local_solution, libMesh::System::current_solution(), libMesh::FEType::default_quadrature_rule(), dim, 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(), mesh, n_vars, libMesh::System::n_vars(), libMesh::System::name(), libMesh::Elem::parent(), libMesh::ErrorEstimator::reduce_error(), libMesh::SERIAL, libMesh::System::solution, std::sqrt(), libMesh::NumericVector< T >::swap(), libMesh::System::update_global_solution(), libMesh::System::variable_name(), libMesh::System::variable_number(), libMesh::DofMap::variable_type(), and libMesh::SystemNorm::weight().

Referenced by main().

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

◆ estimate_errors() [1/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.

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

52 {
53  SystemNorm old_error_norm = this->error_norm;
54 
55  // Sum the error values from each system
56  for (auto s : make_range(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 }
SystemNorm error_norm
When estimating the error in a single system, the error_norm is used to control the scaling and norm ...
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 ...
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:134
template class LIBMESH_EXPORT NumericVector< Number >
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:111

◆ estimate_errors() [2/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.

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

97 {
98  SystemNorm old_error_norm = this->error_norm;
99 
100  // Find the requested error values from each system
101  for (auto s : make_range(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.count(std::make_pair(&sys, v)))
111  continue;
112 
113  // Calculate error in only one variable
114  std::vector<Real> weights(n_vars, 0.0);
115  weights[v] = 1.0;
116  this->error_norm =
117  SystemNorm(std::vector<FEMNormType>(n_vars, old_error_norm.type(v)),
118  weights);
119 
120  const NumericVector<Number> * solution_vector = nullptr;
121  if (solution_vectors &&
122  solution_vectors->find(&sys) != solution_vectors->end())
123  solution_vector = solution_vectors->find(&sys)->second;
124 
125  this->estimate_error
126  (sys, *errors_per_cell[std::make_pair(&sys, v)],
127  solution_vector, estimate_parent_error);
128  }
129  }
130 
131  // Restore our old state before returning
132  this->error_norm = old_error_norm;
133 }
SystemNorm error_norm
When estimating the error in a single system, the error_norm is used to control the scaling and norm ...
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 ...
unsigned int n_vars
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:134
template class LIBMESH_EXPORT NumericVector< Number >

◆ extra_quadrature_order()

void libMesh::ExactErrorEstimator::extra_quadrature_order ( const int  extraorder)
inline

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

The default extraorder is 1, because properly integrating L2 error requires integrating the squares of terms with order p+1, and 2p+2 is 1 higher than what we default to using for reasonable mass matrix integration.

Definition at line 165 of file exact_error_estimator.h.

References _extra_order.

166  { _extra_order = extraorder; }
int _extra_order
Extra order to use for quadrature rule.

◆ find_squared_element_error()

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 402 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().

408 {
409  // The (string) name of this system
410  const std::string & sys_name = system.name();
411  const unsigned int sys_num = system.number();
412 
413  const unsigned int var = system.variable_number(var_name);
414  const unsigned int var_component =
415  system.variable_scalar_number(var, 0);
416 
417  const Parameters & parameters = system.get_equation_systems().parameters;
418 
419  // reinitialize the element-specific data
420  // for the current element
421  fe->reinit (elem);
422 
423  // Get the data we need to compute with
424  const std::vector<Real> & JxW = fe->get_JxW();
425  const std::vector<std::vector<Real>> & phi_values = fe->get_phi();
426  const std::vector<std::vector<RealGradient>> & dphi_values = fe->get_dphi();
427  const std::vector<Point> & q_point = fe->get_xyz();
428 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
429  const std::vector<std::vector<RealTensor>> & d2phi_values = fe->get_d2phi();
430 #endif
431 
432  // The number of shape functions
433  const unsigned int n_sf =
434  cast_int<unsigned int>(Uelem.size());
435 
436  // The number of quadrature points
437  const unsigned int n_qp =
438  cast_int<unsigned int>(JxW.size());
439 
440  Real error_val = 0;
441 
442  // Begin the loop over the Quadrature points.
443  //
444  for (unsigned int qp=0; qp<n_qp; qp++)
445  {
446  // Real u_h = 0.;
447  // RealGradient grad_u_h;
448 
449  Number u_h = 0.;
450 
451  Gradient grad_u_h;
452 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
453  Tensor grad2_u_h;
454 #endif
455 
456  // Compute solution values at the current
457  // quadrature point. This requires a sum
458  // over all the shape functions evaluated
459  // at the quadrature point.
460  for (unsigned int i=0; i<n_sf; i++)
461  {
462  // Values from current solution.
463  u_h += phi_values[i][qp]*Uelem(i);
464  grad_u_h += dphi_values[i][qp]*Uelem(i);
465 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
466  grad2_u_h += d2phi_values[i][qp]*Uelem(i);
467 #endif
468  }
469 
470  // Compute the value of the error at this quadrature point
471  if (error_norm.type(var) == L2 ||
472  error_norm.type(var) == H1 ||
473  error_norm.type(var) == H2)
474  {
475  Number val_error = u_h;
476  if (_exact_value)
477  val_error -= _exact_value(q_point[qp],parameters,sys_name,var_name);
478  else if (_exact_values.size() > sys_num && _exact_values[sys_num])
479  val_error -= _exact_values[sys_num]->
480  component(var_component, q_point[qp], system.time);
481  else if (_equation_systems_fine)
482  val_error -= (*fine_values)(q_point[qp]);
483 
484  // Add the squares of the error to each contribution
485  error_val += JxW[qp]*TensorTools::norm_sq(val_error);
486  }
487 
488  // Compute the value of the error in the gradient at this
489  // quadrature point
490  if (error_norm.type(var) == H1 ||
491  error_norm.type(var) == H1_SEMINORM ||
492  error_norm.type(var) == H2)
493  {
494  Gradient grad_error = grad_u_h;
495  if (_exact_deriv)
496  grad_error -= _exact_deriv(q_point[qp],parameters,sys_name,var_name);
497  else if (_exact_derivs.size() > sys_num && _exact_derivs[sys_num])
498  grad_error -= _exact_derivs[sys_num]->
499  component(var_component, q_point[qp], system.time);
500  else if (_equation_systems_fine)
501  grad_error -= fine_values->gradient(q_point[qp]);
502 
503  error_val += JxW[qp]*grad_error.norm_sq();
504  }
505 
506 
507 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
508  // Compute the value of the error in the hessian at this
509  // quadrature point
510  if ((error_norm.type(var) == H2_SEMINORM ||
511  error_norm.type(var) == H2))
512  {
513  Tensor grad2_error = grad2_u_h;
514  if (_exact_hessian)
515  grad2_error -= _exact_hessian(q_point[qp],parameters,sys_name,var_name);
516  else if (_exact_hessians.size() > sys_num && _exact_hessians[sys_num])
517  grad2_error -= _exact_hessians[sys_num]->
518  component(var_component, q_point[qp], system.time);
519  else if (_equation_systems_fine)
520  grad2_error -= fine_values->hessian(q_point[qp]);
521 
522  error_val += JxW[qp]*grad2_error.norm_sq();
523  }
524 #endif
525 
526  } // end qp loop
527 
528  libmesh_assert_greater_equal (error_val, 0.);
529 
530  return error_val;
531 }
HessianFunctionPointer _exact_hessian
Function pointer to user-provided function which computes the exact hessian of the solution...
SystemNorm error_norm
When estimating the error in a single system, the error_norm is used to control the scaling and norm ...
std::vector< std::unique_ptr< FunctionBase< Tensor > > > _exact_hessians
User-provided functors which compute the exact hessians of the solution for each system.
ValueFunctionPointer _exact_value
Function pointer to user-provided function which computes the exact value of the solution.
std::vector< std::unique_ptr< FunctionBase< Gradient > > > _exact_derivs
User-provided functors which compute the exact derivative of the solution for each system...
FEMNormType type(unsigned int var) const
Definition: system_norm.C:112
NumberVectorValue Gradient
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
GradientFunctionPointer _exact_deriv
Function pointer to user-provided function which computes the exact derivative of the solution...
EquationSystems * _equation_systems_fine
Constant pointer to the EquationSystems object containing a fine grid solution.
NumberTensorValue Tensor
auto norm_sq(const T &a) -> decltype(std::norm(a))
Definition: tensor_tools.h:104
virtual unsigned int size() const override final
Definition: dense_vector.h:104
std::vector< std::unique_ptr< FunctionBase< Number > > > _exact_values
User-provided functors which compute the exact value of the solution for each system.

◆ operator=() [1/2]

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

◆ operator=() [2/2]

ExactErrorEstimator& libMesh::ExactErrorEstimator::operator= ( ExactErrorEstimator &&  )
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.

References TIMPI::Communicator::sum().

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

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 }

◆ type()

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

Implements libMesh::ErrorEstimator.

Definition at line 61 of file exact_error_estimator.C.

References libMesh::EXACT.

62 {
63  return EXACT;
64 }

Member Data Documentation

◆ _equation_systems_fine

EquationSystems* libMesh::ExactErrorEstimator::_equation_systems_fine
private

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

Definition at line 230 of file exact_error_estimator.h.

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

◆ _exact_deriv

GradientFunctionPointer libMesh::ExactErrorEstimator::_exact_deriv
private

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

Definition at line 200 of file exact_error_estimator.h.

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

◆ _exact_derivs

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

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

Definition at line 218 of file exact_error_estimator.h.

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

◆ _exact_hessian

HessianFunctionPointer libMesh::ExactErrorEstimator::_exact_hessian
private

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

Definition at line 206 of file exact_error_estimator.h.

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

◆ _exact_hessians

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

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

Definition at line 224 of file exact_error_estimator.h.

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

◆ _exact_value

ValueFunctionPointer libMesh::ExactErrorEstimator::_exact_value
private

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

Definition at line 194 of file exact_error_estimator.h.

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

◆ _exact_values

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

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

Definition at line 212 of file exact_error_estimator.h.

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

◆ _extra_order

int libMesh::ExactErrorEstimator::_extra_order
private

Extra order to use for quadrature rule.

Definition at line 250 of file exact_error_estimator.h.

Referenced by estimate_error(), and extra_quadrature_order().

◆ 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 158 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::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::PatchRecoveryErrorEstimator::PatchRecoveryErrorEstimator(), libMesh::JumpErrorEstimator::reinit_sides(), and libMesh::UniformRefinementEstimator::UniformRefinementEstimator().


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