ExactErrorEstimator Class Reference

#include <exact_error_estimator.h>

Inheritance diagram for 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_value (Number fptr(const Point &p, const Parameters &Parameters, const std::string &sys_name, const std::string &unknown_name))
void attach_exact_deriv (Gradient fptr(const Point &p, const Parameters &parameters, const std::string &sys_name, const std::string &unknown_name))
void attach_exact_hessian (Tensor fptr(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

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

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


Member Typedef Documentation

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

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

Definition at line 105 of file error_estimator.h.


Constructor & Destructor Documentation

ExactErrorEstimator::ExactErrorEstimator (  )  [inline]

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

Definition at line 70 of file exact_error_estimator.h.

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

00070                         : _exact_value(NULL), 
00071                           _exact_deriv(NULL),
00072                           _exact_hessian(NULL),
00073                           _extra_order(0)
00074   { error_norm = H1; }

ExactErrorEstimator::~ExactErrorEstimator (  )  [inline]

Destructor.

Definition at line 79 of file exact_error_estimator.h.

00079 {}


Member Function Documentation

void ExactErrorEstimator::attach_exact_deriv ( Gradient   fptrconst Point &p,const Parameters &parameters,const std::string &sys_name,const std::string &unknown_name  ) 

Attach function similar to system.h which allows the user to attach an arbitrary function which computes the exact derivative of the solution at any point.

Definition at line 58 of file exact_error_estimator.C.

References _equation_systems_fine, and _exact_deriv.

00062 {
00063   libmesh_assert (fptr != NULL);
00064   _exact_deriv = fptr;
00065 
00066   // If we're not using a fine grid solution
00067   _equation_systems_fine = NULL;
00068 
00069 }

void ExactErrorEstimator::attach_exact_hessian ( Tensor   fptrconst Point &p,const Parameters &parameters,const std::string &sys_name,const std::string &unknown_name  ) 

Attach function similar to system.h which allows the user to attach an arbitrary function which computes the exact second derivatives of the solution at any point.

Definition at line 72 of file exact_error_estimator.C.

References _equation_systems_fine, and _exact_hessian.

00076 {
00077   libmesh_assert (fptr != NULL);
00078   _exact_hessian = fptr;
00079 
00080   // If we're not using a fine grid solution
00081   _equation_systems_fine = NULL;
00082 }

void ExactErrorEstimator::attach_exact_value ( Number   fptrconst Point &p, const Parameters &Parameters, const std::string &sys_name, const std::string &unknown_name  ) 

Attach function similar to system.h which allows the user to attach an arbitrary function which computes the exact value of the solution at any point.

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

References _equation_systems_fine, _exact_deriv, _exact_hessian, and _exact_value.

00085 {
00086   libmesh_assert (es_fine != NULL);
00087   _equation_systems_fine = es_fine;
00088 
00089   // If we're using a fine grid solution, we're not using exact values
00090   _exact_value = NULL;
00091   _exact_deriv = NULL;
00092   _exact_hessian = NULL;
00093 }

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

Definition at line 96 of file exact_error_estimator.C.

References Elem::active(), MeshBase::active_local_elements_begin(), MeshBase::active_local_elements_end(), NumericVector< T >::build(), FEBase::build(), Elem::child(), FEBase::coarsened_dof_values(), System::current_local_solution, System::current_solution(), FEType::default_quadrature_rule(), dim, DofMap::dof_indices(), AutoPtr< Tp >::get(), System::get_dof_map(), System::get_mesh(), DofObject::id(), MeshBase::max_elem_id(), mesh, MeshBase::mesh_dimension(), Elem::n_children(), System::n_vars(), n_vars, System::name(), Elem::parent(), libMeshEnums::SERIAL, System::solution, NumericVector< T >::swap(), System::update_global_solution(), System::variable_name(), System::variable_number(), and DofMap::variable_type().

00106 {
00107   // The current mesh
00108   const MeshBase& mesh = system.get_mesh();
00109 
00110   // The dimensionality of the mesh
00111   const unsigned int dim = mesh.mesh_dimension();
00112   
00113   // The number of variables in the system
00114   const unsigned int n_vars = system.n_vars();
00115   
00116   // The DofMap for this system
00117   const DofMap& dof_map = system.get_dof_map();
00118 
00119   // Resize the error_per_cell vector to be
00120   // the number of elements, initialize it to 0.
00121   error_per_cell.resize (mesh.max_elem_id());
00122   std::fill (error_per_cell.begin(), error_per_cell.end(), 0.);
00123 
00124   // Prepare current_local_solution to localize a non-standard
00125   // solution vector if necessary
00126   if (solution_vector && solution_vector != system.solution.get())
00127     {
00128       NumericVector<Number>* newsol =
00129         const_cast<NumericVector<Number>*>(solution_vector);
00130       System &sys = const_cast<System&>(system);
00131       newsol->swap(*sys.solution);
00132       sys.update();
00133     }
00134   
00135   // Loop over all the variables in the system
00136   for (unsigned int var=0; var<n_vars; var++)
00137     {
00138       // Possibly skip this variable
00139       if (error_norm.weight(var) == 0.0) continue;
00140 
00141       // The (string) name of this variable
00142       const std::string& var_name = system.variable_name(var);
00143       
00144       // The type of finite element to use for this variable
00145       const FEType& fe_type = dof_map.variable_type (var);
00146       
00147       AutoPtr<FEBase> fe (FEBase::build (dim, fe_type));
00148 
00149       // Build an appropriate Gaussian quadrature rule
00150       AutoPtr<QBase> qrule =
00151         fe_type.default_quadrature_rule (dim,
00152                                          _extra_order);
00153 
00154       fe->attach_quadrature_rule (qrule.get());
00155 
00156       // Prepare a global solution and a MeshFunction of the fine system if we need one
00157       AutoPtr<MeshFunction> fine_values;
00158       AutoPtr<NumericVector<Number> > fine_soln = NumericVector<Number>::build();
00159       if (_equation_systems_fine)
00160       {
00161         const System& fine_system = _equation_systems_fine->get_system(system.name());
00162 
00163         std::vector<Number> global_soln;
00164         // FIXME - we're assuming that the fine system solution gets
00165         // used even when a different vector is used for the coarse
00166         // system
00167         fine_system.update_global_solution(global_soln);
00168         fine_soln->init (global_soln.size(), true, SERIAL);
00169         (*fine_soln) = global_soln;
00170 
00171         fine_values = AutoPtr<MeshFunction>
00172           (new MeshFunction(*_equation_systems_fine,
00173                             *fine_soln,
00174                             fine_system.get_dof_map(),
00175                             fine_system.variable_number(var_name)));
00176         fine_values->init();
00177       }
00178 
00179       // Request the data we'll need to compute with
00180       fe->get_JxW();
00181       fe->get_phi();
00182       fe->get_dphi();
00183 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
00184       fe->get_d2phi();
00185 #endif
00186       fe->get_xyz();
00187 
00188       // Iterate over all the active elements in the mesh
00189       // that live on this processor.
00190       MeshBase::const_element_iterator
00191         elem_it  = mesh.active_local_elements_begin();
00192       const MeshBase::const_element_iterator
00193         elem_end = mesh.active_local_elements_end(); 
00194 
00195       for (;elem_it != elem_end; ++elem_it)
00196         {
00197           // e is necessarily an active element on the local processor
00198           const Elem* elem = *elem_it;
00199           const unsigned int e_id = elem->id();
00200 
00201 #ifdef LIBMESH_ENABLE_AMR
00202           // See if the parent of element e has been examined yet;
00203           // if not, we may want to compute the estimator on it
00204           const Elem* parent = elem->parent();
00205 
00206           // We only can compute and only need to compute on
00207           // parents with all active children
00208           bool compute_on_parent = true;
00209           if (!parent || !estimate_parent_error)
00210             compute_on_parent = false;
00211           else
00212             for (unsigned int c=0; c != parent->n_children(); ++c)
00213               if (!parent->child(c)->active())
00214                 compute_on_parent = false;
00215 
00216           if (compute_on_parent &&
00217               !error_per_cell[parent->id()])
00218             {
00219               // Compute a projection onto the parent
00220               DenseVector<Number> Uparent;
00221               FEBase::coarsened_dof_values(*(system.current_local_solution),
00222                                            dof_map, parent, Uparent,
00223                                            var, false);
00224 
00225               error_per_cell[parent->id()] =
00226                 find_squared_element_error(system, var_name, parent, Uparent,
00227                                            fe.get(), fine_values.get());
00228             }
00229 #endif
00230 
00231           // Get the local to global degree of freedom maps
00232           std::vector<unsigned int> dof_indices;
00233           dof_map.dof_indices (elem, dof_indices, var);
00234           DenseVector<Number> Uelem(dof_indices.size());
00235           for (unsigned int i=0; i != dof_indices.size(); ++i)
00236             Uelem(i) = system.current_solution(dof_indices[i]);
00237 
00238           error_per_cell[e_id] =
00239             find_squared_element_error(system, var_name, elem, Uelem,
00240                                        fe.get(), fine_values.get());
00241 
00242         } // End loop over active local elements
00243     } // End loop over variables
00244 
00245 
00246   
00247   // Each processor has now computed the error contribuions
00248   // for its local elements.  We need to sum the vector
00249   // and then take the square-root of each component.  Note
00250   // that we only need to sum if we are running on multiple
00251   // processors, and we only need to take the square-root
00252   // if the value is nonzero.  There will in general be many
00253   // zeros for the inactive elements.
00254 
00255   // First sum the vector of estimated error values
00256   this->reduce_error(error_per_cell);
00257 
00258   // Compute the square-root of each component.
00259   START_LOG("std::sqrt()", "ExactErrorEstimator");
00260   for (unsigned int i=0; i<error_per_cell.size(); i++)
00261     {
00262       
00263       if (error_per_cell[i] != 0.)
00264         {
00265           libmesh_assert (error_per_cell[i] > 0.);
00266           error_per_cell[i] = std::sqrt(error_per_cell[i]);
00267         }
00268       
00269       
00270     }
00271   STOP_LOG("std::sqrt()", "ExactErrorEstimator");
00272 
00273   // If we used a non-standard solution before, now is the time to fix
00274   // the current_local_solution
00275   if (solution_vector && solution_vector != system.solution.get())
00276     {
00277       NumericVector<Number>* newsol =
00278         const_cast<NumericVector<Number>*>(solution_vector);
00279       System &sys = const_cast<System&>(system);
00280       newsol->swap(*sys.solution);
00281       sys.update();
00282     }
00283 }

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

Definition at line 92 of file error_estimator.C.

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

00096 {
00097   SystemNorm old_error_norm = this->error_norm;
00098 
00099   // Find the requested error values from each system
00100   for (unsigned int s = 0; s != equation_systems.n_systems(); ++s)
00101     {
00102       const System &sys = equation_systems.get_system(s);
00103 
00104       unsigned int n_vars = sys.n_vars();
00105 
00106       for (unsigned int v = 0; v != n_vars; ++v)
00107         {
00108           // Only fill in ErrorVectors the user asks for
00109           if (errors_per_cell.find(std::make_pair(&sys, v)) ==
00110               errors_per_cell.end())
00111             continue;
00112 
00113           // Calculate error in only one variable
00114           std::vector<Real> weights(n_vars, 0.0);
00115           weights[v] = 1.0;
00116           this->error_norm =
00117             SystemNorm(std::vector<FEMNormType>(n_vars, old_error_norm.type(0)),
00118                        weights);
00119 
00120           const NumericVector<Number>* solution_vector = NULL;
00121           if (solution_vectors &&
00122               solution_vectors->find(&sys) != solution_vectors->end())
00123             solution_vector = solution_vectors->find(&sys)->second;
00124 
00125           this->estimate_error
00126             (sys, *errors_per_cell[std::make_pair(&sys, v)],
00127              solution_vector, estimate_parent_error);
00128         }
00129     }
00130 
00131   // Restore our old state before returning
00132   this->error_norm = old_error_norm;
00133 }

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

Definition at line 46 of file error_estimator.C.

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

00051 {
00052   SystemNorm old_error_norm = this->error_norm;
00053 
00054   // Sum the error values from each system
00055   for (unsigned int s = 0; s != equation_systems.n_systems(); ++s)
00056     {
00057       ErrorVector system_error_per_cell;
00058       const System &sys = equation_systems.get_system(s);
00059       if (error_norms.find(&sys) == error_norms.end())
00060         this->error_norm = old_error_norm;
00061       else
00062         this->error_norm = error_norms.find(&sys)->second;
00063 
00064       const NumericVector<Number>* solution_vector = NULL;
00065       if (solution_vectors &&
00066           solution_vectors->find(&sys) != solution_vectors->end())
00067         solution_vector = solution_vectors->find(&sys)->second;
00068 
00069       this->estimate_error(sys, system_error_per_cell,
00070                            solution_vector, estimate_parent_error);
00071 
00072       if (s)
00073         {
00074           libmesh_assert(error_per_cell.size() == system_error_per_cell.size());
00075           for (unsigned int i=0; i != error_per_cell.size(); ++i)
00076             error_per_cell[i] += system_error_per_cell[i];
00077         }
00078       else
00079         error_per_cell = system_error_per_cell;
00080     }
00081 
00082   // Restore our old state before returning
00083   this->error_norm = old_error_norm;
00084 }

void ExactErrorEstimator::extra_quadrature_order ( const int  extraorder  )  [inline]

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

Definition at line 127 of file exact_error_estimator.h.

References _extra_order.

00128     { _extra_order = extraorder; }

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

References _equation_systems_fine, _exact_deriv, _exact_hessian, _exact_value, ErrorEstimator::error_norm, FEBase::get_d2phi(), FEBase::get_dphi(), System::get_equation_systems(), FEBase::get_JxW(), FEBase::get_phi(), FEBase::get_xyz(), MeshFunction::gradient(), libMeshEnums::H1, libMeshEnums::H1_SEMINORM, libMeshEnums::H2, libMeshEnums::H2_SEMINORM, MeshFunction::hessian(), libMeshEnums::L2, libmesh_norm(), System::name(), EquationSystems::parameters, FEBase::reinit(), DenseVector< T >::size(), TypeTensor< T >::size_sq(), TypeVector< T >::size_sq(), SystemNorm::type(), and System::variable_number().

00293 {
00294   // The (string) name of this system
00295   const std::string& sys_name = system.name();
00296 
00297   unsigned int var = system.variable_number(var_name);
00298   
00299   const Parameters& parameters = system.get_equation_systems().parameters;
00300 
00301   // reinitialize the element-specific data
00302   // for the current element
00303   fe->reinit (elem);
00304 
00305   // Get the data we need to compute with
00306   const std::vector<Real> &                      JxW          = fe->get_JxW();
00307   const std::vector<std::vector<Real> >&         phi_values   = fe->get_phi();
00308   const std::vector<std::vector<RealGradient> >& dphi_values  = fe->get_dphi();
00309   const std::vector<Point>&                      q_point      = fe->get_xyz();
00310 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
00311   const std::vector<std::vector<RealTensor> >&   d2phi_values = fe->get_d2phi();
00312 #endif
00313 
00314   // The number of shape functions
00315   const unsigned int n_sf = Uelem.size();
00316 
00317   // The number of quadrature points
00318   const unsigned int n_qp = JxW.size();
00319 
00320   Real error_val = 0;
00321 
00322   // Begin the loop over the Quadrature points.
00323   //
00324   for (unsigned int qp=0; qp<n_qp; qp++)
00325     {
00326       // Real u_h = 0.;
00327       // RealGradient grad_u_h;
00328 
00329       Number u_h = 0.;
00330 
00331       Gradient grad_u_h;
00332 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
00333       Tensor grad2_u_h;
00334 #endif
00335 
00336       // Compute solution values at the current
00337       // quadrature point.  This reqiures a sum
00338       // over all the shape functions evaluated
00339       // at the quadrature point.
00340       for (unsigned int i=0; i<n_sf; i++)
00341         {
00342           // Values from current solution.
00343           u_h      += phi_values[i][qp]*Uelem(i);
00344           grad_u_h += dphi_values[i][qp]*Uelem(i);
00345 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
00346           grad2_u_h += d2phi_values[i][qp]*Uelem(i);
00347 #endif
00348         }
00349 
00350       // Compute the value of the error at this quadrature point
00351       if (error_norm.type(var) == L2 ||
00352           error_norm.type(var) == H1 ||
00353           error_norm.type(var) == H2)
00354         {
00355           Number val_error = 0;
00356           if(_exact_value)
00357             val_error = u_h - _exact_value(q_point[qp],parameters,sys_name,var_name);
00358           else if(_equation_systems_fine)
00359             val_error = u_h - (*fine_values)(q_point[qp]);
00360 
00361           // Add the squares of the error to each contribution
00362           error_val += JxW[qp]*libmesh_norm(val_error);
00363         }
00364 
00365       // Compute the value of the error in the gradient at this
00366       // quadrature point
00367       if ((_exact_deriv || _equation_systems_fine) && 
00368           (error_norm.type(var) == H1 ||
00369            error_norm.type(var) == H1_SEMINORM ||
00370            error_norm.type(var) == H2))
00371         {
00372           Gradient grad_error;
00373           if(_exact_deriv)
00374             grad_error = grad_u_h - _exact_deriv(q_point[qp],parameters,sys_name,var_name);
00375           else if(_equation_systems_fine)
00376             grad_error = grad_u_h - fine_values->gradient(q_point[qp]);
00377 
00378           error_val += JxW[qp]*grad_error.size_sq();
00379         }
00380 
00381 
00382 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
00383       // Compute the value of the error in the hessian at this
00384       // quadrature point
00385       if ((_exact_hessian || _equation_systems_fine) && 
00386           (error_norm.type(var) == H2_SEMINORM ||
00387            error_norm.type(var) == H2))
00388         {
00389           Tensor grad2_error;
00390           if(_exact_hessian)
00391             grad2_error = grad2_u_h - _exact_hessian(q_point[qp],parameters,sys_name,var_name);
00392           else if (_equation_systems_fine)
00393             grad2_error = grad2_u_h - fine_values->hessian(q_point[qp]);
00394 
00395           error_val += JxW[qp]*grad2_error.size_sq();
00396         }
00397 #endif
00398 
00399     } // end qp loop
00400 
00401   libmesh_assert (error_val >= 0.);
00402           
00403   return error_val;
00404 }

void ErrorEstimator::reduce_error ( std::vector< float > &  error_per_cell  )  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.

Referenced by PatchRecoveryErrorEstimator::estimate_error(), and JumpErrorEstimator::estimate_error().


Member Data Documentation

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

Definition at line 182 of file exact_error_estimator.h.

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

Gradient(* 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.

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

Tensor(* 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.

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

Number(* 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.

Referenced by attach_reference_solution(), and find_squared_element_error().

Extra order to use for quadrature rule

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

Referenced by KellyErrorEstimator::boundary_side_integration(), DiscontinuityMeasure::boundary_side_integration(), DiscontinuityMeasure::DiscontinuityMeasure(), JumpErrorEstimator::estimate_error(), ErrorEstimator::estimate_errors(), ExactErrorEstimator(), find_squared_element_error(), KellyErrorEstimator::internal_side_integration(), LaplacianErrorEstimator::internal_side_integration(), DiscontinuityMeasure::internal_side_integration(), KellyErrorEstimator::KellyErrorEstimator(), LaplacianErrorEstimator::LaplacianErrorEstimator(), PatchRecoveryErrorEstimator::EstimateError::operator()(), PatchRecoveryErrorEstimator::PatchRecoveryErrorEstimator(), and UniformRefinementEstimator::UniformRefinementEstimator().


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

Site Created By: libMesh Developers
Last modified: November 25 2009 03:44:12.

Hosted By:
SourceForge.net Logo