Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
error_estimator.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2020 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_error_estimator_h
17 #define dealii_error_estimator_h
18 
19 
20 #include <deal.II/base/config.h>
21 
23 #include <deal.II/base/function.h>
24 
26 
28 
29 #include <map>
30 
32 
33 // Forward declarations
34 #ifndef DOXYGEN
35 template <int, int>
36 class Mapping;
37 template <int>
38 class Quadrature;
39 
40 namespace hp
41 {
42  template <int>
43  class QCollection;
44 }
45 #endif
46 
47 
261 template <int dim, int spacedim = dim>
263 {
264 public:
269  enum Strategy
270  {
278  };
279 
340  template <typename InputVector, typename DoFHandlerType>
341  static void
342  estimate(
343  const Mapping<dim, spacedim> &mapping,
344  const DoFHandlerType & dof,
345  const Quadrature<dim - 1> & quadrature,
346  const std::map<types::boundary_id,
348  & neumann_bc,
349  const InputVector & solution,
350  Vector<float> & error,
351  const ComponentMask & component_mask = ComponentMask(),
352  const Function<spacedim> *coefficients = nullptr,
353  const unsigned int n_threads = numbers::invalid_unsigned_int,
356  const Strategy strategy = cell_diameter_over_24);
357 
362  template <typename InputVector, typename DoFHandlerType>
363  static void
364  estimate(
365  const DoFHandlerType & dof,
366  const Quadrature<dim - 1> &quadrature,
367  const std::map<types::boundary_id,
369  & neumann_bc,
370  const InputVector & solution,
371  Vector<float> & error,
372  const ComponentMask & component_mask = ComponentMask(),
373  const Function<spacedim> *coefficients = nullptr,
374  const unsigned int n_threads = numbers::invalid_unsigned_int,
377  const Strategy strategy = cell_diameter_over_24);
378 
392  template <typename InputVector, typename DoFHandlerType>
393  static void
394  estimate(
395  const Mapping<dim, spacedim> &mapping,
396  const DoFHandlerType & dof,
397  const Quadrature<dim - 1> & quadrature,
398  const std::map<types::boundary_id,
400  & neumann_bc,
401  const std::vector<const InputVector *> &solutions,
402  std::vector<Vector<float> *> & errors,
403  const ComponentMask & component_mask = ComponentMask(),
404  const Function<spacedim> * coefficients = nullptr,
405  const unsigned int n_threads = numbers::invalid_unsigned_int,
408  const Strategy strategy = cell_diameter_over_24);
409 
414  template <typename InputVector, typename DoFHandlerType>
415  static void
416  estimate(
417  const DoFHandlerType & dof,
418  const Quadrature<dim - 1> &quadrature,
419  const std::map<types::boundary_id,
421  & neumann_bc,
422  const std::vector<const InputVector *> &solutions,
423  std::vector<Vector<float> *> & errors,
424  const ComponentMask & component_mask = ComponentMask(),
425  const Function<spacedim> * coefficients = nullptr,
426  const unsigned int n_threads = numbers::invalid_unsigned_int,
429  const Strategy strategy = cell_diameter_over_24);
430 
431 
436  template <typename InputVector, typename DoFHandlerType>
437  static void
438  estimate(
439  const Mapping<dim, spacedim> & mapping,
440  const DoFHandlerType & dof,
441  const hp::QCollection<dim - 1> &quadrature,
442  const std::map<types::boundary_id,
444  & neumann_bc,
445  const InputVector & solution,
446  Vector<float> & error,
447  const ComponentMask & component_mask = ComponentMask(),
448  const Function<spacedim> *coefficients = nullptr,
449  const unsigned int n_threads = numbers::invalid_unsigned_int,
452  const Strategy strategy = cell_diameter_over_24);
453 
454 
459  template <typename InputVector, typename DoFHandlerType>
460  static void
461  estimate(
462  const DoFHandlerType & dof,
463  const hp::QCollection<dim - 1> &quadrature,
464  const std::map<types::boundary_id,
466  & neumann_bc,
467  const InputVector & solution,
468  Vector<float> & error,
469  const ComponentMask & component_mask = ComponentMask(),
470  const Function<spacedim> *coefficients = nullptr,
471  const unsigned int n_threads = numbers::invalid_unsigned_int,
474  const Strategy strategy = cell_diameter_over_24);
475 
476 
481  template <typename InputVector, typename DoFHandlerType>
482  static void
483  estimate(
484  const Mapping<dim, spacedim> & mapping,
485  const DoFHandlerType & dof,
486  const hp::QCollection<dim - 1> &quadrature,
487  const std::map<types::boundary_id,
489  & neumann_bc,
490  const std::vector<const InputVector *> &solutions,
491  std::vector<Vector<float> *> & errors,
492  const ComponentMask & component_mask = ComponentMask(),
493  const Function<spacedim> * coefficients = nullptr,
494  const unsigned int n_threads = numbers::invalid_unsigned_int,
497  const Strategy strategy = cell_diameter_over_24);
498 
499 
504  template <typename InputVector, typename DoFHandlerType>
505  static void
506  estimate(
507  const DoFHandlerType & dof,
508  const hp::QCollection<dim - 1> &quadrature,
509  const std::map<types::boundary_id,
511  & neumann_bc,
512  const std::vector<const InputVector *> &solutions,
513  std::vector<Vector<float> *> & errors,
514  const ComponentMask & component_mask = ComponentMask(),
515  const Function<spacedim> * coefficients = nullptr,
516  const unsigned int n_threads = numbers::invalid_unsigned_int,
519  const Strategy strategy = cell_diameter_over_24);
520 
525  "You provided a ComponentMask argument that is invalid. "
526  "Component masks need to be either default constructed "
527  "(in which case they indicate that every component is "
528  "selected) or need to have a length equal to the number "
529  "of vector components of the finite element in use "
530  "by the DoFHandler object. In the latter case, at "
531  "least one component needs to be selected.");
537  "If you do specify the argument for a (possibly "
538  "spatially variable) coefficient function for this function, "
539  "then it needs to refer to a coefficient that is either "
540  "scalar (has one vector component) or has as many vector "
541  "components as there are in the finite element used by "
542  "the DoFHandler argument.");
548  int,
549  int,
550  << "You provided a function map that for boundary indicator "
551  << arg1 << " specifies a function with " << arg2
552  << " vector components. However, the finite "
553  "element in use has "
554  << arg3
555  << " components, and these two numbers need to match.");
560  int,
561  int,
562  << "The number of input vectors, " << arg1
563  << " needs to be equal to the number of output vectors, "
564  << arg2
565  << ". This is not the case in your call of this function.");
570  "You need to specify at least one solution vector as "
571  "input.");
572 };
573 
574 
575 
587 template <int spacedim>
588 class KellyErrorEstimator<1, spacedim>
589 {
590 public:
595  enum Strategy
596  {
604  };
605 
628  template <typename InputVector, typename DoFHandlerType>
629  static void
630  estimate(
631  const Mapping<1, spacedim> &mapping,
632  const DoFHandlerType & dof,
633  const Quadrature<0> & quadrature,
634  const std::map<types::boundary_id,
636  & neumann_bc,
637  const InputVector & solution,
638  Vector<float> & error,
639  const ComponentMask & component_mask = ComponentMask(),
640  const Function<spacedim> *coefficient = nullptr,
641  const unsigned int n_threads = numbers::invalid_unsigned_int,
644  const Strategy strategy = cell_diameter_over_24);
645 
650  template <typename InputVector, typename DoFHandlerType>
651  static void
652  estimate(
653  const DoFHandlerType &dof,
654  const Quadrature<0> & quadrature,
655  const std::map<types::boundary_id,
657  & neumann_bc,
658  const InputVector & solution,
659  Vector<float> & error,
660  const ComponentMask & component_mask = ComponentMask(),
661  const Function<spacedim> *coefficients = nullptr,
662  const unsigned int n_threads = numbers::invalid_unsigned_int,
665  const Strategy strategy = cell_diameter_over_24);
666 
680  template <typename InputVector, typename DoFHandlerType>
681  static void
682  estimate(
683  const Mapping<1, spacedim> &mapping,
684  const DoFHandlerType & dof,
685  const Quadrature<0> & quadrature,
686  const std::map<types::boundary_id,
688  & neumann_bc,
689  const std::vector<const InputVector *> &solutions,
690  std::vector<Vector<float> *> & errors,
691  const ComponentMask & component_mask = ComponentMask(),
692  const Function<spacedim> * coefficients = nullptr,
693  const unsigned int n_threads = numbers::invalid_unsigned_int,
696  const Strategy strategy = cell_diameter_over_24);
697 
702  template <typename InputVector, typename DoFHandlerType>
703  static void
704  estimate(
705  const DoFHandlerType &dof,
706  const Quadrature<0> & quadrature,
707  const std::map<types::boundary_id,
709  & neumann_bc,
710  const std::vector<const InputVector *> &solutions,
711  std::vector<Vector<float> *> & errors,
712  const ComponentMask & component_mask = ComponentMask(),
713  const Function<spacedim> * coefficients = nullptr,
714  const unsigned int n_threads = numbers::invalid_unsigned_int,
717  const Strategy strategy = cell_diameter_over_24);
718 
719 
724  template <typename InputVector, typename DoFHandlerType>
725  static void
726  estimate(
727  const Mapping<1, spacedim> &mapping,
728  const DoFHandlerType & dof,
729  const hp::QCollection<0> & quadrature,
730  const std::map<types::boundary_id,
732  & neumann_bc,
733  const InputVector & solution,
734  Vector<float> & error,
735  const ComponentMask & component_mask = ComponentMask(),
736  const Function<spacedim> *coefficients = nullptr,
737  const unsigned int n_threads = numbers::invalid_unsigned_int,
740  const Strategy strategy = cell_diameter_over_24);
741 
742 
747  template <typename InputVector, typename DoFHandlerType>
748  static void
749  estimate(
750  const DoFHandlerType & dof,
751  const hp::QCollection<0> &quadrature,
752  const std::map<types::boundary_id,
754  & neumann_bc,
755  const InputVector & solution,
756  Vector<float> & error,
757  const ComponentMask & component_mask = ComponentMask(),
758  const Function<spacedim> *coefficients = nullptr,
759  const unsigned int n_threads = numbers::invalid_unsigned_int,
762  const Strategy strategy = cell_diameter_over_24);
763 
764 
769  template <typename InputVector, typename DoFHandlerType>
770  static void
771  estimate(
772  const Mapping<1, spacedim> &mapping,
773  const DoFHandlerType & dof,
774  const hp::QCollection<0> & quadrature,
775  const std::map<types::boundary_id,
777  & neumann_bc,
778  const std::vector<const InputVector *> &solutions,
779  std::vector<Vector<float> *> & errors,
780  const ComponentMask & component_mask = ComponentMask(),
781  const Function<spacedim> * coefficients = nullptr,
782  const unsigned int n_threads = numbers::invalid_unsigned_int,
785  const Strategy strategy = cell_diameter_over_24);
786 
787 
792  template <typename InputVector, typename DoFHandlerType>
793  static void
794  estimate(
795  const DoFHandlerType & dof,
796  const hp::QCollection<0> &quadrature,
797  const std::map<types::boundary_id,
799  & neumann_bc,
800  const std::vector<const InputVector *> &solutions,
801  std::vector<Vector<float> *> & errors,
802  const ComponentMask & component_mask = ComponentMask(),
803  const Function<spacedim> * coefficients = nullptr,
804  const unsigned int n_threads = numbers::invalid_unsigned_int,
807  const Strategy strategy = cell_diameter_over_24);
808 
813  "You provided a ComponentMask argument that is invalid. "
814  "Component masks need to be either default constructed "
815  "(in which case they indicate that every component is "
816  "selected) or need to have a length equal to the number "
817  "of vector components of the finite element in use "
818  "by the DoFHandler object. In the latter case, at "
819  "least one component needs to be selected.");
825  "If you do specify the argument for a (possibly "
826  "spatially variable) coefficient function for this function, "
827  "then it needs to refer to a coefficient that is either "
828  "scalar (has one vector component) or has as many vector "
829  "components as there are in the finite element used by "
830  "the DoFHandler argument.");
836  int,
837  int,
838  << "You provided a function map that for boundary indicator "
839  << arg1 << " specifies a function with " << arg2
840  << " vector components. However, the finite "
841  "element in use has "
842  << arg3
843  << " components, and these two numbers need to match.");
848  int,
849  int,
850  << "The number of input vectors, " << arg1
851  << " needs to be equal to the number of output vectors, "
852  << arg2
853  << ". This is not the case in your call of this function.");
858  "You need to specify at least one solution vector as "
859  "input.");
860 };
861 
862 
863 
865 
866 #endif
DeclExceptionMsg
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:496
KellyErrorEstimator::ExcNoSolutions
static ::ExceptionBase & ExcNoSolutions()
deprecated_function_map.h
KellyErrorEstimator::ExcInvalidComponentMask
static ::ExceptionBase & ExcInvalidComponentMask()
KellyErrorEstimator::cell_diameter_over_24
@ cell_diameter_over_24
Kelly error estimator with the factor .
Definition: error_estimator.h:272
hp::QCollection
Definition: q_collection.h:48
ComponentMask
Definition: component_mask.h:83
KellyErrorEstimator< 1, spacedim >::face_diameter_over_twice_max_degree
@ face_diameter_over_twice_max_degree
Definition: error_estimator.h:601
KellyErrorEstimator::ExcInvalidBoundaryFunction
static ::ExceptionBase & ExcInvalidBoundaryFunction(types::boundary_id arg1, int arg2, int arg3)
KellyErrorEstimator
Definition: error_estimator.h:262
KellyErrorEstimator::cell_diameter
@ cell_diameter
Kelly error estimator with the factor .
Definition: error_estimator.h:277
Mapping
Abstract base class for mapping classes.
Definition: mapping.h:302
types::material_id
unsigned int material_id
Definition: types.h:152
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
KellyErrorEstimator::estimate
static void estimate(const Mapping< dim, spacedim > &mapping, const DoFHandlerType &dof, const Quadrature< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, typename InputVector::value_type > * > &neumann_bc, const InputVector &solution, Vector< float > &error, const ComponentMask &component_mask=ComponentMask(), const Function< spacedim > *coefficients=nullptr, const unsigned int n_threads=numbers::invalid_unsigned_int, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id, const types::material_id material_id=numbers::invalid_material_id, const Strategy strategy=cell_diameter_over_24)
DeclException3
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:564
KellyErrorEstimator::ExcInvalidCoefficient
static ::ExceptionBase & ExcInvalidCoefficient()
component_mask.h
exceptions.h
types::subdomain_id
unsigned int subdomain_id
Definition: types.h:43
numbers::invalid_material_id
const types::material_id invalid_material_id
Definition: types.h:223
unsigned int
function.h
KellyErrorEstimator::ExcIncompatibleNumberOfElements
static ::ExceptionBase & ExcIncompatibleNumberOfElements(int arg1, int arg2)
hp
Definition: hp.h:117
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
KellyErrorEstimator::Strategy
Strategy
Definition: error_estimator.h:269
config.h
KellyErrorEstimator< 1, spacedim >::Strategy
Strategy
Definition: error_estimator.h:595
numbers::invalid_subdomain_id
const types::subdomain_id invalid_subdomain_id
Definition: types.h:285
Function
Definition: function.h:151
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
Quadrature
Definition: quadrature.h:85
KellyErrorEstimator::face_diameter_over_twice_max_degree
@ face_diameter_over_twice_max_degree
Definition: error_estimator.h:275
DeclException2
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:541