DUNE PDELab (2.8)

monomiallocalbasis.hh
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_LOCALFUNCTIONS_MONOMIAL_MONOMIALLOCALBASIS_HH
4#define DUNE_LOCALFUNCTIONS_MONOMIAL_MONOMIALLOCALBASIS_HH
5
6#include <array>
7#include <cassert>
8#include <numeric>
9
11#include <dune/common/math.hh>
12
13#include "../common/localbasis.hh"
14
15namespace Dune
16{
17 namespace MonomImp
18 {
20 template <typename Traits>
21 class EvalAccess {
22 std::vector<typename Traits::RangeType> &out;
23#ifndef NDEBUG
24 unsigned int first_unused_index;
25#endif
26
27 public:
28 EvalAccess(std::vector<typename Traits::RangeType> &out_)
29 : out(out_)
30#ifndef NDEBUG
31 , first_unused_index(0)
32#endif
33 { }
34#ifndef NDEBUG
35 ~EvalAccess() {
36 assert(first_unused_index == out.size());
37 }
38#endif
39 typename Traits::RangeFieldType &operator[](unsigned int index)
40 {
41 assert(index < out.size());
42#ifndef NDEBUG
43 if(first_unused_index <= index)
44 first_unused_index = index+1;
45#endif
46 return out[index][0];
47 }
48 };
49
51 template <typename Traits>
53 std::vector<typename Traits::JacobianType> &out;
54 unsigned int row;
55#ifndef NDEBUG
56 unsigned int first_unused_index;
57#endif
58
59 public:
60 JacobianAccess(std::vector<typename Traits::JacobianType> &out_,
61 unsigned int row_)
62 : out(out_), row(row_)
63#ifndef NDEBUG
64 , first_unused_index(0)
65#endif
66 { }
67#ifndef NDEBUG
69 assert(first_unused_index == out.size());
70 }
71#endif
72 typename Traits::RangeFieldType &operator[](unsigned int index)
73 {
74 assert(index < out.size());
75#ifndef NDEBUG
76 if(first_unused_index <= index)
77 first_unused_index = index+1;
78#endif
79 return out[index][0][row];
80 }
81 };
82
95 template <typename Traits, int c>
96 struct Evaluate
97 {
98 enum {
100 d = Traits::dimDomain - c
101 };
108 template <typename Access>
109 static void eval (
110 const typename Traits::DomainType &in,
113 const std::array<unsigned int, Traits::dimDomain> &derivatives,
116 typename Traits::RangeFieldType prod,
118 int bound,
120 int& index,
122 Access &access)
123 {
124 // start with the highest exponent for this dimension, then work down
125 for (int e = bound; e >= 0; --e)
126 {
127 // the rest of the available exponents, to be used by the other
128 // dimensions
129 int newbound = bound - e;
130 if(e < (int)derivatives[d])
132 eval(in, derivatives, 0, newbound, index, access);
133 else {
134 int coeff = 1;
135 for(int i = e - derivatives[d] + 1; i <= e; ++i)
136 coeff *= i;
137 // call the evaluator for the next dimension
139 eval( // pass the coordinate and the derivatives unchanged
140 in, derivatives,
141 // also pass the product accumulated so far, but also
142 // include the current dimension
143 prod * power(in[d], e-derivatives[d]) * coeff,
144 // pass the number of remaining exponents to the next
145 // dimension
146 newbound,
147 // pass the next index to fill and the output access
148 // wrapper
149 index, access);
150 }
151 }
152 }
153 };
154
159 template <typename Traits>
160 struct Evaluate<Traits, 1>
161 {
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)
169 {
170 if(bound < (int)derivatives[d])
171 prod = 0;
172 else {
173 int coeff = 1;
174 for(int i = bound - derivatives[d] + 1; i <= bound; ++i)
175 coeff *= i;
176 prod *= power(in[d], bound-derivatives[d]) * coeff;
177 }
178 access[index] = prod;
179 ++index;
180 }
181 };
182
183 } //namespace MonomImp
184
198 template<class D, class R, unsigned int d, unsigned int p>
200 {
201 // Helper: Number of shape functions for a k-th order element in dimension dd
202 static constexpr unsigned int size (int dd, int k)
203 {
204 if (dd==0 || k==0)
205 return 1;
206 return size(dd,k-1) + size(dd-1,k);
207 }
208
209 public:
213
215 static constexpr unsigned int size ()
216 {
217 return size(d,p);
218 }
219
221 inline void evaluateFunction (const typename Traits::DomainType& in,
222 std::vector<typename Traits::RangeType>& out) const
223 {
224 out.resize(size());
225 int index = 0;
226 std::array<unsigned int, d> derivatives;
227 std::fill(derivatives.begin(), derivatives.end(), 0);
229 for (unsigned int lp = 0; lp <= p; ++lp)
230 MonomImp::Evaluate<Traits, d>::eval(in, derivatives, 1, lp, index, access);
231 }
232
238 inline void partial(const std::array<unsigned int,d>& order,
239 const typename Traits::DomainType& in,
240 std::vector<typename Traits::RangeType>& out) const
241 {
242 out.resize(size());
243 int index = 0;
245 for (unsigned int lp = 0; lp <= p; ++lp)
246 MonomImp::Evaluate<Traits, d>::eval(in, order, 1, lp, index, access);
247 }
248
250 inline void
251 evaluateJacobian (const typename Traits::DomainType& in, // position
252 std::vector<typename Traits::JacobianType>& out) const // return value
253 {
254 out.resize(size());
255 std::array<unsigned int, d> derivatives;
256 for(unsigned int i = 0; i < d; ++i)
257 derivatives[i] = 0;
258 for(unsigned int i = 0; i < d; ++i)
259 {
260 derivatives[i] = 1;
261 int index = 0;
263 for(unsigned int lp = 0; lp <= p; ++lp)
264 MonomImp::Evaluate<Traits, d>::eval(in, derivatives, 1, lp, index, access);
265 derivatives[i] = 0;
266 }
267 }
268
270 unsigned int order () const
271 {
272 return p;
273 }
274 };
275
276}
277
278#endif // DUNE_LOCALFUNCTIONS_MONOMIAL_MONOMIALLOCALBASIS_HH
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
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 21, 23:30, 2024)