Dune Core Modules (2.7.0)

mimeticall.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_MIMETIC_ALL_HH
4#define DUNE_MIMETIC_ALL_HH
5
6#include <cstddef>
7
11
12#include <dune/geometry/type.hh>
13
14#include "../common/localbasis.hh"
15#include "../common/localkey.hh"
16
17namespace Dune
18{
19 template<class D, class R, int dim>
20 class MimeticLocalBasis
21 {
22 public:
25
26 MimeticLocalBasis (unsigned int variant_)
27 : variant(variant_)
28 {}
29
30 MimeticLocalBasis ()
31 : variant(0)
32 {}
33
34 unsigned int size () const { return variant; }
35
37 inline void evaluateFunction (
38 const typename Traits::DomainType& in,
39 std::vector<typename Traits::RangeType>& out) const
40 {
41 DUNE_THROW(Dune::Exception,"mimetic basis evaluation not available");
42 }
43
45 inline void evaluateJacobian (
46 const typename Traits::DomainType& in,
47 std::vector<typename Traits::JacobianType>& out) const
48 {
49 DUNE_THROW(Dune::Exception,"mimetic basis Jacobian evaluation not available");
50 }
51
53 void partial (const std::array<unsigned int, dim>& /*order*/,
54 const typename Traits::DomainType& /*in*/, // position
55 std::vector<typename Traits::RangeType>& /*out*/) const // return value
56 {
57 DUNE_THROW(Dune::Exception,"mimetic basis partial derivative evaluation not available");
58 }
59
61 unsigned int order () const
62 {
63 DUNE_THROW(Dune::Exception,"mimetic order evaluation not available");
64 }
65
66 private:
67 unsigned int variant;
68 };
69
70 template<class LB>
71 class MimeticLocalInterpolation
72 {
73 public:
74
76 template<typename F, typename C>
77 void interpolate (const F& f, std::vector<C>& out) const {
78 DUNE_THROW(Dune::Exception,"mimetic local interpolation not available");
79 }
80 };
81
86 {
87 public:
88 MimeticLocalCoefficients (unsigned int variant_)
89 : variant(variant_), li(variant_)
90 {
91 for (unsigned int i=0; i<variant; i++)
93 }
94
96 : variant(0), li(0)
97 {}
98
100 std::size_t size () const { return variant; }
101
103 const Dune::LocalKey& localKey (std::size_t i) const {
104 return li[i];
105 }
106
107 private:
108 unsigned int variant;
109 std::vector<Dune::LocalKey> li;
110 };
111}
112
113#endif
Base class for Dune-Exceptions.
Definition: exceptions.hh:94
A dense n x m matrix.
Definition: fmatrix.hh:69
vector space out of a tensor product of fields.
Definition: fvector.hh:96
Describe position of one degree of freedom.
Definition: localkey.hh:21
@ intersectionCodim
Codimension returned by LocalKey::codim() for degrees of freedom attached to an intersection.
Definition: localkey.hh:34
!
Definition: mimeticall.hh:86
const Dune::LocalKey & localKey(std::size_t i) const
map index i to local key
Definition: mimeticall.hh:103
std::size_t size() const
number of coefficients
Definition: mimeticall.hh:100
A few common exception classes.
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Implements a vector constructed from a given type representing a field and a compile-time given size.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
Impl::variant_< T... > variant
Incomplete re-implementation of C++17's std::variant.
Definition: variant.hh:464
Dune namespace.
Definition: alignedallocator.hh:14
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:32
D DomainType
domain type
Definition: localbasis.hh:43
A unique label for each type of element that can occur in a grid.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Oct 13, 22:30, 2024)