FieldMatVec.inc.hpp 14.5 KB
Newer Older
1
2
3
4
5
#pragma once

#include <algorithm>
#include <limits>

6
7
#include <dune/functions/functionspacebases/flatvectorview.hh>

8
9
10
11
12
13
#include <amdis/common/Math.hpp>
#include <amdis/operations/Arithmetic.hpp>
#include <amdis/operations/MaxMin.hpp>

#ifndef DOXYGEN

14
namespace Dune {
15

16
namespace MatVec {
17

18
19
20
21
22
  template <class A>
  auto negate(A const& a)
  {
    return multiplies(a, -1);
  }
23

24
25
26
27
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
  // returns a+b
  template <class A, class B>
  auto plus(A const& a, B const& b)
  {
    using T = std::common_type_t<typename FieldTraits<A>::field_type, typename FieldTraits<B>::field_type>;
    typename MakeMatVec<A,T>::type c{a};

    auto b_ = Dune::Functions::flatVectorView(b);
    auto c_ = Dune::Functions::flatVectorView(c);
    assert(int(b_.size()) == int(c_.size()));
    for(int i = 0; i < int(b_.size()); ++i)
      c_[i] += b_[i];

    return c;
  }

  // returns a-b
  template <class A, class B>
  auto minus(A const& a, B const& b)
  {
    using T = std::common_type_t<typename FieldTraits<A>::field_type, typename FieldTraits<B>::field_type>;
    typename MakeMatVec<A,T>::type c{a};

    auto b_ = Dune::Functions::flatVectorView(b);
    auto c_ = Dune::Functions::flatVectorView(c);
    assert(int(b_.size()) == int(c_.size()));
    for(int i = 0; i < int(b_.size()); ++i)
      c_[i] -= b_[i];

    return c;
  }

56
  template <class A, class B,
57
58
59
60
61
62
63
64
    std::enable_if_t<IsNumber<A>::value && IsNumber<B>::value, int>>
  auto multiplies(A const& a, B const& b)
  {
    return a * b;
  }

  template <class A, class B,
    std::enable_if_t<IsNumber<A>::value != IsNumber<B>::value, int>>
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
  auto multiplies(A const& a, B const& b)
  {
    using T = std::common_type_t<typename FieldTraits<A>::field_type, typename FieldTraits<B>::field_type>;
#if AMDIS_HAS_CXX_CONSTEXPR_IF
    if constexpr(IsNumber<A>::value) {
      typename MakeMatVec<B,T>::type b_{b};
      return b_ *= a;
    } else {
      typename MakeMatVec<A,T>::type a_{a};
      return a_ *= b;
    }
#else
    return Hybrid::ifElse(IsNumber<A>{},
    [&](auto id) { typename MakeMatVec<B,T>::type b_{b}; return id(b_) *= id(a); },
    [&](auto id) { typename MakeMatVec<A,T>::type a_{a}; return id(a_) *= id(b); });
#endif
  }
82

83
84
85
86
87
  template <class T, int N, class S>
  auto multiplies(FieldVector<T,N> const& a, FieldVector<S,N> const& b)
  {
    return a.dot(b);
  }
88

89
90
91
92
93

  template <class Mat, class Vec,
    std::enable_if_t<IsMatrix<Mat>::value && IsVector<Vec>::value, int>>
  auto multiplies(Mat const& mat, Vec const& vec)
  {
Müller, Felix's avatar
Müller, Felix committed
94
    static_assert(int(Mat::cols) == int(Vec::dimension), "");
95
96
97
98
99
100
101
102
103
104
105
    using T = std::common_type_t<typename FieldTraits<Vec>::field_type, typename FieldTraits<Mat>::field_type>;
    FieldVector<T,Mat::rows> y;
    mat.mv(vec, y);
    return y;
  }


  template <class Vec, class Mat,
    std::enable_if_t<IsVector<Vec>::value && IsMatrix<Mat>::value, int>>
  auto multiplies(Vec const& vec, Mat const& mat)
  {
Müller, Felix's avatar
Müller, Felix committed
106
    static_assert(int(Mat::rows) == int(Vec::dimension), "");
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
    using T = std::common_type_t<typename FieldTraits<Vec>::field_type, typename FieldTraits<Mat>::field_type>;
    FieldVector<T,Mat::cols> y;
    mat.mtv(vec, y);
    return y;
  }


  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)
  {
    FieldMatrix<std::common_type_t<T,S>,L,N> C;

    for (int i = 0; i < L; ++i) {
      for (int j = 0; j < N; ++j) {
        C[i][j] = 0;
        for (int k = 0; k < M; ++k)
          C[i][j] += a[i][k]*b[k][j];
      }
    }
    return C;
  }

129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
  template <class T, int SIZE>
  class MatrixView<DiagonalMatrix<T,SIZE>>
  {
    using Matrix = DiagonalMatrix<T,SIZE>;
    using size_type = typename Matrix::size_type;

    struct RowView
    {
      T operator[](size_type c) const
      {
        assert(0 <= c && c < mat_->M());
        return c == r_ ? mat_->diagonal(r_) : T(0);
      }

      Matrix const* mat_;
      size_type r_;
    };

  public:
    MatrixView(Matrix const& mat)
      : mat_(mat)
    {}

    RowView operator[](size_type r) const
    {
      assert(0 <= r && r < mat_.N());
      return {&mat_,r};
    }

    size_type N() const
    {
      return mat_.N();
    }

    size_type M() const
    {
      return mat_.M();
    }

  private:
    Matrix const& mat_;
  };

172
} // end namespace MatVec
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

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

/// Cross-product a 2d-vector = orthogonal vector
template <class T>
FieldVector<T, 2> cross(FieldVector<T, 2> const& a)
{
  return {{ a[1], -a[0] }};
}

/// Cross-product of two vectors (in 3d only)
template <class T>
FieldVector<T, 3> cross(FieldVector<T, 3> const& a, FieldVector<T, 3> const& b)
{
  return {{ a[1]*b[2] - a[2]*b[1],
            a[2]*b[0] - a[0]*b[2],
            a[0]*b[1] - a[1]*b[0] }};
}

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

199
200
201
202
203
204
template <class T, class S, int N>
auto dot(FieldMatrix<T,1,N> const& vec1, FieldMatrix<S,1,N> const& vec2)
{
  return vec1[0].dot(vec2[0]);
}

205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
// ----------------------------------------------------------------------------

namespace Impl
{
  template <class T, int N, class Operation>
  T accumulate(FieldVector<T, N> const& x, T init, Operation op)
  {
    for (int i = 0; i < N; ++i)
      init = op(init, x[i]);
    return init;
  }

  template <class T, int N, class Operation>
  T accumulate(FieldMatrix<T, 1, N> const& x, T init, Operation op)
  {
    for (int i = 0; i < N; ++i)
      init = op(init, x[0][i]);
    return init;
  }

} // end namespace Impl

/// Sum of vector entires.
template <class T, int N>
T sum(FieldVector<T, N> const& x)
{
231
  return Impl::accumulate(x, T(0), AMDiS::Operation::Plus{});
232
233
234
235
236
}

template <class T, int N>
T sum(FieldMatrix<T, 1, N> const& x)
{
237
  return Impl::accumulate(x, T(0), AMDiS::Operation::Plus{});
238
239
240
241
}


/// Dot-product with the vector itself
242
243
244
245
246
247
248
249
template <class T,
  std::enable_if_t<Dune::IsNumber<T>::value, int> >
auto unary_dot(T const& x)
{
  using std::abs;
  return AMDiS::Math::sqr(abs(x));
}

250
251
252
template <class T, int N>
auto unary_dot(FieldVector<T, N> const& x)
{
253
  auto op = [](auto const& a, auto const& b) { using std::abs; return a + AMDiS::Math::sqr(abs(b)); };
254
255
256
257
258
259
  return Impl::accumulate(x, T(0), op);
}

template <class T, int N>
auto unary_dot(FieldMatrix<T, 1, N> const& x)
{
260
  auto op = [](auto const& a, auto const& b) { using std::abs; return a + AMDiS::Math::sqr(abs(b)); };
261
262
263
264
265
266
267
  return Impl::accumulate(x, T(0), op);
}

/// Maximum over all vector entries
template <class T, int N>
auto max(FieldVector<T, N> const& x)
{
268
  return Impl::accumulate(x, std::numeric_limits<T>::lowest(), AMDiS::Operation::Max{});
269
270
271
272
273
}

template <class T, int N>
auto max(FieldMatrix<T, 1, N> const& x)
{
274
  return Impl::accumulate(x, std::numeric_limits<T>::lowest(), AMDiS::Operation::Max{});
275
276
277
278
279
280
}

/// Minimum over all vector entries
template <class T, int N>
auto min(FieldVector<T, N> const& x)
{
281
  return Impl::accumulate(x, std::numeric_limits<T>::max(), AMDiS::Operation::Min{});
282
283
284
285
286
}

template <class T, int N>
auto min(FieldMatrix<T, 1, N> const& x)
{
287
  return Impl::accumulate(x, std::numeric_limits<T>::max(), AMDiS::Operation::Min{});
288
289
290
291
292
293
}

/// Maximum of the absolute values of vector entries
template <class T, int N>
auto abs_max(FieldVector<T, N> const& x)
{
294
  return Impl::accumulate(x, T(0), AMDiS::Operation::AbsMax{});
295
296
297
298
299
}

template <class T, int N>
auto abs_max(FieldMatrix<T, 1, N> const& x)
{
300
  return Impl::accumulate(x, T(0), AMDiS::Operation::AbsMax{});
301
302
303
304
305
306
}

/// Minimum of the absolute values of vector entries
template <class T, int N>
auto abs_min(FieldVector<T, N> const& x)
{
307
  return Impl::accumulate(x, std::numeric_limits<T>::max(), AMDiS::Operation::AbsMin{});
308
309
310
311
312
}

template <class T, int N>
auto abs_min(FieldMatrix<T, 1, N> const& x)
{
313
  return Impl::accumulate(x, std::numeric_limits<T>::max(), AMDiS::Operation::AbsMin{});
314
315
316
317
318
319
320
321
322
323
}

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

/** \ingroup vector_norms
  *  \brief The 1-norm of a vector = sum_i |x_i|
  **/
template <class T, int N>
auto one_norm(FieldVector<T, N> const& x)
{
324
  auto op = [](auto const& a, auto const& b) { using std::abs; return a + abs(b); };
325
326
327
328
329
330
  return Impl::accumulate(x, T(0), op);
}

template <class T, int N>
auto one_norm(FieldMatrix<T, 1, N> const& x)
{
331
  auto op = [](auto const& a, auto const& b) { using std::abs; return a + abs(b); };
332
333
334
335
336
337
  return Impl::accumulate(x, T(0), op);
}

/** \ingroup vector_norms
  *  \brief The euklidean 2-norm of a vector = sqrt( sum_i |x_i|^2 )
  **/
338
339
340
341
342
343
344
345
template <class T,
  std::enable_if_t<Dune::IsNumber<T>::value, int> >
auto two_norm(T const& x)
{
  using std::abs;
  return abs(x);
}

346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
template <class T, int N>
auto two_norm(FieldVector<T, N> const& x)
{
  return std::sqrt(unary_dot(x));
}

template <class T, int N>
auto two_norm(FieldMatrix<T, 1, N> const& x)
{
  return std::sqrt(unary_dot(x));
}

/** \ingroup vector_norms
  *  \brief The p-norm of a vector = ( sum_i |x_i|^p )^(1/p)
  **/
template <int p, class T, int N>
auto p_norm(FieldVector<T, N> const& x)
{
364
  auto op = [](auto const& a, auto const& b) { using std::abs; return a + AMDiS::Math::pow<p>(abs(b)); };
365
366
367
368
369
370
  return std::pow( Impl::accumulate(x, T(0), op), 1.0/p );
}

template <int p, class T, int N>
auto p_norm(FieldMatrix<T, 1, N> const& x)
{
371
  auto op = [](auto const& a, auto const& b) { using std::abs; return a + AMDiS::Math::pow<p>(abs(b)); };
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
  return std::pow( Impl::accumulate(x, T(0), op), 1.0/p );
}

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

template <class T, int N>
auto infty_norm(FieldMatrix<T, 1, N> const& x)
{
  return abs_max(x);
}

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

/// The euklidean distance between two vectors = |lhs-rhs|_2
393
394
395
396
397
398
399
400
template <class T,
  std::enable_if_t<Dune::IsNumber<T>::value, int> >
T distance(T const& lhs, T const& rhs)
{
  using std::abs;
  return abs(lhs - rhs);
}

401
402
403
template <class T, int N>
T distance(FieldVector<T, N> const& lhs, FieldVector<T, N> const& rhs)
{
404
  using std::sqrt;
405
406
  T result = 0;
  for (int i = 0; i < N; ++i)
407
    result += AMDiS::Math::sqr(lhs[i] - rhs[i]);
408
  return sqrt(result);
409
410
411
412
413
414
415
416
}

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

/// Outer product (vec1 * vec2^T)
template <class T, class S, int N, int M, int K>
auto outer(FieldMatrix<T,N,K> const& vec1, FieldMatrix<S,M,K> const& vec2)
{
Praetorius, Simon's avatar
Praetorius, Simon committed
417
  using result_type = FieldMatrix<TYPEOF( std::declval<T>() * std::declval<S>() ), N, M>;
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
  result_type mat;
  for (int i = 0; i < N; ++i)
    for (int j = 0; j < M; ++j)
      mat[i][j] = vec1[i].dot(vec2[j]);
  return mat;
}

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

template <class T>
T det(FieldMatrix<T, 0, 0> const& /*mat*/)
{
  return 0;
}

/// Determinant of a 1x1 matrix
template <class T>
T det(FieldMatrix<T, 1, 1> const& mat)
{
  return mat[0][0];
}

/// Determinant of a 2x2 matrix
template <class T>
T det(FieldMatrix<T, 2, 2> const& mat)
{
  return mat[0][0]*mat[1][1] - mat[0][1]*mat[1][0];
}

/// Determinant of a 3x3 matrix
template <class T>
T det(FieldMatrix<T, 3, 3> const& mat)
{
  return mat[0][0]*mat[1][1]*mat[2][2] + mat[0][1]*mat[1][2]*mat[2][0] + mat[0][2]*mat[1][0]*mat[2][1]
      - (mat[0][2]*mat[1][1]*mat[2][0] + mat[0][1]*mat[1][0]*mat[2][2] + mat[0][0]*mat[1][2]*mat[2][1]);
}

/// Determinant of a NxN matrix
template <class T,  int N>
T det(FieldMatrix<T, N, N> const& mat)
{
  return mat.determinant();
}

/// Return the inverse of the matrix `mat`
template <class T, int N>
auto inv(FieldMatrix<T, N, N> mat)
{
  mat.invert();
  return mat;
}

/// Solve the linear system A*x = b
template <class T, int N>
void solve(FieldMatrix<T, N, N> const& A,  FieldVector<T, N>& x,  FieldVector<T, N> const& b)
{
  A.solve(x, b);
}


/// Gramian determinant = sqrt( det( DT^T * DF ) )
template <class T, int N, int M>
T gramian(FieldMatrix<T,N,M> const& DF)
{
  using std::sqrt;
  return sqrt( det(outer(DF, DF)) );
}

/// Gramian determinant, specialization for 1 column matrices
template <class T, int M>
T gramian(FieldMatrix<T, 1, M> const& DF)
{
  using std::sqrt;
  return sqrt(dot(DF[0], DF[0]));
}

// ----------------------------------------------------------------------------
// some arithmetic operations with FieldMatrix

template <class T, int M, int N>
FieldMatrix<T,N,M> trans(FieldMatrix<T, M, N> const& A)
{
  FieldMatrix<T,N,M> At;
  for (int i = 0; i < M; ++i)
    for (int j = 0; j < N; ++j)
      At[j][i] = A[i][j];

  return At;
}


509
510
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)
511
{
512
  FieldMatrix<std::common_type_t<T1,T2>,M,N> C;
513
514
515
516
517
518
519
520
521
522
523

  for (int m = 0; m < M; ++m) {
    for (int n = 0; n < N; ++n) {
      C[m][n] = 0;
      for (int l = 0; l < L; ++l)
        C[m][n] += A[l][m] * B[n][l];
    }
  }
  return C;
}

524
525
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)
526
{
527
  FieldMatrix<std::common_type_t<T1,T2>,M,N> C;
528
  return multiplies_ABt(A,B,C);
529
530
}

531
532
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)
533
534
535
536
537
538
539
540
541
542
543
{
  for (int m = 0; m < M; ++m) {
    for (int n = 0; n < N; ++n) {
      C[m][n] = 0;
      for (int l = 0; l < L; ++l)
        C[m][n] += A[m][l] * B[n][l];
    }
  }
  return C;
}

544
545
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)
546
547
548
549
550
551
552
553
554
{
  for (int m = 0; m < M; ++m) {
    for (int n = 0; n < N; ++n) {
      C[m][n] = A[m][n] * B.diagonal(n);
    }
  }
  return C;
}

555
556
557
558
559
560
561
562
563
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)
{
  for (int n = 0; n < N; ++n) {
    C[n] = A[0][n] * B.diagonal(n);
  }
  return C;
}

564
565
566
567
568
569
570
571
572
573
574
575
template <class T, int N>
T const& at(FieldMatrix<T,N,1> const& vec, std::size_t i)
{
  return vec[i][0];
}

template <class T, int M>
T const& at(FieldMatrix<T,1,M> const& vec, std::size_t i)
{
  return vec[0][i];
}

576
577
578
579
580
581
template <class T>
T const& at(FieldMatrix<T,1,1> const& vec, std::size_t i)
{
  return vec[0][i];
}

582
583
584
585
586
587
template <class T, int N>
T const& at(FieldVector<T,N> const& vec, std::size_t i)
{
  return vec[i];
}

588
589
590
} // end namespace AMDiS

#endif