FieldMatVec.hpp 12 KB
Newer Older
1
2
#pragma once

3
4
#include <type_traits>

5
#include <dune/common/diagonalmatrix.hh>
6
7
#include <dune/common/fmatrix.hh>
#include <dune/common/fvector.hh>
8
#include <dune/common/typetraits.hh>
9

Praetorius, Simon's avatar
Praetorius, Simon committed
10
11
#include <amdis/common/TypeTraits.hpp>

12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
namespace std
{
  template <class T, int N>
  struct common_type<Dune::FieldVector<T,N>, T>
  {
    using type = T;
  };

  template <class T, int N, int M>
  struct common_type<Dune::FieldMatrix<T,N,M>, T>
  {
    using type = T;
  };
}

27
namespace Dune
28
{
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
  namespace MatVec
  {
    /// Traits to detect fixed size containers like FieldVector and FieldMatrix
    /// @{
    template <class T>
    struct IsMatrix : std::false_type {};

    template <class T, int M, int N>
    struct IsMatrix<FieldMatrix<T,M,N>> : std::true_type {};

    template <class T, int N>
    struct IsMatrix<DiagonalMatrix<T,N>> : std::true_type {};


    template <class T>
    struct IsVector : std::false_type {};

    template <class T, int N>
    struct IsVector<FieldVector<T,N>> : std::true_type {};


    template <class T>
    struct IsMatVec
      : std::integral_constant<bool, IsMatrix<T>::value || IsVector<T>::value> {};
    /// @}

    /// Convert the field_Type of container to type S
    /// @{
    template <class A, class S>
    struct MakeMatVec
    {
      using type = A;
    };

    template <class T, int M, int N, class S>
    struct MakeMatVec<FieldMatrix<T,M,N>,S>
    {
      using type = FieldMatrix<S,M,N>;
    };

    template <class T, int N, class S>
    struct MakeMatVec<DiagonalMatrix<T,N>,S>
    {
      using type = DiagonalMatrix<S,N>;
    };

    template <class T, int N, class S>
    struct MakeMatVec<FieldVector<T,N>,S>
    {
      using type = FieldVector<S,N>;
    };
    /// @}

    /// Convert pseudo-scalar to real scalar type
    /// @{
    template <class T>
85
    decltype(auto) as_scalar(T&& t)
86
    {
Praetorius, Simon's avatar
Praetorius, Simon committed
87
      return FWD(t);
88
89
90
    }

    template <class T>
91
    T as_scalar(FieldVector<T,1> const& t)
92
93
94
95
96
    {
      return t[0];
    }

    template <class T>
97
    T as_scalar(FieldMatrix<T,1,1> const& t)
98
99
100
101
102
    {
      return t[0][0];
    }

    template <class T>
103
    T as_scalar(DiagonalMatrix<T,1> const& t)
104
105
106
107
108
    {
      return t.diagonal(0);
    }
    /// @}

109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
    /// Convert pseudo-vector to real vector type
    /// @{
    template <class T>
    decltype(auto) as_vector(T&& t)
    {
      return FWD(t);
    }

    template <class T, int N>
    FieldVector<T,N> const& as_vector(FieldMatrix<T,1,N> const& t)
    {
      return t[0];
    }

    template <class T, int N>
    FieldVector<T,N>& as_vector(FieldMatrix<T,1,N>& t)
    {
      return t[0];
    }
    /// @}

130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146

    /// Convert pseudo-matrix to real matrix type with proper operator[][]
    /// @{
    template <class T>
    decltype(auto) as_matrix(T&& t)
    {
      return FWD(t);
    }

    template <class Mat>
    class MatrixView;

    template <class T, int N>
    MatrixView<DiagonalMatrix<T,N>> as_matrix(DiagonalMatrix<T,N> const& mat)
    {
      return {mat};
    }
147
    /// @}
148

149
150
151
152
153
154
    // returns -a
    template <class A>
    auto negate(A const& a);

    // returns a+b
    template <class A, class B>
155
    auto plus(A const& a, B const& b);
156
157
158

    // returns a-b
    template <class A, class B>
159
    auto minus(A const& a, B const& b);
160
161
162

    // returns a*b
    template <class A, class B,
163
164
165
166
167
      std::enable_if_t<IsNumber<A>::value && IsNumber<B>::value, int> = 0>
    auto multiplies(A const& a, B const& b);

    template <class A, class B,
      std::enable_if_t<IsNumber<A>::value != IsNumber<B>::value, int> = 0>
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
    auto multiplies(A const& a, B const& b);

    template <class T, int N, class S>
    auto multiplies(FieldVector<T,N> const& a, FieldVector<S,N> const& b);

    template <class Mat, class Vec,
      std::enable_if_t<IsMatrix<Mat>::value && IsVector<Vec>::value, int> = 0>
    auto multiplies(Mat const& mat, Vec const& vec);

    template <class Vec, class Mat,
      std::enable_if_t<IsVector<Vec>::value && IsMatrix<Mat>::value, int> = 0>
    auto multiplies(Vec const& vec, Mat const& mat);

    template <class T, int L, int M, int N, class S>
    auto multiplies(FieldMatrix<T,L,M> const& a, FieldMatrix<S,M,N> const& b);

    // return a/b
    template <class A, class B>
    auto divides(A a, B const& b)
    {
      return a /= b;
    }

  } // end namespace MatVec

  // some arithmetic operations with FieldVector and FieldMatrix

  template <class A,
    std::enable_if_t<MatVec::IsMatVec<A>::value, int> = 0>
  auto operator-(A const& a)
  {
199
    return MatVec::negate(MatVec::as_scalar(a));
200
  }
201

202
203
204
205
  template <class A, class B,
    std::enable_if_t<MatVec::IsMatVec<A>::value || MatVec::IsMatVec<B>::value, int> = 0>
  auto operator+(A const& a, B const& b)
  {
206
    return MatVec::plus(MatVec::as_scalar(a), MatVec::as_scalar(b));
207
  }
208

209
210
211
212
  template <class A, class B,
    std::enable_if_t<MatVec::IsMatVec<A>::value || MatVec::IsMatVec<B>::value, int> = 0>
  auto operator-(A const& a, B const& b)
  {
213
    return MatVec::minus(MatVec::as_scalar(a), MatVec::as_scalar(b));
214
  }
215

216
217
218
219
  template <class A, class B,
    std::enable_if_t<MatVec::IsMatVec<A>::value || MatVec::IsMatVec<B>::value, int> = 0>
  auto operator*(A const& a, B const& b)
  {
220
    return MatVec::multiplies(MatVec::as_scalar(a), MatVec::as_scalar(b));
221
  }
222

223
224
225
226
  template <class A, class B,
    std::enable_if_t<MatVec::IsMatVec<A>::value || MatVec::IsMatVec<B>::value, int> = 0>
  auto operator/(A const& a, B const& b)
  {
227
    return MatVec::divides(MatVec::as_scalar(a), MatVec::as_scalar(b));
228
  }
229
230
231
232
233

  // ----------------------------------------------------------------------------

  /// Cross-product a 2d-vector = orthogonal vector
  template <class T>
234
  FieldVector<T, 2> cross(FieldVector<T, 2> const& a);
235

236
  /// Cross-product of two 3d-vectors
237
  template <class T>
238
  FieldVector<T, 3> cross(FieldVector<T, 3> const& a, FieldVector<T, 3> const& b);
239
240
241

  /// Dot product (vec1^T * vec2)
  template <class T, class S, int N>
242
  auto dot(FieldVector<T,N> const& vec1, FieldVector<S,N> const& vec2);
243

244
245
246
  template <class T, class S, int N>
  auto dot(FieldMatrix<T,1,N> const& vec1, FieldMatrix<S,1,N> const& vec2);

247
248
249
250
  // ----------------------------------------------------------------------------

  /// Sum of vector entires.
  template <class T, int N>
251
  T sum(FieldVector<T, N> const& x);
252

253
  template <class T, int N>
254
  T sum(FieldMatrix<T, 1, N> const& x);
255

256
257

  /// Dot-product with the vector itself
258
259
260
261
  template <class T,
    std::enable_if_t<Dune::IsNumber<T>::value, int> = 0>
  auto unary_dot(T const& x);

262
  template <class T, int N>
263
  auto unary_dot(FieldVector<T, N> const& x);
264

265
  template <class T, int N>
266
  auto unary_dot(FieldMatrix<T, 1, N> const& x);
267

268
269
  /// Maximum over all vector entries
  template <class T, int N>
270
  auto max(FieldVector<T, N> const& x);
271

272
  template <class T, int N>
273
  auto max(FieldMatrix<T, 1, N> const& x);
274

275
276
  /// Minimum over all vector entries
  template <class T, int N>
277
  auto min(FieldVector<T, N> const& x);
278

279
  template <class T, int N>
280
  auto min(FieldMatrix<T, 1, N> const& x);
281

282
283
  /// Maximum of the absolute values of vector entries
  template <class T, int N>
284
  auto abs_max(FieldVector<T, N> const& x);
285

286
  template <class T, int N>
287
  auto abs_max(FieldMatrix<T, 1, N> const& x);
288

289
290
  /// Minimum of the absolute values of vector entries
  template <class T, int N>
291
  auto abs_min(FieldVector<T, N> const& x);
292

293
  template <class T, int N>
294
  auto abs_min(FieldMatrix<T, 1, N> const& x);
295

296
297
298
299
300
301
  // ----------------------------------------------------------------------------

  /** \ingroup vector_norms
   *  \brief The 1-norm of a vector = sum_i |x_i|
   **/
  template <class T, int N>
302
  auto one_norm(FieldVector<T, N> const& x);
303

304
  template <class T, int N>
305
  auto one_norm(FieldMatrix<T, 1, N> const& x);
306

307
308
309
  /** \ingroup vector_norms
   *  \brief The euklidean 2-norm of a vector = sqrt( sum_i |x_i|^2 )
   **/
310
311
312
313
  template <class T,
    std::enable_if_t<Dune::IsNumber<T>::value, int> = 0>
  auto two_norm(T const& x);

314
  template <class T, int N>
315
  auto two_norm(FieldVector<T, N> const& x);
316

317
  template <class T, int N>
318
  auto two_norm(FieldMatrix<T, 1, N> const& x);
319

320
321
322
323
  /** \ingroup vector_norms
   *  \brief The p-norm of a vector = ( sum_i |x_i|^p )^(1/p)
   **/
  template <int p, class T, int N>
324
  auto p_norm(FieldVector<T, N> const& x);
325

326
  template <int p, class T, int N>
327
  auto p_norm(FieldMatrix<T, 1, N> const& x);
328

329
330
331
332
  /** \ingroup vector_norms
   *  \brief The infty-norm of a vector = max_i |x_i| = alias for \ref abs_max
   **/
  template <class T, int N>
333
  auto infty_norm(FieldVector<T, N> const& x);
334

335
  template <class T, int N>
336
  auto infty_norm(FieldMatrix<T, 1, N> const& x);
337

338
339
340
  // ----------------------------------------------------------------------------

  /// The euklidean distance between two vectors = |lhs-rhs|_2
341
342
343
344
  template <class T,
    std::enable_if_t<Dune::IsNumber<T>::value, int> = 0>
  T distance(T const& lhs, T const& rhs);

345
  template <class T, int N>
346
  T distance(FieldVector<T, N> const& lhs, FieldVector<T, N> const& rhs);
347
348
349
350
351

  // ----------------------------------------------------------------------------

  /// Outer product (vec1 * vec2^T)
  template <class T, class S, int N, int M, int K>
352
  auto outer(FieldMatrix<T,N,K> const& vec1, FieldMatrix<S,M,K> const& vec2);
353
354
355
356

  // ----------------------------------------------------------------------------

  template <class T>
357
  T det(FieldMatrix<T, 0, 0> const& /*mat*/);
358
359
360

  /// Determinant of a 1x1 matrix
  template <class T>
361
  T det(FieldMatrix<T, 1, 1> const& mat);
362
363
364

  /// Determinant of a 2x2 matrix
  template <class T>
365
  T det(FieldMatrix<T, 2, 2> const& mat);
366
367
368

  /// Determinant of a 3x3 matrix
  template <class T>
369
  T det(FieldMatrix<T, 3, 3> const& mat);
370
371
372

  /// Determinant of a NxN matrix
  template <class T,  int N>
373
374
  T det(FieldMatrix<T, N, N> const& mat);

375
376
377

  /// Return the inverse of the matrix `mat`
  template <class T, int N>
378
  auto inv(FieldMatrix<T, N, N> mat);
379
380
381

  /// Solve the linear system A*x = b
  template <class T, int N>
382
  void solve(FieldMatrix<T, N, N> const& A,  FieldVector<T, N>& x,  FieldVector<T, N> const& b);
383
384
385
386


  /// Gramian determinant = sqrt( det( DT^T * DF ) )
  template <class T, int N, int M>
387
  T gramian(FieldMatrix<T,N,M> const& DF);
388
389
390

  /// Gramian determinant, specialization for 1 column matrices
  template <class T, int M>
391
  T gramian(FieldMatrix<T, 1, M> const& DF);
392

393
394
395
  // ----------------------------------------------------------------------------
  // some arithmetic operations with FieldMatrix

396
  template <class T, int M, int N>
397
  FieldMatrix<T,N,M> trans(FieldMatrix<T, M, N> const& A);
398

399
400
401
402
403
404
  template <class T, int N>
  DiagonalMatrix<T,N> const& trans(DiagonalMatrix<T,N> const& A)
  {
    return A;
  }

405
  // -----------------------------------------------------------------------------
406

407
408
  template <class T1, class T2, int M, int N, int L>
  FieldMatrix<std::common_type_t<T1,T2>,M,N> multiplies_AtB(FieldMatrix<T1, L, M> const& A,  FieldMatrix<T2, N, L> const& B);
409

410
411
  template <class T1, class T2, int M, int N, int L>
  FieldMatrix<std::common_type_t<T1,T2>,M,N> multiplies_ABt(FieldMatrix<T1, M, L> const& A,  FieldMatrix<T2, N, L> const& B);
412

413
414
  template <class T1, class T2, class T3, int M, int N, int L>
  FieldMatrix<T3,M,N>& multiplies_ABt(FieldMatrix<T1, M, L> const& A,  FieldMatrix<T2, N, L> const& B, FieldMatrix<T3,M,N>& C);
415

416
417
418
419
420
  template <class T1, class T2, class T3, int M, int N>
  FieldMatrix<T3,M,N>& multiplies_ABt(FieldMatrix<T1, M, N> const& A,  DiagonalMatrix<T2, N> const& B, FieldMatrix<T3,M,N>& C);

  template <class T1, class T2, class T3, int N>
  FieldVector<T3,N>& multiplies_ABt(FieldMatrix<T1, 1, N> const& A,  DiagonalMatrix<T2, N> const& B, FieldVector<T3,N>& C);
421

422
  // -----------------------------------------------------------------------------
423
424
425
426
427
428
429

  template <class T, int N>
  T const& at(FieldMatrix<T,N,1> const& vec, std::size_t i);

  template <class T, int M>
  T const& at(FieldMatrix<T,1,M> const& vec, std::size_t i);

430
  // necessary specialization to resolve ambiguous calls
431
432
433
  template <class T>
  T const& at(FieldMatrix<T,1,1> const& vec, std::size_t i);

434
435
436
  template <class T, int N>
  T const& at(FieldVector<T,N> const& vec, std::size_t i);

437
438
439
440
441
442
443
} // end namespace Dune

namespace AMDiS
{
  using Dune::FieldMatrix;
  using Dune::FieldVector;
}
444
445

#include "FieldMatVec.inc.hpp"