5#ifndef DUNE_LOCALFUNCTIONS_MONOMIAL_MONOMIALLOCALBASIS_HH
6#define DUNE_LOCALFUNCTIONS_MONOMIAL_MONOMIALLOCALBASIS_HH
15#include "../common/localbasis.hh"
22 template <
typename Traits>
24 std::vector<typename Traits::RangeType> &out;
26 unsigned int first_unused_index;
30 EvalAccess(std::vector<typename Traits::RangeType> &out_)
33 , first_unused_index(0)
38 assert(first_unused_index == out.size());
41 typename Traits::RangeFieldType &operator[](
unsigned int index)
43 assert(index < out.size());
45 if(first_unused_index <= index)
46 first_unused_index = index+1;
53 template <
typename Traits>
55 std::vector<typename Traits::JacobianType> &out;
58 unsigned int first_unused_index;
64 : out(out_), row(row_)
66 , first_unused_index(0)
71 assert(first_unused_index == out.size());
74 typename Traits::RangeFieldType &operator[](
unsigned int index)
76 assert(index < out.size());
78 if(first_unused_index <= index)
79 first_unused_index = index+1;
81 return out[index][0][row];
97 template <
typename Traits,
int c>
101 constexpr static int d = Traits::dimDomain - c;
109 template <
typename Access>
111 const typename Traits::DomainType &in,
114 const std::array<unsigned int, Traits::dimDomain> &derivatives,
117 typename Traits::RangeFieldType prod,
126 for (
int e = bound; e >= 0; --e)
130 int newbound = bound - e;
131 if(e < (
int)derivatives[
d])
133 eval(in, derivatives, 0, newbound, index, access);
136 for(
int i = e - derivatives[
d] + 1; i <= e; ++i)
144 prod *
power(in[
d], e-derivatives[
d]) * coeff,
160 template <
typename Traits>
163 constexpr static int d = Traits::dimDomain-1;
165 template <
typename Access>
166 static void eval (
const typename Traits::DomainType &in,
167 const std::array<unsigned int, Traits::dimDomain> &derivatives,
168 typename Traits::RangeFieldType prod,
169 int bound,
int& index, Access &access)
171 if(bound < (
int)derivatives[
d])
175 for(
int i = bound - derivatives[
d] + 1; i <= bound; ++i)
177 prod *=
power(in[
d], bound-derivatives[
d]) * coeff;
179 access[index] = prod;
199 template<
class D,
class R,
unsigned int d,
unsigned int p>
203 static constexpr unsigned int size (
int dd,
int k)
216 static constexpr unsigned int size ()
223 std::vector<typename Traits::RangeType>& out)
const
227 std::array<unsigned int, d> derivatives;
228 std::fill(derivatives.begin(), derivatives.end(), 0);
230 for (
unsigned int lp = 0; lp <= p; ++lp)
241 std::vector<typename Traits::RangeType>& out)
const
246 for (
unsigned int lp = 0; lp <= p; ++lp)
253 std::vector<typename Traits::JacobianType>& out)
const
256 std::array<unsigned int, d> derivatives;
257 for(
unsigned int i = 0; i < d; ++i)
259 for(
unsigned int i = 0; i < d; ++i)
264 for(
unsigned int lp = 0; lp <= p; ++lp)
A dense n x m matrix.
Definition: fmatrix.hh:117
vector space out of a tensor product of fields.
Definition: fvector.hh:95
Access output vector of evaluateFunction() and evaluate()
Definition: monomiallocalbasis.hh:23
Access output vector of evaluateJacobian()
Definition: monomiallocalbasis.hh:54
Constant shape function.
Definition: monomiallocalbasis.hh:201
unsigned int order() const
Polynomial order of the shape functions.
Definition: monomiallocalbasis.hh:271
void partial(const std::array< unsigned int, d > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate partial derivatives of any order of all shape functions.
Definition: monomiallocalbasis.hh:239
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: monomiallocalbasis.hh:252
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: monomiallocalbasis.hh:222
LocalBasisTraits< D, d, Dune::FieldVector< D, d >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, d > > Traits
export type traits for function signature
Definition: monomiallocalbasis.hh:213
static constexpr unsigned int size()
Number of shape functions.
Definition: monomiallocalbasis.hh:216
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Some useful basic math stuff.
Dune namespace.
Definition: alignedallocator.hh:13
constexpr Base power(Base m, Exponent p)
Power method for integer exponents.
Definition: math.hh:75
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:34
D DomainType
domain type
Definition: localbasis.hh:42
static void eval(const typename Traits::DomainType &in, const std::array< unsigned int, Traits::dimDomain > &derivatives, typename Traits::RangeFieldType prod, int bound, int &index, Access &access)
Definition: monomiallocalbasis.hh:166
Definition: monomiallocalbasis.hh:99
static constexpr int d
The next dimension to try for factors.
Definition: monomiallocalbasis.hh:101
static void eval(const typename Traits::DomainType &in, const std::array< unsigned int, Traits::dimDomain > &derivatives, typename Traits::RangeFieldType prod, int bound, int &index, Access &access)
Definition: monomiallocalbasis.hh:110