DUNE-FEM (unstable)

interpolation.hh
1#ifndef DUNE_FEM_SPACE_LOCALFINITEELEMENT_INTERPOLATION_HH
2#define DUNE_FEM_SPACE_LOCALFINITEELEMENT_INTERPOLATION_HH
3
4#include <cstddef>
5
6#include <utility>
7#include <vector>
8#include <optional>
9
11
12#include <dune/fem/misc/threads/threadsafevalue.hh>
13
14#include <dune/fem/function/common/localcontribution.hh>
15#include <dune/fem/space/basisfunctionset/transformed.hh>
16#include <dune/fem/space/basisfunctionset/vectorial.hh>
17#include <dune/fem/space/combinedspace/interpolation.hh>
18#include <dune/fem/space/common/localinterpolation.hh>
19
20namespace Dune
21{
22
23 namespace Fem
24 {
25
26 namespace Impl
27 {
28 // LocalFunctionWrapper
29 // --------------------
30 template< class LocalFunction, class BasisFunctionSet >
31 struct LocalFunctionWrapper
32 {
33 struct Traits
34 {
35 typedef typename LocalFunction::DomainType DomainType;
36 typedef typename LocalFunction::RangeType RangeType;
37 };
38 typedef typename LocalFunction::DomainType DomainType;
39 typedef typename LocalFunction::RangeType RangeType;
40
41 LocalFunctionWrapper ( const LocalFunction &lf, const BasisFunctionSet &bset ) : lf_( lf ) {}
42
43 template< class Arg >
44 void evaluate ( const Arg &x, typename Traits::RangeType &y ) const
45 {
46 lf_.evaluate( x, y );
47 }
48 template< class Arg >
49 typename Traits::RangeType operator()(const Arg &x) const
50 {
51 typename Traits::RangeType y;
52 evaluate(x,y);
53 return y;
54 }
55
56 private:
57 const LocalFunction &lf_;
58 };
59 // LocalFunctionWrapper
60 // --------------------
61 template< class LocalFunction, class Entity, class ShapeFunctionSet, class Transformation >
62 struct LocalFunctionWrapper< LocalFunction, TransformedBasisFunctionSet< Entity, ShapeFunctionSet, Transformation > >
63 {
64 typedef TransformedBasisFunctionSet< Entity, ShapeFunctionSet, Transformation > BasisFunctionSetType;
65
66 struct Traits
67 {
68 typedef typename LocalFunction::DomainType DomainType;
69 typedef typename LocalFunction::RangeType RangeType;
70 };
71 typedef typename LocalFunction::DomainType DomainType;
72 typedef typename LocalFunction::RangeType RangeType;
73
74 LocalFunctionWrapper ( const LocalFunction &lf, const BasisFunctionSetType &bset )
75 : lf_( lf ), bset_( bset )
76 {
77 bset_.geometry();
78 }
79
80 template< class Arg >
81 void evaluate ( const Arg &x, typename Traits::RangeType &y ) const
82 {
83 typename Traits::RangeType help;
84 lf_.evaluate( x, help );
85 // geometry needs to come from local function. For some reasons the
86 // geometry from the basis set was invalid. Something to investigate.
87 typename Transformation::InverseTransformationType transf( lf_.geometry(), x );
88 y = transf.apply( help );
89 }
90 template< class Arg >
91 typename Traits::RangeType operator() ( const Arg &x ) const
92 {
93 typename Traits::RangeType y;
94 evaluate(x,y);
95 return y;
96 }
97
98 private:
99 const LocalFunction &lf_;
100 const BasisFunctionSetType &bset_;
101 };
102
103 } // namespace Impl
104
105
106 // LocalFiniteElementInterpolation
107 // -------------------------------
108
109 template< class Space, class LocalInterpolation, bool scalarSFS >
110 class LocalFiniteElementInterpolation;
111
112 // a vector valued shape basis function set taken from
113 // dune-localfunction - the vector value 'blow up' is not yet supported
114 // for this case
115 template< class Space, class LocalInterpolation >
116 class LocalFiniteElementInterpolation<Space,LocalInterpolation,false>
117 {
118 typedef LocalFiniteElementInterpolation< Space, LocalInterpolation,false > ThisType;
119
120 public:
121 typedef typename Space::BasisFunctionSetType BasisFunctionSetType;
122 typedef LocalInterpolation LocalInterpolationType;
123
124 private:
125 typedef typename BasisFunctionSetType::FunctionSpaceType FunctionSpaceType;
126
127 typedef typename FunctionSpaceType::RangeType RangeType;
128 typedef typename FunctionSpaceType::RangeFieldType RangeFieldType;
129 static const int dimRange = FunctionSpaceType::dimRange;
130
131 typedef std::size_t size_type;
132
133 template< class LocalFunction >
134 using LocalFunctionWrapper = Impl::LocalFunctionWrapper< LocalFunction, BasisFunctionSetType >;
135
136 public:
137 LocalFiniteElementInterpolation()
138 : basisFunctionSet_()
139 , localInterpolation_()
140 {
141 }
142
143 explicit LocalFiniteElementInterpolation ( const BasisFunctionSetType &basisFunctionSet,
144 const LocalInterpolationType &localInterpolation = LocalInterpolationType() )
145 : basisFunctionSet_( basisFunctionSet ),
146 localInterpolation_( localInterpolation )
147 {}
148
149 ThisType& operator=( const ThisType& other )
150 {
151 basisFunctionSet_ = other.basisFunctionSet_;
152 localInterpolation_.emplace( other.localInterpolation() );
153 return *this;
154 }
155
156 void unbind()
157 {
158 // basisFunctionSet_ = BasisFunctionSetType();
159 }
160
161 template< class LocalFunction, class Dof >
162 void operator() ( const LocalFunction &localFunction, std::vector< Dof > &localDofVector ) const
163 {
164 LocalFunctionWrapper< LocalFunction > wrapper( localFunction, basisFunctionSet() );
165 localInterpolation().interpolate( wrapper, localDofVector );
166 }
167
168 template< class LocalFunction, class Dof, class Allocator >
169 void operator() ( const LocalFunction &localFunction, Dune::DynamicVector< Dof, Allocator > &localDofVector ) const
170 {
171 (*this)(localFunction,localDofVector.container() );
172 }
173
174 template< class LocalFunction, class DiscreteFunction, template< class > class Assembly >
175 void operator() ( const LocalFunction &localFunction, LocalContribution< DiscreteFunction, Assembly > &localContribution ) const
176 {
177 (*this)(localFunction,localContribution.localDofVector());
178 }
179
180 BasisFunctionSetType basisFunctionSet () const { return basisFunctionSet_; }
181 const LocalInterpolationType &localInterpolation () const
182 {
183 assert( localInterpolation_.has_value() );
184 return *localInterpolation_;
185 }
186
187 private:
188 BasisFunctionSetType basisFunctionSet_;
189 std::optional< LocalInterpolationType > localInterpolation_;
190 };
191
192 namespace Impl
193 {
194 template <int dimRange>
195 struct RangeConverter
196 {
197 RangeConverter ( std::size_t range ) : range_( range ) {}
198
199 template< class T >
200 FieldVector< T, 1 > operator() ( const FieldVector< T, dimRange > &in ) const
201 {
202 return in[ range_ ];
203 }
204
205 template< class T, int j >
206 FieldMatrix< T, 1, j > operator() ( const FieldMatrix< T, dimRange, j > &in ) const
207 {
208 // return in[ range_ ]; // implicit conversion fails
209 FieldMatrix<T,1,j> ret;
210 ret[0] = in[range_];
211 return ret;
212 }
213
214 protected:
215 std::size_t range_;
216 };
217 template <class DofVector, class DofAlignment>
218 struct SubDofVectorWrapper
219 : public SubDofVector< DofVector, DofAlignment >
220 {
221 typedef SubDofVector< DofVector, DofAlignment > BaseType;
222
223 SubDofVectorWrapper( DofVector& dofs, int coordinate, const DofAlignment &dofAlignment )
224 : BaseType( dofs, coordinate, dofAlignment )
225 {}
226
228 void clear() {}
229 void resize( const unsigned int) {}
230 };
231
232 }
233
234 // a scalar valued shape basis function set taken from
235 // dune-localfunction - for this the vector value 'blow up' is supported
236 // as with other spaces
237 template< class Space, class LocalInterpolation >
238 class LocalFiniteElementInterpolation<Space,LocalInterpolation,true>
239 {
240 typedef LocalFiniteElementInterpolation< Space, LocalInterpolation,true > ThisType;
241
242 public:
243 typedef typename Space::BasisFunctionSetType BasisFunctionSetType;
244 typedef LocalInterpolation LocalInterpolationType;
245
246 private:
247 typedef typename BasisFunctionSetType::FunctionSpaceType FunctionSpaceType;
248 // typedef typename Space::FunctionSpaceType FunctionSpaceType;
249
250 typedef typename FunctionSpaceType::RangeType RangeType;
251 typedef typename FunctionSpaceType::RangeFieldType RangeFieldType;
252 static const int dimRange = FunctionSpaceType::dimRange;
253 static const int dimR = Space::FunctionSpaceType::dimRange;
254
255 typedef std::size_t size_type;
256
257 typedef VerticalDofAlignment< BasisFunctionSetType, RangeType> DofAlignmentType;
258
259 public:
260 LocalFiniteElementInterpolation ()
261 : basisFunctionSet_(),
262 localInterpolation_(),
263 dofAlignment_()
264 {}
265
266 explicit LocalFiniteElementInterpolation ( const BasisFunctionSetType &basisFunctionSet,
267 const LocalInterpolationType &localInterpolation = LocalInterpolationType() )
268 : basisFunctionSet_( basisFunctionSet ),
269 localInterpolation_( localInterpolation ),
270 dofAlignment_( basisFunctionSet_ )
271 {}
272
273 LocalFiniteElementInterpolation ( const ThisType& other )
274 : basisFunctionSet_( other.basisFunctionSet_ ),
275 localInterpolation_( other.localInterpolation_ ),
276 dofAlignment_( basisFunctionSet_ )
277 {}
278
279 ThisType& operator=( const ThisType& other )
280 {
281 basisFunctionSet_ = other.basisFunctionSet_;
282 localInterpolation_.emplace( other.localInterpolation() );
283 return *this;
284 }
285
286 void unbind()
287 {
288 // basisFunctionSet_ = BasisFunctionSetType();
289 }
290
291 template< class LocalFunction, class Vector>
292 void operator() ( const LocalFunction &localFunction, Vector &localDofVector ) const
293 {
294 // clear dofs before something is adedd
295 // localDofVector.clear(); // does not exist on DynVector so use 'fill' instead
296 std::fill(localDofVector.begin(),localDofVector.end(),0);
297 for( std::size_t i = 0; i < dimR; ++i )
298 {
299 Impl::SubDofVectorWrapper< Vector, DofAlignmentType > subLdv( localDofVector, i, dofAlignment_ );
300 (*this)(
301 localFunctionConverter( localFunction, Impl::RangeConverter<dimR>(i) ),
302 subLdv, PriorityTag<42>()
303 );
304 }
305 }
306
307 template< class LocalFunction, class Dof, class Allocator >
308 void operator() ( const LocalFunction &localFunction, Dune::DynamicVector< Dof, Allocator > &localDofVector ) const
309 {
310 (*this)(localFunction,localDofVector.container() );
311 }
312
313 template< class LocalFunction, class DiscreteFunction, template< class > class Assembly >
314 void operator() ( const LocalFunction &localFunction, LocalContribution< DiscreteFunction, Assembly > &localContribution ) const
315 {
316 (*this)(localFunction,localContribution.localDofVector());
317 }
318
319 BasisFunctionSetType basisFunctionSet () const { return basisFunctionSet_; }
320 const LocalInterpolationType &localInterpolation () const
321 {
322 assert( localInterpolation_.has_value() );
323 return *localInterpolation_;
324 }
325
326 private:
327 template< class LocalFunction, class Vector>
328 auto operator() ( const LocalFunction &localFunction, Vector &localDofVector, PriorityTag<1> ) const
329 -> void_t< decltype(
330 std::declval<LocalInterpolationType>().interpolate(
331 localFunction, localDofVector)) >
332 {
333 localInterpolation().interpolate( localFunction, localDofVector );
334 }
335 template< class LocalFunction, class Vector>
336 void operator() ( const LocalFunction &localFunction, Vector &localDofVector, PriorityTag<0> ) const
337 {
338 std::vector<typename Vector::value_type> tmp(basisFunctionSet_.size());
339 localInterpolation().interpolate( localFunction, tmp);
340 for (unsigned int i=0;i<tmp.size();++i)
341 localDofVector[i] = tmp[i];
342 }
343 BasisFunctionSetType basisFunctionSet_;
344 std::optional< LocalInterpolationType > localInterpolation_;
345 DofAlignmentType dofAlignment_;
346 };
347
348 template <class DiscreteFunctionSpace >
349 class LocalFEInterpolationWrapper
350 : public LocalInterpolationWrapper< DiscreteFunctionSpace >
351 {
352 typedef LocalInterpolationWrapper< DiscreteFunctionSpace > BaseType;
353 protected:
354 using BaseType :: interpolation_;
355 typedef std::vector< typename DiscreteFunctionSpace::RangeFieldType >
356 TemporarayDofVectorType;
358
359 using BaseType::interpolation;
360 public:
361 LocalFEInterpolationWrapper( const DiscreteFunctionSpace& space )
362 : BaseType( space )
363 {}
364
365 using BaseType :: bind;
366 using BaseType :: unbind;
367
368 template< class LocalFunction, class Dof >
369 void operator() ( const LocalFunction &localFunction, std::vector< Dof > &localDofVector ) const
370 {
371 // interpolation only works for std::vector storage
372 interpolation()( localFunction, localDofVector );
373 }
374
375 template< class LocalFunction, class Dof, class Allocator >
376 void operator() ( const LocalFunction &localFunction, Dune::DynamicVector< Dof, Allocator > &localDofVector ) const
377 {
378 // DynamicVector internally stores a std::vector
379 (*this)(localFunction, localDofVector.container() );
380 }
381
382
383 template< class LocalFunction, class DiscreteFunction, template< class > class Assembly >
384 void operator() ( const LocalFunction &localFunction, LocalContribution< DiscreteFunction, Assembly > &localContribution ) const
385 {
386 // LocalContribution internally stores a std::vector
387 (*this)(localFunction, localContribution.localDofVector());
388 }
389
391 template< class LocalFunction, class LocalDofVector >
392 void operator () ( const LocalFunction &localFunction, LocalDofVector &dofs ) const
393 {
394 const int size = dofs.size();
395 TemporarayDofVectorType& tmpDofs = *tmpDofs_;
396 tmpDofs.resize( size );
397 (*this)(localFunction, tmpDofs );
398
399 // copy to dof vector
400 for( int i=0; i<size; ++i )
401 {
402 dofs[ i ] = tmpDofs[ i ];
403 }
404 }
405 };
406
407 } // namespace Fem
408
409} // namespace Dune
410
411#endif // #ifndef DUNE_FEM_SPACE_LOCALFINITEELEMENT_INTERPOLATION_HH
discrete function space
Construct a vector with a dynamic size.
Definition: dynvector.hh:59
FunctionSpaceTraits::RangeType RangeType
Type of range vector (using type of range field) has a Dune::FieldVector type interface.
Definition: functionspaceinterface.hh:71
@ dimRange
dimension of range vector space
Definition: functionspaceinterface.hh:48
FunctionSpaceTraits::RangeFieldType RangeFieldType
Intrinsic type used for values in the range field (usually a double)
Definition: functionspaceinterface.hh:63
A vector valued function space.
Definition: functionspace.hh:60
FunctionSpaceType::DomainType DomainType
type of domain vectors, i.e., type of coordinates
Definition: localfunction.hh:108
FunctionSpaceType::RangeType RangeType
type of range vectors, i.e., type of function values
Definition: localfunction.hh:110
Implements a vector constructed from a given type representing a field and a compile-time given size.
typename Impl::voider< Types... >::type void_t
Is void for all valid input types. The workhorse for C++11 SFINAE-techniques.
Definition: typetraits.hh:40
Dune namespace.
Definition: alignedallocator.hh:13
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 21, 23:30, 2024)