3 #ifndef DUNE_LOCALFUNCTIONS_MONOMIAL_MONOMIALLOCALBASIS_HH 4 #define DUNE_LOCALFUNCTIONS_MONOMIAL_MONOMIALLOCALBASIS_HH 10 #include <dune/common/fmatrix.hh> 12 #include "../common/localbasis.hh" 20 template<
int d,
int k>
52 template <
typename Traits>
54 std::vector<typename Traits::RangeType> &out;
56 unsigned int first_unused_index;
60 EvalAccess(std::vector<typename Traits::RangeType> &out_)
63 , first_unused_index(0)
68 assert(first_unused_index == out.size());
71 typename Traits::RangeFieldType &
operator[](
unsigned int index)
73 assert(index < out.size());
75 if(first_unused_index <= index)
76 first_unused_index = index+1;
83 template <
typename Traits>
85 std::vector<typename Traits::JacobianType> &out;
88 unsigned int first_unused_index;
94 : out(out_), row(row_)
96 , first_unused_index(0)
101 assert(first_unused_index == out.size());
104 typename Traits::RangeFieldType &
operator[](
unsigned int index)
106 assert(index < out.size());
108 if(first_unused_index <= index)
109 first_unused_index = index+1;
111 return out[index][0][row];
127 template <
typename Traits,
int c>
132 d = Traits::dimDomain - c
140 template <
typename Access>
142 const typename Traits::DomainType &in,
145 const std::array<int, Traits::dimDomain> &derivatives,
148 typename Traits::RangeFieldType prod,
157 for (
int e = bound; e >= 0; --e)
161 int newbound = bound - e;
162 if(e < derivatives[d])
164 eval(in, derivatives, 0, newbound, index, access);
167 for(
int i = e - derivatives[d] + 1; i <= e; ++i)
175 prod *
ipow(in[d], e-derivatives[d]) * coeff,
191 template <
typename Traits>
194 enum { d = Traits::dimDomain-1 };
196 template <
typename Access>
197 static void eval (
const typename Traits::DomainType &in,
198 const std::array<int, Traits::dimDomain> &derivatives,
199 typename Traits::RangeFieldType prod,
200 int bound,
int& index, Access &access)
202 if(bound < derivatives[d])
206 for(
int i = bound - derivatives[d] + 1; i <= bound; ++i)
208 prod *=
ipow(in[d], bound-derivatives[d]) * coeff;
210 access[index] = prod;
230 template<
class D,
class R,
unsigned int d,
unsigned int p>
246 std::vector<typename Traits::RangeType>& out)
const 250 std::array<int, d> derivatives;
251 std::fill(derivatives.begin(), derivatives.end(), 0);
253 for (
unsigned int lp = 0; lp <= p; ++lp)
262 inline void partial(
const std::array<unsigned int,d>& order,
264 std::vector<typename Traits::RangeType>& out)
const 268 std::array<int, d> derivatives;
269 std::copy(order.begin(), order.end(), derivatives.begin());
271 for (
unsigned int lp = 0; lp <= p; ++lp)
278 std::vector<typename Traits::JacobianType>& out)
const 281 std::array<int, d> derivatives;
282 for(
unsigned int i = 0; i < d; ++i)
284 for(
unsigned int i = 0; i < d; ++i)
289 for(
unsigned int lp = 0; lp <= p; ++lp)
304 #endif // DUNE_LOCALFUNCTIONS_MONOMIAL_MONOMIALLOCALBASIS_HH Definition: monomiallocalbasis.hh:128
Traits::RangeFieldType & operator[](unsigned int index)
Definition: monomiallocalbasis.hh:71
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:262
Definition: monomiallocalbasis.hh:21
~EvalAccess()
Definition: monomiallocalbasis.hh:67
T ipow(T base, int exp)
Definition: monomiallocalbasis.hh:38
Constant shape function.
Definition: monomiallocalbasis.hh:231
Access output vector of evaluateJacobian()
Definition: monomiallocalbasis.hh:84
static void eval(const typename Traits::DomainType &in, const std::array< int, Traits::dimDomain > &derivatives, typename Traits::RangeFieldType prod, int bound, int &index, Access &access)
Definition: monomiallocalbasis.hh:197
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:31
JacobianAccess(std::vector< typename Traits::JacobianType > &out_, unsigned int row_)
Definition: monomiallocalbasis.hh:92
D DomainType
domain type
Definition: localbasis.hh:43
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:236
static void eval(const typename Traits::DomainType &in, const std::array< int, Traits::dimDomain > &derivatives, typename Traits::RangeFieldType prod, int bound, int &index, Access &access)
Definition: monomiallocalbasis.hh:141
Traits::RangeFieldType & operator[](unsigned int index)
Definition: monomiallocalbasis.hh:104
unsigned int order() const
Polynomial order of the shape functions.
Definition: monomiallocalbasis.hh:296
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:15
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: monomiallocalbasis.hh:277
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: monomiallocalbasis.hh:245
Access output vector of evaluateFunction() and evaluate()
Definition: monomiallocalbasis.hh:53
EvalAccess(std::vector< typename Traits::RangeType > &out_)
Definition: monomiallocalbasis.hh:60
unsigned int size() const
number of shape functions
Definition: monomiallocalbasis.hh:239
Definition: monomiallocalbasis.hh:22
~JacobianAccess()
Definition: monomiallocalbasis.hh:100