Reference documentation for deal.II version 9.1.0-pre
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.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 
22 #include <deal.II/base/exceptions.h>
23 #include <deal.II/base/function.h>
24 
25 #include <deal.II/dofs/deprecated_function_map.h>
26 
27 #include <deal.II/fe/component_mask.h>
28 
29 #include <map>
30 
31 DEAL_II_NAMESPACE_OPEN
32 
33 
34 template <int, int>
35 class Mapping;
36 template <int>
37 class Quadrature;
38 
39 namespace hp
40 {
41  template <int>
42  class QCollection;
43 }
44 
45 
46 
258 template <int dim, int spacedim = dim>
260 {
261 public:
266  enum Strategy
267  {
269  cell_diameter_over_24 = 0,
274  cell_diameter
275  };
276 
337  template <typename InputVector, typename DoFHandlerType>
338  static void
339  estimate(
340  const Mapping<dim, spacedim> &mapping,
341  const DoFHandlerType & dof,
342  const Quadrature<dim - 1> & quadrature,
343  const std::map<types::boundary_id,
345  & neumann_bc,
346  const InputVector & solution,
347  Vector<float> & error,
348  const ComponentMask & component_mask = ComponentMask(),
349  const Function<spacedim> *coefficients = nullptr,
350  const unsigned int n_threads = numbers::invalid_unsigned_int,
353  const Strategy strategy = cell_diameter_over_24);
354 
359  template <typename InputVector, typename DoFHandlerType>
360  static void
361  estimate(
362  const DoFHandlerType & dof,
363  const Quadrature<dim - 1> &quadrature,
364  const std::map<types::boundary_id,
366  & neumann_bc,
367  const InputVector & solution,
368  Vector<float> & error,
369  const ComponentMask & component_mask = ComponentMask(),
370  const Function<spacedim> *coefficients = nullptr,
371  const unsigned int n_threads = numbers::invalid_unsigned_int,
374  const Strategy strategy = cell_diameter_over_24);
375 
389  template <typename InputVector, typename DoFHandlerType>
390  static void
391  estimate(
392  const Mapping<dim, spacedim> &mapping,
393  const DoFHandlerType & dof,
394  const Quadrature<dim - 1> & quadrature,
395  const std::map<types::boundary_id,
397  & neumann_bc,
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 
411  template <typename InputVector, typename DoFHandlerType>
412  static void
413  estimate(
414  const DoFHandlerType & dof,
415  const Quadrature<dim - 1> &quadrature,
416  const std::map<types::boundary_id,
418  & neumann_bc,
419  const std::vector<const InputVector *> &solutions,
420  std::vector<Vector<float> *> & errors,
421  const ComponentMask & component_mask = ComponentMask(),
422  const Function<spacedim> * coefficients = 0,
423  const unsigned int n_threads = numbers::invalid_unsigned_int,
426  const Strategy strategy = cell_diameter_over_24);
427 
428 
433  template <typename InputVector, typename DoFHandlerType>
434  static void
435  estimate(
436  const Mapping<dim, spacedim> & mapping,
437  const DoFHandlerType & dof,
438  const hp::QCollection<dim - 1> &quadrature,
439  const std::map<types::boundary_id,
441  & neumann_bc,
442  const InputVector & solution,
443  Vector<float> & error,
444  const ComponentMask & component_mask = ComponentMask(),
445  const Function<spacedim> *coefficients = 0,
446  const unsigned int n_threads = numbers::invalid_unsigned_int,
449  const Strategy strategy = cell_diameter_over_24);
450 
451 
456  template <typename InputVector, typename DoFHandlerType>
457  static void
458  estimate(
459  const DoFHandlerType & dof,
460  const hp::QCollection<dim - 1> &quadrature,
461  const std::map<types::boundary_id,
463  & neumann_bc,
464  const InputVector & solution,
465  Vector<float> & error,
466  const ComponentMask & component_mask = ComponentMask(),
467  const Function<spacedim> *coefficients = nullptr,
468  const unsigned int n_threads = numbers::invalid_unsigned_int,
471  const Strategy strategy = cell_diameter_over_24);
472 
473 
478  template <typename InputVector, typename DoFHandlerType>
479  static void
480  estimate(
481  const Mapping<dim, spacedim> & mapping,
482  const DoFHandlerType & dof,
483  const hp::QCollection<dim - 1> &quadrature,
484  const std::map<types::boundary_id,
486  & neumann_bc,
487  const std::vector<const InputVector *> &solutions,
488  std::vector<Vector<float> *> & errors,
489  const ComponentMask & component_mask = ComponentMask(),
490  const Function<spacedim> * coefficients = 0,
491  const unsigned int n_threads = numbers::invalid_unsigned_int,
494  const Strategy strategy = cell_diameter_over_24);
495 
496 
501  template <typename InputVector, typename DoFHandlerType>
502  static void
503  estimate(
504  const DoFHandlerType & dof,
505  const hp::QCollection<dim - 1> &quadrature,
506  const std::map<types::boundary_id,
508  & neumann_bc,
509  const std::vector<const InputVector *> &solutions,
510  std::vector<Vector<float> *> & errors,
511  const ComponentMask & component_mask = ComponentMask(),
512  const Function<spacedim> * coefficients = nullptr,
513  const unsigned int n_threads = numbers::invalid_unsigned_int,
516  const Strategy strategy = cell_diameter_over_24);
517 
521  DeclExceptionMsg(ExcInvalidComponentMask,
522  "You provided a ComponentMask argument that is invalid. "
523  "Component masks need to be either default constructed "
524  "(in which case they indicate that every component is "
525  "selected) or need to have a length equal to the number "
526  "of vector components of the finite element in use "
527  "by the DoFHandler object. In the latter case, at "
528  "least one component needs to be selected.");
533  ExcInvalidCoefficient,
534  "If you do specify the argument for a (possibly "
535  "spatially variable) coefficient function for this function, "
536  "then it needs to refer to a coefficient that is either "
537  "scalar (has one vector component) or has as many vector "
538  "components as there are in the finite element used by "
539  "the DoFHandler argument.");
543  DeclException3(ExcInvalidBoundaryFunction,
545  int,
546  int,
547  << "You provided a function map that for boundary indicator "
548  << arg1 << " specifies a function with " << arg2
549  << " vector components. However, the finite "
550  "element in use has "
551  << arg2
552  << " components, and these two numbers need to match.");
556  DeclException2(ExcIncompatibleNumberOfElements,
557  int,
558  int,
559  << "The number of input vectors, " << arg1
560  << " needs to be equal to the number of output vectors, "
561  << arg2
562  << ". This is not the case in your call of this function.");
566  DeclExceptionMsg(ExcNoSolutions,
567  "You need to specify at least one solution vector as "
568  "input.");
569 };
570 
571 
572 
584 template <int spacedim>
585 class KellyErrorEstimator<1, spacedim>
586 {
587 public:
592  enum Strategy
593  {
595  cell_diameter_over_24 = 0,
600  cell_diameter
601  };
602 
625  template <typename InputVector, typename DoFHandlerType>
626  static void
627  estimate(
628  const Mapping<1, spacedim> &mapping,
629  const DoFHandlerType & dof,
630  const Quadrature<0> & quadrature,
631  const std::map<types::boundary_id,
633  & neumann_bc,
634  const InputVector & solution,
635  Vector<float> & error,
636  const ComponentMask & component_mask = ComponentMask(),
637  const Function<spacedim> *coefficient = nullptr,
638  const unsigned int n_threads = numbers::invalid_unsigned_int,
641  const Strategy strategy = cell_diameter_over_24);
642 
647  template <typename InputVector, typename DoFHandlerType>
648  static void
649  estimate(
650  const DoFHandlerType &dof,
651  const Quadrature<0> & quadrature,
652  const std::map<types::boundary_id,
654  & neumann_bc,
655  const InputVector & solution,
656  Vector<float> & error,
657  const ComponentMask & component_mask = ComponentMask(),
658  const Function<spacedim> *coefficients = nullptr,
659  const unsigned int n_threads = numbers::invalid_unsigned_int,
662  const Strategy strategy = cell_diameter_over_24);
663 
677  template <typename InputVector, typename DoFHandlerType>
678  static void
679  estimate(
680  const Mapping<1, spacedim> &mapping,
681  const DoFHandlerType & dof,
682  const Quadrature<0> & quadrature,
683  const std::map<types::boundary_id,
685  & neumann_bc,
686  const std::vector<const InputVector *> &solutions,
687  std::vector<Vector<float> *> & errors,
688  const ComponentMask & component_mask = ComponentMask(),
689  const Function<spacedim> * coefficients = 0,
690  const unsigned int n_threads = numbers::invalid_unsigned_int,
693  const Strategy strategy = cell_diameter_over_24);
694 
699  template <typename InputVector, typename DoFHandlerType>
700  static void
701  estimate(
702  const DoFHandlerType &dof,
703  const Quadrature<0> & quadrature,
704  const std::map<types::boundary_id,
706  & neumann_bc,
707  const std::vector<const InputVector *> &solutions,
708  std::vector<Vector<float> *> & errors,
709  const ComponentMask & component_mask = ComponentMask(),
710  const Function<spacedim> * coefficients = 0,
711  const unsigned int n_threads = numbers::invalid_unsigned_int,
714  const Strategy strategy = cell_diameter_over_24);
715 
716 
721  template <typename InputVector, typename DoFHandlerType>
722  static void
723  estimate(
724  const Mapping<1, spacedim> &mapping,
725  const DoFHandlerType & dof,
726  const hp::QCollection<0> & quadrature,
727  const std::map<types::boundary_id,
729  & neumann_bc,
730  const InputVector & solution,
731  Vector<float> & error,
732  const ComponentMask & component_mask = ComponentMask(),
733  const Function<spacedim> *coefficients = 0,
734  const unsigned int n_threads = numbers::invalid_unsigned_int,
737  const Strategy strategy = cell_diameter_over_24);
738 
739 
744  template <typename InputVector, typename DoFHandlerType>
745  static void
746  estimate(
747  const DoFHandlerType & dof,
748  const hp::QCollection<0> &quadrature,
749  const std::map<types::boundary_id,
751  & neumann_bc,
752  const InputVector & solution,
753  Vector<float> & error,
754  const ComponentMask & component_mask = ComponentMask(),
755  const Function<spacedim> *coefficients = 0,
756  const unsigned int n_threads = numbers::invalid_unsigned_int,
759  const Strategy strategy = cell_diameter_over_24);
760 
761 
766  template <typename InputVector, typename DoFHandlerType>
767  static void
768  estimate(
769  const Mapping<1, spacedim> &mapping,
770  const DoFHandlerType & dof,
771  const hp::QCollection<0> & quadrature,
772  const std::map<types::boundary_id,
774  & neumann_bc,
775  const std::vector<const InputVector *> &solutions,
776  std::vector<Vector<float> *> & errors,
777  const ComponentMask & component_mask = ComponentMask(),
778  const Function<spacedim> * coefficients = 0,
779  const unsigned int n_threads = numbers::invalid_unsigned_int,
782  const Strategy strategy = cell_diameter_over_24);
783 
784 
789  template <typename InputVector, typename DoFHandlerType>
790  static void
791  estimate(
792  const DoFHandlerType & dof,
793  const hp::QCollection<0> &quadrature,
794  const std::map<types::boundary_id,
796  & neumann_bc,
797  const std::vector<const InputVector *> &solutions,
798  std::vector<Vector<float> *> & errors,
799  const ComponentMask & component_mask = ComponentMask(),
800  const Function<spacedim> * coefficients = 0,
801  const unsigned int n_threads = numbers::invalid_unsigned_int,
804  const Strategy strategy = cell_diameter_over_24);
805 
809  DeclExceptionMsg(ExcInvalidComponentMask,
810  "You provided a ComponentMask argument that is invalid. "
811  "Component masks need to be either default constructed "
812  "(in which case they indicate that every component is "
813  "selected) or need to have a length equal to the number "
814  "of vector components of the finite element in use "
815  "by the DoFHandler object. In the latter case, at "
816  "least one component needs to be selected.");
821  ExcInvalidCoefficient,
822  "If you do specify the argument for a (possibly "
823  "spatially variable) coefficient function for this function, "
824  "then it needs to refer to a coefficient that is either "
825  "scalar (has one vector component) or has as many vector "
826  "components as there are in the finite element used by "
827  "the DoFHandler argument.");
831  DeclException3(ExcInvalidBoundaryFunction,
833  int,
834  int,
835  << "You provided a function map that for boundary indicator "
836  << arg1 << " specifies a function with " << arg2
837  << " vector components. However, the finite "
838  "element in use has "
839  << arg3
840  << " components, and these two numbers need to match.");
844  DeclException2(ExcIncompatibleNumberOfElements,
845  int,
846  int,
847  << "The number of input vectors, " << arg1
848  << " needs to be equal to the number of output vectors, "
849  << arg2
850  << ". This is not the case in your call of this function.");
854  DeclExceptionMsg(ExcNoSolutions,
855  "You need to specify at least one solution vector as "
856  "input.");
857 };
858 
859 
860 
861 DEAL_II_NAMESPACE_CLOSE
862 
863 #endif
static const unsigned int invalid_unsigned_int
Definition: types.h:173
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:420
const types::subdomain_id invalid_subdomain_id
Definition: types.h:255
unsigned int material_id
Definition: types.h:134
unsigned int subdomain_id
Definition: types.h:43
Abstract base class for mapping classes.
Definition: dof_tools.h:57
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:397
Definition: hp.h:102
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:432
const types::material_id invalid_material_id
Definition: types.h:196
unsigned int boundary_id
Definition: types.h:111