KASKADE 7 development version
concepts.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-2023 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
71{
72public:
74 typedef unspecified Scalar;
75
77 using GridView = unspecified;
78
80 typedef unspecified Grid;
81
83 typedef unspecified IndexSet;
84
86 typedef unspecified ShapeFunctionSet;
87
91 typedef unspecified GlobalIndexRange;
92
96 typedef unspecified SortedIndexRange;
97
100
102
107
110
112 GridView const& gridView() const;
113
120 ShapeFunctionSet const& shapefunctions(Cell const& cell) const;
121
123 ShapeFunctionSet const& shapefunctions(size_t index) const;
124
125
127 size_t size() const;
128
132 int maxOrder() const;
133
139
145
152
158 GlobalIndexRange globalIndices(size_t index) const;
159
160
166 Combiner combiner(Cell const& cell, size_t index) const;
167
169 void update();
170};
171
172
191{
192public:
200
203 ConverterConcept(Cell<Grid> const& cell);
204
211 void moveTo(Cell<Grid> const& cell);
212
214 void setLocalPosition(LocalPosition<Grid> const& x);
215
218
220 VariationalArg<RT,dim,k> global(VariationalArg<RT,dim,k> const& sf, int deriv=1) const;
221
225 template <class Scalar>
226 VariationalArg<Scalar,dim,1> global(std::tuple<Dune::FieldVector<Scalar,1>,
228 Tensor<Scalar,1,dim,dim>> const& sf) const;
229
232};
233
234
271{
272public:
287 template <int n, int m>
288 void rightTransform(DynamicMatrix<Dune::FieldMatrix<Real,n,m>>& A) const;
289
296 template <int n, int m>
297 void rightTransform(std::vector<VariationalArg<Real,n,m>>& v) const;
298
304 template <int n, int m>
305 void leftPseudoInverse(DynamicMatrix<Dune::FieldMatrix<Real,n,m>>& A) const;
306
308 operator Dune::BCRSMatrix<Dune::FieldMatrix<Real,1,1>>() const;
309};
310
311//---------------------------------------------------------------------
312
324{
325public:
330 typedef unspecified Scalar;
331
335 static int const components = unspecified;
336
341
347 template <class Cell>
348 int order(Cell const& cell) const;
349
354 template <class Cell>
355 ValueType value(Cell const& cell,
356 Dune::FieldVector<typename Cell::Geometry::ctype,
357 Cell::dimension> const& localCoordinate) const;
358};
359
360//---------------------------------------------------------------------
361
375{
379 template<class Cell>
380 double operator()(Cell const& c) const;
384 double scalingFactor() const;
388 double offset() const;
389};
390
391//---------------------------------------------------------------------
392
405{
406public:
416 template <class Cell, class Sequence>
417 void operator()(Cell const& cell,
419 Sequence& values) const;
420};
421
422//---------------------------------------------------------------------
423
437{
438 template <int row, int col>
439 class D2
440 {
441 static bool const present;
442 static bool const lumped;
443 };
444};
445
446//---------------------------------------------------------------------
447
476template <class X, class Y>
477class MatrixRepresentedOperator: public Dune::LinearOperator<X,Y>
478{
479public:
481 typedef X domain_type;
483 typedef Y range_type;
484
486 template <class MatrixType>
487 MatrixType get() const;
488
490
499 template <class Vector>
500 void rangeToVector(Y const& y, Vector& coord) const;
501
503
512 template <class Vector>
513 void vectorToDomain(Vector const& coord, X& x) const;
514};
515
516
526{
527public:
532
533 template <class Cell>
534 int integrationOrder(Cell const& /* ci */, int shapeFunctionOrder) const;
535
547 template <class Cell, class Index, class Values>
548 void operator()(Cell const& cell, Index idx,
549 LocalPosition xi, double weight,
550 Values const& x);
551
557 void join(Collector& c);
558};
559
573{
574public:
578 using Space = unspecified;
579
583 Space const& space() const;
584
588 Dune::FieldVector<unspecified,unspecified> value(typename Space::Evaluator const& evaluator) const;
589};
590
602{
603public:
609 Dune::FieldVector<unspecified,unspecified> value(Cell const& cell, Position const& xloc) const;
610};
611
612
613
621{
622public:
626 using Scalar = unspecified;
627
631 static int const dim = unspecified;
632
636 using Tensor = unspecified;
637
641 using VoigtTensor = unspecified;
642
648
652 Scalar d0() const;
653
657 Scalar d1(Tensor const& e1) const;
658
664 Scalar d2(Tensor const& e1, Tensor const& e2) const;
665
674 Tensor stress() const;
675
695};
696
704{
705public:
709 using Scalar = unspecified;
710
714 static int const dim = unspecified;
715
720
735
739 Scalar d0() const;
740
744 Scalar d1(Invariants const& i1) const;
745
749 Scalar d2(Invariants const& i1, Invariants const& i2) const;
750};
751
752// ------------------------------------------------------------------------------------------------
753
763{
764 template <int row> struct D1 { static bool const assemble = unspecified; };
765 template <int row, int col> struct D2 { static bool const assemble = unspecified; };
766};
Defines an interface for integrating values computed over a grid.
Definition: concepts.hh:526
Collector(Collector const &)
Copy constructor.
void join(Collector &c)
Merges the result of two concurrently filled collectors.
int integrationOrder(Cell const &, int shapeFunctionOrder) const
void operator()(Cell const &cell, Index idx, LocalPosition xi, double weight, Values const &x)
Collects values computed on the grid.
Algebraic combiner of global degrees of freedom to shape function coefficients.
Definition: concepts.hh:271
void rightTransform(std::vector< VariationalArg< Real, n, m > > &v) const
In-place computation of row vectors .
void rightTransform(DynamicMatrix< Dune::FieldMatrix< Real, n, m > > &A) const
In-place computation of .
void leftPseudoInverse(DynamicMatrix< Dune::FieldMatrix< Real, n, m > > &A) const
In-place computation of . Given a set of shape function coefficients, this computes the best fit of g...
Geometrical transformation of shape function values.
Definition: concepts.hh:191
void setLocalPosition(LocalPosition< Grid > const &x)
Sets the argument of the transformation .
Dune::FieldMatrix< RT, k, 1 > local(Dune::FieldMatrix< RT, k, 1 >) const
Transforms global values to local shape function values.
VariationalArg< RT, dim, k > global(VariationalArg< RT, dim, k > const &sf, int deriv=1) const
Applies the transformation to shape function value and derivatives up to order deriv.
Dune::FieldMatrix< RT, k, 1 > global(Dune::FieldMatrix< RT, k, 1 > const &sf) const
Applies the transformation to shape function value.
ConverterConcept(Cell< Grid > const &cell)
VariationalArg< Scalar, dim, 1 > global(std::tuple< Dune::FieldVector< Scalar, 1 >, Dune::FieldMatrix< Scalar, 1, dim >, Tensor< Scalar, 1, dim, dim > > const &sf) const
Applies the transformation to shape function value, gradient, and Hessian, returning a filled Variat...
void moveTo(Cell< Grid > const &cell)
Defines the cell on which the converter lives.
ConverterConcept()
Create a converter.
Function is the interface for evaluatable functions .
Definition: concepts.hh:324
int order(Cell const &cell) const
static int const components
Definition: concepts.hh:335
unspecified Scalar
Definition: concepts.hh:330
Dune::FieldVector< Scalar, components > ValueType
Definition: concepts.hh:340
ValueType value(Cell const &cell, Dune::FieldVector< typename Cell::Geometry::ctype, Cell::dimension > const &localCoordinate) const
A function that supports efficient evaluation via an evaluator.
Definition: concepts.hh:573
unspecified Space
The type of finite element space the evaluators of which are required for evaluating this FunctionVie...
Definition: concepts.hh:578
Space const & space() const
Provides a reference to the finite element space required for evaluating this FunctionView.
Dune::FieldVector< unspecified, unspecified > value(typename Space::Evaluator const &evaluator) const
Evaluates the function view.
Defines the block structure of defect systems to be solved by the hierarchical error estimator.
Definition: concepts.hh:437
A hyperelastic material law.
Definition: concepts.hh:621
Scalar d0() const
Evaluates the stored energy density .
unspecified VoigtTensor
The type of stress and strain tensors (usually Dune::FieldMatrix<Scalar,dim*(dim+1)/2,...
Definition: concepts.hh:641
unspecified Tensor
The type of stress and strain tensors (usually Dune::FieldMatrix<Scalar,dim,dim>).
Definition: concepts.hh:636
VoigtTensor strainToStressMatrix() const
Returns the tangent stiffness tensor mapping strain tensor variations to stress tensor (2nd Piola-Ki...
Scalar d1(Tensor const &e1) const
Evaluates the first directional derivative .
Tensor stress() const
Returns the 2nd Piola-Kirchhoff stress tensor corresponding to the current strain .
static int const dim
The spatial dimension (usually 2 or 3).
Definition: concepts.hh:631
Scalar d2(Tensor const &e1, Tensor const &e2) const
Evaluates the second directional derivative .
void setLinearizationPoint(Tensor const &e)
Defines a new evaluation/linearization point.
unspecified Scalar
The scalar field type (usually double).
Definition: concepts.hh:626
A hyperelastic material law formulated in terms of the invariants of the Cauchy-Green strain tensor ...
Definition: concepts.hh:704
Scalar d1(Invariants const &i1) const
Evaluates the first directional derivative .
static int const dim
The spatial dimension (usually 2 or 3).
Definition: concepts.hh:714
unspecified Scalar
The scalar field type (usually double).
Definition: concepts.hh:709
void InvariantsMaterialConcept(Invariants const &i,...)
Constructs the material at shifted invariants .
Scalar d0() const
Evaluates the stored energy density .
Scalar d2(Invariants const &i1, Invariants const &i2) const
Evaluates the second directional derivative .
Management of degrees of freedom.
Definition: concepts.hh:71
unspecified Grid
The grid type of the grid view.
Definition: concepts.hh:80
GlobalIndexRange globalIndices(size_t index) const
The global indices associated with the cell with given index. Returns a const sequence containing th...
CombinerConcept Combiner
The type of the algebraic combiner associated with this cell.
Definition: concepts.hh:106
unspecified GridView
The grid view on which the finite element functions are defined.
Definition: concepts.hh:77
size_t size() const
Returns the number of global degrees of freedom managed.
Kaskade::Cell< Grid > Cell
The type of codimension 0 entities in the grid.
Definition: concepts.hh:109
SortedIndexRange sortedIndices(Cell const &cell) const
Returns an immutable sequence of (global index, local index) pairs sorted in ascending global index o...
unspecified IndexSet
The entity indexing of the grid view.
Definition: concepts.hh:83
void update()
Recomputes the node management, e.g. after grid modifications.
ConverterConcept Converter
The type of the geometrical transformation .
Definition: concepts.hh:99
GridView const & gridView() const
Returns the grid view used.
unspecified SortedIndexRange
A standard sequence with value type std::pair<IndexSet::IndexType,int>.
Definition: concepts.hh:96
ShapeFunctionSet const & shapefunctions(Cell const &cell) const
Returns the shape function set living on the given cell.
unspecified Scalar
The type of scalars used for shape function values.
Definition: concepts.hh:74
int maxOrder() const
Returns the maximal polynomial order of shape functions encountered in any cell.
GlobalIndexRange globalIndices(Cell const &cell) const
The global indices associated with this cell. Returns a const sequence containing the global indices...
unspecified ShapeFunctionSet
The type of shape function sets living on the cells.
Definition: concepts.hh:86
unspecified GlobalIndexRange
A standard sequence with value type IndexSet::IndexType.
Definition: concepts.hh:91
SortedIndexRange sortedIndices(size_t n) const
Returns an immutable sequence of (global index, local index) pairs sorted in ascending global index o...
Combiner combiner(Cell const &cell, size_t index) const
Returns the algebraic combinator.
ShapeFunctionSet const & shapefunctions(size_t index) const
Returns the shape function set living on the cell with given index.
Defines a linear operator with matrix representation .
Definition: concepts.hh:478
void vectorToDomain(Vector const &coord, X &x) const
Get from coordinate vector .
void rangeToVector(Y const &y, Vector &coord) const
Get coordinate vector from .
MatrixType get() const
Get matrix-representation of the operator.
A callable that allows to implement weighted (scaled) norms.
Definition: concepts.hh:405
void operator()(Cell const &cell, Dune::FieldVector< typename Cell::Geometry::ctype, Cell::dimension > const &localPos, Sequence &values) const
A function that supports efficient evaluation via cell and local coordinate.
Definition: concepts.hh:602
Dune::FieldVector< unspecified, unspecified > value(Cell const &cell, Position const &xloc) const
Evaluates the function view.
typename GridView::template Codim< 0 >::Entity Cell
The type of cells (entities of codimension 0) in the grid view.
Definition: gridBasics.hh:35
Dune::FieldVector< typename GridView::ctype, GridView::dimension > LocalPosition
The type of local positions within cells of the grid view.
Definition: gridBasics.hh:59
Dune::FieldVector< Scalar, dim > Vector
static bool const assemble
Definition: concepts.hh:764
static bool const assemble
Definition: concepts.hh:765
A block filter specifying which subblocks of a stiffness matrix or right hand side shall be assembled...
Definition: concepts.hh:763
A Weighing is a class that allows to compute weighted averages of values associated to cells of a gri...
Definition: concepts.hh:375
double scalingFactor() const
Defines the global scaling factor .
double operator()(Cell const &c) const
This defines a weight for each cell of the grid.
double offset() const
Defines the global offset .