Dune Core Modules (2.4.2)

fvector.hh
Go to the documentation of this file.
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_FVECTOR_HH
4#define DUNE_FVECTOR_HH
5
6#include <cmath>
7#include <cstddef>
8#include <cstdlib>
9#include <complex>
10#include <cstring>
11#include <utility>
12#include <initializer_list>
13#include <algorithm>
14
16
17#include "typetraits.hh"
18#include "exceptions.hh"
19#include "array.hh"
20
21#include "ftraits.hh"
22#include "densevector.hh"
23#include "unused.hh"
24
25namespace Dune {
26
36 template< class K, int SIZE > class FieldVector;
37 template< class K, int SIZE >
38 struct DenseMatVecTraits< FieldVector<K,SIZE> >
39 {
40 typedef FieldVector<K,SIZE> derived_type;
41 typedef Dune::array<K,SIZE> container_type;
42 typedef K value_type;
43 typedef typename container_type::size_type size_type;
44 };
45
46 template< class K, int SIZE >
47 struct FieldTraits< FieldVector<K,SIZE> >
48 {
49 typedef typename FieldTraits<K>::field_type field_type;
50 typedef typename FieldTraits<K>::real_type real_type;
51 };
52
61 template<typename C, int SIZE>
63 {
64 enum {
69 value = true
70 };
71 };
72
73 template<typename T, int SIZE>
74 struct IsFieldVectorSizeCorrect<FieldVector<T,SIZE>,SIZE>
75 {
76 enum {value = true};
77 };
78
79 template<typename T, int SIZE, int SIZE1>
80 struct IsFieldVectorSizeCorrect<FieldVector<T,SIZE1>,SIZE>
81 {
82 enum {value = false};
83 };
84
85
91 template< class K, int SIZE >
93 public DenseVector< FieldVector<K,SIZE> >
94 {
95 Dune::array<K,SIZE> _data;
97 public:
99 enum {
101 dimension = SIZE
102 };
103
104 typedef typename Base::size_type size_type;
105 typedef typename Base::value_type value_type;
106
109 : _data{}
110 {}
111
113 explicit FieldVector (const K& t)
114 {
115 fill(t);
116 }
117
119 FieldVector (const FieldVector & x) : Base(), _data(x._data)
120 {}
121
122 FieldVector (std::initializer_list<K> const &l)
123 {
124 assert(l.size() == dimension);// Actually, this is not needed any more!
125 std::copy_n(l.begin(), std::min(static_cast<std::size_t>(dimension),
126 l.size()),
127 _data.begin());
128 }
129
141 template<class C>
142 FieldVector (const DenseVector<C> & x, typename Dune::enable_if<IsFieldVectorSizeCorrect<C,SIZE>::value>::type* dummy=0 )
143 {
145 // do a run-time size check, for the case that x is not a FieldVector
146 assert(x.size() == SIZE); // Actually this is not needed any more!
147 std::copy_n(x.begin(), std::min(static_cast<std::size_t>(SIZE),x.size()), _data.begin());
148 }
149
151 template<class K1, int SIZE1>
153 {
154 static_assert(SIZE1 == SIZE, "FieldVector in constructor has wrong size");
155 for (size_type i = 0; i<SIZE; i++)
156 _data[i] = x[i];
157 }
158 using Base::operator=;
159
160 DUNE_CONSTEXPR size_type size () const { return vec_size(); }
161
162 // make this thing a vector
163 DUNE_CONSTEXPR size_type vec_size () const { return SIZE; }
164 K & vec_access(size_type i) { return _data[i]; }
165 const K & vec_access(size_type i) const { return _data[i]; }
166 private:
167 void fill(const K& t)
168 {
169 for (int i=0; i<SIZE; i++) _data[i]=t;
170 }
171 };
172
184 template<class K, int SIZE>
185 inline std::istream &operator>> ( std::istream &in,
187 {
189 for( typename FieldVector<K, SIZE>::size_type i = 0; i < SIZE; ++i )
190 in >> w[ i ];
191 if(in)
192 v = w;
193 return in;
194 }
195
196#ifndef DOXYGEN
197 template< class K >
198 struct DenseMatVecTraits< FieldVector<K,1> >
199 {
200 typedef FieldVector<K,1> derived_type;
201 typedef K container_type;
202 typedef K value_type;
203 typedef size_t size_type;
204 };
205
208 template<class K>
209 class FieldVector<K, 1> :
210 public DenseVector< FieldVector<K,1> >
211 {
212 K _data;
213 typedef DenseVector< FieldVector<K,1> > Base;
214 public:
216 enum {
218 dimension = 1
219 };
220
221 typedef typename Base::size_type size_type;
222
223 //===== construction
224
226 FieldVector ()
227 : _data()
228 {}
229
231 template<typename T,
232 typename EnableIf = typename std::enable_if<
233 std::is_convertible<T, K>::value &&
234 ! std::is_same<K, DenseVector<typename FieldTraits<T>::field_type>
235 >::value
236 >::type
237 >
238 FieldVector (const T& k) : _data(k) {}
239
241 template<class C>
242 FieldVector (const DenseVector<C> & x)
243 {
244 static_assert(((bool)IsFieldVectorSizeCorrect<C,1>::value), "FieldVectors do not match in dimension!");
245 assert(x.size() == 1);
246 _data = x[0];
247 }
248
250 FieldVector ( const FieldVector &other )
251 : Base(), _data( other._data )
252 {}
253
255 template<typename T,
256 typename EnableIf = typename std::enable_if<
257 std::is_convertible<T, K>::value &&
258 ! std::is_same<K, DenseVector<typename FieldTraits<T>::field_type>
259 >::value
260 >::type
261 >
262 inline FieldVector& operator= (const T& k)
263 {
264 _data = k;
265 return *this;
266 }
267
268 DUNE_CONSTEXPR size_type size () const { return vec_size(); }
269
270 //===== forward methods to container
271 DUNE_CONSTEXPR size_type vec_size () const { return 1; }
272 K & vec_access(size_type i)
273 {
275 assert(i == 0);
276 return _data;
277 }
278 const K & vec_access(size_type i) const
279 {
281 assert(i == 0);
282 return _data;
283 }
284
285 //===== conversion operator
286
288 operator K& () { return _data; }
289
291 operator const K& () const { return _data; }
292 };
293
294 /* ----- FV / FV ----- */
295 /* mostly not necessary as these operations are already covered via the cast operator */
296
298 template<class K>
299 inline bool operator> (const FieldVector<K,1>& a, const FieldVector<K,1>& b)
300 {
301 return a[0]>b[0];
302 }
303
305 template<class K>
306 inline bool operator>= (const FieldVector<K,1>& a, const FieldVector<K,1>& b)
307 {
308 return a[0]>=b[0];
309 }
310
312 template<class K>
313 inline bool operator< (const FieldVector<K,1>& a, const FieldVector<K,1>& b)
314 {
315 return a[0]<b[0];
316 }
317
319 template<class K>
320 inline bool operator<= (const FieldVector<K,1>& a, const FieldVector<K,1>& b)
321 {
322 return a[0]<=b[0];
323 }
324
325 /* ----- FV / scalar ----- */
326
328 template<class K>
329 inline FieldVector<K,1> operator+ (const FieldVector<K,1>& a, const K b)
330 {
331 return a[0]+b;
332 }
333
335 template<class K>
336 inline FieldVector<K,1> operator- (const FieldVector<K,1>& a, const K b)
337 {
338 return a[0]-b;
339 }
340
342 template<class K>
343 inline FieldVector<K,1> operator* (const FieldVector<K,1>& a, const K b)
344 {
345 return a[0]*b;
346 }
347
349 template<class K>
350 inline FieldVector<K,1> operator/ (const FieldVector<K,1>& a, const K b)
351 {
352 return a[0]/b;
353 }
354
356 template<class K>
357 inline bool operator> (const FieldVector<K,1>& a, const K b)
358 {
359 return a[0]>b;
360 }
361
363 template<class K>
364 inline bool operator>= (const FieldVector<K,1>& a, const K b)
365 {
366 return a[0]>=b;
367 }
368
370 template<class K>
371 inline bool operator< (const FieldVector<K,1>& a, const K b)
372 {
373 return a[0]<b;
374 }
375
377 template<class K>
378 inline bool operator<= (const FieldVector<K,1>& a, const K b)
379 {
380 return a[0]<=b;
381 }
382
384 template<class K>
385 inline bool operator== (const FieldVector<K,1>& a, const K b)
386 {
387 return a[0]==b;
388 }
389
391 template<class K>
392 inline bool operator!= (const FieldVector<K,1>& a, const K b)
393 {
394 return a[0]!=b;
395 }
396
397 /* ----- scalar / FV ------ */
398
400 template<class K>
401 inline FieldVector<K,1> operator+ (const K a, const FieldVector<K,1>& b)
402 {
403 return a+b[0];
404 }
405
407 template<class K>
408 inline FieldVector<K,1> operator- (const K a, const FieldVector<K,1>& b)
409 {
410 return a-b[0];
411 }
412
414 template<class K>
415 inline FieldVector<K,1> operator* (const K a, const FieldVector<K,1>& b)
416 {
417 return a*b[0];
418 }
419
421 template<class K>
422 inline FieldVector<K,1> operator/ (const K a, const FieldVector<K,1>& b)
423 {
424 return a/b[0];
425 }
426
428 template<class K>
429 inline bool operator> (const K a, const FieldVector<K,1>& b)
430 {
431 return a>b[0];
432 }
433
435 template<class K>
436 inline bool operator>= (const K a, const FieldVector<K,1>& b)
437 {
438 return a>=b[0];
439 }
440
442 template<class K>
443 inline bool operator< (const K a, const FieldVector<K,1>& b)
444 {
445 return a<b[0];
446 }
447
449 template<class K>
450 inline bool operator<= (const K a, const FieldVector<K,1>& b)
451 {
452 return a<=b[0];
453 }
454
456 template<class K>
457 inline bool operator== (const K a, const FieldVector<K,1>& b)
458 {
459 return a==b[0];
460 }
461
463 template<class K>
464 inline bool operator!= (const K a, const FieldVector<K,1>& b)
465 {
466 return a!=b[0];
467 }
468#endif
469
472} // end namespace
473
474#endif
Fallback implementation of the std::array class (a static array)
Interface for a class of dense vectors over a given field.
Definition: densevector.hh:234
Traits::value_type value_type
export the type representing the field
Definition: densevector.hh:256
Iterator begin()
begin iterator
Definition: densevector.hh:307
size_type size() const
size method
Definition: densevector.hh:296
derived_type & operator=(const value_type &k)
Assignment operator for scalar.
Definition: densevector.hh:275
Traits::size_type size_type
The type used for the index access and size operation.
Definition: densevector.hh:265
vector space out of a tensor product of fields.
Definition: fvector.hh:94
FieldVector(const K &t)
Constructor making vector with identical coordinates.
Definition: fvector.hh:113
FieldVector(const FieldVector &x)
Constructor making vector with identical coordinates.
Definition: fvector.hh:119
@ dimension
The size of this vector.
Definition: fvector.hh:101
FieldVector(const FieldVector< K1, SIZE1 > &x)
Constructor making vector with identical coordinates.
Definition: fvector.hh:152
FieldVector()
Constructor making default-initialized vector.
Definition: fvector.hh:108
FieldVector(const DenseVector< C > &x, typename Dune::enable_if< IsFieldVectorSizeCorrect< C, SIZE >::value >::type *dummy=0)
Copy constructor from a second vector of possibly different type.
Definition: fvector.hh:142
Definition of the DUNE_CONSTEXPR macro.
#define DUNE_CONSTEXPR
Set method or expression constexpr if supported by the compiler.
Definition: constexpr.hh:21
Implements the dense vector interface, with an exchangeable storage class.
Type traits to determine the type of reals (when working with complex numbers)
std::istream & operator>>(std::istream &is, tuple< T1 > &t)
Read a tuple.
Definition: tuples.hh:203
EnableIfInterOperable< T1, T2, bool >::type operator<(const RandomAccessIteratorFacade< T1, V1, R1, D > &lhs, const RandomAccessIteratorFacade< T2, V2, R2, D > &rhs)
Comparison operator.
Definition: iteratorfacades.hh:626
EnableIfInterOperable< T1, T2, bool >::type operator>(const RandomAccessIteratorFacade< T1, V1, R1, D > &lhs, const RandomAccessIteratorFacade< T2, V2, R2, D > &rhs)
Comparison operator.
Definition: iteratorfacades.hh:672
EnableIfInterOperable< T1, T2, bool >::type operator<=(const RandomAccessIteratorFacade< T1, V1, R1, D > &lhs, const RandomAccessIteratorFacade< T2, V2, R2, D > &rhs)
Comparison operator.
Definition: iteratorfacades.hh:649
EnableIfInterOperable< T1, T2, bool >::type operator==(const ForwardIteratorFacade< T1, V1, R1, D > &lhs, const ForwardIteratorFacade< T2, V2, R2, D > &rhs)
Checks for equality.
Definition: iteratorfacades.hh:230
EnableIfInterOperable< T1, T2, bool >::type operator>=(const RandomAccessIteratorFacade< T1, V1, R1, D > &lhs, const RandomAccessIteratorFacade< T2, V2, R2, D > &rhs)
Comparison operator.
Definition: iteratorfacades.hh:694
EnableIfInterOperable< T1, T2, bool >::type operator!=(const ForwardIteratorFacade< T1, V1, R1, D > &lhs, const ForwardIteratorFacade< T2, V2, R2, D > &rhs)
Checks for inequality.
Definition: iteratorfacades.hh:252
Dune namespace.
Definition: alignment.hh:10
TMP to check the size of a DenseVectors statically, if possible.
Definition: fvector.hh:63
@ value
Definition: fvector.hh:69
Traits for type conversions and type information.
Definition of the DUNE_UNUSED macro for the case that config.h is not available.
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentional unused function parameters with.
Definition: unused.hh:18
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)