libMesh::KellyErrorEstimator Class Reference

#include <kelly_error_estimator.h>

Inheritance diagram for libMesh::KellyErrorEstimator:

List of all members.

Public Types

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

Public Member Functions

 KellyErrorEstimator ()
 ~KellyErrorEstimator ()
void attach_flux_bc_function (std::pair< bool, Real > fptr(const System &system, const Point &p, const std::string &var_name))
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

bool scale_by_n_flux_faces
SystemNorm error_norm

Protected Member Functions

virtual void initialize (const System &system, ErrorVector &error_per_cell, bool estimate_parent_error)
virtual void internal_side_integration ()
virtual bool boundary_side_integration ()
void reinit_sides ()
float coarse_n_flux_faces_increment ()
void reduce_error (std::vector< float > &error_per_cell) const

Protected Attributes

const Systemmy_system
std::pair< bool, Real >(* _bc_function )(const System &system, const Point &p, const std::string &var_name)
bool integrate_boundary_sides
const Elemfine_elem
const Elemcoarse_elem
Real fine_error
Real coarse_error
unsigned int fine_side
unsigned int var
DenseVector< NumberUfine
DenseVector< NumberUcoarse
AutoPtr< FEBasefe_fine
AutoPtr< FEBasefe_coarse

Detailed Description

This class implements the Kelly error indicator which is based on the flux jumps between elements. See the JumpErrorEstimator class for most user APIs

Full BibteX reference:

@Article{Kelly83error,
  author = {D.~W.~Kelly and J.~P.~Gago and O.~C.~Zienkiewicz and I.~Babuska},
  title  = {{A posteriori error analysis and adaptive
             processes in the finite element method: Part I Error analysis}},
  journal = {Int. J. Num. Meth. Engng.},
  volume  = {19},
  pages   = {1593--1619},
  year    = {1983}
}
Author:
Benjamin S. Kirk, 2003.

Definition at line 61 of file kelly_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 107 of file error_estimator.h.


Constructor & Destructor Documentation

libMesh::KellyErrorEstimator::KellyErrorEstimator (  )  [inline]

Constructor. Responsible for initializing the _bc_function function pointer to NULL. Defaults to H1 seminorm; changes to system norm are ignored.

Definition at line 70 of file kelly_error_estimator.h.

References libMesh::ErrorEstimator::error_norm, and libMeshEnums::H1_SEMINORM.

00070                         : _bc_function(NULL)
00071   { error_norm = H1_SEMINORM; }

libMesh::KellyErrorEstimator::~KellyErrorEstimator (  )  [inline]

Destructor.

Definition at line 76 of file kelly_error_estimator.h.

00076 {}


Member Function Documentation

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

Register a user function to use in computing the flux BCs. The return value is std::pair<bool, Real>

Definition at line 169 of file kelly_error_estimator.C.

References _bc_function, and libMesh::JumpErrorEstimator::integrate_boundary_sides.

00172 {
00173   _bc_function = fptr;
00174 
00175 // We may be turning boundary side integration on or off
00176   if (fptr)
00177     integrate_boundary_sides = true;
00178   else
00179     integrate_boundary_sides = false;
00180 }

bool libMesh::KellyErrorEstimator::boundary_side_integration (  )  [protected, virtual]

The function which calculates a normal derivative jump based error term on a boundary side. Returns true if the flux bc function is in fact defined on the current side.

Reimplemented from libMesh::JumpErrorEstimator.

Definition at line 100 of file kelly_error_estimator.C.

References _bc_function, libMesh::TypeVector< T >::add_scaled(), libMesh::ErrorEstimator::error_norm, libMesh::JumpErrorEstimator::fe_fine, libMesh::JumpErrorEstimator::fine_elem, libMesh::JumpErrorEstimator::fine_error, libMesh::Elem::hmax(), my_system, libMesh::TensorTools::norm_sq(), libMesh::Real, libMesh::DenseVector< T >::size(), libMesh::JumpErrorEstimator::Ufine, libMesh::JumpErrorEstimator::var, libMesh::System::variable_name(), and libMesh::SystemNorm::weight().

00101 {
00102   const std::string &var_name = my_system->variable_name(var);
00103   const unsigned int n_fine_dofs = Ufine.size();
00104 
00105   std::vector<std::vector<RealGradient> > dphi_fine = fe_fine->get_dphi();
00106   std::vector<Point> face_normals = fe_fine->get_normals();
00107   std::vector<Real> JxW_face = fe_fine->get_JxW();
00108   std::vector<Point> qface_point = fe_fine->get_xyz();
00109 
00110   // The reinitialization also recomputes the locations of
00111   // the quadrature points on the side.  By checking if the
00112   // first quadrature point on the side is on a flux boundary
00113   // for a particular variable, we will determine if the whole
00114   // element is on a flux boundary (assuming quadrature points
00115   // are strictly contained in the side).
00116   if (this->_bc_function(*my_system, qface_point[0], var_name).first)
00117     {
00118       const Real h = fine_elem->hmax();
00119 
00120       // The number of quadrature points
00121       const unsigned int n_qp = fe_fine->n_quadrature_points();
00122 
00123       // The error contribution from this face
00124       Real error = 1.e-30;
00125 
00126       // loop over the integration points on the face.
00127       for (unsigned int qp=0; qp<n_qp; qp++)
00128         {
00129           // Value of the imposed flux BC at this quadrature point.
00130           const std::pair<bool,Real> flux_bc =
00131             this->_bc_function(*my_system, qface_point[qp], var_name);
00132 
00133           // Be sure the BC function still thinks we're on the
00134           // flux boundary.
00135           libmesh_assert_equal_to (flux_bc.first, true);
00136 
00137           // The solution gradient from each element
00138           Gradient grad_fine;
00139 
00140           // Compute the solution gradient on element e
00141           for (unsigned int i=0; i != n_fine_dofs; i++)
00142             grad_fine.add_scaled (dphi_fine[i][qp], Ufine(i));
00143 
00144           // The difference between the desired BC and the approximate solution.
00145           const Number jump = flux_bc.second - grad_fine*face_normals[qp];
00146 
00147           // The flux jump squared.  If using complex numbers,
00148           // TensorTools::norm_sq(z) returns |z|^2, where |z| is the modulus of z.
00149           const Real jump2 = TensorTools::norm_sq(jump);
00150 
00151           // Integrate the error on the face.  The error is
00152           // scaled by an additional power of h, where h is
00153           // the maximum side length for the element.  This
00154           // arises in the definition of the indicator.
00155           error += JxW_face[qp]*jump2;
00156 
00157         } // End quadrature point loop
00158 
00159       fine_error = error*h*error_norm.weight(var);
00160 
00161       return true;
00162     } // end if side on flux boundary
00163   return false;
00164 }

float libMesh::JumpErrorEstimator::coarse_n_flux_faces_increment (  )  [protected, inherited]

A utility function to correctly increase n_flux_faces for the coarse element

Definition at line 445 of file jump_error_estimator.C.

References libMesh::JumpErrorEstimator::coarse_elem, libMesh::Elem::dim(), libMesh::JumpErrorEstimator::fine_elem, and libMesh::Elem::level().

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

00446 {
00447   // Keep track of the number of internal flux sides found on each
00448   // element
00449   unsigned int dim = coarse_elem->dim();
00450 
00451   const unsigned int divisor =
00452     1 << (dim-1)*(fine_elem->level() - coarse_elem->level());
00453 
00454   // With a difference of n levels between fine and coarse elements,
00455   // we compute a fractional flux face for the coarse element by adding:
00456   // 1/2^n in 2D
00457   // 1/4^n in 3D
00458   // each time.  This code will get hit 2^n times in 2D and 4^n
00459   // times in 3D so that the final flux face count for the coarse
00460   // element will be an integer value.
00461 
00462   return 1.0f / static_cast<float>(divisor);
00463 }

void libMesh::JumpErrorEstimator::estimate_error ( const System system,
ErrorVector error_per_cell,
const NumericVector< Number > *  solution_vector = NULL,
bool  estimate_parent_error = false 
) [virtual, inherited]

This function uses the derived class's jump error estimate formula to estimate the error on each cell. The estimated error is output in the vector error_per_cell

Implements libMesh::ErrorEstimator.

Definition at line 54 of file jump_error_estimator.C.

References libMesh::Elem::active(), libMesh::MeshBase::active_local_elements_begin(), libMesh::MeshBase::active_local_elements_end(), libMesh::JumpErrorEstimator::boundary_side_integration(), libMesh::FEGenericBase< Real >::build(), libMesh::Elem::child(), libMesh::JumpErrorEstimator::coarse_elem, libMesh::JumpErrorEstimator::coarse_error, libMesh::JumpErrorEstimator::coarse_n_flux_faces_increment(), libMesh::FEGenericBase< Real >::coarsened_dof_values(), libMesh::System::current_solution(), libMesh::FEType::default_quadrature_order(), libMesh::DofMap::dof_indices(), libMesh::ErrorEstimator::error_norm, libMesh::ErrorVectorReal, libMesh::JumpErrorEstimator::fe_coarse, libMesh::JumpErrorEstimator::fe_fine, libMesh::JumpErrorEstimator::fine_elem, libMesh::JumpErrorEstimator::fine_error, libMesh::JumpErrorEstimator::fine_side, libMesh::AutoPtr< Tp >::get(), libMesh::System::get_dof_map(), libMesh::System::get_mesh(), libMesh::DofObject::id(), libMesh::JumpErrorEstimator::initialize(), libMesh::JumpErrorEstimator::integrate_boundary_sides, libMesh::JumpErrorEstimator::internal_side_integration(), libMesh::Elem::level(), libMesh::MeshBase::max_elem_id(), mesh, libMesh::MeshBase::mesh_dimension(), libMesh::Elem::n_children(), libMesh::Elem::n_neighbors(), libMesh::System::n_vars(), n_vars, libMesh::Elem::neighbor(), libMesh::Elem::parent(), libMesh::ErrorEstimator::reduce_error(), libMesh::JumpErrorEstimator::reinit_sides(), libMesh::DenseVector< T >::resize(), libMesh::JumpErrorEstimator::scale_by_n_flux_faces, libMesh::System::solution, libMesh::NumericVector< T >::swap(), libMesh::JumpErrorEstimator::Ucoarse, libMesh::JumpErrorEstimator::Ufine, libMesh::JumpErrorEstimator::var, libMesh::DofMap::variable_type(), and libMesh::SystemNorm::weight().

00058 {
00059   START_LOG("estimate_error()", "JumpErrorEstimator");
00060   /*
00061 
00062   Conventions for assigning the direction of the normal:
00063 
00064   - e & f are global element ids
00065 
00066   Case (1.) Elements are at the same level, e<f
00067             Compute the flux jump on the face and
00068             add it as a contribution to error_per_cell[e]
00069             and error_per_cell[f]
00070 
00071                    ----------------------
00072                   |           |          |
00073                   |           |    f     |
00074                   |           |          |
00075                   |    e      |---> n    |
00076                   |           |          |
00077                   |           |          |
00078                    ----------------------
00079 
00080 
00081    Case (2.) The neighbor is at a higher level.
00082              Compute the flux jump on e's face and
00083              add it as a contribution to error_per_cell[e]
00084              and error_per_cell[f]
00085 
00086                    ----------------------
00087                   |     |     |          |
00088                   |     |  e  |---> n    |
00089                   |     |     |          |
00090                   |-----------|    f     |
00091                   |     |     |          |
00092                   |     |     |          |
00093                   |     |     |          |
00094                    ----------------------
00095   */
00096 
00097   // The current mesh
00098   const MeshBase& mesh = system.get_mesh();
00099 
00100   // The dimensionality of the mesh
00101   const unsigned int dim = mesh.mesh_dimension();
00102 
00103   // The number of variables in the system
00104   const unsigned int n_vars = system.n_vars();
00105 
00106   // The DofMap for this system
00107   const DofMap& dof_map = system.get_dof_map();
00108 
00109   // Resize the error_per_cell vector to be
00110   // the number of elements, initialize it to 0.
00111   error_per_cell.resize (mesh.max_elem_id());
00112   std::fill (error_per_cell.begin(), error_per_cell.end(), 0.);
00113 
00114   // Declare a vector of floats which is as long as
00115   // error_per_cell above, and fill with zeros.  This vector will be
00116   // used to keep track of the number of edges (faces) on each active
00117   // element which are either:
00118   // 1) an internal edge
00119   // 2) an edge on a Neumann boundary for which a boundary condition
00120   //    function has been specified.
00121   // The error estimator can be scaled by the number of flux edges (faces)
00122   // which the element actually has to obtain a more uniform measure
00123   // of the error.  Use floats instead of ints since in case 2 (above)
00124   // f gets 1/2 of a flux face contribution from each of his
00125   // neighbors
00126   std::vector<float> n_flux_faces (error_per_cell.size());
00127 
00128   // Prepare current_local_solution to localize a non-standard
00129   // solution vector if necessary
00130   if (solution_vector && solution_vector != system.solution.get())
00131     {
00132       NumericVector<Number>* newsol =
00133         const_cast<NumericVector<Number>*>(solution_vector);
00134       System &sys = const_cast<System&>(system);
00135       newsol->swap(*sys.solution);
00136       sys.update();
00137     }
00138 
00139   // Loop over all the variables in the system
00140   for (var=0; var<n_vars; var++)
00141     {
00142       // Possibly skip this variable
00143       if (error_norm.weight(var) == 0.0) continue;
00144 
00145       // The type of finite element to use for this variable
00146       const FEType& fe_type = dof_map.variable_type (var);
00147 
00148       // Finite element objects for the same face from
00149       // different sides
00150       fe_fine = FEBase::build (dim, fe_type);
00151       fe_coarse = FEBase::build (dim, fe_type);
00152 
00153       // Build an appropriate Gaussian quadrature rule
00154       QGauss qrule (dim-1, fe_type.default_quadrature_order());
00155 
00156       // Tell the finite element for the fine element about the quadrature
00157       // rule.  The finite element for the coarse element need not know about it
00158       fe_fine->attach_quadrature_rule (&qrule);
00159 
00160       // By convention we will always do the integration
00161       // on the face of element e.  We'll need its Jacobian values and
00162       // physical point locations, at least
00163       fe_fine->get_JxW();
00164       fe_fine->get_xyz();
00165 
00166       // Our derived classes may want to do some initialization here
00167       this->initialize(system, error_per_cell, estimate_parent_error);
00168 
00169       // The global DOF indices for elements e & f
00170       std::vector<dof_id_type> dof_indices_fine;
00171       std::vector<dof_id_type> dof_indices_coarse;
00172 
00173 
00174 
00175       // Iterate over all the active elements in the mesh
00176       // that live on this processor.
00177       MeshBase::const_element_iterator       elem_it  = mesh.active_local_elements_begin();
00178       const MeshBase::const_element_iterator elem_end = mesh.active_local_elements_end();
00179 
00180       for (; elem_it != elem_end; ++elem_it)
00181         {
00182           // e is necessarily an active element on the local processor
00183           const Elem* e = *elem_it;
00184           const dof_id_type e_id = e->id();
00185 
00186 #ifdef LIBMESH_ENABLE_AMR
00187           // See if the parent of element e has been examined yet;
00188           // if not, we may want to compute the estimator on it
00189           const Elem* parent = e->parent();
00190 
00191           // We only can compute and only need to compute on
00192           // parents with all active children
00193           bool compute_on_parent = true;
00194           if (!parent || !estimate_parent_error)
00195             compute_on_parent = false;
00196           else
00197             for (unsigned int c=0; c != parent->n_children(); ++c)
00198               if (!parent->child(c)->active())
00199                 compute_on_parent = false;
00200 
00201           if (compute_on_parent &&
00202               !error_per_cell[parent->id()])
00203             {
00204               // Compute a projection onto the parent
00205               DenseVector<Number> Uparent;
00206               FEBase::coarsened_dof_values(*(system.solution),
00207                                            dof_map, parent, Uparent,
00208                                            var, false);
00209 
00210               // Loop over the neighbors of the parent
00211               for (unsigned int n_p=0; n_p<parent->n_neighbors(); n_p++)
00212                 {
00213                   if (parent->neighbor(n_p) != NULL) // parent has a neighbor here
00214                     {
00215                       // Find the active neighbors in this direction
00216                       std::vector<const Elem*> active_neighbors;
00217                       parent->neighbor(n_p)->
00218                         active_family_tree_by_neighbor(active_neighbors,
00219                                                        parent);
00220                       // Compute the flux to each active neighbor
00221                       for (unsigned int a=0;
00222                            a != active_neighbors.size(); ++a)
00223                         {
00224                           const Elem *f = active_neighbors[a];
00225                       // FIXME - what about when f->level <
00226                       // parent->level()??
00227                           if (f->level() >= parent->level())
00228                             {
00229                               fine_elem = f;
00230                               coarse_elem = parent;
00231                               Ucoarse = Uparent;
00232 
00233                               dof_map.dof_indices (fine_elem, dof_indices_fine, var);
00234                               const unsigned int n_dofs_fine =
00235                                 libmesh_cast_int<unsigned int>(dof_indices_fine.size());
00236                               Ufine.resize(n_dofs_fine);
00237 
00238                               for (unsigned int i=0; i<n_dofs_fine; i++)
00239                                 Ufine(i) = system.current_solution(dof_indices_fine[i]);
00240                               this->reinit_sides();
00241                               this->internal_side_integration();
00242 
00243                               error_per_cell[fine_elem->id()] +=
00244                                 static_cast<ErrorVectorReal>(fine_error);
00245                               error_per_cell[coarse_elem->id()] += 
00246                                 static_cast<ErrorVectorReal>(coarse_error);
00247 
00248                               // Keep track of the number of internal flux
00249                               // sides found on each element
00250                               n_flux_faces[fine_elem->id()]++;
00251                               n_flux_faces[coarse_elem->id()] += this->coarse_n_flux_faces_increment();
00252                             }
00253                         }
00254                     }
00255                   else if (integrate_boundary_sides)
00256                     {
00257                       fine_elem = parent;
00258                       Ufine = Uparent;
00259 
00260                       // Reinitialize shape functions on the fine element side
00261                       fe_fine->reinit (fine_elem, fine_side);
00262 
00263                       if (this->boundary_side_integration())
00264                         {
00265                           error_per_cell[fine_elem->id()] +=
00266                             static_cast<ErrorVectorReal>(fine_error);
00267                           n_flux_faces[fine_elem->id()]++;
00268                         }
00269                     }
00270                 }
00271             }
00272 #endif // #ifdef LIBMESH_ENABLE_AMR
00273 
00274           // If we do any more flux integration, e will be the fine element
00275           fine_elem = e;
00276 
00277           // Loop over the neighbors of element e
00278           for (unsigned int n_e=0; n_e<e->n_neighbors(); n_e++)
00279             {
00280               fine_side = n_e;
00281 
00282               if (e->neighbor(n_e) != NULL) // e is not on the boundary
00283                 {
00284                   const Elem* f           = e->neighbor(n_e);
00285                   const dof_id_type f_id = f->id();
00286 
00287                   // Compute flux jumps if we are in case 1 or case 2.
00288                   if ((f->active() && (f->level() == e->level()) && (e_id < f_id))
00289                       || (f->level() < e->level()))
00290                     {
00291                       // f is now the coarse element
00292                       coarse_elem = f;
00293 
00294                       // Get the DOF indices for the two elements
00295                       dof_map.dof_indices (fine_elem, dof_indices_fine, var);
00296                       dof_map.dof_indices (coarse_elem, dof_indices_coarse, var);
00297 
00298                       // The number of DOFS on each element
00299                       const unsigned int n_dofs_fine = 
00300                         libmesh_cast_int<unsigned int>(dof_indices_fine.size());
00301                       const unsigned int n_dofs_coarse =
00302                         libmesh_cast_int<unsigned int>(dof_indices_coarse.size());
00303                       Ufine.resize(n_dofs_fine);
00304                       Ucoarse.resize(n_dofs_coarse);
00305 
00306                       // The local solutions on each element
00307                       for (unsigned int i=0; i<n_dofs_fine; i++)
00308                         Ufine(i) = system.current_solution(dof_indices_fine[i]);
00309                       for (unsigned int i=0; i<n_dofs_coarse; i++)
00310                         Ucoarse(i) = system.current_solution(dof_indices_coarse[i]);
00311 
00312                       this->reinit_sides();
00313                       this->internal_side_integration();
00314 
00315                       error_per_cell[fine_elem->id()] +=
00316                         static_cast<ErrorVectorReal>(fine_error);
00317                       error_per_cell[coarse_elem->id()] +=
00318                         static_cast<ErrorVectorReal>(coarse_error);
00319 
00320                       // Keep track of the number of internal flux
00321                       // sides found on each element
00322                       n_flux_faces[fine_elem->id()]++;
00323                       n_flux_faces[coarse_elem->id()] += this->coarse_n_flux_faces_increment();
00324                     } // end if (case1 || case2)
00325                 } // if (e->neigbor(n_e) != NULL)
00326 
00327               // Otherwise, e is on the boundary.  If it happens to
00328               // be on a Dirichlet boundary, we need not do anything.
00329               // On the other hand, if e is on a Neumann (flux) boundary
00330               // with grad(u).n = g, we need to compute the additional residual
00331               // (h * \int |g - grad(u_h).n|^2 dS)^(1/2).
00332               // We can only do this with some knowledge of the boundary
00333               // conditions, i.e. the user must have attached an appropriate
00334               // BC function.
00335               else
00336                 {
00337                   if (integrate_boundary_sides)
00338                     {
00339                       // Reinitialize shape functions on the fine element side
00340                       fe_fine->reinit (fine_elem, fine_side);
00341 
00342                       // Get the DOF indices
00343                       dof_map.dof_indices (fine_elem, dof_indices_fine, var);
00344 
00345                       // The number of DOFS on each element
00346                       const unsigned int n_dofs_fine =
00347                         libmesh_cast_int<unsigned int>(dof_indices_fine.size());
00348                       Ufine.resize(n_dofs_fine);
00349 
00350                       for (unsigned int i=0; i<n_dofs_fine; i++)
00351                         Ufine(i) = system.current_solution(dof_indices_fine[i]);
00352 
00353                       if (this->boundary_side_integration())
00354                         {
00355                           error_per_cell[fine_elem->id()] +=
00356                             static_cast<ErrorVectorReal>(fine_error);
00357                           n_flux_faces[fine_elem->id()]++;
00358                         }
00359                     } // end if _bc_function != NULL
00360                 } // end if (e->neighbor(n_e) == NULL)
00361             } // end loop over neighbors
00362         } // End loop over active local elements
00363     } // End loop over variables
00364 
00365 
00366 
00367   // Each processor has now computed the error contribuions
00368   // for its local elements.  We need to sum the vector
00369   // and then take the square-root of each component.  Note
00370   // that we only need to sum if we are running on multiple
00371   // processors, and we only need to take the square-root
00372   // if the value is nonzero.  There will in general be many
00373   // zeros for the inactive elements.
00374 
00375   // First sum the vector of estimated error values
00376   this->reduce_error(error_per_cell);
00377 
00378   // Compute the square-root of each component.
00379   for (std::size_t i=0; i<error_per_cell.size(); i++)
00380     if (error_per_cell[i] != 0.)
00381       error_per_cell[i] = std::sqrt(error_per_cell[i]);
00382 
00383 
00384   if (this->scale_by_n_flux_faces)
00385     {
00386       // Sum the vector of flux face counts
00387       this->reduce_error(n_flux_faces);
00388 
00389       // Sanity check: Make sure the number of flux faces is
00390       // always an integer value
00391 #ifdef DEBUG
00392       for (unsigned int i=0; i<n_flux_faces.size(); ++i)
00393         libmesh_assert_equal_to (n_flux_faces[i], static_cast<float>(static_cast<unsigned int>(n_flux_faces[i])) );
00394 #endif
00395 
00396       // Scale the error by the number of flux faces for each element
00397       for (unsigned int i=0; i<n_flux_faces.size(); ++i)
00398         {
00399           if (n_flux_faces[i] == 0.0) // inactive or non-local element
00400             continue;
00401 
00402           //libMesh::out << "Element " << i << " has " << n_flux_faces[i] << " flux faces." << std::endl;
00403           error_per_cell[i] /= static_cast<ErrorVectorReal>(n_flux_faces[i]);
00404         }
00405     }
00406 
00407   // If we used a non-standard solution before, now is the time to fix
00408   // the current_local_solution
00409   if (solution_vector && solution_vector != system.solution.get())
00410     {
00411       NumericVector<Number>* newsol =
00412         const_cast<NumericVector<Number>*>(solution_vector);
00413       System &sys = const_cast<System&>(system);
00414       newsol->swap(*sys.solution);
00415       sys.update();
00416     }
00417 
00418   STOP_LOG("estimate_error()", "JumpErrorEstimator");
00419 }

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 93 of file error_estimator.C.

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

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

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 47 of file error_estimator.C.

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

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

void libMesh::KellyErrorEstimator::initialize ( const System system,
ErrorVector error_per_cell,
bool  estimate_parent_error 
) [protected, virtual]

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

Reimplemented from libMesh::JumpErrorEstimator.

Definition at line 43 of file kelly_error_estimator.C.

References libMesh::JumpErrorEstimator::fe_coarse, libMesh::JumpErrorEstimator::fe_fine, and my_system.

00046 {
00047   // Hang onto the system - we may need it for variable names later.
00048   my_system = &system;
00049 
00050   // We'll need gradients and normal vectors for flux jump computation
00051   fe_fine->get_dphi();
00052   fe_fine->get_normals();
00053   fe_coarse->get_dphi();
00054 }

void libMesh::KellyErrorEstimator::internal_side_integration (  )  [protected, virtual]

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

Implements libMesh::JumpErrorEstimator.

Definition at line 59 of file kelly_error_estimator.C.

References libMesh::TypeVector< T >::add_scaled(), libMesh::JumpErrorEstimator::coarse_elem, libMesh::JumpErrorEstimator::coarse_error, libMesh::ErrorEstimator::error_norm, libMesh::JumpErrorEstimator::fe_coarse, libMesh::JumpErrorEstimator::fe_fine, libMesh::JumpErrorEstimator::fine_elem, libMesh::JumpErrorEstimator::fine_error, libMesh::Elem::hmax(), libMesh::TensorTools::norm_sq(), libMesh::Real, libMesh::DenseVector< T >::size(), libMesh::JumpErrorEstimator::Ucoarse, libMesh::JumpErrorEstimator::Ufine, libMesh::JumpErrorEstimator::var, and libMesh::SystemNorm::weight().

00060 {
00061   Real error = 1.e-30;
00062   unsigned int n_qp = fe_fine->n_quadrature_points();
00063   unsigned int n_fine_dofs = Ufine.size();
00064   unsigned int n_coarse_dofs = Ucoarse.size();
00065 
00066   std::vector<std::vector<RealGradient> > dphi_coarse = fe_coarse->get_dphi();
00067   std::vector<std::vector<RealGradient> > dphi_fine = fe_fine->get_dphi();
00068   std::vector<Point> face_normals = fe_fine->get_normals();
00069   std::vector<Real> JxW_face = fe_fine->get_JxW();
00070 
00071   for (unsigned int qp=0; qp != n_qp; ++qp)
00072     {
00073       // Calculate solution gradients on fine and coarse elements
00074       // at this quadrature point
00075       Gradient grad_fine, grad_coarse;
00076       for (unsigned int i=0; i != n_coarse_dofs; ++i)
00077         grad_coarse.add_scaled (dphi_coarse[i][qp], Ucoarse(i));
00078 
00079       for (unsigned int i=0; i != n_fine_dofs; ++i)
00080         grad_fine.add_scaled (dphi_fine[i][qp], Ufine(i));
00081 
00082       // Find the jump in the normal derivative
00083       // at this quadrature point
00084       const Number jump = (grad_fine - grad_coarse)*face_normals[qp];
00085       const Real jump2 = TensorTools::norm_sq(jump);
00086 
00087       // Accumulate the jump integral
00088       error += JxW_face[qp] * jump2;
00089     }
00090 
00091   // Add the h-weighted jump integral to each error term
00092   fine_error =
00093     error * fine_elem->hmax() * error_norm.weight(var);
00094   coarse_error =
00095     error * coarse_elem->hmax() * error_norm.weight(var);
00096 }

void libMesh::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 libMesh::WeightedPatchRecoveryErrorEstimator::estimate_error(), libMesh::PatchRecoveryErrorEstimator::estimate_error(), and libMesh::JumpErrorEstimator::estimate_error().

void libMesh::JumpErrorEstimator::reinit_sides (  )  [protected, inherited]

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

Definition at line 424 of file jump_error_estimator.C.

References libMesh::JumpErrorEstimator::coarse_elem, libMesh::Elem::dim(), libMesh::JumpErrorEstimator::fe_coarse, libMesh::JumpErrorEstimator::fe_fine, libMesh::JumpErrorEstimator::fine_elem, libMesh::JumpErrorEstimator::fine_side, and libMesh::FEInterface::inverse_map().

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

00425 {
00426   // The master quadrature point locations on the coarse element
00427   std::vector<Point> qp_coarse;
00428 
00429   // Reinitialize shape functions on the fine element side
00430   fe_fine->reinit (fine_elem, fine_side);
00431 
00432   // Get the physical locations of the fine element quadrature points
00433   std::vector<Point> qface_point = fe_fine->get_xyz();
00434 
00435   // Find their locations on the coarse element
00436   FEInterface::inverse_map (coarse_elem->dim(), fe_coarse->get_fe_type(),
00437                             coarse_elem, qface_point, qp_coarse);
00438 
00439   // Calculate the coarse element shape functions at those locations
00440   fe_coarse->reinit (coarse_elem, &qp_coarse);
00441 }


Member Data Documentation

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

Pointer to function that returns BC information.

Referenced by attach_flux_bc_function(), and boundary_side_integration().

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

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

unsigned int libMesh::JumpErrorEstimator::fine_side [protected, inherited]

Which side of the fine element is this?

Definition at line 142 of file jump_error_estimator.h.

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

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

Definition at line 127 of file jump_error_estimator.h.

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

A pointer to the current System

Definition at line 112 of file kelly_error_estimator.h.

Referenced by boundary_side_integration(), and initialize().

This boolean flag allows you to scale the error indicator result for each element by the number of "flux faces" the element actually has. This tends to weight more evenly cells which are on the boundaries and thus have fewer contributions to their flux. The value is initialized to false, simply set it to true if you want to use the feature.

Definition at line 88 of file jump_error_estimator.h.

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


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

Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:27 UTC

Hosted By:
SourceForge.net Logo