dune-localfunctions  2.9.0
mimeticall.hh
Go to the documentation of this file.
1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 // SPDX-FileCopyrightInfo: Copyright (C) DUNE Project contributors, see file LICENSE.md in module root
4 // SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
5 #ifndef DUNE_MIMETIC_ALL_HH
6 #define DUNE_MIMETIC_ALL_HH
7 
8 #include <cstddef>
9 
10 #include <dune/common/exceptions.hh>
11 #include <dune/common/fvector.hh>
12 #include <dune/common/fmatrix.hh>
13 
14 #include <dune/geometry/type.hh>
15 
16 #include "../common/localbasis.hh"
17 #include "../common/localkey.hh"
18 
19 namespace Dune
20 {
21  template<class D, class R, int dim>
23  {
24  public:
26  R,1,Dune::FieldVector<R,1>, Dune::FieldMatrix<R,1,dim> > Traits;
27 
28  MimeticLocalBasis (unsigned int variant_)
29  : variant(variant_)
30  {}
31 
33  : variant(0)
34  {}
35 
36  unsigned int size () const { return variant; }
37 
39  inline void evaluateFunction (
40  const typename Traits::DomainType& in,
41  std::vector<typename Traits::RangeType>& out) const
42  {
43  DUNE_THROW(Dune::Exception,"mimetic basis evaluation not available");
44  }
45 
47  inline void evaluateJacobian (
48  const typename Traits::DomainType& in,
49  std::vector<typename Traits::JacobianType>& out) const
50  {
51  DUNE_THROW(Dune::Exception,"mimetic basis Jacobian evaluation not available");
52  }
53 
55  void partial (const std::array<unsigned int, dim>& /*order*/,
56  const typename Traits::DomainType& /*in*/, // position
57  std::vector<typename Traits::RangeType>& /*out*/) const // return value
58  {
59  DUNE_THROW(Dune::Exception,"mimetic basis partial derivative evaluation not available");
60  }
61 
63  unsigned int order () const
64  {
65  DUNE_THROW(Dune::Exception,"mimetic order evaluation not available");
66  }
67 
68  private:
69  unsigned int variant;
70  };
71 
72  template<class LB>
74  {
75  public:
76 
78  template<typename F, typename C>
79  void interpolate (const F& f, std::vector<C>& out) const {
80  DUNE_THROW(Dune::Exception,"mimetic local interpolation not available");
81  }
82  };
83 
88  {
89  public:
90  MimeticLocalCoefficients (unsigned int variant_)
91  : variant(variant_), li(variant_)
92  {
93  for (unsigned int i=0; i<variant; i++)
95  }
96 
98  : variant(0), li(0)
99  {}
100 
102  std::size_t size () const { return variant; }
103 
105  const Dune::LocalKey& localKey (std::size_t i) const {
106  return li[i];
107  }
108 
109  private:
110  unsigned int variant;
111  std::vector<Dune::LocalKey> li;
112  };
113 }
114 
115 #endif
Definition: bdfmcube.hh:18
Type traits for LocalBasisVirtualInterface.
Definition: common/localbasis.hh:34
D DomainType
domain type
Definition: common/localbasis.hh:42
Describe position of one degree of freedom.
Definition: localkey.hh:23
@ intersectionCodim
Codimension returned by LocalKey::codim() for degrees of freedom attached to an intersection.
Definition: localkey.hh:36
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
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
const Dune::LocalKey & localKey(std::size_t i) const
map index i to local key
Definition: mimeticall.hh:105