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 ), geometry_( lf_.entity().geometry() ), bset_( bset ) {}
76
77 template< class Arg >
78 void evaluate ( const Arg &x, typename Traits::RangeType &y ) const
79 {
80 typename Traits::RangeType help;
81 lf_.evaluate( x, help );
82 typename Transformation::InverseTransformationType transf( geometry_, x );
83 y = transf.apply( help );
84 }
85 template< class Arg >
86 typename Traits::RangeType operator() ( const Arg &x ) const
87 {
88 typename Traits::RangeType y;
89 evaluate(x,y);
90 return y;
91 }
92
93 private:
94 const LocalFunction &lf_;
95 const typename Entity::Geometry geometry_;
96 const BasisFunctionSetType &bset_;
97 };
98
99 } // namespace Impl
100
101
102 // LocalFiniteElementInterpolation
103 // -------------------------------
104
105 template< class Space, class LocalInterpolation, bool scalarSFS >
106 class LocalFiniteElementInterpolation;
107
108 // a vector valued shape basis function set taken from
109 // dune-localfunction - the vector value 'blow up' is not yet supported
110 // for this case
111 template< class Space, class LocalInterpolation >
112 class LocalFiniteElementInterpolation<Space,LocalInterpolation,false>
113 {
114 typedef LocalFiniteElementInterpolation< Space, LocalInterpolation,false > ThisType;
115
116 public:
117 typedef typename Space::BasisFunctionSetType BasisFunctionSetType;
118 typedef LocalInterpolation LocalInterpolationType;
119
120 private:
121 typedef typename BasisFunctionSetType::FunctionSpaceType FunctionSpaceType;
122
123 typedef typename FunctionSpaceType::RangeType RangeType;
124 typedef typename FunctionSpaceType::RangeFieldType RangeFieldType;
125 static const int dimRange = FunctionSpaceType::dimRange;
126
127 typedef std::size_t size_type;
128
129 template< class LocalFunction >
130 using LocalFunctionWrapper = Impl::LocalFunctionWrapper< LocalFunction, BasisFunctionSetType >;
131
132 public:
133 LocalFiniteElementInterpolation()
134 : basisFunctionSet_()
135 , localInterpolation_()
136 {
137 }
138
139 explicit LocalFiniteElementInterpolation ( const BasisFunctionSetType &basisFunctionSet,
140 const LocalInterpolationType &localInterpolation = LocalInterpolationType() )
141 : basisFunctionSet_( basisFunctionSet ),
142 localInterpolation_( localInterpolation )
143 {}
144
145 ThisType& operator=( const ThisType& other )
146 {
147 basisFunctionSet_ = other.basisFunctionSet_;
148 localInterpolation_.emplace( other.localInterpolation() );
149 return *this;
150 }
151
152 void unbind()
153 {
154 // basisFunctionSet_ = BasisFunctionSetType();
155 }
156
157 template< class LocalFunction, class Dof >
158 void operator() ( const LocalFunction &localFunction, std::vector< Dof > &localDofVector ) const
159 {
160 LocalFunctionWrapper< LocalFunction > wrapper( localFunction, basisFunctionSet() );
161 localInterpolation().interpolate( wrapper, localDofVector );
162 }
163
164 template< class LocalFunction, class Dof, class Allocator >
165 void operator() ( const LocalFunction &localFunction, Dune::DynamicVector< Dof, Allocator > &localDofVector ) const
166 {
167 (*this)(localFunction,localDofVector.container() );
168 }
169
170 template< class LocalFunction, class DiscreteFunction, template< class > class Assembly >
171 void operator() ( const LocalFunction &localFunction, LocalContribution< DiscreteFunction, Assembly > &localContribution ) const
172 {
173 (*this)(localFunction,localContribution.localDofVector());
174 }
175
176 BasisFunctionSetType basisFunctionSet () const { return basisFunctionSet_; }
177 const LocalInterpolationType &localInterpolation () const
178 {
179 assert( localInterpolation_.has_value() );
180 return *localInterpolation_;
181 }
182
183 private:
184 BasisFunctionSetType basisFunctionSet_;
185 std::optional< LocalInterpolationType > localInterpolation_;
186 };
187
188 namespace Impl
189 {
190 template <int dimRange>
191 struct RangeConverter
192 {
193 RangeConverter ( std::size_t range ) : range_( range ) {}
194
195 template< class T >
196 FieldVector< T, 1 > operator() ( const FieldVector< T, dimRange > &in ) const
197 {
198 return in[ range_ ];
199 }
200
201 template< class T, int j >
202 FieldMatrix< T, 1, j > operator() ( const FieldMatrix< T, dimRange, j > &in ) const
203 {
204 // return in[ range_ ]; // implicit conversion fails
205 FieldMatrix<T,1,j> ret;
206 ret[0] = in[range_];
207 return ret;
208 }
209
210 protected:
211 std::size_t range_;
212 };
213 template <class DofVector, class DofAlignment>
214 struct SubDofVectorWrapper
215 : public SubDofVector< DofVector, DofAlignment >
216 {
217 typedef SubDofVector< DofVector, DofAlignment > BaseType;
218
219 SubDofVectorWrapper( DofVector& dofs, int coordinate, const DofAlignment &dofAlignment )
220 : BaseType( dofs, coordinate, dofAlignment )
221 {}
222
224 void clear() {}
225 void resize( const unsigned int) {}
226 };
227
228 }
229
230 // a scalar valued shape basis function set taken from
231 // dune-localfunction - for this the vector value 'blow up' is supported
232 // as with other spaces
233 template< class Space, class LocalInterpolation >
234 class LocalFiniteElementInterpolation<Space,LocalInterpolation,true>
235 {
236 typedef LocalFiniteElementInterpolation< Space, LocalInterpolation,true > ThisType;
237
238 public:
239 typedef typename Space::BasisFunctionSetType BasisFunctionSetType;
240 typedef LocalInterpolation LocalInterpolationType;
241
242 private:
243 typedef typename BasisFunctionSetType::FunctionSpaceType FunctionSpaceType;
244 // typedef typename Space::FunctionSpaceType FunctionSpaceType;
245
246 typedef typename FunctionSpaceType::RangeType RangeType;
247 typedef typename FunctionSpaceType::RangeFieldType RangeFieldType;
248 static const int dimRange = FunctionSpaceType::dimRange;
249 static const int dimR = Space::FunctionSpaceType::dimRange;
250
251 typedef std::size_t size_type;
252
253 typedef VerticalDofAlignment< BasisFunctionSetType, RangeType> DofAlignmentType;
254
255 public:
256 LocalFiniteElementInterpolation ()
257 : basisFunctionSet_(),
258 localInterpolation_(),
259 dofAlignment_()
260 {}
261
262 explicit LocalFiniteElementInterpolation ( const BasisFunctionSetType &basisFunctionSet,
263 const LocalInterpolationType &localInterpolation = LocalInterpolationType() )
264 : basisFunctionSet_( basisFunctionSet ),
265 localInterpolation_( localInterpolation ),
266 dofAlignment_( basisFunctionSet_ )
267 {}
268
269 LocalFiniteElementInterpolation ( const ThisType& other )
270 : basisFunctionSet_( other.basisFunctionSet_ ),
271 localInterpolation_( other.localInterpolation_ ),
272 dofAlignment_( basisFunctionSet_ )
273 {}
274
275 ThisType& operator=( const ThisType& other )
276 {
277 basisFunctionSet_ = other.basisFunctionSet_;
278 localInterpolation_.emplace( other.localInterpolation() );
279 return *this;
280 }
281
282 void unbind()
283 {
284 // basisFunctionSet_ = BasisFunctionSetType();
285 }
286
287 template< class LocalFunction, class Vector>
288 void operator() ( const LocalFunction &localFunction, Vector &localDofVector ) const
289 {
290 // clear dofs before something is adedd
291 // localDofVector.clear(); // does not exist on DynVector so use 'fill' instead
292 std::fill(localDofVector.begin(),localDofVector.end(),0);
293 for( std::size_t i = 0; i < dimR; ++i )
294 {
295 Impl::SubDofVectorWrapper< Vector, DofAlignmentType > subLdv( localDofVector, i, dofAlignment_ );
296 (*this)(
297 localFunctionConverter( localFunction, Impl::RangeConverter<dimR>(i) ),
298 subLdv, PriorityTag<42>()
299 );
300 }
301 }
302
303 template< class LocalFunction, class Dof, class Allocator >
304 void operator() ( const LocalFunction &localFunction, Dune::DynamicVector< Dof, Allocator > &localDofVector ) const
305 {
306 (*this)(localFunction,localDofVector.container() );
307 }
308
309 template< class LocalFunction, class DiscreteFunction, template< class > class Assembly >
310 void operator() ( const LocalFunction &localFunction, LocalContribution< DiscreteFunction, Assembly > &localContribution ) const
311 {
312 (*this)(localFunction,localContribution.localDofVector());
313 }
314
315 BasisFunctionSetType basisFunctionSet () const { return basisFunctionSet_; }
316 const LocalInterpolationType &localInterpolation () const
317 {
318 assert( localInterpolation_.has_value() );
319 return *localInterpolation_;
320 }
321
322 private:
323 template< class LocalFunction, class Vector>
324 auto operator() ( const LocalFunction &localFunction, Vector &localDofVector, PriorityTag<1> ) const
325 -> void_t< decltype(
326 std::declval<LocalInterpolationType>().interpolate(
327 localFunction, localDofVector)) >
328 {
329 localInterpolation().interpolate( localFunction, localDofVector );
330 }
331 template< class LocalFunction, class Vector>
332 void operator() ( const LocalFunction &localFunction, Vector &localDofVector, PriorityTag<0> ) const
333 {
334 std::vector<typename Vector::value_type> tmp(basisFunctionSet_.size());
335 localInterpolation().interpolate( localFunction, tmp);
336 for (unsigned int i=0;i<tmp.size();++i)
337 localDofVector[i] = tmp[i];
338 }
339 BasisFunctionSetType basisFunctionSet_;
340 std::optional< LocalInterpolationType > localInterpolation_;
341 DofAlignmentType dofAlignment_;
342 };
343
344 template <class DiscreteFunctionSpace >
345 class LocalFEInterpolationWrapper
346 : public LocalInterpolationWrapper< DiscreteFunctionSpace >
347 {
348 typedef LocalInterpolationWrapper< DiscreteFunctionSpace > BaseType;
349 protected:
350 using BaseType :: interpolation_;
351 typedef std::vector< typename DiscreteFunctionSpace::RangeFieldType >
352 TemporarayDofVectorType;
354
355 using BaseType::interpolation;
356 public:
357 LocalFEInterpolationWrapper( const DiscreteFunctionSpace& space )
358 : BaseType( space )
359 {}
360
361 using BaseType :: bind;
362 using BaseType :: unbind;
363
364 template< class LocalFunction, class Dof >
365 void operator() ( const LocalFunction &localFunction, std::vector< Dof > &localDofVector ) const
366 {
367 // interpolation only works for std::vector storage
368 interpolation()( localFunction, localDofVector );
369 }
370
371 template< class LocalFunction, class Dof, class Allocator >
372 void operator() ( const LocalFunction &localFunction, Dune::DynamicVector< Dof, Allocator > &localDofVector ) const
373 {
374 // DynamicVector internally stores a std::vector
375 (*this)(localFunction, localDofVector.container() );
376 }
377
378
379 template< class LocalFunction, class DiscreteFunction, template< class > class Assembly >
380 void operator() ( const LocalFunction &localFunction, LocalContribution< DiscreteFunction, Assembly > &localContribution ) const
381 {
382 // LocalContribution internally stores a std::vector
383 (*this)(localFunction, localContribution.localDofVector());
384 }
385
387 template< class LocalFunction, class LocalDofVector >
388 void operator () ( const LocalFunction &localFunction, LocalDofVector &dofs ) const
389 {
390 const int size = dofs.size();
391 TemporarayDofVectorType& tmpDofs = *tmpDofs_;
392 tmpDofs.resize( size );
393 (*this)(localFunction, tmpDofs );
394
395 // copy to dof vector
396 for( int i=0; i<size; ++i )
397 {
398 dofs[ i ] = tmpDofs[ i ];
399 }
400 }
401 };
402
403 } // namespace Fem
404
405} // namespace Dune
406
407#endif // #ifndef DUNE_FEM_SPACE_LOCALFINITEELEMENT_INTERPOLATION_HH
discrete function space
Construct a vector with a dynamic size.
Definition: dynvector.hh:59
GridImp::template Codim< cd >::Geometry Geometry
The corresponding geometry type.
Definition: entity.hh:100
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:105
FunctionSpaceType::RangeType RangeType
type of range vectors, i.e., type of function values
Definition: localfunction.hh:107
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 (Jul 27, 22:29, 2024)