libMesh::PatchRecoveryErrorEstimator Class Reference

#include <patch_recovery_error_estimator.h>

Inheritance diagram for libMesh::PatchRecoveryErrorEstimator:

Classes

class  EstimateError
 

Public Types

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

Public Member Functions

 PatchRecoveryErrorEstimator ()
 
 ~PatchRecoveryErrorEstimator ()
 
virtual void estimate_error (const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=NULL, bool estimate_parent_error=false)
 
void set_patch_reuse (bool)
 
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

unsigned int target_patch_size
 
Patch::PMF patch_growth_strategy
 
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
 

Static Protected Member Functions

static std::vector< Realspecpoly (const unsigned int dim, const Order order, const Point p, const unsigned int matsize)
 

Protected Attributes

bool patch_reuse
 

Friends

class EstimateError
 

Detailed Description

This class implements the Patch Recovery error indicator.

Author
Varis Carey, Benjamin S. Kirk, 2004.

Definition at line 47 of file patch_recovery_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 109 of file error_estimator.h.

Constructor & Destructor Documentation

libMesh::PatchRecoveryErrorEstimator::PatchRecoveryErrorEstimator ( )
inline

Constructor. Defaults to H1 seminorm. All Hilbert norms and seminorms should be supported now. W1,p and W2,p norms would be natural to support if any contributors make the effort.

Definition at line 56 of file patch_recovery_error_estimator.h.

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

libMesh::PatchRecoveryErrorEstimator::~PatchRecoveryErrorEstimator ( )
inline

Destructor.

Definition at line 66 of file patch_recovery_error_estimator.h.

66 {}

Member Function Documentation

void libMesh::PatchRecoveryErrorEstimator::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 Patch Recovery error estimate to estimate the error on each cell. The estimated error is output in the vector error_per_cell

Implements libMesh::ErrorEstimator.

Reimplemented in libMesh::WeightedPatchRecoveryErrorEstimator.

Definition at line 129 of file patch_recovery_error_estimator.C.

References libMesh::MeshBase::active_local_elements_begin(), libMesh::MeshBase::active_local_elements_end(), libMesh::ParallelObject::comm(), EstimateError, libMesh::System::get_mesh(), libMesh::MeshBase::max_elem_id(), mesh, libMesh::Threads::parallel_for(), libMesh::ErrorEstimator::reduce_error(), libMesh::System::solution, libMesh::START_LOG(), libMesh::STOP_LOG(), and libMesh::NumericVector< T >::swap().

133 {
134  START_LOG("estimate_error()", "PatchRecoveryErrorEstimator");
135 
136  // The current mesh
137  const MeshBase& mesh = system.get_mesh();
138 
139  // Resize the error_per_cell vector to be
140  // the number of elements, initialize it to 0.
141  error_per_cell.resize (mesh.max_elem_id());
142  std::fill (error_per_cell.begin(), error_per_cell.end(), 0.);
143 
144  // Prepare current_local_solution to localize a non-standard
145  // solution vector if necessary
146  if (solution_vector && solution_vector != system.solution.get())
147  {
148  NumericVector<Number>* newsol =
149  const_cast<NumericVector<Number>*>(solution_vector);
150  System &sys = const_cast<System&>(system);
151  newsol->swap(*sys.solution);
152  sys.update();
153  }
154 
155  //------------------------------------------------------------
156  // Iterate over all the active elements in the mesh
157  // that live on this processor.
158  Threads::parallel_for (ConstElemRange(mesh.active_local_elements_begin(),
159  mesh.active_local_elements_end(),
160  200),
161  EstimateError(system,
162  *this,
163  error_per_cell)
164  );
165 
166  // Each processor has now computed the error contributions
167  // for its local elements, and error_per_cell contains 0 for all the
168  // non-local elements. Summing the vector will provide the true
169  // value for each element, local or remote
170  this->reduce_error(error_per_cell, system.comm());
171 
172  // If we used a non-standard solution before, now is the time to fix
173  // the current_local_solution
174  if (solution_vector && solution_vector != system.solution.get())
175  {
176  NumericVector<Number>* newsol =
177  const_cast<NumericVector<Number>*>(solution_vector);
178  System &sys = const_cast<System&>(system);
179  newsol->swap(*sys.solution);
180  sys.update();
181  }
182 
183  STOP_LOG("estimate_error()", "PatchRecoveryErrorEstimator");
184 }
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 
)
virtualinherited

This virtual function can be redefined in derived classes, but by default computes the sum of the error_per_cell for each system in the equation_systems.

Currently this function ignores the error_norm member variable, and uses the function argument error_norms instead.

This function is named estimate_errors instead of estimate_error because otherwise C++ can get confused.

Reimplemented in libMesh::UniformRefinementEstimator.

Definition at line 48 of file error_estimator.C.

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

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

This virtual function can be redefined in derived classes, but by default it calls estimate_error repeatedly to calculate the requested error vectors.

Currently this function ignores the error_norm.weight() values because it calculates each variable's error individually, unscaled.

The user selects which errors get computed by filling a map with error vectors: If errors_per_cell[&system][v] exists, it will be filled with the error values in variable v of system

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

Reimplemented in libMesh::UniformRefinementEstimator.

Definition at line 94 of file error_estimator.C.

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

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

This method takes the local error contributions in error_per_cell from each processor and combines them to get the global error vector.

Definition at line 33 of file error_estimator.C.

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

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

35 {
36  // This function must be run on all processors at once
37  // parallel_object_only();
38 
39  // Each processor has now computed the error contribuions
40  // for its local elements. We may need to sum the vector to
41  // recover the error for each element.
42 
43  comm.sum(error_per_cell);
44 }
void libMesh::PatchRecoveryErrorEstimator::set_patch_reuse ( bool  patch_reuse_flag)

Definition at line 47 of file patch_recovery_error_estimator.C.

References patch_reuse.

48  {
49  patch_reuse = patch_reuse_flag;
50  }
std::vector< Real > libMesh::PatchRecoveryErrorEstimator::specpoly ( const unsigned int  dim,
const Order  order,
const Point  p,
const unsigned int  matsize 
)
staticprotected

Returns the spectral polynomial basis function values at a point x,y,z

Definition at line 54 of file patch_recovery_error_estimator.C.

References libMesh::Real.

Referenced by libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), and libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()().

58 {
59  std::vector<Real> psi;
60  psi.reserve(matsize);
61  Real x = p(0), y=0., z=0.;
62  unsigned int npows = order+1;
63  std::vector<Real> xpow(npows,1.), ypow, zpow;
64  for (unsigned int i=1; i != npows; ++i)
65  xpow[i] = xpow[i-1] * x;
66  if (dim > 1)
67  {
68  y = p(1);
69  ypow.resize(npows,1.);
70  for (unsigned int i=1; i != npows; ++i)
71  ypow[i] = ypow[i-1] * y;
72  }
73  if (dim > 2)
74  {
75  z = p(2);
76  zpow.resize(npows,1.);
77  for (unsigned int i=1; i != npows; ++i)
78  zpow[i] = zpow[i-1] * z;
79  }
80 
81  // builds psi vector of form 1 x y z x^2 xy xz y^2 yz z^2 etc..
82  // I haven't added 1D support here
83  for (unsigned int poly_deg=0; poly_deg <= static_cast<unsigned int>(order) ; poly_deg++)
84  { // loop over all polynomials of total degreee = poly_deg
85 
86  switch (dim)
87  {
88  // 3D spectral polynomial basis functions
89  case 3:
90  {
91  for (int xexp=poly_deg; xexp >= 0; xexp--) // use an int for xexp since we -- it
92  for (int yexp=poly_deg-xexp; yexp >= 0; yexp--)
93  {
94  int zexp = poly_deg - xexp - yexp;
95  psi.push_back(xpow[xexp]*ypow[yexp]*zpow[zexp]);
96  }
97  break;
98  }
99 
100  // 2D spectral polynomial basis functions
101  case 2:
102  {
103  for (int xexp=poly_deg; xexp >= 0; xexp--) // use an int for xexp since we -- it
104  {
105  int yexp = poly_deg - xexp;
106  psi.push_back(xpow[xexp]*ypow[yexp]);
107  }
108  break;
109  }
110 
111  // 1D spectral polynomial basis functions
112  case 1:
113  {
114  int xexp = poly_deg;
115  psi.push_back(xpow[xexp]);
116  break;
117  }
118 
119  default:
120  libmesh_error();
121  }
122  }
123 
124  return psi;
125 }

Friends And Related Function Documentation

friend class EstimateError
friend

Definition at line 140 of file patch_recovery_error_estimator.h.

Referenced by estimate_error().

Member Data Documentation

SystemNorm libMesh::ErrorEstimator::error_norm
inherited

When estimating the error in a single system, the error_norm is used to control the scaling and norm choice for each variable. Not all estimators will support all norm choices. The default scaling is for all variables to be weighted equally. The default norm choice depends on the error estimator.

Part of this functionality was supported via component_scale and sobolev_order in older libMesh versions, and a small part was supported via component_mask in even older versions. Hopefully the encapsulation here will allow us to avoid changing this API again.

Definition at line 141 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(), libMesh::ExactErrorEstimator::ExactErrorEstimator(), libMesh::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()(), PatchRecoveryErrorEstimator(), and libMesh::UniformRefinementEstimator::UniformRefinementEstimator().

Patch::PMF libMesh::PatchRecoveryErrorEstimator::patch_growth_strategy

The PatchErrorEstimator will use this pointer to a Patch member function when growing patches. The default strategy used is Patch::add_local_face_neighbors. Patch::add_local_point_neighbors may be more reliable but slower.

Definition at line 92 of file patch_recovery_error_estimator.h.

Referenced by libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), and libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()().

unsigned int libMesh::PatchRecoveryErrorEstimator::target_patch_size

The PatchErrorEstimator will build patches of at least this many elements to perform estimates

Definition at line 84 of file patch_recovery_error_estimator.h.

Referenced by libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), and libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()().


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

Site Created By: libMesh Developers
Last modified: February 07 2014 16:58:01 UTC

Hosted By:
SourceForge.net Logo