KASKADE 7 development version
variationalfunctional.hh
Go to the documentation of this file.
1/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
2/* */
3/* This file is part of the library KASKADE 7 */
4/* https://www.zib.de/research/projects/kaskade7-finite-element-toolbox */
5/* */
6/* Copyright (C) 2002-2021 Zuse Institute Berlin */
7/* */
8/* KASKADE 7 is distributed under the terms of the ZIB Academic License. */
9/* see $KASKADE/academic.txt */
10/* */
11/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
12
18class unspecified;
19
64{
65public:
69 typedef unspecified Scalar;
70
76 typedef unspecified AnsatzVars;
77
78
84 typedef unspecified TestVars;
85
99 typedef unspecified OriginVars;
100
104 static ProblemType const type;
105
110 typedef typename Vars::Grid::template Codim<0>::Entity Entity;
111
143 {
144
149 void moveTo(Cell const&)
150 {
151 }
152
157 template<class Position, class Evaluators>
158 void evaluateAt(Position const& x, Evaluators const& evaluators)
159 {
160 }
161
165 Scalar d0() const;
166
172 template <int row, int dim>
173 Dune::FieldVector<Scalar,TestVars::template Components<row>::m> d1(VariationalArg<Scalar,dim> const& arg) const;
174
190 template <int row, int col, int dim>
191 Dune::FieldMatrix<Scalar,TestVars::template Components<row>::m,AnsatzVars::template Components<row>::m>
192 d2(VariationalArg<Scalar,dim> const& argT, VariationalArg<Scalar,dim> const& argA) const;
193 };
194
206 {
212 void moveTo(typename Vars::GridView::IntersectionIterator const&);
213
219 template<class Position, class Evaluators>
220 void evaluateAt(Position const& xi, Evaluators const& evaluators);
221
225 Scalar d0() const;
226
233 template <int row, int dim>
234 Dune::FieldVector<Scalar,TestVars::template Components<row>::m> d1(VariationalArg<Scalar,dim> const& argT) const;
235
246 template <int row, int col, int dim>
247 Dune::FieldMatrix<Scalar,TestVars::template Components<row>::m,AnsatzVars::template Components<row>::m>
248 d2(VariationalArg<Scalar,dim> const& argT, VariationalArg<Scalar,dim> const& argA) const;
249 };
250
260 {
261 template <class FaceIterator>
262 void moveTo(FaceIterator const&);
263
271 template <class Position, class Evaluators>
272 void evaluateAt(Position const& x, Evaluators const& evaluators, Evaluators const& neighbourEvaluators);
273
277 Scalar d0() const;
278
286 template <int row, int dim>
287 Dune::FieldVector<Scalar,TestVars::template Components<row>::m> d1(VariationalArg<Scalar,dim> const& argT) const;
288
301 template <int row, int col, int dim>
303 d2(VariationalArg<Scalar,dim> const& argT, VariationalArg<Scalar,dim> const& argA, bool centerCell) const;
304 };
305
315 template <int row>
316 struct D1
317 {
324 static bool const present;
325
329 static int const derivatives;
330 };
331
341 template <int row, int col>
342 struct D2
343 {
357 static bool const present;
358
374 static bool const symmetric;
375
395 static const bool makePositive;
396
403 static bool const lumped;
404
410 static int const derivatives;
411 };
412
420 DomainCache createDomainCache(int flags) const {}
421
430
441
442
452 int integrationOrder(typename Variables::Grid::template Codim<0>::Entity const& cell,
453 int shapeFunctionOrder, bool boundary) const;
454
466 bool considerFace(Face<typename Variables::GridView> const& face) const;
467
468};
469
470// ------------------------------------------------------------------------------------------------
471// ------------------------------------------------------------------------------------------------
472
509{
510public:
514 typedef unspecified Scalar;
515
521 typedef unspecified AnsatzVars;
522
523
529 typedef unspecified TestVars;
530
534 typedef unspecified OriginVars;
535
540 static ProblemType const type;
541
546 typedef typename Vars::Grid::template Codim<0>::Entity Entity;
547
555 struct DomainCache : public EvalCacheBase
556 {
557
567 DomainCache(Functional const& f, typename Vars::VariableSet const& u, int flags) {};
568
573 template<class Entity>
574 void moveTo(Entity const&)
575 {
576 }
577
582 template<class Position, class Evaluators>
583 void evaluateAt(Position const& x, Evaluators const& evaluators)
584 {
585 }
586
590 Scalar d0() const;
591
597 template <int row, int dim>
598 Dune::FieldVector<Scalar,Variables::template Components<row>::m> d1(VariationalArg<Scalar,dim> const& arg) const;
599
609 template <int row, int col, int dim>
610 Dune::FieldMatrix<Scalar,Variables::template Components<row>::m,Variables::template Components<col>::m>
611 d2(VariationalArg<Scalar,dim> const& argT, VariationalArg<Scalar,dim> const& argA) const;
612 };
613
624 struct BoundaryCache : public EvalCacheBase
625 {
629 typedef typename AnsatzVars::Grid::LeafIntersectionIterator FaceIterator;
630
640 BoundaryCache(Functional const& f, typename Vars::VariableSet const& u, int const flags=7) {};
641
645 template<class FaceIterator>
646 void moveTo(FaceIterator const&)
647 {
648 }
649
653 template<class Position, class Evaluators>
654 void evaluateAt(Position const& x, Evaluators const& evaluators)
655 {
656 }
657
661 Scalar d0() const;
662
669 template <int row, int dim>
670 Dune::FieldVector<Scalar,Variables::template Components<row>::m> d1(VariationalArg<Scalar,dim> const& arg) const;
671
682 template <int row, int col, int dim>
683 Dune::FieldMatrix<Scalar,Variables::template Components<row>::m,Variables::template Components<col>::m>
684 d2(VariationalArg<Scalar,dim> const& argT, VariationalArg<Scalar,dim> const& argA) const;
685 };
686
687
695 template <int row>
696 struct D1
697 {
702 static bool const present;
703 };
704
712 template <int row, int col>
713 struct D2
714 {
723 static bool const present;
724
732 static bool const symmetric;
733
739 static bool const lumped;
740 };
741
749 int integrationOrder(typename Variables::Grid::template Codim<0>::Entity const& cell,
750 int shapeFunctionOrder, bool boundary) const;
751
763 bool considerFace(Face<typename Variables::GridView> const& face) const;
764};
765
766// ------------------------------------------------------------------------------------------------
767// ------------------------------------------------------------------------------------------------
768
805{
806public:
810 typedef unspecified Scalar;
811
817 typedef unspecified AnsatzVars;
818
819
824 typedef unspecified TestVars;
825
829 typedef unspecified OriginVars;
830
835 static ProblemType const type;
836
842
850 struct DomainCache : public EvalCacheBase
851 {
852
865 DomainCache(Functional const& f, typename OriginVars::VariableSet const& u, int flags) {};
866
871 void moveTo(Cell const&)
872 {
873 }
874
879 template<class Position, class Evaluators>
880 void evaluateAt(Position const& x, Evaluators const& evaluators)
881 {
882 }
883
888 template <int row, int dim>
890 d1(VariationalArg<Scalar,dim> const& arg) const;
891
901 template <int row, int col, int dim>
903 Variables::template components<col>>
904 d2(VariationalArg<Scalar,dim> const& argT, VariationalArg<Scalar,dim> const& argA) const;
905
914 template <int row, int col, int dim>
916 Variables::template components<col>>
917 b2(VariationalArg<Scalar,dim> const& argT, VariationalArg<Scalar,dim> const& argA) const;
918 };
919
930 struct BoundaryCache : public EvalCacheBase
931 {
935 typedef typename AnsatzVars::Grid::LeafIntersectionIterator FaceIterator;
936
946 BoundaryCache(Functional const& f, typename OriginVars::VariableSet const& u, int const flags=7) {};
947
951 template<class FaceIterator>
952 void moveTo(FaceIterator const&)
953 {
954 }
955
959 template<class Position, class Evaluators>
960 void evaluateAt(Position const& x, Evaluators const& evaluators)
961 {
962 }
963
966 template <int row, int dim>
967 Dune::FieldVector<Scalar,Variables::template Components<row>::m> d1(VariationalArg<Scalar,dim> const& arg) const;
968
974 template <int row, int col, int dim>
976 Variables::template Components<col>>
977 d2(VariationalArg<Scalar,dim> const& argT, VariationalArg<Scalar,dim> const& argA) const;
978
984 template <int row, int col, int dim>
986 Variables::template Components<col>>
987 b2(VariationalArg<Scalar,dim> const& argT, VariationalArg<Scalar,dim> const& argA) const;
988 };
989
990
998 template <int row>
999 struct D1
1000 {
1005 static bool const present;
1006 };
1007
1015 template <int row, int col>
1016 struct D2
1017 {
1026 static bool const present;
1027
1035 static bool const symmetric;
1036
1042 static bool const lumped;
1043 };
1044
1052 template <int row, int col>
1053 struct B2
1054 {
1059 static bool const present;
1060
1068 static bool const symmetric;
1069
1075 static bool const lumped;
1076
1083 static bool const constant;
1084 };
1085
1089 double time() const;
1090
1097
1108 void temporalEvaluationRange(double t0, double t1);
1109
1117 int integrationOrder(typename Variables::Grid::template Codim<0>::Entity const& cell,
1118 int shapeFunctionOrder, bool boundary) const;
1119};
1120
1121// ------------------------------------------------------------------------------------------------
1122// ------------------------------------------------------------------------------------------------
1123
1134{
1138 template <int row, int col> class D2;
1139};
Documentation of the concept of a nonlinear variational functional The nonlinear variational function...
static ProblemType const type
The type of problem, either VariationalFunctional or WeakFormulation.
unspecified AnsatzVars
A description of the ansatz variables occuring in the variational problem. This should be an instance...
int integrationOrder(typename Variables::Grid::template Codim< 0 >::Entity const &cell, int shapeFunctionOrder, bool boundary) const
Vars::Grid::template Codim< 0 >::Entity Entity
bool considerFace(Face< typename Variables::GridView > const &face) const
unspecified Scalar
The scalar type to be used for this variational problem.
unspecified OriginVars
A description of the type of the point of linearization.
Documentation of the concept of a parabolic parabolic equation.
unspecified Scalar
The scalar type to be used for this variational problem.
unspecified OriginVars
A description of the type of current evaluation state.
unspecified AnsatzVars
A description of the ansatz variables occuring in the weak formulation. This should be an instance of...
int integrationOrder(typename Variables::Grid::template Codim< 0 >::Entity const &cell, int shapeFunctionOrder, bool boundary) const
void temporalEvaluationRange(double t0, double t1)
Notifies the PDE in which time interval evaluations will be done.
static ProblemType const type
The type of problem, either VariationalFunctional (for symmetric problems) or WeakFormulation.
double time() const
Reports the current simulated time for which the evaluations are done.
Kaskade::Cell< typename AnsatzVars::GridView > Cell
ParabolicEquation & setTime(double t)
Sets a new time for evaluating PDE coefficients.
Documentation of the concept of a quadratic variational functionalThe variational functional concept ...
DomainCache createDomainCache(int flags) const
Creates a DomainCache.
Vars::Grid::template Codim< 0 >::Entity Entity
int integrationOrder(typename Variables::Grid::template Codim< 0 >::Entity const &cell, int shapeFunctionOrder, bool boundary) const
This method defines a suitable quadrature order to be used on the given cell.
BoundaryCache createBoundaryCache(int flags) const
Creates a BoundaryCache.
static ProblemType const type
The type of problem, either VariationalFunctional or WeakFormulation.
bool considerFace(Face< typename Variables::GridView > const &face) const
unspecified Scalar
The scalar type to be used for this variational problem.
unspecified OriginVars
A description of the variables defining the linearization point.
InnerBoundaryCache createInnerBoundaryCache(int flags) const
Creates an InnerBoundaryCache.
unspecified AnsatzVars
A description of the ansatz variables occuring in the variational problem.
unspecified TestVars
A description of the test variables occuring in the "variational" problem.
typename boost::fusion::result_of::as_vector< typename boost::fusion::result_of::transform< Spaces, GetEvaluatorTypes >::type >::type Evaluators
the heterogeneous sequence type of Evaluators for the given spaces
ProblemType
A type for describing the nature of a PDE problem.
typename GridView::template Codim< 0 >::Entity Cell
The type of cells (entities of codimension 0) in the grid view.
Definition: gridBasics.hh:35
Concept for providing block information to hierarchical error estimator.
Dune::FieldVector< Scalar, Variables::template Components< row >::m > d1(VariationalArg< Scalar, dim > const &arg) const
void evaluateAt(Position const &x, Evaluators const &evaluators)
BoundaryCache(Functional const &f, typename Vars::VariableSet const &u, int const flags=7)
Dune::FieldMatrix< Scalar, Variables::template Components< row >::m, Variables::template Components< col >::m > d2(VariationalArg< Scalar, dim > const &argT, VariationalArg< Scalar, dim > const &argA) const
AnsatzVars::Grid::LeafIntersectionIterator FaceIterator
Type of the Face Information.
This evaluation cache concept defines the interface that is accessed by the assembler.
void evaluateAt(Position const &x, Evaluators const &evaluators)
DomainCache(Functional const &f, typename Vars::VariableSet const &u, int flags)
Dune::FieldMatrix< Scalar, Variables::template Components< row >::m, Variables::template Components< col >::m > d2(VariationalArg< Scalar, dim > const &argT, VariationalArg< Scalar, dim > const &argA) const
Dune::FieldVector< Scalar, Variables::template Components< row >::m > d1(VariationalArg< Scalar, dim > const &arg) const
static bool const present
Presence flag. Is true if that block is statically present, i.e. .
Defining boundary conditions.
Dune::FieldMatrix< Scalar, Variables::template components< row >, Variables::template Components< col > > b2(VariationalArg< Scalar, dim > const &argT, VariationalArg< Scalar, dim > const &argA) const
void evaluateAt(Position const &x, Evaluators const &evaluators)
Dune::FieldVector< Scalar, Variables::template Components< row >::m > d1(VariationalArg< Scalar, dim > const &arg) const
AnsatzVars::Grid::LeafIntersectionIterator FaceIterator
Type of the Face Information.
Dune::FieldMatrix< Scalar, Variables::template components< row >, Variables::template Components< col > > d2(VariationalArg< Scalar, dim > const &argT, VariationalArg< Scalar, dim > const &argA) const
BoundaryCache(Functional const &f, typename OriginVars::VariableSet const &u, int const flags=7)
This evaluation cache concept defines the interface that is accessed by the assembler.
void evaluateAt(Position const &x, Evaluators const &evaluators)
Dune::FieldMatrix< Scalar, Variables::template components< row >, Variables::template components< col > > b2(VariationalArg< Scalar, dim > const &argT, VariationalArg< Scalar, dim > const &argA) const
Dune::FieldMatrix< Scalar, Variables::template components< row >, Variables::template components< col > > d2(VariationalArg< Scalar, dim > const &argT, VariationalArg< Scalar, dim > const &argA) const
DomainCache(Functional const &f, typename OriginVars::VariableSet const &u, int flags)
Dune::FieldVector< Scalar, Variables::template components< row > > d1(VariationalArg< Scalar, dim > const &arg) const
Evaluation of boundary conditions.
void moveTo(typename Vars::GridView::IntersectionIterator const &)
Construct all data that is constant in an entity.
void evaluateAt(Position const &xi, Evaluators const &evaluators)
Construct all data that is constant at one point.
Dune::FieldVector< Scalar, TestVars::template Components< row >::m > d1(VariationalArg< Scalar, dim > const &argT) const
Dune::FieldMatrix< Scalar, TestVars::template Components< row >::m, AnsatzVars::template Components< row >::m > d2(VariationalArg< Scalar, dim > const &argT, VariationalArg< Scalar, dim > const &argA) const
Provides static information about the right hand side blocks.
static bool const present
Provides static information about the subvector blocks of the right hand side.
static int const derivatives
Determines the highest derivatives that need to be available for right hand side evaluation.
Provides static information about the submatrix blocks of the stiffness block matrix.
static bool const present
Specifies the presence of a submatrix block.
static int const derivatives
Determines the highest derivatives that need to be available for evaluation of this block.
static bool const lumped
Specifies whether only the diagonal of the subblock shall be assembled. Is true if only the diagonal ...
static const bool makePositive
If this flag is true (and symmetric==true), the assembler will enforce positive semidefiniteness of a...
static bool const symmetric
Specifies whether the subblock is symmetric.
This evaluation cache concept defines the interface that is accessed by the assembler.
void moveTo(Cell const &)
Construct all data that is constant in an entity. Default implementation in EvalCacheBase: does nothi...
void evaluateAt(Position const &x, Evaluators const &evaluators)
Construct all data that is constant at one point. Default implementation in EvalCacheBase: assert(fal...
Dune::FieldMatrix< Scalar, TestVars::template Components< row >::m, AnsatzVars::template Components< row >::m > d2(VariationalArg< Scalar, dim > const &argT, VariationalArg< Scalar, dim > const &argA) const
Dune::FieldVector< Scalar, TestVars::template Components< row >::m > d1(VariationalArg< Scalar, dim > const &arg) const
Evaluates jump contributions at interior faces, suitable for discontinuous Galerkin methods (optional...
Scalar d0() const
Returns the value .
Dune::FieldMatrix< Scalar, TestVars::::template Components< row >::m, AnsatzVars::::template Components< row >::m > d2(VariationalArg< Scalar, dim > const &argT, VariationalArg< Scalar, dim > const &argA, bool centerCell) const
Dune::FieldVector< Scalar, TestVars::template Components< row >::m > d1(VariationalArg< Scalar, dim > const &argT) const
void moveTo(FaceIterator const &)
void evaluateAt(Position const &x, Evaluators const &evaluators, Evaluators const &neighbourEvaluators)
Prepare the cache for evaluation at given quadrature point.