Reference documentation for deal.II version 9.0.0
error_estimator.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2018 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 at
12 // the top level of the deal.II distribution.
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 #include <deal.II/base/exceptions.h>
22 #include <deal.II/base/function.h>
23 #include <deal.II/dofs/function_map.h>
24 #include <deal.II/fe/component_mask.h>
25 #include <map>
26 
27 DEAL_II_NAMESPACE_OPEN
28 
29 
30 template <int, int> class Mapping;
31 template <int> class Quadrature;
32 
33 namespace hp
34 {
35  template <int> class QCollection;
36 }
37 
38 
39 
250 template <int dim, int spacedim=dim>
252 {
253 public:
258  enum Strategy
259  {
266  };
267 
328  template <typename InputVector, typename DoFHandlerType>
329  static void estimate
330  (const Mapping<dim, spacedim> &mapping,
331  const DoFHandlerType &dof,
332  const Quadrature<dim-1> &quadrature,
334  const InputVector &solution,
335  Vector<float> &error,
336  const ComponentMask &component_mask = ComponentMask(),
337  const Function<spacedim> *coefficients = nullptr,
338  const unsigned int n_threads = numbers::invalid_unsigned_int,
341  const Strategy strategy = cell_diameter_over_24);
342 
347  template <typename InputVector, typename DoFHandlerType>
348  static void estimate
349  (const DoFHandlerType &dof,
350  const Quadrature<dim-1> &quadrature,
352  const InputVector &solution,
353  Vector<float> &error,
354  const ComponentMask &component_mask = ComponentMask(),
355  const Function<spacedim> *coefficients = nullptr,
356  const unsigned int n_threads = numbers::invalid_unsigned_int,
359  const Strategy strategy = cell_diameter_over_24);
360 
374  template <typename InputVector, typename DoFHandlerType>
375  static void estimate
376  (const Mapping<dim, spacedim> &mapping,
377  const DoFHandlerType &dof,
378  const Quadrature<dim-1> &quadrature,
380  const std::vector<const InputVector *> &solutions,
381  std::vector<Vector<float>*> &errors,
382  const ComponentMask &component_mask = ComponentMask(),
383  const Function<spacedim> *coefficients = 0,
384  const unsigned int n_threads = numbers::invalid_unsigned_int,
387  const Strategy strategy = cell_diameter_over_24);
388 
393  template <typename InputVector, typename DoFHandlerType>
394  static void estimate
395  (const DoFHandlerType &dof,
396  const Quadrature<dim-1> &quadrature,
398  const std::vector<const InputVector *> &solutions,
399  std::vector<Vector<float>*> &errors,
400  const ComponentMask &component_mask = ComponentMask(),
401  const Function<spacedim> *coefficients = 0,
402  const unsigned int n_threads = numbers::invalid_unsigned_int,
405  const Strategy strategy = cell_diameter_over_24);
406 
407 
412  template <typename InputVector, typename DoFHandlerType>
413  static void estimate
414  (const Mapping<dim, spacedim> &mapping,
415  const DoFHandlerType &dof,
416  const hp::QCollection<dim-1> &quadrature,
418  const InputVector &solution,
419  Vector<float> &error,
420  const ComponentMask &component_mask = ComponentMask(),
421  const Function<spacedim> *coefficients = 0,
422  const unsigned int n_threads = numbers::invalid_unsigned_int,
425  const Strategy strategy = cell_diameter_over_24);
426 
427 
432  template <typename InputVector, typename DoFHandlerType>
433  static void estimate
434  (const DoFHandlerType &dof,
435  const hp::QCollection<dim-1> &quadrature,
437  const InputVector &solution,
438  Vector<float> &error,
439  const ComponentMask &component_mask = ComponentMask(),
440  const Function<spacedim> *coefficients = nullptr,
441  const unsigned int n_threads = numbers::invalid_unsigned_int,
444  const Strategy strategy = cell_diameter_over_24);
445 
446 
451  template <typename InputVector, typename DoFHandlerType>
452  static void estimate
453  (const Mapping<dim, spacedim> &mapping,
454  const DoFHandlerType &dof,
455  const hp::QCollection<dim-1> &quadrature,
457  const std::vector<const InputVector *> &solutions,
458  std::vector<Vector<float>*> &errors,
459  const ComponentMask &component_mask = ComponentMask(),
460  const Function<spacedim> *coefficients = 0,
461  const unsigned int n_threads = numbers::invalid_unsigned_int,
464  const Strategy strategy = cell_diameter_over_24);
465 
466 
471  template <typename InputVector, typename DoFHandlerType>
472  static void estimate
473  (const DoFHandlerType &dof,
474  const hp::QCollection<dim-1> &quadrature,
476  const std::vector<const InputVector *> &solutions,
477  std::vector<Vector<float>*> &errors,
478  const ComponentMask &component_mask = ComponentMask(),
479  const Function<spacedim> *coefficients = nullptr,
480  const unsigned int n_threads = numbers::invalid_unsigned_int,
483  const Strategy strategy = cell_diameter_over_24);
484 
489  "You provided a ComponentMask argument that is invalid. "
490  "Component masks need to be either default constructed "
491  "(in which case they indicate that every component is "
492  "selected) or need to have a length equal to the number "
493  "of vector components of the finite element in use "
494  "by the DoFHandler object. In the latter case, at "
495  "least one component needs to be selected.");
500  "If you do specify the argument for a (possibly "
501  "spatially variable) coefficient function for this function, "
502  "then it needs to refer to a coefficient that is either "
503  "scalar (has one vector component) or has as many vector "
504  "components as there are in the finite element used by "
505  "the DoFHandler argument.");
511  int,
512  int,
513  << "You provided a function map that for boundary indicator "
514  << arg1 << " specifies a function with "
515  << arg2 << " vector components. However, the finite "
516  "element in use has "
517  << arg2 << " components, and these two numbers need to match.");
522  int, int,
523  << "The number of input vectors, " << arg1
524  << " needs to be equal to the number of output vectors, "
525  << arg2
526  << ". This is not the case in your call of this function.");
531  "You need to specify at least one solution vector as "
532  "input.");
533 };
534 
535 
536 
548 template <int spacedim>
549 class KellyErrorEstimator<1,spacedim>
550 {
551 public:
556  enum Strategy
557  {
564  };
565 
588  template <typename InputVector, typename DoFHandlerType>
589  static void estimate
590  (const Mapping<1,spacedim> &mapping,
591  const DoFHandlerType &dof,
592  const Quadrature<0> &quadrature,
594  const InputVector &solution,
595  Vector<float> &error,
596  const ComponentMask &component_mask = ComponentMask(),
597  const Function<spacedim> *coefficient = nullptr,
598  const unsigned int n_threads = numbers::invalid_unsigned_int,
601  const Strategy strategy = cell_diameter_over_24);
602 
607  template <typename InputVector, typename DoFHandlerType>
608  static void estimate
609  (const DoFHandlerType &dof,
610  const Quadrature<0> &quadrature,
612  const InputVector &solution,
613  Vector<float> &error,
614  const ComponentMask &component_mask = ComponentMask(),
615  const Function<spacedim> *coefficients = nullptr,
616  const unsigned int n_threads = numbers::invalid_unsigned_int,
619  const Strategy strategy = cell_diameter_over_24);
620 
634  template <typename InputVector, typename DoFHandlerType>
635  static void estimate
636  (const Mapping<1,spacedim> &mapping,
637  const DoFHandlerType &dof,
638  const Quadrature<0> &quadrature,
640  const std::vector<const InputVector *> &solutions,
641  std::vector<Vector<float>*> &errors,
642  const ComponentMask &component_mask = ComponentMask(),
643  const Function<spacedim> *coefficients = 0,
644  const unsigned int n_threads = numbers::invalid_unsigned_int,
647  const Strategy strategy = cell_diameter_over_24);
648 
653  template <typename InputVector, typename DoFHandlerType>
654  static void estimate
655  (const DoFHandlerType &dof,
656  const Quadrature<0> &quadrature,
658  const std::vector<const InputVector *> &solutions,
659  std::vector<Vector<float>*> &errors,
660  const ComponentMask &component_mask = ComponentMask(),
661  const Function<spacedim> *coefficients = 0,
662  const unsigned int n_threads = numbers::invalid_unsigned_int,
665  const Strategy strategy = cell_diameter_over_24);
666 
667 
672  template <typename InputVector, typename DoFHandlerType>
673  static void estimate
674  (const Mapping<1,spacedim> &mapping,
675  const DoFHandlerType &dof,
676  const hp::QCollection<0> &quadrature,
678  const InputVector &solution,
679  Vector<float> &error,
680  const ComponentMask &component_mask = ComponentMask(),
681  const Function<spacedim> *coefficients = 0,
682  const unsigned int n_threads = numbers::invalid_unsigned_int,
685  const Strategy strategy = cell_diameter_over_24);
686 
687 
692  template <typename InputVector, typename DoFHandlerType>
693  static void estimate
694  (const DoFHandlerType &dof,
695  const hp::QCollection<0> &quadrature,
697  const InputVector &solution,
698  Vector<float> &error,
699  const ComponentMask &component_mask = ComponentMask(),
700  const Function<spacedim> *coefficients = 0,
701  const unsigned int n_threads = numbers::invalid_unsigned_int,
704  const Strategy strategy = cell_diameter_over_24);
705 
706 
711  template <typename InputVector, typename DoFHandlerType>
712  static void estimate
713  (const Mapping<1,spacedim> &mapping,
714  const DoFHandlerType &dof,
715  const hp::QCollection<0> &quadrature,
717  const std::vector<const InputVector *> &solutions,
718  std::vector<Vector<float>*> &errors,
719  const ComponentMask &component_mask = ComponentMask(),
720  const Function<spacedim> *coefficients = 0,
721  const unsigned int n_threads = numbers::invalid_unsigned_int,
724  const Strategy strategy = cell_diameter_over_24);
725 
726 
731  template <typename InputVector, typename DoFHandlerType>
732  static void estimate
733  (const DoFHandlerType &dof,
734  const hp::QCollection<0> &quadrature,
736  const std::vector<const InputVector *> &solutions,
737  std::vector<Vector<float>*> &errors,
738  const ComponentMask &component_mask = ComponentMask(),
739  const Function<spacedim> *coefficients = 0,
740  const unsigned int n_threads = numbers::invalid_unsigned_int,
743  const Strategy strategy = cell_diameter_over_24);
744 
749  "You provided a ComponentMask argument that is invalid. "
750  "Component masks need to be either default constructed "
751  "(in which case they indicate that every component is "
752  "selected) or need to have a length equal to the number "
753  "of vector components of the finite element in use "
754  "by the DoFHandler object. In the latter case, at "
755  "least one component needs to be selected.");
760  "If you do specify the argument for a (possibly "
761  "spatially variable) coefficient function for this function, "
762  "then it needs to refer to a coefficient that is either "
763  "scalar (has one vector component) or has as many vector "
764  "components as there are in the finite element used by "
765  "the DoFHandler argument.");
771  int,
772  int,
773  << "You provided a function map that for boundary indicator "
774  << arg1 << " specifies a function with "
775  << arg2 << " vector components. However, the finite "
776  "element in use has "
777  << arg3 << " components, and these two numbers need to match.");
782  int, int,
783  << "The number of input vectors, " << arg1
784  << " needs to be equal to the number of output vectors, "
785  << arg2
786  << ". This is not the case in your call of this function.");
791  "You need to specify at least one solution vector as "
792  "input.");
793 };
794 
795 
796 
797 DEAL_II_NAMESPACE_CLOSE
798 
799 #endif
static ::ExceptionBase & ExcInvalidBoundaryFunction(types::boundary_id arg1, int arg2, int arg3)
static const unsigned int invalid_unsigned_int
Definition: types.h:173
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:358
static ::ExceptionBase & ExcInvalidCoefficient()
the boundary residual estimator with the factor .
const types::subdomain_id invalid_subdomain_id
Definition: types.h:248
static ::ExceptionBase & ExcInvalidComponentMask()
unsigned int boundary_id
Definition: types.h:110
unsigned int material_id
Definition: types.h:133
Abstract base class for mapping classes.
Definition: dof_tools.h:46
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:335
static void estimate(const Mapping< dim, spacedim > &mapping, const DoFHandlerType &dof, const Quadrature< dim-1 > &quadrature, const typename FunctionMap< spacedim, typename InputVector::value_type >::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)
std::map< types::boundary_id, const Function< dim, Number > * > type
Definition: function_map.h:81
static ::ExceptionBase & ExcIncompatibleNumberOfElements(int arg1, int arg2)
unsigned int subdomain_id
Definition: types.h:42
Definition: hp.h:102
the boundary residual estimator with the factor .
Kelly error estimator with the factor .
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:370
Kelly error estimator with the factor .
const types::material_id invalid_material_id
Definition: types.h:194
static ::ExceptionBase & ExcNoSolutions()