exact_error_estimator.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 // C++ includes
20 #include <algorithm> // for std::fill
21 #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC
22 #include <cmath> // for sqrt
23 
24 
25 // Local Includes
26 #include "libmesh/libmesh_common.h"
28 #include "libmesh/dof_map.h"
30 #include "libmesh/error_vector.h"
31 #include "libmesh/fe_base.h"
33 #include "libmesh/elem.h"
34 #include "libmesh/mesh_base.h"
35 #include "libmesh/mesh_function.h"
36 #include "libmesh/numeric_vector.h"
37 #include "libmesh/quadrature.h"
38 #include "libmesh/system.h"
39 #include "libmesh/tensor_tools.h"
41 #include "libmesh/enum_norm_type.h"
42 
43 namespace libMesh
44 {
45 
46 //-----------------------------------------------------------------
47 // ErrorEstimator implementations
50  _exact_value(nullptr),
51  _exact_deriv(nullptr),
52  _exact_hessian(nullptr),
53  _equation_systems_fine(nullptr),
54  _extra_order(0)
55 {
56  error_norm = H1;
57 }
58 
59 
61 {
62  return EXACT;
63 }
64 
65 
66 void ExactErrorEstimator::attach_exact_value (ValueFunctionPointer fptr)
67 {
68  libmesh_assert(fptr);
69  _exact_value = fptr;
70 
71  // We're not using a fine grid solution
72  _equation_systems_fine = nullptr;
73 
74  // We're not using user-provided functors
75  this->clear_functors();
76 }
77 
78 
80 {
81  // Automatically delete any previous _exact_values entries, then add a new
82  // entry for each system.
83  _exact_values.clear();
84 
85  for (auto ptr : f)
86  _exact_values.emplace_back(ptr ? ptr->clone() : nullptr);
87 }
88 
89 
90 void ExactErrorEstimator::attach_exact_value (unsigned int sys_num,
92 {
93  if (_exact_values.size() <= sys_num)
94  _exact_values.resize(sys_num+1);
95 
96  if (f)
97  _exact_values[sys_num] = f->clone();
98 }
99 
100 
101 void ExactErrorEstimator::attach_exact_deriv (GradientFunctionPointer gptr)
102 {
103  libmesh_assert(gptr);
104  _exact_deriv = gptr;
105 
106  // We're not using a fine grid solution
107  _equation_systems_fine = nullptr;
108 
109  // We're not using user-provided functors
110  this->clear_functors();
111 }
112 
113 
115 {
116  // Automatically delete any previous _exact_derivs entries, then add a new
117  // entry for each system.
118  _exact_derivs.clear();
119 
120  for (auto ptr : g)
121  _exact_derivs.emplace_back(ptr ? ptr->clone() : nullptr);
122 }
123 
124 
125 void ExactErrorEstimator::attach_exact_deriv (unsigned int sys_num,
127 {
128  if (_exact_derivs.size() <= sys_num)
129  _exact_derivs.resize(sys_num+1);
130 
131  if (g)
132  _exact_derivs[sys_num] = g->clone();
133 }
134 
135 
136 
137 
138 void ExactErrorEstimator::attach_exact_hessian (HessianFunctionPointer hptr)
139 {
140  libmesh_assert(hptr);
141  _exact_hessian = hptr;
142 
143  // We're not using a fine grid solution
144  _equation_systems_fine = nullptr;
145 
146  // We're not using user-provided functors
147  this->clear_functors();
148 }
149 
150 
152 {
153  // Automatically delete any previous _exact_hessians entries, then add a new
154  // entry for each system.
155  _exact_hessians.clear();
156 
157  for (auto ptr : h)
158  _exact_hessians.emplace_back(ptr ? ptr->clone() : nullptr);
159 }
160 
161 
162 void ExactErrorEstimator::attach_exact_hessian (unsigned int sys_num,
164 {
165  if (_exact_hessians.size() <= sys_num)
166  _exact_hessians.resize(sys_num+1);
167 
168  if (h)
169  _exact_hessians[sys_num] = h->clone();
170 }
171 
172 
174 {
175  libmesh_assert(es_fine);
176  _equation_systems_fine = es_fine;
177 
178  // If we're using a fine grid solution, we're not using exact value
179  // function pointers or functors.
180  _exact_value = nullptr;
181  _exact_deriv = nullptr;
182  _exact_hessian = nullptr;
183 
184  this->clear_functors();
185 }
186 
188  ErrorVector & error_per_cell,
189  const NumericVector<Number> * solution_vector,
190  bool estimate_parent_error)
191 {
192  // Ignore the fact that this variable is unused when !LIBMESH_ENABLE_AMR
193  libmesh_ignore(estimate_parent_error);
194 
195  // The current mesh
196  const MeshBase & mesh = system.get_mesh();
197 
198  // The dimensionality of the mesh
199  const unsigned int dim = mesh.mesh_dimension();
200 
201  // The number of variables in the system
202  const unsigned int n_vars = system.n_vars();
203 
204  // The DofMap for this system
205  const DofMap & dof_map = system.get_dof_map();
206 
207  // Resize the error_per_cell vector to be
208  // the number of elements, initialize it to 0.
209  error_per_cell.resize (mesh.max_elem_id());
210  std::fill (error_per_cell.begin(), error_per_cell.end(), 0.);
211 
212  // Prepare current_local_solution to localize a non-standard
213  // solution vector if necessary
214  if (solution_vector && solution_vector != system.solution.get())
215  {
216  NumericVector<Number> * newsol =
217  const_cast<NumericVector<Number> *>(solution_vector);
218  System & sys = const_cast<System &>(system);
219  newsol->swap(*sys.solution);
220  sys.update();
221  }
222 
223  // Loop over all the variables in the system
224  for (unsigned int var=0; var<n_vars; var++)
225  {
226  // Possibly skip this variable
227  if (error_norm.weight(var) == 0.0) continue;
228 
229  // The (string) name of this variable
230  const std::string & var_name = system.variable_name(var);
231 
232  // The type of finite element to use for this variable
233  const FEType & fe_type = dof_map.variable_type (var);
234 
235  std::unique_ptr<FEBase> fe (FEBase::build (dim, fe_type));
236 
237  // Build an appropriate Gaussian quadrature rule
238  std::unique_ptr<QBase> qrule =
239  fe_type.default_quadrature_rule (dim,
240  _extra_order);
241 
242  fe->attach_quadrature_rule (qrule.get());
243 
244  // Prepare a global solution and a MeshFunction of the fine system if we need one
245  std::unique_ptr<MeshFunction> fine_values;
246  std::unique_ptr<NumericVector<Number>> fine_soln = NumericVector<Number>::build(system.comm());
248  {
249  const System & fine_system = _equation_systems_fine->get_system(system.name());
250 
251  std::vector<Number> global_soln;
252  // FIXME - we're assuming that the fine system solution gets
253  // used even when a different vector is used for the coarse
254  // system
255  fine_system.update_global_solution(global_soln);
256  fine_soln->init
257  (cast_int<numeric_index_type>(global_soln.size()), true,
258  SERIAL);
259  (*fine_soln) = global_soln;
260 
261  fine_values = std::unique_ptr<MeshFunction>
263  *fine_soln,
264  fine_system.get_dof_map(),
265  fine_system.variable_number(var_name)));
266  fine_values->init();
267  } else {
268  // Initialize functors if we're using them
269  for (auto & ev : _exact_values)
270  if (ev)
271  ev->init();
272 
273  for (auto & ed : _exact_derivs)
274  if (ed)
275  ed->init();
276 
277  for (auto & eh : _exact_hessians)
278  if (eh)
279  eh->init();
280  }
281 
282  // Request the data we'll need to compute with
283  fe->get_JxW();
284  fe->get_phi();
285  fe->get_dphi();
286 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
287  fe->get_d2phi();
288 #endif
289  fe->get_xyz();
290 
291 #ifdef LIBMESH_ENABLE_AMR
292  // If we compute on parent elements, we'll want to do so only
293  // once on each, so we need to keep track of which we've done.
294  std::vector<bool> computed_var_on_parent;
295 
296  if (estimate_parent_error)
297  computed_var_on_parent.resize(error_per_cell.size(), false);
298 #endif
299 
300  // TODO: this ought to be threaded (and using subordinate
301  // MeshFunction objects in each thread rather than a single
302  // master)
303 
304  // Iterate over all the active elements in the mesh
305  // that live on this processor.
306  for (const auto & elem : mesh.active_local_element_ptr_range())
307  {
308  const dof_id_type e_id = elem->id();
309 
310 #ifdef LIBMESH_ENABLE_AMR
311  // See if the parent of element e has been examined yet;
312  // if not, we may want to compute the estimator on it
313  const Elem * parent = elem->parent();
314 
315  // We only can compute and only need to compute on
316  // parents with all active children
317  bool compute_on_parent = true;
318  if (!parent || !estimate_parent_error)
319  compute_on_parent = false;
320  else
321  for (auto & child : parent->child_ref_range())
322  if (!child.active())
323  compute_on_parent = false;
324 
325  if (compute_on_parent &&
326  !computed_var_on_parent[parent->id()])
327  {
328  computed_var_on_parent[parent->id()] = true;
329 
330  // Compute a projection onto the parent
331  DenseVector<Number> Uparent;
333  dof_map, parent, Uparent,
334  var, false);
335 
336  error_per_cell[parent->id()] +=
337  static_cast<ErrorVectorReal>
338  (find_squared_element_error(system, var_name,
339  parent, Uparent,
340  fe.get(),
341  fine_values.get()));
342  }
343 #endif
344 
345  // Get the local to global degree of freedom maps
346  std::vector<dof_id_type> dof_indices;
347  dof_map.dof_indices (elem, dof_indices, var);
348  const unsigned int n_dofs =
349  cast_int<unsigned int>(dof_indices.size());
350  DenseVector<Number> Uelem(n_dofs);
351  for (unsigned int i=0; i != n_dofs; ++i)
352  Uelem(i) = system.current_solution(dof_indices[i]);
353 
354  error_per_cell[e_id] +=
355  static_cast<ErrorVectorReal>
356  (find_squared_element_error(system, var_name, elem,
357  Uelem, fe.get(),
358  fine_values.get()));
359 
360  } // End loop over active local elements
361  } // End loop over variables
362 
363 
364 
365  // Each processor has now computed the error contributions
366  // for its local elements. We need to sum the vector
367  // and then take the square-root of each component. Note
368  // that we only need to sum if we are running on multiple
369  // processors, and we only need to take the square-root
370  // if the value is nonzero. There will in general be many
371  // zeros for the inactive elements.
372 
373  // First sum the vector of estimated error values
374  this->reduce_error(error_per_cell, system.comm());
375 
376  // Compute the square-root of each component.
377  {
378  LOG_SCOPE("std::sqrt()", "ExactErrorEstimator");
379  for (auto & val : error_per_cell)
380  if (val != 0.)
381  {
382  libmesh_assert_greater (val, 0.);
383  val = std::sqrt(val);
384  }
385  }
386 
387  // If we used a non-standard solution before, now is the time to fix
388  // the current_local_solution
389  if (solution_vector && solution_vector != system.solution.get())
390  {
391  NumericVector<Number> * newsol =
392  const_cast<NumericVector<Number> *>(solution_vector);
393  System & sys = const_cast<System &>(system);
394  newsol->swap(*sys.solution);
395  sys.update();
396  }
397 }
398 
399 
400 
402  const std::string & var_name,
403  const Elem * elem,
404  const DenseVector<Number> & Uelem,
405  FEBase * fe,
406  MeshFunction * fine_values) const
407 {
408  // The (string) name of this system
409  const std::string & sys_name = system.name();
410  const unsigned int sys_num = system.number();
411 
412  const unsigned int var = system.variable_number(var_name);
413  const unsigned int var_component =
414  system.variable_scalar_number(var, 0);
415 
416  const Parameters & parameters = system.get_equation_systems().parameters;
417 
418  // reinitialize the element-specific data
419  // for the current element
420  fe->reinit (elem);
421 
422  // Get the data we need to compute with
423  const std::vector<Real> & JxW = fe->get_JxW();
424  const std::vector<std::vector<Real>> & phi_values = fe->get_phi();
425  const std::vector<std::vector<RealGradient>> & dphi_values = fe->get_dphi();
426  const std::vector<Point> & q_point = fe->get_xyz();
427 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
428  const std::vector<std::vector<RealTensor>> & d2phi_values = fe->get_d2phi();
429 #endif
430 
431  // The number of shape functions
432  const unsigned int n_sf =
433  cast_int<unsigned int>(Uelem.size());
434 
435  // The number of quadrature points
436  const unsigned int n_qp =
437  cast_int<unsigned int>(JxW.size());
438 
439  Real error_val = 0;
440 
441  // Begin the loop over the Quadrature points.
442  //
443  for (unsigned int qp=0; qp<n_qp; qp++)
444  {
445  // Real u_h = 0.;
446  // RealGradient grad_u_h;
447 
448  Number u_h = 0.;
449 
450  Gradient grad_u_h;
451 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
452  Tensor grad2_u_h;
453 #endif
454 
455  // Compute solution values at the current
456  // quadrature point. This requires a sum
457  // over all the shape functions evaluated
458  // at the quadrature point.
459  for (unsigned int i=0; i<n_sf; i++)
460  {
461  // Values from current solution.
462  u_h += phi_values[i][qp]*Uelem(i);
463  grad_u_h += dphi_values[i][qp]*Uelem(i);
464 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
465  grad2_u_h += d2phi_values[i][qp]*Uelem(i);
466 #endif
467  }
468 
469  // Compute the value of the error at this quadrature point
470  if (error_norm.type(var) == L2 ||
471  error_norm.type(var) == H1 ||
472  error_norm.type(var) == H2)
473  {
474  Number val_error = u_h;
475  if (_exact_value)
476  val_error -= _exact_value(q_point[qp],parameters,sys_name,var_name);
477  else if (_exact_values.size() > sys_num && _exact_values[sys_num])
478  val_error -= _exact_values[sys_num]->
479  component(var_component, q_point[qp], system.time);
480  else if (_equation_systems_fine)
481  val_error -= (*fine_values)(q_point[qp]);
482 
483  // Add the squares of the error to each contribution
484  error_val += JxW[qp]*TensorTools::norm_sq(val_error);
485  }
486 
487  // Compute the value of the error in the gradient at this
488  // quadrature point
489  if (error_norm.type(var) == H1 ||
490  error_norm.type(var) == H1_SEMINORM ||
491  error_norm.type(var) == H2)
492  {
493  Gradient grad_error = grad_u_h;
494  if (_exact_deriv)
495  grad_error -= _exact_deriv(q_point[qp],parameters,sys_name,var_name);
496  else if (_exact_derivs.size() > sys_num && _exact_derivs[sys_num])
497  grad_error -= _exact_derivs[sys_num]->
498  component(var_component, q_point[qp], system.time);
499  else if (_equation_systems_fine)
500  grad_error -= fine_values->gradient(q_point[qp]);
501 
502  error_val += JxW[qp]*grad_error.norm_sq();
503  }
504 
505 
506 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
507  // Compute the value of the error in the hessian at this
508  // quadrature point
509  if ((error_norm.type(var) == H2_SEMINORM ||
510  error_norm.type(var) == H2))
511  {
512  Tensor grad2_error = grad2_u_h;
513  if (_exact_hessian)
514  grad2_error -= _exact_hessian(q_point[qp],parameters,sys_name,var_name);
515  else if (_exact_hessians.size() > sys_num && _exact_hessians[sys_num])
516  grad2_error -= _exact_hessians[sys_num]->
517  component(var_component, q_point[qp], system.time);
518  else if (_equation_systems_fine)
519  grad2_error -= fine_values->hessian(q_point[qp]);
520 
521  error_val += JxW[qp]*grad2_error.norm_sq();
522  }
523 #endif
524 
525  } // end qp loop
526 
527  libmesh_assert_greater_equal (error_val, 0.);
528 
529  return error_val;
530 }
531 
532 
533 
535 {
536  _exact_values.clear();
537  _exact_derivs.clear();
538  _exact_hessians.clear();
539 }
540 
541 
542 
543 } // namespace libMesh
Manages the family, order, etc. parameters for a given FE.
Definition: fe_type.h:179
virtual unsigned int size() const override
Definition: dense_vector.h:92
HessianFunctionPointer _exact_hessian
Manages multiples systems of equations.
void attach_exact_hessian(unsigned int sys_num, FunctionBase< Tensor > *h)
const Elem * parent() const
Definition: elem.h:2479
unsigned int variable_scalar_number(const std::string &var, unsigned int component) const
Definition: system.h:2164
void update_global_solution(std::vector< Number > &global_soln) const
Definition: system.C:642
const std::vector< Point > & get_xyz() const
Definition: fe_abstract.h:238
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Definition: dof_map.C:1930
std::vector< std::unique_ptr< FunctionBase< Tensor > > > _exact_hessians
void attach_exact_hessians(std::vector< FunctionBase< Tensor > *> h)
const T_sys & get_system(const std::string &name) const
const FEType & variable_type(const unsigned int c) const
Definition: dof_map.h:1792
const EquationSystems & get_equation_systems() const
Definition: system.h:712
The base class for all geometric element types.
Definition: elem.h:100
MeshBase & mesh
Gradient gradient(const Point &p, const Real time=0.)
const Parallel::Communicator & comm() const
const unsigned int n_vars
Definition: tecplot_io.C:69
void attach_exact_value(unsigned int sys_num, FunctionBase< Number > *f)
const MeshBase & get_mesh() const
Definition: system.h:2033
DIE A HORRIBLE DEATH HERE typedef float ErrorVectorReal
std::vector< std::unique_ptr< FunctionBase< Gradient > > > _exact_derivs
Base class for Mesh.
Definition: mesh_base.h:77
SimpleRange< ChildRefIter > child_ref_range()
Definition: elem.h:1779
Number current_solution(const dof_id_type global_dof_number) const
Definition: system.C:194
virtual void estimate_error(const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=nullptr, bool estimate_parent_error=false) override
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:791
virtual ErrorEstimatorType type() const override
Manages the degrees of freedom (DOFs) in a simulation.
Definition: dof_map.h:176
std::unique_ptr< QBase > default_quadrature_rule(const unsigned int dim, const int extraorder=0) const
Definition: fe_type.C:31
const std::vector< std::vector< OutputGradient > > & get_dphi() const
Definition: fe_base.h:215
void libmesh_ignore(const Args &...)
unsigned int number() const
Definition: system.h:2025
dof_id_type id() const
Definition: dof_object.h:655
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:92
unsigned short int variable_number(const std::string &var) const
Definition: system.C:1243
const std::vector< Real > & get_JxW() const
Definition: fe_abstract.h:245
FEMNormType type(unsigned int var) const
Definition: system_norm.C:110
std::unique_ptr< NumericVector< Number > > solution
Definition: system.h:1523
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
virtual void reinit(const Elem *elem, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr)=0
const std::string & variable_name(const unsigned int i) const
Definition: system.h:2153
void reduce_error(std::vector< ErrorVectorReal > &error_per_cell, const Parallel::Communicator &comm) const
const std::vector< std::vector< OutputTensor > > & get_d2phi() const
Definition: fe_base.h:289
void attach_exact_derivs(std::vector< FunctionBase< Gradient > *> g)
static std::unique_ptr< NumericVector< T > > build(const Parallel::Communicator &comm, const SolverPackage solver_package=libMesh::default_solver_package())
Real norm_sq() const
Definition: type_vector.h:943
Real weight(unsigned int var) const
Definition: system_norm.C:132
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
GradientFunctionPointer _exact_deriv
virtual void swap(NumericVector< T > &v)
Real find_squared_element_error(const System &system, const std::string &var_name, const Elem *elem, const DenseVector< Number > &Uelem, FEBase *fe, MeshFunction *fine_values) const
std::unique_ptr< NumericVector< Number > > current_local_solution
Definition: system.h:1535
virtual std::unique_ptr< FunctionBase< Output > > clone() const =0
void attach_exact_values(std::vector< FunctionBase< Number > *> f)
const std::string & name() const
Definition: system.h:2017
void attach_reference_solution(EquationSystems *es_fine)
void attach_exact_deriv(unsigned int sys_num, FunctionBase< Gradient > *g)
unsigned int n_vars() const
Definition: system.h:2105
const DofMap & get_dof_map() const
Definition: system.h:2049
Tensor hessian(const Point &p, const Real time=0.)
std::vector< std::unique_ptr< FunctionBase< Number > > > _exact_values
Real norm_sq() const
Definition: type_tensor.h:1299
uint8_t dof_id_type
Definition: id_types.h:64
const std::vector< std::vector< OutputShape > > & get_phi() const
Definition: fe_base.h:207