5#ifndef DUNE_MIMETIC_ALL_HH
6#define DUNE_MIMETIC_ALL_HH
10#include <dune/common/exceptions.hh>
11#include <dune/common/fvector.hh>
12#include <dune/common/fmatrix.hh>
14#include <dune/geometry/type.hh>
16#include "../common/localbasis.hh"
17#include "../common/localkey.hh"
21 template<
class D,
class R,
int dim>
26 R,1,Dune::FieldVector<R,1>, Dune::FieldMatrix<R,1,dim> >
Traits;
36 unsigned int size ()
const {
return variant; }
41 std::vector<typename Traits::RangeType>& out)
const
43 DUNE_THROW(Dune::Exception,
"mimetic basis evaluation not available");
49 std::vector<typename Traits::JacobianType>& out)
const
51 DUNE_THROW(Dune::Exception,
"mimetic basis Jacobian evaluation not available");
55 void partial (
const std::array<unsigned int, dim>& ,
57 std::vector<typename Traits::RangeType>& )
const
59 DUNE_THROW(Dune::Exception,
"mimetic basis partial derivative evaluation not available");
65 DUNE_THROW(Dune::Exception,
"mimetic order evaluation not available");
78 template<
typename F,
typename C>
80 DUNE_THROW(Dune::Exception,
"mimetic local interpolation not available");
91 : variant(variant_), li(variant_)
93 for (
unsigned int i=0; i<variant; i++)
102 std::size_t
size ()
const {
return variant; }
110 unsigned int variant;
111 std::vector<Dune::LocalKey> li;
Definition bdfmcube.hh:18
Type traits for LocalBasisVirtualInterface.
Definition common/localbasis.hh:35
D DomainType
domain type
Definition common/localbasis.hh:43
Describe position of one degree of freedom.
Definition localkey.hh:24
@ intersectionCodim
Codimension returned by LocalKey::codim() for degrees of freedom attached to an intersection.
Definition localkey.hh:37
Definition mimeticall.hh:23
MimeticLocalBasis(unsigned int variant_)
Definition mimeticall.hh:28
Dune::LocalBasisTraits< D, dim, Dune::FieldVector< D, dim >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, dim > > Traits
Definition mimeticall.hh:26
MimeticLocalBasis()
Definition mimeticall.hh:32
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition mimeticall.hh:47
unsigned int order() const
Polynomial order of the shape functions.
Definition mimeticall.hh:63
void partial(const std::array< unsigned int, dim > &, const typename Traits::DomainType &, std::vector< typename Traits::RangeType > &) const
Evaluate partial derivatives of all shape functions.
Definition mimeticall.hh:55
unsigned int size() const
Definition mimeticall.hh:36
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition mimeticall.hh:39
Definition mimeticall.hh:74
void interpolate(const F &f, std::vector< C > &out) const
Local interpolation of a function.
Definition mimeticall.hh:79
!
Definition mimeticall.hh:88
const Dune::LocalKey & localKey(std::size_t i) const
map index i to local key
Definition mimeticall.hh:105
MimeticLocalCoefficients(unsigned int variant_)
Definition mimeticall.hh:90
std::size_t size() const
number of coefficients
Definition mimeticall.hh:102
MimeticLocalCoefficients()
Definition mimeticall.hh:97