3 #ifndef DUNE_MIMETIC_ALL_HH
4 #define DUNE_MIMETIC_ALL_HH
8 #include <dune/common/exceptions.hh>
9 #include <dune/common/fvector.hh>
10 #include <dune/common/fmatrix.hh>
12 #include <dune/geometry/type.hh>
14 #include "../common/localbasis.hh"
15 #include "../common/localkey.hh"
19 template<
class D,
class R,
int dim>
24 R,1,Dune::FieldVector<R,1>, Dune::FieldMatrix<R,1,dim> >
Traits;
34 unsigned int size ()
const {
return variant; }
39 std::vector<typename Traits::RangeType>& out)
const
41 DUNE_THROW(Dune::Exception,
"mimetic basis evaluation not available");
47 std::vector<typename Traits::JacobianType>& out)
const
49 DUNE_THROW(Dune::Exception,
"mimetic basis Jacobian evaluation not available");
55 DUNE_THROW(Dune::Exception,
"mimetic order evaluation not available");
68 template<
typename F,
typename C>
70 DUNE_THROW(Dune::Exception,
"mimetic local interpolation not available");
81 : variant(variant_), li(variant_)
83 for (
unsigned int i=0; i<variant; i++)
92 std::size_t
size ()
const {
return variant; }
100 unsigned int variant;
101 std::vector<Dune::LocalKey> li;
MimeticLocalCoefficients()
Definition: mimeticall.hh:87
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: mimeticall.hh:45
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:14
Definition: mimeticall.hh:20
unsigned int size() const
Definition: mimeticall.hh:34
Definition: mimeticall.hh:63
MimeticLocalBasis(unsigned int variant_)
Definition: mimeticall.hh:26
Dune::LocalBasisTraits< D, dim, Dune::FieldVector< D, dim >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, dim > > Traits
Definition: mimeticall.hh:24
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: mimeticall.hh:37
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:37
MimeticLocalBasis()
Definition: mimeticall.hh:30
const Dune::LocalKey & localKey(std::size_t i) const
map index i to local key
Definition: mimeticall.hh:95
unsigned int order() const
Polynomial order of the shape functions.
Definition: mimeticall.hh:53
MimeticLocalCoefficients(unsigned int variant_)
Definition: mimeticall.hh:80
D DomainType
domain type
Definition: localbasis.hh:49
std::size_t size() const
number of coefficients
Definition: mimeticall.hh:92
Codimension returned by LocalKey::codim() for degrees of freedom attached to an intersection.
Definition: localkey.hh:35
!
Definition: mimeticall.hh:77
Describe position of one degree of freedom.
Definition: localkey.hh:21
void interpolate(const F &f, std::vector< C > &out) const
Local interpolation of a function.
Definition: mimeticall.hh:69