3#ifndef DUNE_LOCALFUNCTIONS_MONOMIAL_MONOMIALLOCALBASIS_HH
4#define DUNE_LOCALFUNCTIONS_MONOMIAL_MONOMIALLOCALBASIS_HH
13#include "../common/localbasis.hh"
20 template <
typename Traits>
22 std::vector<typename Traits::RangeType> &out;
24 unsigned int first_unused_index;
28 EvalAccess(std::vector<typename Traits::RangeType> &out_)
31 , first_unused_index(0)
36 assert(first_unused_index == out.size());
39 typename Traits::RangeFieldType &operator[](
unsigned int index)
41 assert(index < out.size());
43 if(first_unused_index <= index)
44 first_unused_index = index+1;
51 template <
typename Traits>
53 std::vector<typename Traits::JacobianType> &out;
56 unsigned int first_unused_index;
62 : out(out_), row(row_)
64 , first_unused_index(0)
69 assert(first_unused_index == out.size());
72 typename Traits::RangeFieldType &operator[](
unsigned int index)
74 assert(index < out.size());
76 if(first_unused_index <= index)
77 first_unused_index = index+1;
79 return out[index][0][row];
95 template <
typename Traits,
int c>
100 d = Traits::dimDomain - c
108 template <
typename Access>
110 const typename Traits::DomainType &in,
113 const std::array<unsigned int, Traits::dimDomain> &derivatives,
116 typename Traits::RangeFieldType prod,
125 for (
int e = bound; e >= 0; --e)
129 int newbound = bound - e;
130 if(e < (
int)derivatives[
d])
132 eval(in, derivatives, 0, newbound, index, access);
135 for(
int i = e - derivatives[
d] + 1; i <= e; ++i)
143 prod *
power(in[
d], e-derivatives[
d]) * coeff,
159 template <
typename Traits>
162 enum {
d = Traits::dimDomain-1 };
164 template <
typename Access>
165 static void eval (
const typename Traits::DomainType &in,
166 const std::array<unsigned int, Traits::dimDomain> &derivatives,
167 typename Traits::RangeFieldType prod,
168 int bound,
int& index, Access &access)
170 if(bound < (
int)derivatives[
d])
174 for(
int i = bound - derivatives[
d] + 1; i <= bound; ++i)
176 prod *=
power(in[
d], bound-derivatives[
d]) * coeff;
178 access[index] = prod;
198 template<
class D,
class R,
unsigned int d,
unsigned int p>
202 static constexpr unsigned int size (
int dd,
int k)
215 static constexpr unsigned int size ()
222 std::vector<typename Traits::RangeType>& out)
const
226 std::array<unsigned int, d> derivatives;
227 std::fill(derivatives.begin(), derivatives.end(), 0);
229 for (
unsigned int lp = 0; lp <= p; ++lp)
240 std::vector<typename Traits::RangeType>& out)
const
245 for (
unsigned int lp = 0; lp <= p; ++lp)
252 std::vector<typename Traits::JacobianType>& out)
const
255 std::array<unsigned int, d> derivatives;
256 for(
unsigned int i = 0; i < d; ++i)
258 for(
unsigned int i = 0; i < d; ++i)
263 for(
unsigned int lp = 0; lp <= p; ++lp)
A dense n x m matrix.
Definition: fmatrix.hh:69
vector space out of a tensor product of fields.
Definition: fvector.hh:95
Access output vector of evaluateFunction() and evaluate()
Definition: monomiallocalbasis.hh:21
Access output vector of evaluateJacobian()
Definition: monomiallocalbasis.hh:52
Constant shape function.
Definition: monomiallocalbasis.hh:200
unsigned int order() const
Polynomial order of the shape functions.
Definition: monomiallocalbasis.hh:270
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:238
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: monomiallocalbasis.hh:251
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: monomiallocalbasis.hh:221
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:212
static constexpr unsigned int size()
Number of shape functions.
Definition: monomiallocalbasis.hh:215
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:11
constexpr Mantissa power(Mantissa m, Exponent p)
Power method for integer exponents.
Definition: math.hh:73
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:32
D DomainType
domain type
Definition: localbasis.hh:43
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:165
Definition: monomiallocalbasis.hh:97
@ d
The next dimension to try for factors.
Definition: monomiallocalbasis.hh:100
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:109