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

3
#include <dune/common/diagonalmatrix.hh>
4
5
6
#include <dune/common/fmatrix.hh>
#include <dune/common/fvector.hh>

7
#include <amdis/common/ConceptsBase.hpp>
8
#include <amdis/common/ScalarTypes.hpp>
9
10
11

namespace AMDiS
{
Praetorius, Simon's avatar
Praetorius, Simon committed
12
13
14
15
  using Dune::FieldVector;
  using Dune::FieldMatrix;

  // some arithmetic operations with FieldVector
16
17
18

  template <class T, int N, class S,
    REQUIRES(Concepts::Arithmetic<S>) >
19
  FieldVector<T,N> operator*(FieldVector<T,N> v, S factor);
20
21
22

  template <class S, class T, int N,
    REQUIRES(Concepts::Arithmetic<S>) >
23
  FieldVector<T,N> operator*(S factor, FieldVector<T,N> v);
24
25
26

  template <class T, int N, class S,
    REQUIRES(Concepts::Arithmetic<S>) >
27
  FieldVector<T,N> operator/(FieldVector<T,N> v, S factor);
28

29
  template <class T>
30
  FieldVector<T,1> operator*(FieldVector<T,1> const& v, FieldVector<T,1> const& w);
31

32
  template <class T, int N>
33
  FieldVector<T,N> operator*(FieldVector<T,1> const& factor, FieldVector<T,N> v);
34

35
  template <class T, int N>
36
  FieldVector<T,N> operator*(FieldVector<T,N> v, FieldVector<T,1> const& factor);
37

38
39
40
41
  // ----------------------------------------------------------------------------

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

44
  /// Cross-product of two 3d-vectors
45
  template <class T>
46
  FieldVector<T, 3> cross(FieldVector<T, 3> const& a, FieldVector<T, 3> const& b);
47
48
49

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

52
53
  template <class T, int N, int M,
    REQUIRES( N!=1 && M!=1 )>
54
  auto operator*(FieldVector<T,N> const& v, FieldVector<T,M> const& w);
55

56
57
58
59
  // ----------------------------------------------------------------------------

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

62
  template <class T, int N>
63
  T sum(FieldMatrix<T, 1, N> const& x);
64

65
66
67

  /// Dot-product with the vector itself
  template <class T, int N>
68
  auto unary_dot(FieldVector<T, N> const& x);
69

70
  template <class T, int N>
71
  auto unary_dot(FieldMatrix<T, 1, N> const& x);
72

73
74
  /// Maximum over all vector entries
  template <class T, int N>
75
  auto max(FieldVector<T, N> const& x);
76

77
  template <class T, int N>
78
  auto max(FieldMatrix<T, 1, N> const& x);
79

80
81
  /// Minimum over all vector entries
  template <class T, int N>
82
  auto min(FieldVector<T, N> const& x);
83

84
  template <class T, int N>
85
  auto min(FieldMatrix<T, 1, N> const& x);
86

87
88
  /// Maximum of the absolute values of vector entries
  template <class T, int N>
89
  auto abs_max(FieldVector<T, N> const& x);
90

91
  template <class T, int N>
92
  auto abs_max(FieldMatrix<T, 1, N> const& x);
93

94
95
  /// Minimum of the absolute values of vector entries
  template <class T, int N>
96
  auto abs_min(FieldVector<T, N> const& x);
97

98
  template <class T, int N>
99
  auto abs_min(FieldMatrix<T, 1, N> const& x);
100

101
102
103
104
105
106
  // ----------------------------------------------------------------------------

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

109
  template <class T, int N>
110
  auto one_norm(FieldMatrix<T, 1, N> const& x);
111

112
113
114
115
  /** \ingroup vector_norms
   *  \brief The euklidean 2-norm of a vector = sqrt( sum_i |x_i|^2 )
   **/
  template <class T, int N>
116
  auto two_norm(FieldVector<T, N> const& x);
117

118
  template <class T, int N>
119
  auto two_norm(FieldMatrix<T, 1, N> const& x);
120

121
122
123
124
  /** \ingroup vector_norms
   *  \brief The p-norm of a vector = ( sum_i |x_i|^p )^(1/p)
   **/
  template <int p, class T, int N>
125
  auto p_norm(FieldVector<T, N> const& x);
126

127
  template <int p, class T, int N>
128
  auto p_norm(FieldMatrix<T, 1, N> const& x);
129

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

136
  template <class T, int N>
137
  auto infty_norm(FieldMatrix<T, 1, N> const& x);
138

139
140
141
142
  // ----------------------------------------------------------------------------

  /// The euklidean distance between two vectors = |lhs-rhs|_2
  template <class T, int N>
143
  T distance(FieldVector<T, N> const& lhs, FieldVector<T, N> const& rhs);
144
145
146
147
148

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

  /// Outer product (vec1 * vec2^T)
  template <class T, class S, int N, int M, int K>
149
  auto outer(FieldMatrix<T,N,K> const& vec1, FieldMatrix<S,M,K> const& vec2);
150
151
152
153

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

  template <class T>
154
  T det(FieldMatrix<T, 0, 0> const& /*mat*/);
155
156
157

  /// Determinant of a 1x1 matrix
  template <class T>
158
  T det(FieldMatrix<T, 1, 1> const& mat);
159
160
161

  /// Determinant of a 2x2 matrix
  template <class T>
162
  T det(FieldMatrix<T, 2, 2> const& mat);
163
164
165

  /// Determinant of a 3x3 matrix
  template <class T>
166
  T det(FieldMatrix<T, 3, 3> const& mat);
167
168
169

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

172
173
174

  /// Return the inverse of the matrix `mat`
  template <class T, int N>
175
  auto inv(FieldMatrix<T, N, N> mat);
176
177
178

  /// Solve the linear system A*x = b
  template <class T, int N>
179
  void solve(FieldMatrix<T, N, N> const& A,  FieldVector<T, N>& x,  FieldVector<T, N> const& b);
180
181
182
183


  /// Gramian determinant = sqrt( det( DT^T * DF ) )
  template <class T, int N, int M>
184
  T gramian(FieldMatrix<T,N,M> const& DF);
185
186
187

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

190
191
192
  // ----------------------------------------------------------------------------
  // some arithmetic operations with FieldMatrix

193
  template <class T, int M, int N>
194
  FieldMatrix<T,N,M> trans(FieldMatrix<T, M, N> const& A);
195
196


197
198
  template <class T, int M, int N, class S,
    REQUIRES(Concepts::Arithmetic<S>) >
199
  FieldMatrix<T,M,N> operator*(S scalar, FieldMatrix<T, M, N> A);
200

201
202
  template <class T, int M, int N, class S,
    REQUIRES(Concepts::Arithmetic<S>) >
203
  FieldMatrix<T,M,N> operator*(FieldMatrix<T, M, N> A, S scalar);
204

205
206
  template <class T, int M, int N, class S,
    REQUIRES(Concepts::Arithmetic<S>) >
207
208
  FieldMatrix<T,M,N> operator/(FieldMatrix<T, M, N> A, S scalar);

209
210

  template <class T, int M, int N>
211
  FieldMatrix<T,M,N> operator+(FieldMatrix<T, M, N> A, FieldMatrix<T, M, N> const& B);
212
213

  template <class T, int M, int N>
214
  FieldMatrix<T,M,N> operator-(FieldMatrix<T, M, N> A, FieldMatrix<T, M, N> const& B);
215
216
217


  template <class T, int N, int M>
218
  FieldVector<T,N> operator*(FieldMatrix<T,N,M> const& mat, FieldVector<T,M> const& vec);
219

220
221
222
223
224
225
  template <class T, int N, int M>
  FieldMatrix<T,N,M> operator*(FieldMatrix<T,N,M> mat, FieldVector<T,1> const& scalar);

  template <class T, int N>
  FieldMatrix<T,N,1> operator*(FieldMatrix<T,N,1> mat, FieldVector<T,1> const& scalar);

226
227


228
  template <class T, int M, int N, int L>
229
  FieldMatrix<T,M,N> multiplies(FieldMatrix<T, M, L> const& A,  FieldMatrix<T, L, N> const& B);
230
231
232
233



  template <class T, int M, int N, int L>
234
  FieldMatrix<T,M,N> multiplies_AtB(FieldMatrix<T, L, M> const& A,  FieldMatrix<T, N, L> const& B);
235
236

  template <class T, int M, int N, int L>
237
  FieldMatrix<T,M,N> multiplies_ABt(FieldMatrix<T, M, L> const& A,  FieldMatrix<T, N, L> const& B);
238
239

  template <class T, int M, int N, int L>
240
  FieldMatrix<T,M,N>& multiplies_ABt(FieldMatrix<T, M, L> const& A,  FieldMatrix<T, N, L> const& B, FieldMatrix<T,M,N>& C);
241

242
  template <class T, int M, int N>
243
  FieldMatrix<T,M,N>& multiplies_ABt(FieldMatrix<T, M, N> const& A,  Dune::DiagonalMatrix<T, N> const& B, FieldMatrix<T,M,N>& C);
244

245
} // end namespace AMDiS
246
247

#include "FieldMatVec.inc.hpp"