libMesh::WeightedPatchRecoveryErrorEstimator Class Reference

#include <weighted_patch_recovery_error_estimator.h>

Inheritance diagram for libMesh::WeightedPatchRecoveryErrorEstimator:

Classes

class  EstimateError
 

Public Types

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

Public Member Functions

 WeightedPatchRecoveryErrorEstimator ()=default
 
 WeightedPatchRecoveryErrorEstimator (const WeightedPatchRecoveryErrorEstimator &)=default
 
 WeightedPatchRecoveryErrorEstimator (WeightedPatchRecoveryErrorEstimator &&)=default
 
WeightedPatchRecoveryErrorEstimatoroperator= (const WeightedPatchRecoveryErrorEstimator &)=default
 
WeightedPatchRecoveryErrorEstimatoroperator= (WeightedPatchRecoveryErrorEstimator &&)=default
 
virtual ~WeightedPatchRecoveryErrorEstimator ()=default
 
virtual void estimate_error (const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=nullptr, bool estimate_parent_error=false) override
 
virtual ErrorEstimatorType type () const override
 
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=nullptr, 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=nullptr, bool estimate_parent_error=false)
 

Public Attributes

std::vector< FEMFunctionBase< Number > * > weight_functions
 
unsigned int target_patch_size
 
Patch::PMF patch_growth_strategy
 
SystemNorm error_norm
 

Protected Member Functions

void reduce_error (std::vector< ErrorVectorReal > &error_per_cell, const Parallel::Communicator &comm) 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
Vikram Garg
Date
2012

Definition at line 48 of file weighted_patch_recovery_error_estimator.h.

Member Typedef Documentation

◆ ErrorMap

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

Constructor & Destructor Documentation

◆ WeightedPatchRecoveryErrorEstimator() [1/3]

libMesh::WeightedPatchRecoveryErrorEstimator::WeightedPatchRecoveryErrorEstimator ( )
default

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.

◆ WeightedPatchRecoveryErrorEstimator() [2/3]

libMesh::WeightedPatchRecoveryErrorEstimator::WeightedPatchRecoveryErrorEstimator ( const WeightedPatchRecoveryErrorEstimator )
default

Copy/move ctor, copy/move assignment operator, and destructor are all explicitly defaulted for this class.

◆ WeightedPatchRecoveryErrorEstimator() [3/3]

libMesh::WeightedPatchRecoveryErrorEstimator::WeightedPatchRecoveryErrorEstimator ( WeightedPatchRecoveryErrorEstimator &&  )
default

◆ ~WeightedPatchRecoveryErrorEstimator()

virtual libMesh::WeightedPatchRecoveryErrorEstimator::~WeightedPatchRecoveryErrorEstimator ( )
virtualdefault

Member Function Documentation

◆ estimate_error()

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

This function uses the Patch Recovery error estimate to estimate the error on each cell. The estimated error is output in the vector error_per_cell

Reimplemented from libMesh::PatchRecoveryErrorEstimator.

Definition at line 55 of file weighted_patch_recovery_error_estimator.C.

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

59 {
60  LOG_SCOPE("estimate_error()", "WeightedPatchRecoveryErrorEstimator");
61 
62  // The current mesh
63  const MeshBase & mesh = system.get_mesh();
64 
65  // Resize the error_per_cell vector to be
66  // the number of elements, initialize it to 0.
67  error_per_cell.resize (mesh.max_elem_id());
68  std::fill (error_per_cell.begin(), error_per_cell.end(), 0.);
69 
70  // Prepare current_local_solution to localize a non-standard
71  // solution vector if necessary
72  if (solution_vector && solution_vector != system.solution.get())
73  {
74  NumericVector<Number> * newsol =
75  const_cast<NumericVector<Number> *>(solution_vector);
76  System & sys = const_cast<System &>(system);
77  newsol->swap(*sys.solution);
78  sys.update();
79  }
80 
81  //------------------------------------------------------------
82  // Iterate over all the active elements in the mesh
83  // that live on this processor.
84  Threads::parallel_for (ConstElemRange(mesh.active_local_elements_begin(),
85  mesh.active_local_elements_end(),
86  200),
87  EstimateError(system,
88  *this,
89  error_per_cell)
90  );
91 
92  // Each processor has now computed the error contributions
93  // for its local elements, and error_per_cell contains 0 for all the
94  // non-local elements. Summing the vector will provide the true
95  // value for each element, local or remote
96  this->reduce_error(error_per_cell, system.comm());
97 
98  // If we used a non-standard solution before, now is the time to fix
99  // the current_local_solution
100  if (solution_vector && solution_vector != system.solution.get())
101  {
102  NumericVector<Number> * newsol =
103  const_cast<NumericVector<Number> *>(solution_vector);
104  System & sys = const_cast<System &>(system);
105  newsol->swap(*sys.solution);
106  sys.update();
107  }
108 }
void parallel_for(const Range &range, const Body &body)
Definition: threads_none.h:73
MeshBase & mesh
StoredRange< MeshBase::const_element_iterator, const Elem * > ConstElemRange
Definition: elem_range.h:34
void reduce_error(std::vector< ErrorVectorReal > &error_per_cell, const Parallel::Communicator &comm) const

◆ estimate_errors() [1/2]

void libMesh::ErrorEstimator::estimate_errors ( const EquationSystems equation_systems,
ErrorVector error_per_cell,
const std::map< const System *, SystemNorm > &  error_norms,
const std::map< const System *, const NumericVector< Number > *> *  solution_vectors = nullptr,
bool  estimate_parent_error = false 
)
virtualinherited

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

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

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

Reimplemented in libMesh::UniformRefinementEstimator.

Definition at line 47 of file error_estimator.C.

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

52 {
53  SystemNorm old_error_norm = this->error_norm;
54 
55  // Sum the error values from each system
56  for (unsigned int s = 0; s != equation_systems.n_systems(); ++s)
57  {
58  ErrorVector system_error_per_cell;
59  const System & sys = equation_systems.get_system(s);
60  if (error_norms.find(&sys) == error_norms.end())
61  this->error_norm = old_error_norm;
62  else
63  this->error_norm = error_norms.find(&sys)->second;
64 
65  const NumericVector<Number> * solution_vector = nullptr;
66  if (solution_vectors &&
67  solution_vectors->find(&sys) != solution_vectors->end())
68  solution_vector = solution_vectors->find(&sys)->second;
69 
70  this->estimate_error(sys, system_error_per_cell,
71  solution_vector, estimate_parent_error);
72 
73  if (s)
74  {
75  libmesh_assert_equal_to (error_per_cell.size(), system_error_per_cell.size());
76  for (auto i : index_range(error_per_cell))
77  error_per_cell[i] += system_error_per_cell[i];
78  }
79  else
80  error_per_cell = system_error_per_cell;
81  }
82 
83  // Restore our old state before returning
84  this->error_norm = old_error_norm;
85 }
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Definition: int_range.h:104
virtual void estimate_error(const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=nullptr, bool estimate_parent_error=false)=0

◆ estimate_errors() [2/2]

void libMesh::ErrorEstimator::estimate_errors ( const EquationSystems equation_systems,
ErrorMap errors_per_cell,
const std::map< const System *, const NumericVector< Number > *> *  solution_vectors = nullptr,
bool  estimate_parent_error = false 
)
virtualinherited

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

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

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

◆ operator=() [1/2]

WeightedPatchRecoveryErrorEstimator& libMesh::WeightedPatchRecoveryErrorEstimator::operator= ( const WeightedPatchRecoveryErrorEstimator )
default

◆ operator=() [2/2]

WeightedPatchRecoveryErrorEstimator& libMesh::WeightedPatchRecoveryErrorEstimator::operator= ( WeightedPatchRecoveryErrorEstimator &&  )
default

◆ reduce_error()

void libMesh::ErrorEstimator::reduce_error ( std::vector< ErrorVectorReal > &  error_per_cell,
const Parallel::Communicator comm 
) const
protectedinherited

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

Definition at line 32 of file error_estimator.C.

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

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

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

◆ set_patch_reuse()

void libMesh::PatchRecoveryErrorEstimator::set_patch_reuse ( bool  patch_reuse_flag)
inherited

◆ specpoly()

std::vector< Real > libMesh::PatchRecoveryErrorEstimator::specpoly ( const unsigned int  dim,
const Order  order,
const Point  p,
const unsigned int  matsize 
)
staticprotectedinherited
Returns
The spectral polynomial basis function values at a point (x,y,z).

Definition at line 76 of file patch_recovery_error_estimator.C.

References libMesh::Real.

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

80 {
81  std::vector<Real> psi;
82  psi.reserve(matsize);
83  unsigned int npows = order+1;
84  std::vector<Real> xpow(npows,1.), ypow, zpow;
85  {
86  Real x = p(0);
87  for (auto i : IntRange<int>(1, npows))
88  xpow[i] = xpow[i-1] * x;
89  }
90  if (dim > 1)
91  {
92  Real y = p(1);
93  ypow.resize(npows,1.);
94  for (auto i : IntRange<int>(1, npows))
95  ypow[i] = ypow[i-1] * y;
96  }
97  if (dim > 2)
98  {
99  Real z = p(2);
100  zpow.resize(npows,1.);
101  for (auto i : IntRange<int>(1, npows))
102  zpow[i] = zpow[i-1] * z;
103  }
104 
105  // builds psi vector of form 1 x y z x^2 xy xz y^2 yz z^2 etc..
106  // I haven't added 1D support here
107  for (unsigned int poly_deg=0; poly_deg <= static_cast<unsigned int>(order) ; poly_deg++)
108  { // loop over all polynomials of total degree = poly_deg
109 
110  switch (dim)
111  {
112  // 3D spectral polynomial basis functions
113  case 3:
114  {
115  for (int xexp=poly_deg; xexp >= 0; xexp--) // use an int for xexp since we -- it
116  for (int yexp=poly_deg-xexp; yexp >= 0; yexp--)
117  {
118  int zexp = poly_deg - xexp - yexp;
119  psi.push_back(xpow[xexp]*ypow[yexp]*zpow[zexp]);
120  }
121  break;
122  }
123 
124  // 2D spectral polynomial basis functions
125  case 2:
126  {
127  for (int xexp=poly_deg; xexp >= 0; xexp--) // use an int for xexp since we -- it
128  {
129  int yexp = poly_deg - xexp;
130  psi.push_back(xpow[xexp]*ypow[yexp]);
131  }
132  break;
133  }
134 
135  // 1D spectral polynomial basis functions
136  case 1:
137  {
138  int xexp = poly_deg;
139  psi.push_back(xpow[xexp]);
140  break;
141  }
142 
143  default:
144  libmesh_error_msg("Invalid dimension dim " << dim);
145  }
146  }
147 
148  return psi;
149 }
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ type()

ErrorEstimatorType libMesh::WeightedPatchRecoveryErrorEstimator::type ( ) const
overridevirtual

Friends And Related Function Documentation

◆ EstimateError

friend class EstimateError
friend

Definition at line 119 of file weighted_patch_recovery_error_estimator.h.

Referenced by estimate_error().

Member Data Documentation

◆ error_norm

SystemNorm libMesh::ErrorEstimator::error_norm
inherited

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

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

Definition at line 161 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::ExactErrorEstimator::estimate_error(), libMesh::ErrorEstimator::estimate_errors(), libMesh::ExactErrorEstimator::ExactErrorEstimator(), libMesh::ExactErrorEstimator::find_squared_element_error(), libMesh::LaplacianErrorEstimator::init_context(), libMesh::DiscontinuityMeasure::init_context(), libMesh::KellyErrorEstimator::init_context(), libMesh::LaplacianErrorEstimator::internal_side_integration(), libMesh::DiscontinuityMeasure::internal_side_integration(), libMesh::KellyErrorEstimator::internal_side_integration(), libMesh::KellyErrorEstimator::KellyErrorEstimator(), libMesh::LaplacianErrorEstimator::LaplacianErrorEstimator(), libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::PatchRecoveryErrorEstimator::PatchRecoveryErrorEstimator(), libMesh::JumpErrorEstimator::reinit_sides(), and libMesh::UniformRefinementEstimator::UniformRefinementEstimator().

◆ patch_growth_strategy

Patch::PMF libMesh::PatchRecoveryErrorEstimator::patch_growth_strategy
inherited

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 99 of file patch_recovery_error_estimator.h.

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

◆ patch_reuse

◆ target_patch_size

unsigned int libMesh::PatchRecoveryErrorEstimator::target_patch_size
inherited

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

Definition at line 91 of file patch_recovery_error_estimator.h.

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

◆ weight_functions

std::vector<FEMFunctionBase<Number> *> libMesh::WeightedPatchRecoveryErrorEstimator::weight_functions

Vector of fem function base pointers, the user will fill this in with pointers to the appropriate weight functions.

Definition at line 84 of file weighted_patch_recovery_error_estimator.h.

Referenced by libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()().


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