libMesh::JumpErrorEstimator Class Referenceabstract

#include <jump_error_estimator.h>

Inheritance diagram for libMesh::JumpErrorEstimator:

Public Types

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

Public Member Functions

 JumpErrorEstimator ()
 
virtual ~JumpErrorEstimator ()
 
virtual void estimate_error (const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=libmesh_nullptr, bool estimate_parent_error=false) libmesh_override
 
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=libmesh_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=libmesh_nullptr, bool estimate_parent_error=false)
 
virtual ErrorEstimatorType type () const =0
 

Public Attributes

bool scale_by_n_flux_faces
 
SystemNorm error_norm
 

Protected Member Functions

void reinit_sides ()
 
float coarse_n_flux_faces_increment ()
 
virtual void init_context (FEMContext &c)
 
virtual void internal_side_integration ()=0
 
virtual bool boundary_side_integration ()
 
void reduce_error (std::vector< ErrorVectorReal > &error_per_cell, const Parallel::Communicator &comm LIBMESH_CAN_DEFAULT_TO_COMMWORLD) const
 

Protected Attributes

bool integrate_boundary_sides
 
UniquePtr< FEMContextfine_context
 
UniquePtr< FEMContextcoarse_context
 
Real fine_error
 
Real coarse_error
 
unsigned int var
 

Detailed Description

This abstract base class implements utility functions for error estimators which are based on integrated jumps between elements.

Author
Roy H. Stogner
Date
2006

Definition at line 48 of file jump_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 112 of file error_estimator.h.

Constructor & Destructor Documentation

libMesh::JumpErrorEstimator::JumpErrorEstimator ( )
inline
virtual libMesh::JumpErrorEstimator::~JumpErrorEstimator ( )
inlinevirtual

Destructor.

Definition at line 67 of file jump_error_estimator.h.

References estimate_error(), and libmesh_nullptr.

67 {}

Member Function Documentation

virtual bool libMesh::JumpErrorEstimator::boundary_side_integration ( )
inlineprotectedvirtual

The function, to be implemented by derived classes, which calculates an error term on a boundary side Returns true if the flux bc function is in fact defined on the current side

Reimplemented in libMesh::KellyErrorEstimator, and libMesh::DiscontinuityMeasure.

Definition at line 120 of file jump_error_estimator.h.

Referenced by estimate_error().

120 { return false; }
float libMesh::JumpErrorEstimator::coarse_n_flux_faces_increment ( )
protected

A utility function to correctly increase n_flux_faces for the coarse element

Definition at line 459 of file jump_error_estimator.C.

References coarse_context, and fine_context.

Referenced by estimate_error().

460 {
461  // Keep track of the number of internal flux sides found on each
462  // element
463  unsigned int dim = coarse_context->get_elem().dim();
464 
465  const unsigned int divisor =
466  1 << (dim-1)*(fine_context->get_elem().level() -
467  coarse_context->get_elem().level());
468 
469  // With a difference of n levels between fine and coarse elements,
470  // we compute a fractional flux face for the coarse element by adding:
471  // 1/2^n in 2D
472  // 1/4^n in 3D
473  // each time. This code will get hit 2^n times in 2D and 4^n
474  // times in 3D so that the final flux face count for the coarse
475  // element will be an integer value.
476 
477  return 1.0f / static_cast<float>(divisor);
478 }
UniquePtr< FEMContext > fine_context
UniquePtr< FEMContext > coarse_context
void libMesh::JumpErrorEstimator::estimate_error ( const System system,
ErrorVector error_per_cell,
const NumericVector< Number > *  solution_vector = libmesh_nullptr,
bool  estimate_parent_error = false 
)
virtual

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

Conventions for assigning the direction of the normal:

  • e & f are global element ids

Case (1.) Elements are at the same level, e<f Compute the flux jump on the face and add it as a contribution to error_per_cell[e] and error_per_cell[f]


| | |

f
e —> n

Case (2.) The neighbor is at a higher level. Compute the flux jump on e's face and add it as a contribution to error_per_cell[e] and error_per_cell[f]


| | | | | | e |—> n | | | | | |--------—| f |


Implements libMesh::ErrorEstimator.

Definition at line 53 of file jump_error_estimator.C.

References libMesh::Elem::active(), libMesh::MeshBase::active_local_elements_begin(), libMesh::MeshBase::active_local_elements_end(), boundary_side_integration(), libMesh::Elem::child_ptr(), coarse_context, coarse_error, coarse_n_flux_faces_increment(), libMesh::FEGenericBase< OutputType >::coarsened_dof_values(), libMesh::ParallelObject::comm(), libMesh::ErrorEstimator::error_norm, libMesh::ErrorVectorReal, libMesh::FEType::family, fine_context, fine_error, libMesh::System::get_dof_map(), libMesh::FEAbstract::get_fe_type(), libMesh::System::get_mesh(), libMesh::FEAbstract::get_xyz(), libMesh::DofObject::id(), init_context(), integrate_boundary_sides, internal_side_integration(), libMesh::Elem::level(), libmesh_nullptr, libMesh::MeshBase::max_elem_id(), mesh, libMesh::Elem::n_children(), libMesh::Elem::n_neighbors(), n_vars, libMesh::System::n_vars(), libMesh::Elem::neighbor_ptr(), libMesh::Elem::parent(), libMesh::ErrorEstimator::reduce_error(), reinit_sides(), libMesh::SCALAR, scale_by_n_flux_faces, libMesh::DenseVector< T >::size(), libMesh::System::solution, libMesh::NumericVector< T >::swap(), libMesh::sys, var, and libMesh::SystemNorm::weight().

Referenced by ~JumpErrorEstimator().

57 {
58  LOG_SCOPE("estimate_error()", "JumpErrorEstimator");
59 
96  // The current mesh
97  const MeshBase & mesh = system.get_mesh();
98 
99  // The number of variables in the system
100  const unsigned int n_vars = system.n_vars();
101 
102  // The DofMap for this system
103  const DofMap & dof_map = system.get_dof_map();
104 
105  // Resize the error_per_cell vector to be
106  // the number of elements, initialize it to 0.
107  error_per_cell.resize (mesh.max_elem_id());
108  std::fill (error_per_cell.begin(), error_per_cell.end(), 0.);
109 
110  // Declare a vector of floats which is as long as
111  // error_per_cell above, and fill with zeros. This vector will be
112  // used to keep track of the number of edges (faces) on each active
113  // element which are either:
114  // 1) an internal edge
115  // 2) an edge on a Neumann boundary for which a boundary condition
116  // function has been specified.
117  // The error estimator can be scaled by the number of flux edges (faces)
118  // which the element actually has to obtain a more uniform measure
119  // of the error. Use floats instead of ints since in case 2 (above)
120  // f gets 1/2 of a flux face contribution from each of his
121  // neighbors
122  std::vector<float> n_flux_faces;
124  n_flux_faces.resize(error_per_cell.size(), 0);
125 
126  // Prepare current_local_solution to localize a non-standard
127  // solution vector if necessary
128  if (solution_vector && solution_vector != system.solution.get())
129  {
130  NumericVector<Number> * newsol =
131  const_cast<NumericVector<Number> *>(solution_vector);
132  System & sys = const_cast<System &>(system);
133  newsol->swap(*sys.solution);
134  sys.update();
135  }
136 
137  fine_context.reset(new FEMContext(system));
138  coarse_context.reset(new FEMContext(system));
139 
140  // Loop over all the variables we've been requested to find jumps in, to
141  // pre-request
142  for (var=0; var<n_vars; var++)
143  {
144  // Possibly skip this variable
145  if (error_norm.weight(var) == 0.0) continue;
146 
147  // FIXME: Need to generalize this to vector-valued elements. [PB]
148  FEBase * side_fe = libmesh_nullptr;
149 
150  const std::set<unsigned char> & elem_dims =
151  fine_context->elem_dimensions();
152 
153  for (std::set<unsigned char>::const_iterator dim_it =
154  elem_dims.begin(); dim_it != elem_dims.end(); ++dim_it)
155  {
156  const unsigned char dim = *dim_it;
157 
158  fine_context->get_side_fe( var, side_fe, dim );
159 
160  libmesh_assert_not_equal_to(side_fe->get_fe_type().family, SCALAR);
161 
162  side_fe->get_xyz();
163  }
164  }
165 
166  this->init_context(*fine_context);
168 
169  // Iterate over all the active elements in the mesh
170  // that live on this processor.
171  MeshBase::const_element_iterator elem_it = mesh.active_local_elements_begin();
172  const MeshBase::const_element_iterator elem_end = mesh.active_local_elements_end();
173 
174  for (; elem_it != elem_end; ++elem_it)
175  {
176  // e is necessarily an active element on the local processor
177  const Elem * e = *elem_it;
178  const dof_id_type e_id = e->id();
179 
180 #ifdef LIBMESH_ENABLE_AMR
181  // See if the parent of element e has been examined yet;
182  // if not, we may want to compute the estimator on it
183  const Elem * parent = e->parent();
184 
185  // We only can compute and only need to compute on
186  // parents with all active children
187  bool compute_on_parent = true;
188  if (!parent || !estimate_parent_error)
189  compute_on_parent = false;
190  else
191  for (unsigned int c=0; c != parent->n_children(); ++c)
192  if (!parent->child_ptr(c)->active())
193  compute_on_parent = false;
194 
195  if (compute_on_parent &&
196  !error_per_cell[parent->id()])
197  {
198  // Compute a projection onto the parent
199  DenseVector<Number> Uparent;
201  (*(system.solution), dof_map, parent, Uparent, false);
202 
203  // Loop over the neighbors of the parent
204  for (unsigned int n_p=0; n_p<parent->n_neighbors(); n_p++)
205  {
206  if (parent->neighbor_ptr(n_p) != libmesh_nullptr) // parent has a neighbor here
207  {
208  // Find the active neighbors in this direction
209  std::vector<const Elem *> active_neighbors;
210  parent->neighbor_ptr(n_p)->
211  active_family_tree_by_neighbor(active_neighbors,
212  parent);
213  // Compute the flux to each active neighbor
214  for (unsigned int a=0;
215  a != active_neighbors.size(); ++a)
216  {
217  const Elem * f = active_neighbors[a];
218  // FIXME - what about when f->level <
219  // parent->level()??
220  if (f->level() >= parent->level())
221  {
222  fine_context->pre_fe_reinit(system, f);
223  coarse_context->pre_fe_reinit(system, parent);
224  libmesh_assert_equal_to
225  (coarse_context->get_elem_solution().size(),
226  Uparent.size());
227  coarse_context->get_elem_solution() = Uparent;
228 
229  this->reinit_sides();
230 
231  // Loop over all significant variables in the system
232  for (var=0; var<n_vars; var++)
233  if (error_norm.weight(var) != 0.0)
234  {
236 
237  error_per_cell[fine_context->get_elem().id()] +=
238  static_cast<ErrorVectorReal>(fine_error);
239  error_per_cell[coarse_context->get_elem().id()] +=
240  static_cast<ErrorVectorReal>(coarse_error);
241  }
242 
243  // Keep track of the number of internal flux
244  // sides found on each element
246  {
247  n_flux_faces[fine_context->get_elem().id()]++;
248  n_flux_faces[coarse_context->get_elem().id()] +=
250  }
251  }
252  }
253  }
254  else if (integrate_boundary_sides)
255  {
256  fine_context->pre_fe_reinit(system, parent);
257  libmesh_assert_equal_to
258  (fine_context->get_elem_solution().size(),
259  Uparent.size());
260  fine_context->get_elem_solution() = Uparent;
261  fine_context->side = n_p;
262  fine_context->side_fe_reinit();
263 
264  // If we find a boundary flux for any variable,
265  // let's just count it as a flux face for all
266  // variables. Otherwise we'd need to keep track of
267  // a separate n_flux_faces and error_per_cell for
268  // every single var.
269  bool found_boundary_flux = false;
270 
271  for (var=0; var<n_vars; var++)
272  if (error_norm.weight(var) != 0.0)
273  {
274  if (this->boundary_side_integration())
275  {
276  error_per_cell[fine_context->get_elem().id()] +=
277  static_cast<ErrorVectorReal>(fine_error);
278  found_boundary_flux = true;
279  }
280  }
281 
282  if (scale_by_n_flux_faces && found_boundary_flux)
283  n_flux_faces[fine_context->get_elem().id()]++;
284  }
285  }
286  }
287 #endif // #ifdef LIBMESH_ENABLE_AMR
288 
289  // If we do any more flux integration, e will be the fine element
290  fine_context->pre_fe_reinit(system, e);
291 
292  // Loop over the neighbors of element e
293  for (unsigned int n_e=0; n_e<e->n_neighbors(); n_e++)
294  {
295  if ((e->neighbor_ptr(n_e) != libmesh_nullptr) ||
297  {
298  fine_context->side = n_e;
299  fine_context->side_fe_reinit();
300  }
301 
302  if (e->neighbor_ptr(n_e) != libmesh_nullptr) // e is not on the boundary
303  {
304  const Elem * f = e->neighbor_ptr(n_e);
305  const dof_id_type f_id = f->id();
306 
307  // Compute flux jumps if we are in case 1 or case 2.
308  if ((f->active() && (f->level() == e->level()) && (e_id < f_id))
309  || (f->level() < e->level()))
310  {
311  // f is now the coarse element
312  coarse_context->pre_fe_reinit(system, f);
313 
314  this->reinit_sides();
315 
316  // Loop over all significant variables in the system
317  for (var=0; var<n_vars; var++)
318  if (error_norm.weight(var) != 0.0)
319  {
321 
322  error_per_cell[fine_context->get_elem().id()] +=
323  static_cast<ErrorVectorReal>(fine_error);
324  error_per_cell[coarse_context->get_elem().id()] +=
325  static_cast<ErrorVectorReal>(coarse_error);
326  }
327 
328  // Keep track of the number of internal flux
329  // sides found on each element
331  {
332  n_flux_faces[fine_context->get_elem().id()]++;
333  n_flux_faces[coarse_context->get_elem().id()] +=
335  }
336  } // end if (case1 || case2)
337  } // if (e->neigbor(n_e) != libmesh_nullptr)
338 
339  // Otherwise, e is on the boundary. If it happens to
340  // be on a Dirichlet boundary, we need not do anything.
341  // On the other hand, if e is on a Neumann (flux) boundary
342  // with grad(u).n = g, we need to compute the additional residual
343  // (h * \int |g - grad(u_h).n|^2 dS)^(1/2).
344  // We can only do this with some knowledge of the boundary
345  // conditions, i.e. the user must have attached an appropriate
346  // BC function.
347  else if (integrate_boundary_sides)
348  {
349  bool found_boundary_flux = false;
350 
351  for (var=0; var<n_vars; var++)
352  if (error_norm.weight(var) != 0.0)
353  if (this->boundary_side_integration())
354  {
355  error_per_cell[fine_context->get_elem().id()] +=
356  static_cast<ErrorVectorReal>(fine_error);
357  found_boundary_flux = true;
358  }
359 
360  if (scale_by_n_flux_faces && found_boundary_flux)
361  n_flux_faces[fine_context->get_elem().id()]++;
362  } // end if (e->neighbor_ptr(n_e) == libmesh_nullptr)
363  } // end loop over neighbors
364  } // End loop over active local elements
365 
366 
367  // Each processor has now computed the error contribuions
368  // for its local elements. We need to sum the vector
369  // and then take the square-root of each component. Note
370  // that we only need to sum if we are running on multiple
371  // processors, and we only need to take the square-root
372  // if the value is nonzero. There will in general be many
373  // zeros for the inactive elements.
374 
375  // First sum the vector of estimated error values
376  this->reduce_error(error_per_cell, system.comm());
377 
378  // Compute the square-root of each component.
379  for (std::size_t i=0; i<error_per_cell.size(); i++)
380  if (error_per_cell[i] != 0.)
381  error_per_cell[i] = std::sqrt(error_per_cell[i]);
382 
383 
384  if (this->scale_by_n_flux_faces)
385  {
386  // Sum the vector of flux face counts
387  this->reduce_error(n_flux_faces, system.comm());
388 
389  // Sanity check: Make sure the number of flux faces is
390  // always an integer value
391 #ifdef DEBUG
392  for (std::size_t i=0; i<n_flux_faces.size(); ++i)
393  libmesh_assert_equal_to (n_flux_faces[i], static_cast<float>(static_cast<unsigned int>(n_flux_faces[i])) );
394 #endif
395 
396  // Scale the error by the number of flux faces for each element
397  for (std::size_t i=0; i<n_flux_faces.size(); ++i)
398  {
399  if (n_flux_faces[i] == 0.0) // inactive or non-local element
400  continue;
401 
402  //libMesh::out << "Element " << i << " has " << n_flux_faces[i] << " flux faces." << std::endl;
403  error_per_cell[i] /= static_cast<ErrorVectorReal>(n_flux_faces[i]);
404  }
405  }
406 
407  // If we used a non-standard solution before, now is the time to fix
408  // the current_local_solution
409  if (solution_vector && solution_vector != system.solution.get())
410  {
411  NumericVector<Number> * newsol =
412  const_cast<NumericVector<Number> *>(solution_vector);
413  System & sys = const_cast<System &>(system);
414  newsol->swap(*sys.solution);
415  sys.update();
416  }
417 }
ImplicitSystem & sys
Real weight(unsigned int var) const
Definition: system_norm.h:292
MeshBase & mesh
const class libmesh_nullptr_t libmesh_nullptr
virtual void internal_side_integration()=0
UniquePtr< FEMContext > fine_context
const unsigned int n_vars
Definition: tecplot_io.C:68
DIE A HORRIBLE DEATH HERE typedef float ErrorVectorReal
static void coarsened_dof_values(const NumericVector< Number > &global_vector, const DofMap &dof_map, const Elem *coarse_elem, DenseVector< Number > &coarse_dofs, const unsigned int var, const bool use_old_dof_indices=false)
Definition: fe_base.C:804
FEGenericBase< Real > FEBase
void reduce_error(std::vector< ErrorVectorReal > &error_per_cell, const Parallel::Communicator &comm LIBMESH_CAN_DEFAULT_TO_COMMWORLD) const
virtual void update()
Definition: system.C:420
virtual void get(const std::vector< numeric_index_type > &index, T *values) const
UniquePtr< FEMContext > coarse_context
virtual void init_context(FEMContext &c)
uint8_t dof_id_type
Definition: id_types.h:64
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 = libmesh_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 48 of file error_estimator.C.

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

Referenced by libMesh::ErrorEstimator::~ErrorEstimator().

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 = libmesh_nullptr;
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 (std::size_t 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 }
ImplicitSystem & sys
const class libmesh_nullptr_t libmesh_nullptr
virtual void estimate_error(const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=libmesh_nullptr, bool estimate_parent_error=false)=0
void libMesh::ErrorEstimator::estimate_errors ( const EquationSystems equation_systems,
ErrorMap errors_per_cell,
const std::map< const System *, const NumericVector< Number > * > *  solution_vectors = libmesh_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 94 of file error_estimator.C.

References libMesh::ErrorEstimator::error_norm, libMesh::ErrorEstimator::estimate_error(), libMesh::EquationSystems::get_system(), libmesh_nullptr, libMesh::EquationSystems::n_systems(), n_vars, libMesh::System::n_vars(), libMesh::sys, 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 = libmesh_nullptr;
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 }
ImplicitSystem & sys
const class libmesh_nullptr_t libmesh_nullptr
const unsigned int n_vars
Definition: tecplot_io.C:68
virtual void estimate_error(const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=libmesh_nullptr, bool estimate_parent_error=false)=0
void libMesh::JumpErrorEstimator::init_context ( FEMContext c)
protectedvirtual

An initialization function, to give derived classes a chance to request specific data from the FE objects

Reimplemented in libMesh::KellyErrorEstimator, libMesh::DiscontinuityMeasure, and libMesh::LaplacianErrorEstimator.

Definition at line 47 of file jump_error_estimator.C.

Referenced by estimate_error().

48 {
49 }
virtual void libMesh::JumpErrorEstimator::internal_side_integration ( )
protectedpure virtual

The function, to be implemented by derived classes, which calculates an error term on an internal side

Implemented in libMesh::KellyErrorEstimator, libMesh::DiscontinuityMeasure, and libMesh::LaplacianErrorEstimator.

Referenced by estimate_error().

void libMesh::ErrorEstimator::reduce_error ( std::vector< ErrorVectorReal > &  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::PatchRecoveryErrorEstimator::estimate_error(), libMesh::AdjointRefinementEstimator::estimate_error(), and libMesh::ExactErrorEstimator::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::JumpErrorEstimator::reinit_sides ( )
protected

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

Definition at line 422 of file jump_error_estimator.C.

References coarse_context, libMesh::ErrorEstimator::error_norm, fine_context, libMesh::FEAbstract::get_fe_type(), libMesh::FEInterface::inverse_map(), libmesh_nullptr, n_vars, libMesh::FEAbstract::reinit(), and libMesh::SystemNorm::weight().

Referenced by estimate_error().

423 {
424  fine_context->side_fe_reinit();
425 
426  unsigned int dim = fine_context->get_elem().dim();
427  libmesh_assert_equal_to(dim, coarse_context->get_elem().dim());
428 
429  FEBase * fe_fine = libmesh_nullptr;
430  fine_context->get_side_fe( 0, fe_fine, dim );
431 
432  // Get the physical locations of the fine element quadrature points
433  std::vector<Point> qface_point = fe_fine->get_xyz();
434 
435  // Find the master quadrature point locations on the coarse element
436  FEBase * fe_coarse = libmesh_nullptr;
437  coarse_context->get_side_fe( 0, fe_coarse, dim );
438 
439  std::vector<Point> qp_coarse;
440 
442  (coarse_context->get_elem().dim(), fe_coarse->get_fe_type(),
443  &coarse_context->get_elem(), qface_point, qp_coarse);
444 
445  // The number of variables in the system
446  const unsigned int n_vars = fine_context->n_vars();
447 
448  // Calculate all coarse element shape functions at those locations
449  for (unsigned int v=0; v<n_vars; v++)
450  if (error_norm.weight(v) != 0.0)
451  {
452  coarse_context->get_side_fe( v, fe_coarse, dim );
453  fe_coarse->reinit (&coarse_context->get_elem(), &qp_coarse);
454  }
455 }
Real weight(unsigned int var) const
Definition: system_norm.h:292
const class libmesh_nullptr_t libmesh_nullptr
UniquePtr< FEMContext > fine_context
const unsigned int n_vars
Definition: tecplot_io.C:68
static Point inverse_map(const unsigned int dim, const FEType &fe_t, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
Definition: fe_interface.C:556
FEGenericBase< Real > FEBase
UniquePtr< FEMContext > coarse_context

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 149 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(), 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::PatchRecoveryErrorEstimator::PatchRecoveryErrorEstimator(), reinit_sides(), and libMesh::UniformRefinementEstimator::UniformRefinementEstimator().

bool libMesh::JumpErrorEstimator::integrate_boundary_sides
protected

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

Definition at line 126 of file jump_error_estimator.h.

Referenced by libMesh::DiscontinuityMeasure::attach_essential_bc_function(), libMesh::KellyErrorEstimator::attach_flux_bc_function(), and estimate_error().

bool libMesh::JumpErrorEstimator::scale_by_n_flux_faces

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 89 of file jump_error_estimator.h.

Referenced by estimate_error().


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