version 0.4.1
fetraits.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2021-2024 The Ikarus Developers mueller@ibb.uni-stuttgart.de
2// SPDX-License-Identifier: LGPL-3.0-or-later
3
9#pragma once
10
13
14namespace Ikarus {
15
24template <typename BH, typename FER = FERequirements<>, bool useEigenRef = false, bool useFlat = true>
26{
28 using BasisHandler = BH;
29
31 using FERequirementType = FER;
32
34 static constexpr bool useFlatBasis = useFlat;
35
38
41
43 using Basis = std::conditional_t<useFlat, FlatBasis, UntouchedBasis>;
44
46 using LocalView = typename Basis::LocalView;
47
49 using GridView = typename Basis::GridView;
50
52 using Element = typename LocalView::Element;
53
55 using Geometry = typename Element::Geometry;
56
58 using GlobalIndex = typename LocalView::MultiIndex;
59
61 using ctype = double;
62
64 static constexpr int worlddim = Geometry::coorddimension;
65
67 static constexpr int mydim = Element::mydimension;
68
70 static constexpr int dimension = Element::dimension;
71
73 using GlobalCoordinates = Eigen::Matrix<ctype, worlddim, 1>;
74
76 using ParameterSpaceType = Eigen::Matrix<ctype, mydim, 1>;
77
79 template <typename ScalarType = ctype>
80 using VectorType =
81 std::conditional_t<useEigenRef, Eigen::Ref<Eigen::VectorX<ScalarType>>, Eigen::VectorX<ScalarType>&>;
82
84 template <typename ScalarType = ctype>
85 using MatrixType =
86 std::conditional_t<useEigenRef, Eigen::Ref<Eigen::MatrixX<ScalarType>>, Eigen::MatrixX<ScalarType>&>;
87};
88
89} // namespace Ikarus
Several concepts.
Definition of the LinearElastic class for finite element mechanics computations.
Definition: simpleassemblers.hh:22
Traits for handling finite elements.
Definition: fetraits.hh:26
FER FERequirementType
Type of the requirements for the finite element.
Definition: fetraits.hh:31
std::conditional_t< useEigenRef, Eigen::Ref< Eigen::VectorX< ScalarType > >, Eigen::VectorX< ScalarType > & > VectorType
Type of the internal forces.
Definition: fetraits.hh:81
typename Basis::LocalView LocalView
Type of the local view.
Definition: fetraits.hh:46
static constexpr int dimension
Dimension of the grid.
Definition: fetraits.hh:70
BH BasisHandler
Type of the basis of the finite element.
Definition: fetraits.hh:28
typename BasisHandler::UntouchedBasis UntouchedBasis
Type of the untouched basis.
Definition: fetraits.hh:40
typename Basis::GridView GridView
Type of the grid view.
Definition: fetraits.hh:49
Eigen::Matrix< ctype, mydim, 1 > ParameterSpaceType
Type of the ParameterSpace coordinate.
Definition: fetraits.hh:76
static constexpr bool useFlatBasis
A bool to indicate if the provided basishandler should hand out the flat basis.
Definition: fetraits.hh:34
std::conditional_t< useFlat, FlatBasis, UntouchedBasis > Basis
Type of the basis version.
Definition: fetraits.hh:43
typename BasisHandler::FlatBasis FlatBasis
Type of the flat basis.
Definition: fetraits.hh:37
static constexpr int worlddim
Dimension of the world space.
Definition: fetraits.hh:64
static constexpr int mydim
Dimension of the geometry.
Definition: fetraits.hh:67
double ctype
Type used for coordinates.
Definition: fetraits.hh:61
typename LocalView::Element Element
Type of the grid element.
Definition: fetraits.hh:52
typename Element::Geometry Geometry
Type of the element geometry.
Definition: fetraits.hh:55
Eigen::Matrix< ctype, worlddim, 1 > GlobalCoordinates
Type of the coordinate.
Definition: fetraits.hh:73
typename LocalView::MultiIndex GlobalIndex
Type of the global index.
Definition: fetraits.hh:58
std::conditional_t< useEigenRef, Eigen::Ref< Eigen::MatrixX< ScalarType > >, Eigen::MatrixX< ScalarType > & > MatrixType
Type of the stiffness matrix.
Definition: fetraits.hh:86
decltype(Dune::Functions::DefaultGlobalBasis(std::declval< PreBasis >())) UntouchedBasis
The type of the untouched basis.
Definition: utils/basis.hh:35
decltype(Dune::Functions::DefaultGlobalBasis(Ikarus::flatPreBasis(std::declval< PreBasis >()))) FlatBasis
The type of the flattened basis.
Definition: utils/basis.hh:37