deal.II version 9.7.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\}}\)
Loading...
Searching...
No Matches
error_estimator.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 1998 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_error_estimator_h
16#define dealii_error_estimator_h
17
18
19#include <deal.II/base/config.h>
20
23
25
27
28#include <map>
29
31
32// Forward declarations
33#ifndef DOXYGEN
34template <int dim, int spacedim>
36class DoFHandler;
37template <int, int>
38class Mapping;
39template <int>
40class Quadrature;
41
42namespace hp
43{
44 template <int>
45 class QCollection;
46
47 template <int, int>
49} // namespace hp
50#endif
51
52
264template <int dim, int spacedim = dim>
266{
267public:
282
343 template <typename Number>
344 static void
347 const DoFHandler<dim, spacedim> &dof,
348 const Quadrature<dim - 1> &quadrature,
349 const std::map<types::boundary_id, const Function<spacedim, Number> *>
350 &neumann_bc,
351 const ReadVector<Number> &solution,
352 Vector<float> &error,
353 const ComponentMask &component_mask = {},
354 const Function<spacedim> *coefficients = nullptr,
355 const unsigned int n_threads = numbers::invalid_unsigned_int,
358 const Strategy strategy = cell_diameter_over_24);
359
364 template <typename Number>
365 static void
367 const DoFHandler<dim, spacedim> &dof,
368 const Quadrature<dim - 1> &quadrature,
369 const std::map<types::boundary_id, const Function<spacedim, Number> *>
370 &neumann_bc,
371 const ReadVector<Number> &solution,
372 Vector<float> &error,
373 const ComponentMask &component_mask = {},
374 const Function<spacedim> *coefficients = nullptr,
375 const unsigned int n_threads = numbers::invalid_unsigned_int,
378 const Strategy strategy = cell_diameter_over_24);
379
393 template <typename Number>
394 static void
397 const DoFHandler<dim, spacedim> &dof,
398 const Quadrature<dim - 1> &quadrature,
399 const std::map<types::boundary_id, const Function<spacedim, Number> *>
400 &neumann_bc,
401 const ArrayView<const ReadVector<Number> *> &solutions,
402 ArrayView<Vector<float> *> &errors,
403 const ComponentMask &component_mask = {},
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 Number>
415 static void
417 const DoFHandler<dim, spacedim> &dof,
418 const Quadrature<dim - 1> &quadrature,
419 const std::map<types::boundary_id, const Function<spacedim, Number> *>
420 &neumann_bc,
421 const ArrayView<const ReadVector<Number> *> &solutions,
422 ArrayView<Vector<float> *> &errors,
423 const ComponentMask &component_mask = {},
424 const Function<spacedim> *coefficients = nullptr,
425 const unsigned int n_threads = numbers::invalid_unsigned_int,
428 const Strategy strategy = cell_diameter_over_24);
429
430
435 template <typename Number>
436 static void
439 const DoFHandler<dim, spacedim> &dof,
440 const hp::QCollection<dim - 1> &quadrature,
441 const std::map<types::boundary_id, const Function<spacedim, Number> *>
442 &neumann_bc,
443 const ReadVector<Number> &solution,
444 Vector<float> &error,
445 const ComponentMask &component_mask = {},
446 const Function<spacedim> *coefficients = nullptr,
447 const unsigned int n_threads = numbers::invalid_unsigned_int,
450 const Strategy strategy = cell_diameter_over_24);
451
452
460 template <typename Number>
461 DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
462 "Use the version of this function which takes a hp::MappingCollection instead.")
463 static void estimate(
464 const Mapping<dim, spacedim> &mapping,
465 const DoFHandler<dim, spacedim> &dof,
466 const hp::QCollection<dim - 1> &quadrature,
467 const std::map<types::boundary_id, const Function<spacedim, Number> *>
468 &neumann_bc,
469 const ReadVector<Number> &solution,
470 Vector<float> &error,
471 const ComponentMask &component_mask = {},
472 const Function<spacedim> *coefficients = nullptr,
473 const unsigned int n_threads = numbers::invalid_unsigned_int,
476 const Strategy strategy = cell_diameter_over_24);
477
478
483 template <typename Number>
484 static void
486 const DoFHandler<dim, spacedim> &dof,
487 const hp::QCollection<dim - 1> &quadrature,
488 const std::map<types::boundary_id, const Function<spacedim, Number> *>
489 &neumann_bc,
490 const ReadVector<Number> &solution,
491 Vector<float> &error,
492 const ComponentMask &component_mask = {},
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 Number>
505 static void
508 const DoFHandler<dim, spacedim> &dof,
509 const hp::QCollection<dim - 1> &quadrature,
510 const std::map<types::boundary_id, const Function<spacedim, Number> *>
511 &neumann_bc,
512 const ArrayView<const ReadVector<Number> *> &solutions,
513 ArrayView<Vector<float> *> &errors,
514 const ComponentMask &component_mask = {},
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
521
529 template <typename Number>
530 DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
531 "Use the version of this function which takes a hp::MappingCollection instead.")
532 static void estimate(
533 const Mapping<dim, spacedim> &mapping,
534 const DoFHandler<dim, spacedim> &dof,
535 const hp::QCollection<dim - 1> &quadrature,
536 const std::map<types::boundary_id, const Function<spacedim, Number> *>
537 &neumann_bc,
538 const ArrayView<const ReadVector<Number> *> &solutions,
539 ArrayView<Vector<float> *> &errors,
540 const ComponentMask &component_mask = {},
541 const Function<spacedim> *coefficients = nullptr,
542 const unsigned int n_threads = numbers::invalid_unsigned_int,
545 const Strategy strategy = cell_diameter_over_24);
546
547
552 template <typename Number>
553 static void
555 const DoFHandler<dim, spacedim> &dof,
556 const hp::QCollection<dim - 1> &quadrature,
557 const std::map<types::boundary_id, const Function<spacedim, Number> *>
558 &neumann_bc,
559 const ArrayView<const ReadVector<Number> *> &solutions,
560 ArrayView<Vector<float> *> &errors,
561 const ComponentMask &component_mask = {},
562 const Function<spacedim> *coefficients = nullptr,
563 const unsigned int n_threads = numbers::invalid_unsigned_int,
566 const Strategy strategy = cell_diameter_over_24);
567
572 "You provided a ComponentMask argument that is invalid. "
573 "Component masks need to be either default constructed "
574 "(in which case they indicate that every component is "
575 "selected) or need to have a length equal to the number "
576 "of vector components of the finite element in use "
577 "by the DoFHandler object. In the latter case, at "
578 "least one component needs to be selected.");
584 "If you do specify the argument for a (possibly "
585 "spatially variable) coefficient function for this function, "
586 "then it needs to refer to a coefficient that is either "
587 "scalar (has one vector component) or has as many vector "
588 "components as there are in the finite element used by "
589 "the DoFHandler argument.");
595 int,
596 int,
597 << "You provided a function map that for boundary indicator "
598 << arg1 << " specifies a function with " << arg2
599 << " vector components. However, the finite "
600 "element in use has "
601 << arg3
602 << " components, and these two numbers need to match.");
607 int,
608 int,
609 << "The number of input vectors, " << arg1
610 << " needs to be equal to the number of output vectors, "
611 << arg2
612 << ". This is not the case in your call of this function.");
617 "You need to specify at least one solution vector as "
618 "input.");
619};
620
621
622
632template <int spacedim>
633class KellyErrorEstimator<1, spacedim>
634{
635public:
650
651
652
675 template <typename Number>
676 static void
677 estimate(
679 const DoFHandler<1, spacedim> &dof,
680 const Quadrature<0> &quadrature,
681 const std::map<types::boundary_id, const Function<spacedim, Number> *>
682 &neumann_bc,
683 const ReadVector<Number> &solution,
684 Vector<float> &error,
685 const ComponentMask &component_mask = {},
686 const Function<spacedim> *coefficient = nullptr,
687 const unsigned int n_threads = numbers::invalid_unsigned_int,
690 const Strategy strategy = cell_diameter_over_24);
691
692
693
698 template <typename Number>
699 static void
700 estimate(
701 const DoFHandler<1, spacedim> &dof,
702 const Quadrature<0> &quadrature,
703 const std::map<types::boundary_id, const Function<spacedim, Number> *>
704 &neumann_bc,
705 const ReadVector<Number> &solution,
706 Vector<float> &error,
707 const ComponentMask &component_mask = {},
708 const Function<spacedim> *coefficients = nullptr,
709 const unsigned int n_threads = numbers::invalid_unsigned_int,
712 const Strategy strategy = cell_diameter_over_24);
713
714
715
729 template <typename Number>
730 static void
731 estimate(
732 const Mapping<1, spacedim> &mapping,
733 const DoFHandler<1, spacedim> &dof,
734 const Quadrature<0> &quadrature,
735 const std::map<types::boundary_id, const Function<spacedim, Number> *>
736 &neumann_bc,
737 const ArrayView<const ReadVector<Number> *> &solutions,
738 ArrayView<Vector<float> *> &errors,
739 const ComponentMask &component_mask = {},
740 const Function<spacedim> *coefficients = nullptr,
741 const unsigned int n_threads = numbers::invalid_unsigned_int,
744 const Strategy strategy = cell_diameter_over_24);
745
746
747
752 template <typename Number>
753 static void
754 estimate(
755 const DoFHandler<1, spacedim> &dof,
756 const Quadrature<0> &quadrature,
757 const std::map<types::boundary_id, const Function<spacedim, Number> *>
758 &neumann_bc,
759 const ArrayView<const ReadVector<Number> *> &solutions,
760 ArrayView<Vector<float> *> &errors,
761 const ComponentMask &component_mask = {},
762 const Function<spacedim> *coefficients = nullptr,
763 const unsigned int n_threads = numbers::invalid_unsigned_int,
766 const Strategy strategy = cell_diameter_over_24);
767
768
769
774 template <typename Number>
775 static void
776 estimate(
777 const hp::MappingCollection<1, spacedim> &mapping,
778 const DoFHandler<1, spacedim> &dof,
779 const hp::QCollection<0> &quadrature,
780 const std::map<types::boundary_id, const Function<spacedim, Number> *>
781 &neumann_bc,
782 const ReadVector<Number> &solution,
783 Vector<float> &error,
784 const ComponentMask &component_mask = {},
785 const Function<spacedim> *coefficients = nullptr,
786 const unsigned int n_threads = numbers::invalid_unsigned_int,
789 const Strategy strategy = cell_diameter_over_24);
790
791
792
800 template <typename Number>
801 DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
802 "Use the version of this function which takes a hp::MappingCollection instead.")
803 static void estimate(
804 const Mapping<1, spacedim> &mapping,
805 const DoFHandler<1, spacedim> &dof,
806 const hp::QCollection<0> &quadrature,
807 const std::map<types::boundary_id, const Function<spacedim, Number> *>
808 &neumann_bc,
809 const ReadVector<Number> &solution,
810 Vector<float> &error,
811 const ComponentMask &component_mask = {},
812 const Function<spacedim> *coefficients = nullptr,
813 const unsigned int n_threads = numbers::invalid_unsigned_int,
816 const Strategy strategy = cell_diameter_over_24);
817
818
819
824 template <typename Number>
825 static void
826 estimate(
827 const DoFHandler<1, spacedim> &dof,
828 const hp::QCollection<0> &quadrature,
829 const std::map<types::boundary_id, const Function<spacedim, Number> *>
830 &neumann_bc,
831 const ReadVector<Number> &solution,
832 Vector<float> &error,
833 const ComponentMask &component_mask = {},
834 const Function<spacedim> *coefficients = nullptr,
835 const unsigned int n_threads = numbers::invalid_unsigned_int,
838 const Strategy strategy = cell_diameter_over_24);
839
840
841
846 template <typename Number>
847 static void
848 estimate(
849 const hp::MappingCollection<1, spacedim> &mapping,
850 const DoFHandler<1, spacedim> &dof,
851 const hp::QCollection<0> &quadrature,
852 const std::map<types::boundary_id, const Function<spacedim, Number> *>
853 &neumann_bc,
854 const ArrayView<const ReadVector<Number> *> &solutions,
855 ArrayView<Vector<float> *> &errors,
856 const ComponentMask &component_mask = {},
857 const Function<spacedim> *coefficients = nullptr,
858 const unsigned int n_threads = numbers::invalid_unsigned_int,
861 const Strategy strategy = cell_diameter_over_24);
862
863
864
872 template <typename Number>
873 DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
874 "Use the version of this function which takes a hp::MappingCollection instead.")
875 static void estimate(
876 const Mapping<1, spacedim> &mapping,
877 const DoFHandler<1, spacedim> &dof,
878 const hp::QCollection<0> &quadrature,
879 const std::map<types::boundary_id, const Function<spacedim, Number> *>
880 &neumann_bc,
881 const ArrayView<const ReadVector<Number> *> &solutions,
882 ArrayView<Vector<float> *> &errors,
883 const ComponentMask &component_mask = {},
884 const Function<spacedim> *coefficients = nullptr,
885 const unsigned int n_threads = numbers::invalid_unsigned_int,
888 const Strategy strategy = cell_diameter_over_24);
889
890
891
896 template <typename Number>
897 static void
898 estimate(
899 const DoFHandler<1, spacedim> &dof,
900 const hp::QCollection<0> &quadrature,
901 const std::map<types::boundary_id, const Function<spacedim, Number> *>
902 &neumann_bc,
903 const ArrayView<const ReadVector<Number> *> &solutions,
904 ArrayView<Vector<float> *> &errors,
905 const ComponentMask &component_mask = {},
906 const Function<spacedim> *coefficients = nullptr,
907 const unsigned int n_threads = numbers::invalid_unsigned_int,
910 const Strategy strategy = cell_diameter_over_24);
911
912
913
918 "You provided a ComponentMask argument that is invalid. "
919 "Component masks need to be either default constructed "
920 "(in which case they indicate that every component is "
921 "selected) or need to have a length equal to the number "
922 "of vector components of the finite element in use "
923 "by the DoFHandler object. In the latter case, at "
924 "least one component needs to be selected.");
925
931 "If you do specify the argument for a (possibly "
932 "spatially variable) coefficient function for this function, "
933 "then it needs to refer to a coefficient that is either "
934 "scalar (has one vector component) or has as many vector "
935 "components as there are in the finite element used by "
936 "the DoFHandler argument.");
937
943 int,
944 int,
945 << "You provided a function map that for boundary indicator "
946 << arg1 << " specifies a function with " << arg2
947 << " vector components. However, the finite "
948 "element in use has "
949 << arg3
950 << " components, and these two numbers need to match.");
951
956 int,
957 int,
958 << "The number of input vectors, " << arg1
959 << " needs to be equal to the number of output vectors, "
960 << arg2
961 << ". This is not the case in your call of this function.");
962
967 "You need to specify at least one solution vector as "
968 "input.");
969};
970
971
973
974#endif
@ cell_diameter
Kelly error estimator with the factor .
@ cell_diameter_over_24
Kelly error estimator with the factor .
static void estimate(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, Number > * > &neumann_bc, const ReadVector< Number > &solution, Vector< float > &error, const ComponentMask &component_mask={}, 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)
static void estimate(const DoFHandler< dim, spacedim > &dof, const Quadrature< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, Number > * > &neumann_bc, const ArrayView< const ReadVector< Number > * > &solutions, ArrayView< Vector< float > * > &errors, const ComponentMask &component_mask={}, 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)
static void estimate(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, Number > * > &neumann_bc, const ArrayView< const ReadVector< Number > * > &solutions, ArrayView< Vector< float > * > &errors, const ComponentMask &component_mask={}, 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)
static void estimate(const DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, Number > * > &neumann_bc, const ArrayView< const ReadVector< Number > * > &solutions, ArrayView< Vector< float > * > &errors, const ComponentMask &component_mask={}, 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)
static void estimate(const hp::MappingCollection< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, Number > * > &neumann_bc, const ReadVector< Number > &solution, Vector< float > &error, const ComponentMask &component_mask={}, 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)
@ cell_diameter_over_24
Kelly error estimator with the factor .
@ cell_diameter
Kelly error estimator with the factor .
static void estimate(const DoFHandler< dim, spacedim > &dof, const Quadrature< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, Number > * > &neumann_bc, const ReadVector< Number > &solution, Vector< float > &error, const ComponentMask &component_mask={}, 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)
static void estimate(const DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, Number > * > &neumann_bc, const ReadVector< Number > &solution, Vector< float > &error, const ComponentMask &component_mask={}, 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)
static void estimate(const hp::MappingCollection< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const hp::QCollection< dim - 1 > &quadrature, const std::map< types::boundary_id, const Function< spacedim, Number > * > &neumann_bc, const ArrayView< const ReadVector< Number > * > &solutions, ArrayView< Vector< float > * > &errors, const ComponentMask &component_mask={}, 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)
Abstract base class for mapping classes.
Definition mapping.h:320
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:248
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
static ::ExceptionBase & ExcInvalidBoundaryFunction(types::boundary_id arg1, int arg2, int arg3)
static ::ExceptionBase & ExcInvalidBoundaryFunction(types::boundary_id arg1, int arg2, int arg3)
static ::ExceptionBase & ExcIncompatibleNumberOfElements(int arg1, int arg2)
static ::ExceptionBase & ExcNoSolutions()
#define DeclException2(Exception2, type1, type2, outsequence)
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcInvalidComponentMask()
#define DeclException3(Exception3, type1, type2, type3, outsequence)
static ::ExceptionBase & ExcInvalidCoefficient()
static ::ExceptionBase & ExcInvalidCoefficient()
static ::ExceptionBase & ExcNoSolutions()
static ::ExceptionBase & ExcIncompatibleNumberOfElements(int arg1, int arg2)
MappingQ< dim, spacedim > StaticMappingQ1< dim, spacedim >::mapping
Definition mapping_q1.h:104
Definition hp.h:117
constexpr unsigned int invalid_unsigned_int
Definition types.h:238
constexpr types::material_id invalid_material_id
Definition types.h:294
constexpr types::subdomain_id invalid_subdomain_id
Definition types.h:381
STL namespace.
Definition types.h:32
unsigned int material_id
Definition types.h:184
unsigned int boundary_id
Definition types.h:161
unsigned int subdomain_id
Definition types.h:52