$extrastylesheet
libMesh::ExactErrorEstimator Class Reference

#include <exact_error_estimator.h>

Inheritance diagram for libMesh::ExactErrorEstimator:

List of all members.

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=NULL, bool estimate_parent_error=false)
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=NULL, 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=NULL, bool estimate_parent_error=false)

Public Attributes

SystemNorm error_norm

Protected Member Functions

void reduce_error (std::vector< float > &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, 2006.

Definition at line 66 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 110 of file error_estimator.h.


Constructor & Destructor Documentation

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

Definition at line 74 of file exact_error_estimator.h.

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

Destructor.

Definition at line 86 of file exact_error_estimator.h.

{}

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 _exact_derivs, and libMesh::FunctionBase< Output >::clone().

{
  if (_exact_derivs.size() <= sys_num)
    _exact_derivs.resize(sys_num+1, NULL);

  if (g)
    _exact_derivs[sys_num] = g->clone().release();
}
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.

Definition at line 92 of file exact_error_estimator.C.

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

{
  libmesh_assert(gptr);
  _exact_deriv = gptr;

  // We're not using a fine grid solution
  _equation_systems_fine = NULL;

  // We're not using user-provided functors
  this->clear_functors();
}

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.

{
  // Clear out any previous _exact_derivs entries, then add a new
  // entry for each system.
  for (unsigned int i=0; i != _exact_derivs.size(); ++i)
    delete (_exact_derivs[i]);

  _exact_derivs.clear();
  _exact_derivs.resize(g.size(), NULL);

  // We use clone() to get non-sliced copies of FunctionBase
  // subclasses, but we don't currently put the resulting UniquePtrs
  // into an STL container.
  for (unsigned int i=0; i != g.size(); ++i)
    if (g[i])
      _exact_derivs[i] = g[i]->clone().release();
}
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, and libMesh::FunctionBase< Output >::clone().

{
  if (_exact_hessians.size() <= sys_num)
    _exact_hessians.resize(sys_num+1, NULL);

  if (h)
    _exact_hessians[sys_num] = h->clone().release();
}
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.

Definition at line 140 of file exact_error_estimator.C.

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

{
  libmesh_assert(hptr);
  _exact_hessian = hptr;

  // We're not using a fine grid solution
  _equation_systems_fine = NULL;

  // We're not using user-provided functors
  this->clear_functors();
}

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.

{
  // Clear out any previous _exact_hessians entries, then add a new
  // entry for each system.
  for (unsigned int i=0; i != _exact_hessians.size(); ++i)
    delete (_exact_hessians[i]);

  _exact_hessians.clear();
  _exact_hessians.resize(h.size(), NULL);

  // We use clone() to get non-sliced copies of FunctionBase
  // subclasses, but we don't currently put the resulting UniquePtrs
  // into an STL container.
  for (unsigned int i=0; i != h.size(); ++i)
    if (h[i])
      _exact_hessians[i] = h[i]->clone().release();
}
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 _exact_values, and libMesh::FunctionBase< Output >::clone().

{
  if (_exact_values.size() <= sys_num)
    _exact_values.resize(sys_num+1, NULL);

  if (f)
    _exact_values[sys_num] = f->clone().release();
}
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.

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.

{
  // Clear out any previous _exact_values entries, then add a new
  // entry for each system.
  for (unsigned int i=0; i != _exact_values.size(); ++i)
    delete (_exact_values[i]);

  _exact_values.clear();
  _exact_values.resize(f.size(), NULL);

  // We use clone() to get non-sliced copies of FunctionBase
  // subclasses, but we don't currently put the resulting UniquePtrs
  // into an STL container.
  for (unsigned int i=0; i != f.size(); ++i)
    if (f[i])
      _exact_values[i] = f[i]->clone().release();
}

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(), and libMesh::libmesh_assert().

{
  libmesh_assert(es_fine);
  _equation_systems_fine = es_fine;

  // If we're using a fine grid solution, we're not using exact value
  // function pointers or functors.
  _exact_value = NULL;
  _exact_deriv = NULL;
  _exact_hessian = NULL;

  this->clear_functors();
}

Helper method for cleanup

Definition at line 562 of file exact_error_estimator.C.

References _exact_derivs, _exact_hessians, and _exact_values.

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

{
  // delete will clean up any cloned functors and no-op on any NULL
  // pointers

  for (unsigned int i=0; i != _exact_values.size(); ++i)
    delete (_exact_values[i]);

  for (unsigned int i=0; i != _exact_derivs.size(); ++i)
    delete (_exact_derivs[i]);

  for (unsigned int i=0; i != _exact_hessians.size(); ++i)
    delete (_exact_hessians[i]);
}
void libMesh::ExactErrorEstimator::estimate_error ( const System system,
ErrorVector error_per_cell,
const NumericVector< Number > *  solution_vector = NULL,
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 201 of file exact_error_estimator.C.

References libMesh::Elem::active(), libMesh::MeshBase::active_local_elements_begin(), libMesh::MeshBase::active_local_elements_end(), libMesh::NumericVector< T >::build(), libMesh::FEGenericBase< OutputType >::build(), libMesh::Elem::child(), 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::ErrorVectorReal, libMesh::System::get_dof_map(), libMesh::System::get_mesh(), libMesh::DofObject::id(), 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::SERIAL, libMesh::System::solution, libMesh::START_LOG(), libMesh::NumericVector< T >::swap(), libMesh::sys, libMesh::System::update_global_solution(), libMesh::System::variable_name(), libMesh::System::variable_number(), and libMesh::DofMap::variable_type().

{
  // The current mesh
  const MeshBase& mesh = system.get_mesh();

  // The dimensionality of the mesh
  const unsigned int dim = mesh.mesh_dimension();

  // The number of variables in the system
  const unsigned int n_vars = system.n_vars();

  // The DofMap for this system
  const DofMap& dof_map = system.get_dof_map();

  // Resize the error_per_cell vector to be
  // the number of elements, initialize it to 0.
  error_per_cell.resize (mesh.max_elem_id());
  std::fill (error_per_cell.begin(), error_per_cell.end(), 0.);

  // Prepare current_local_solution to localize a non-standard
  // solution vector if necessary
  if (solution_vector && solution_vector != system.solution.get())
    {
      NumericVector<Number>* newsol =
        const_cast<NumericVector<Number>*>(solution_vector);
      System &sys = const_cast<System&>(system);
      newsol->swap(*sys.solution);
      sys.update();
    }

  // Loop over all the variables in the system
  for (unsigned int var=0; var<n_vars; var++)
    {
      // Possibly skip this variable
      if (error_norm.weight(var) == 0.0) continue;

      // The (string) name of this variable
      const std::string& var_name = system.variable_name(var);

      // The type of finite element to use for this variable
      const FEType& fe_type = dof_map.variable_type (var);

      UniquePtr<FEBase> fe (FEBase::build (dim, fe_type));

      // Build an appropriate Gaussian quadrature rule
      UniquePtr<QBase> qrule =
        fe_type.default_quadrature_rule (dim,
                                         _extra_order);

      fe->attach_quadrature_rule (qrule.get());

      // Prepare a global solution and a MeshFunction of the fine system if we need one
      UniquePtr<MeshFunction> fine_values;
      UniquePtr<NumericVector<Number> > fine_soln = NumericVector<Number>::build(system.comm());
      if (_equation_systems_fine)
        {
          const System& fine_system = _equation_systems_fine->get_system(system.name());

          std::vector<Number> global_soln;
          // FIXME - we're assuming that the fine system solution gets
          // used even when a different vector is used for the coarse
          // system
          fine_system.update_global_solution(global_soln);
          fine_soln->init
            (cast_int<numeric_index_type>(global_soln.size()), true,
             SERIAL);
          (*fine_soln) = global_soln;

          fine_values = UniquePtr<MeshFunction>
            (new MeshFunction(*_equation_systems_fine,
                              *fine_soln,
                              fine_system.get_dof_map(),
                              fine_system.variable_number(var_name)));
          fine_values->init();
        } else {
        // Initialize functors if we're using them
        for (unsigned int i=0; i != _exact_values.size(); ++i)
          if (_exact_values[i])
            _exact_values[i]->init();

        for (unsigned int i=0; i != _exact_derivs.size(); ++i)
          if (_exact_derivs[i])
            _exact_derivs[i]->init();

        for (unsigned int i=0; i != _exact_hessians.size(); ++i)
          if (_exact_hessians[i])
            _exact_hessians[i]->init();
      }

      // Request the data we'll need to compute with
      fe->get_JxW();
      fe->get_phi();
      fe->get_dphi();
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
      fe->get_d2phi();
#endif
      fe->get_xyz();

#ifdef LIBMESH_ENABLE_AMR
      // If we compute on parent elements, we'll want to do so only
      // once on each, so we need to keep track of which we've done.
      std::vector<bool> computed_var_on_parent;

      if (estimate_parent_error)
        computed_var_on_parent.resize(error_per_cell.size(), false);
#endif

      // TODO: this ought to be threaded (and using subordinate
      // MeshFunction objects in each thread rather than a single
      // master)

      // Iterate over all the active elements in the mesh
      // that live on this processor.
      MeshBase::const_element_iterator
        elem_it  = mesh.active_local_elements_begin();
      const MeshBase::const_element_iterator
        elem_end = mesh.active_local_elements_end();

      for (;elem_it != elem_end; ++elem_it)
        {
          // e is necessarily an active element on the local processor
          const Elem* elem = *elem_it;
          const dof_id_type e_id = elem->id();

#ifdef LIBMESH_ENABLE_AMR
          // See if the parent of element e has been examined yet;
          // if not, we may want to compute the estimator on it
          const Elem* parent = elem->parent();

          // We only can compute and only need to compute on
          // parents with all active children
          bool compute_on_parent = true;
          if (!parent || !estimate_parent_error)
            compute_on_parent = false;
          else
            for (unsigned int c=0; c != parent->n_children(); ++c)
              if (!parent->child(c)->active())
                compute_on_parent = false;

          if (compute_on_parent &&
              !computed_var_on_parent[parent->id()])
            {
              computed_var_on_parent[parent->id()] = true;

              // Compute a projection onto the parent
              DenseVector<Number> Uparent;
              FEBase::coarsened_dof_values(*(system.current_local_solution),
                                           dof_map, parent, Uparent,
                                           var, false);

              error_per_cell[parent->id()] +=
                static_cast<ErrorVectorReal>
                (find_squared_element_error(system, var_name,
                                            parent, Uparent,
                                            fe.get(),
                                            fine_values.get()));
            }
#endif

          // Get the local to global degree of freedom maps
          std::vector<dof_id_type> dof_indices;
          dof_map.dof_indices (elem, dof_indices, var);
          const unsigned int n_dofs =
            cast_int<unsigned int>(dof_indices.size());
          DenseVector<Number> Uelem(n_dofs);
          for (unsigned int i=0; i != n_dofs; ++i)
            Uelem(i) = system.current_solution(dof_indices[i]);

          error_per_cell[e_id] +=
            static_cast<ErrorVectorReal>
            (find_squared_element_error(system, var_name, elem,
                                        Uelem, fe.get(),
                                        fine_values.get()));

        } // End loop over active local elements
    } // End loop over variables



  // Each processor has now computed the error contribuions
  // for its local elements.  We need to sum the vector
  // and then take the square-root of each component.  Note
  // that we only need to sum if we are running on multiple
  // processors, and we only need to take the square-root
  // if the value is nonzero.  There will in general be many
  // zeros for the inactive elements.

  // First sum the vector of estimated error values
  this->reduce_error(error_per_cell, system.comm());

  // Compute the square-root of each component.
  START_LOG("std::sqrt()", "ExactErrorEstimator");
  for (dof_id_type i=0; i<error_per_cell.size(); i++)
    {

      if (error_per_cell[i] != 0.)
        {
          libmesh_assert_greater (error_per_cell[i], 0.);
          error_per_cell[i] = std::sqrt(error_per_cell[i]);
        }


    }
  STOP_LOG("std::sqrt()", "ExactErrorEstimator");

  // If we used a non-standard solution before, now is the time to fix
  // the current_local_solution
  if (solution_vector && solution_vector != system.solution.get())
    {
      NumericVector<Number>* newsol =
        const_cast<NumericVector<Number>*>(solution_vector);
      System &sys = const_cast<System&>(system);
      newsol->swap(*sys.solution);
      sys.update();
    }
}
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 = NULL,
bool  estimate_parent_error = false 
) [virtual, inherited]

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::EquationSystems::n_systems(), and libMesh::sys.

{
  SystemNorm old_error_norm = this->error_norm;

  // Sum the error values from each system
  for (unsigned int s = 0; s != equation_systems.n_systems(); ++s)
    {
      ErrorVector system_error_per_cell;
      const System &sys = equation_systems.get_system(s);
      if (error_norms.find(&sys) == error_norms.end())
        this->error_norm = old_error_norm;
      else
        this->error_norm = error_norms.find(&sys)->second;

      const NumericVector<Number>* solution_vector = NULL;
      if (solution_vectors &&
          solution_vectors->find(&sys) != solution_vectors->end())
        solution_vector = solution_vectors->find(&sys)->second;

      this->estimate_error(sys, system_error_per_cell,
                           solution_vector, estimate_parent_error);

      if (s)
        {
          libmesh_assert_equal_to (error_per_cell.size(), system_error_per_cell.size());
          for (unsigned int i=0; i != error_per_cell.size(); ++i)
            error_per_cell[i] += system_error_per_cell[i];
        }
      else
        error_per_cell = system_error_per_cell;
    }

  // Restore our old state before returning
  this->error_norm = old_error_norm;
}
void libMesh::ErrorEstimator::estimate_errors ( const EquationSystems equation_systems,
ErrorMap errors_per_cell,
const std::map< const System *, const NumericVector< Number > * > *  solution_vectors = NULL,
bool  estimate_parent_error = false 
) [virtual, inherited]

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::EquationSystems::n_systems(), n_vars, libMesh::System::n_vars(), libMesh::sys, and libMesh::SystemNorm::type().

{
  SystemNorm old_error_norm = this->error_norm;

  // Find the requested error values from each system
  for (unsigned int s = 0; s != equation_systems.n_systems(); ++s)
    {
      const System &sys = equation_systems.get_system(s);

      unsigned int n_vars = sys.n_vars();

      for (unsigned int v = 0; v != n_vars; ++v)
        {
          // Only fill in ErrorVectors the user asks for
          if (errors_per_cell.find(std::make_pair(&sys, v)) ==
              errors_per_cell.end())
            continue;

          // Calculate error in only one variable
          std::vector<Real> weights(n_vars, 0.0);
          weights[v] = 1.0;
          this->error_norm =
            SystemNorm(std::vector<FEMNormType>(n_vars, old_error_norm.type(v)),
                       weights);

          const NumericVector<Number>* solution_vector = NULL;
          if (solution_vectors &&
              solution_vectors->find(&sys) != solution_vectors->end())
            solution_vector = solution_vectors->find(&sys)->second;

          this->estimate_error
            (sys, *errors_per_cell[std::make_pair(&sys, v)],
             solution_vector, estimate_parent_error);
        }
    }

  // Restore our old state before returning
  this->error_norm = old_error_norm;
}
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 166 of file exact_error_estimator.h.

References _extra_order.

  { _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 429 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::System::number(), libMesh::EquationSystems::parameters, libMesh::Real, libMesh::FEAbstract::reinit(), libMesh::DenseVector< T >::size(), libMesh::TypeVector< T >::size_sq(), libMesh::TypeTensor< T >::size_sq(), libMesh::System::time, libMesh::SystemNorm::type(), libMesh::System::variable_number(), and libMesh::System::variable_scalar_number().

{
  // The (string) name of this system
  const std::string& sys_name = system.name();
  const unsigned int sys_num = system.number();

  const unsigned int var = system.variable_number(var_name);
  const unsigned int var_component =
    system.variable_scalar_number(var, 0);

  const Parameters& parameters = system.get_equation_systems().parameters;

  // reinitialize the element-specific data
  // for the current element
  fe->reinit (elem);

  // Get the data we need to compute with
  const std::vector<Real> &                      JxW          = fe->get_JxW();
  const std::vector<std::vector<Real> >&         phi_values   = fe->get_phi();
  const std::vector<std::vector<RealGradient> >& dphi_values  = fe->get_dphi();
  const std::vector<Point>&                      q_point      = fe->get_xyz();
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
  const std::vector<std::vector<RealTensor> >&   d2phi_values = fe->get_d2phi();
#endif

  // The number of shape functions
  const unsigned int n_sf =
    cast_int<unsigned int>(Uelem.size());

  // The number of quadrature points
  const unsigned int n_qp =
    cast_int<unsigned int>(JxW.size());

  Real error_val = 0;

  // Begin the loop over the Quadrature points.
  //
  for (unsigned int qp=0; qp<n_qp; qp++)
    {
      // Real u_h = 0.;
      // RealGradient grad_u_h;

      Number u_h = 0.;

      Gradient grad_u_h;
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
      Tensor grad2_u_h;
#endif

      // Compute solution values at the current
      // quadrature point.  This reqiures a sum
      // over all the shape functions evaluated
      // at the quadrature point.
      for (unsigned int i=0; i<n_sf; i++)
        {
          // Values from current solution.
          u_h      += phi_values[i][qp]*Uelem(i);
          grad_u_h += dphi_values[i][qp]*Uelem(i);
#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
          grad2_u_h += d2phi_values[i][qp]*Uelem(i);
#endif
        }

      // Compute the value of the error at this quadrature point
      if (error_norm.type(var) == L2 ||
          error_norm.type(var) == H1 ||
          error_norm.type(var) == H2)
        {
          Number val_error = u_h;
          if (_exact_value)
            val_error -= _exact_value(q_point[qp],parameters,sys_name,var_name);
          else if (_exact_values.size() > sys_num && _exact_values[sys_num])
            val_error -= _exact_values[sys_num]->
              component(var_component, q_point[qp], system.time);
          else if (_equation_systems_fine)
            val_error -= (*fine_values)(q_point[qp]);

          // Add the squares of the error to each contribution
          error_val += JxW[qp]*TensorTools::norm_sq(val_error);
        }

      // Compute the value of the error in the gradient at this
      // quadrature point
      if (error_norm.type(var) == H1 ||
          error_norm.type(var) == H1_SEMINORM ||
          error_norm.type(var) == H2)
        {
          Gradient grad_error = grad_u_h;
          if(_exact_deriv)
            grad_error -= _exact_deriv(q_point[qp],parameters,sys_name,var_name);
          else if (_exact_derivs.size() > sys_num && _exact_derivs[sys_num])
            grad_error -= _exact_derivs[sys_num]->
              component(var_component, q_point[qp], system.time);
          else if(_equation_systems_fine)
            grad_error -= fine_values->gradient(q_point[qp]);

          error_val += JxW[qp]*grad_error.size_sq();
        }


#ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
      // Compute the value of the error in the hessian at this
      // quadrature point
      if ((error_norm.type(var) == H2_SEMINORM ||
           error_norm.type(var) == H2))
        {
          Tensor grad2_error = grad2_u_h;
          if(_exact_hessian)
            grad2_error -= _exact_hessian(q_point[qp],parameters,sys_name,var_name);
          else if (_exact_hessians.size() > sys_num && _exact_hessians[sys_num])
            grad2_error -= _exact_hessians[sys_num]->
              component(var_component, q_point[qp], system.time);
          else if (_equation_systems_fine)
            grad2_error -= fine_values->hessian(q_point[qp]);

          error_val += JxW[qp]*grad2_error.size_sq();
        }
#endif

    } // end qp loop

  libmesh_assert_greater_equal (error_val, 0.);

  return error_val;
}
void libMesh::ErrorEstimator::reduce_error ( std::vector< float > &  error_per_cell,
const Parallel::Communicator &comm  LIBMESH_CAN_DEFAULT_TO_COMMWORLD 
) const [protected, inherited]

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(), and libMesh::AdjointRefinementEstimator::estimate_error().

{
  // This function must be run on all processors at once
  // parallel_object_only();

  // Each processor has now computed the error contribuions
  // for its local elements.  We may need to sum the vector to
  // recover the error for each element.

  comm.sum(error_per_cell);
}

Member Data Documentation

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

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

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

Definition at line 227 of file exact_error_estimator.h.

Referenced by attach_exact_deriv(), attach_exact_derivs(), clear_functors(), 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 212 of file exact_error_estimator.h.

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

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

Definition at line 233 of file exact_error_estimator.h.

Referenced by attach_exact_hessian(), attach_exact_hessians(), clear_functors(), 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 194 of file exact_error_estimator.h.

Referenced by attach_reference_solution(), and find_squared_element_error().

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

Definition at line 221 of file exact_error_estimator.h.

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

Extra order to use for quadrature rule

Definition at line 259 of file exact_error_estimator.h.

Referenced by extra_quadrature_order().

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 142 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(), libMesh::ErrorEstimator::estimate_errors(), ExactErrorEstimator(), find_squared_element_error(), 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(), and libMesh::UniformRefinementEstimator::UniformRefinementEstimator().


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