dune-localfunctions 2.8.0
refinedp0localbasis.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#ifndef DUNE_REFINED_P0_LOCALBASIS_HH
4#define DUNE_REFINED_P0_LOCALBASIS_HH
5
6#include <numeric>
7
8#include <dune/common/fvector.hh>
9#include <dune/common/fmatrix.hh>
10
13
14namespace Dune
15{
16
35 template<class D, class R, int dim>
37 : public RefinedSimplexLocalBasis<D,dim>
38 {
39 // 2 to the k-th power
40 enum {N = 1<<dim};
41 public:
43 typedef LocalBasisTraits<D,dim,Dune::FieldVector<D,dim>,R,1,Dune::FieldVector<R,1>, Dune::FieldMatrix<R,1,dim> > Traits;
44
46 unsigned int size () const
47 {
48 return N;
49 }
50
52 inline void evaluateFunction (const typename Traits::DomainType& in,
53 std::vector<typename Traits::RangeType>& out) const
54 {
55 int subElement = this->getSubElement(in);
56 out.resize(N);
57 for(int i=0; i<N; ++i)
58 out[i] = (i==subElement) ? 1 : 0;
59 }
60
61 inline void
62 evaluateJacobian (const typename Traits::DomainType& in, // position
63 std::vector<typename Traits::JacobianType>& out) const // return value
64 {
65 out.resize(N);
66 for(int i=0; i<N; ++i)
67 out[i][0] = 0;
68 }
69
71 void partial (const std::array<unsigned int, dim>& order,
72 const typename Traits::DomainType& in, // position
73 std::vector<typename Traits::RangeType>& out) const // return value
74 {
75 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
76 if (totalOrder == 0) {
77 evaluateFunction(in, out);
78 } else {
79 out.resize(size());
80 for (std::size_t i = 0; i < size(); ++i)
81 out[i] = 0;
82 }
83 }
84
89 unsigned int order () const
90 {
91 return 0;
92 }
93
94 };
95
96}
97#endif
Contains a base class for LocalBasis classes based on uniform refinement.
Definition: bdfmcube.hh:16
Type traits for LocalBasisVirtualInterface.
Definition: common/localbasis.hh:32
D DomainType
domain type
Definition: common/localbasis.hh:43
Definition: refinedsimplexlocalbasis.hh:18
Uniformly refined constant shape functions on a unit simplex in R^dim.
Definition: refinedp0localbasis.hh:38
void partial(const std::array< unsigned int, dim > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate partial derivatives of all shape functions.
Definition: refinedp0localbasis.hh:71
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: refinedp0localbasis.hh:52
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Definition: refinedp0localbasis.hh:62
unsigned int order() const
Polynomial order of the shape functions.
Definition: refinedp0localbasis.hh:89
LocalBasisTraits< D, dim, Dune::FieldVector< D, dim >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, dim > > Traits
export type traits for function signature
Definition: refinedp0localbasis.hh:43
unsigned int size() const
number of shape functions
Definition: refinedp0localbasis.hh:46