OGS
LinAlg.h
Go to the documentation of this file.
1
11#pragma once
12
13#include <cassert>
14
15#include "BaseLib/Error.h"
16#include "LinAlgEnums.h"
17
18#ifdef USE_PETSC
19#include <petscsystypes.h>
20#endif
21
22namespace MathLib
23{
24
35namespace LinAlg
36{
37
38// Matrix or Vector
39
41template <typename MatrixOrVector>
43{
44 y = x;
45}
46
48template <typename MatrixOrVector>
49void scale(MatrixOrVector& x, double const a)
50{
51 x *= a;
52}
53
55template <typename MatrixOrVector>
56void aypx(MatrixOrVector& y, double const a, MatrixOrVector const& x)
57{
58 y = a * y + x;
59}
60
62template <typename MatrixOrVector>
63void axpy(MatrixOrVector& y, double const a, MatrixOrVector const& x)
64{
65 y += a * x;
66}
67
69template <typename MatrixOrVector>
70void axpby(MatrixOrVector& y, double const a, double const b,
71 MatrixOrVector const& x)
72{
73 y = a * x + b * y;
74}
75
77template <typename MatrixOrVector>
79 MatrixOrVector const& y);
80
82template <typename MatrixOrVector>
83double norm1(MatrixOrVector const& x);
84
86template <typename MatrixOrVector>
87double norm2(MatrixOrVector const& x);
88
90template <typename MatrixOrVector>
91double normMax(MatrixOrVector const& x);
92
93template <typename MatrixOrVector>
95{
96 switch (type)
97 {
99 return norm1(x);
101 return norm2(x);
103 return normMax(x);
104 default:
105 OGS_FATAL("Invalid norm type given.");
106 }
107}
108
109template <typename Matrix>
110void finalizeAssembly(Matrix& /*A*/);
111
112// Matrix and Vector
113
120template <typename Matrix, typename Vector>
121void matMult(Matrix const& A, Vector const& x, Vector& y)
122{
123 assert(&x != &y);
124 y = A * x;
125}
126
133template <typename Matrix, typename Vector>
134void matMultAdd(Matrix const& A, Vector const& v1, Vector const& v2, Vector& v3)
135{
136 assert(&v1 != &v3);
137 v3 = v2 + A * v1;
138}
139
140} // namespace LinAlg
141} // namespace MathLib
142
143// Global PETScMatrix/PETScVector //////////////////////////////////////////
144#ifdef USE_PETSC
145
146namespace MathLib
147{
148
149class PETScMatrix;
150class PETScVector;
151
152namespace LinAlg
153{
154
155// Vector
156
159void setLocalAccessibleVector(PETScVector const& x);
160
161void set(PETScVector& x, PetscScalar const a);
162
163void copy(PETScVector const& x, PETScVector& y);
164
165void scale(PETScVector& x, PetscScalar const a);
166
167// y = a*y + X
168void aypx(PETScVector& y, PetscScalar const a, PETScVector const& x);
169
170// y = a*x + y
171void axpy(PETScVector& y, PetscScalar const a, PETScVector const& x);
172
173// y = a*x + b*y
174void axpby(PETScVector& y, PetscScalar const a, PetscScalar const b,
175 PETScVector const& x);
176
177// Matrix
178
179void copy(PETScMatrix const& A, PETScMatrix& B);
180
181// A = a*A
182void scale(PETScMatrix& A, PetscScalar const a);
183
184// Y = a*Y + X
185void aypx(PETScMatrix& Y, PetscScalar const a, PETScMatrix const& X);
186
187// Y = a*X + Y
188void axpy(PETScMatrix& Y, PetscScalar const a, PETScMatrix const& X);
189
190// Matrix and Vector
191
192// v3 = A*v1 + v2
193void matMult(PETScMatrix const& A, PETScVector const& x, PETScVector& y);
194
195// y = A*x
196void matMultAdd(PETScMatrix const& A, PETScVector const& v1,
197 PETScVector const& v2, PETScVector& v3);
198
199// new_A = A^T * A
200// new_b = A^T * b
201void linearSysNormalize(PETScMatrix const& A, PETScMatrix& new_A,
202 PETScVector const& b, PETScVector& new_b);
203
204void finalizeAssembly(PETScMatrix& A);
205void finalizeAssembly(PETScVector& x);
206
207} // namespace LinAlg
208} // namespace MathLib
209
210// Sparse global EigenMatrix/EigenVector
211// //////////////////////////////////////////
212#else
213
214namespace MathLib
215{
216
217class EigenMatrix;
218class EigenVector;
219
220namespace LinAlg
221{
222
223// Vector
224
233void setLocalAccessibleVector(EigenVector const& x);
234
235void set(EigenVector& x, double const a);
236
237void copy(EigenVector const& x, EigenVector& y);
238
239void scale(EigenVector& x, double const a);
240
241// y = a*y + X
242void aypx(EigenVector& y, double const a, EigenVector const& x);
243
244// y = a*x + y
245void axpy(EigenVector& y, double const a, EigenVector const& x);
246
247// y = a*x + b*y
248void axpby(EigenVector& y, double const a, double const b,
249 EigenVector const& x);
250
251// Matrix
252
253void copy(EigenMatrix const& A, EigenMatrix& B);
254
255// A = a*A
256void scale(EigenMatrix& A, double const a);
257
258// Y = a*Y + X
259void aypx(EigenMatrix& Y, double const a, EigenMatrix const& X);
260
261// Y = a*X + Y
262void axpy(EigenMatrix& Y, double const a, EigenMatrix const& X);
263
264// Matrix and Vector
265
266// y = A*x
267void matMult(EigenMatrix const& A, EigenVector const& x, EigenVector& y);
268
269// v3 = A*v1 + v2
270void matMultAdd(EigenMatrix const& A, EigenVector const& v1,
271 EigenVector const& v2, EigenVector& v3);
272// new_A = A^T * A
273// new_b = A^T * b
274void linearSysNormalize(EigenMatrix const& A, EigenMatrix& new_A,
275 EigenVector const& b, EigenVector& new_b);
276void finalizeAssembly(EigenMatrix& x);
277void finalizeAssembly(EigenVector& A);
278
279} // namespace LinAlg
280
281} // namespace MathLib
282
283#endif
284
285namespace MathLib::LinAlg
286{
287
296template <typename VectorType>
297double computeRelativeNorm(VectorType const& x, VectorType const& y,
299{
301 {
302 OGS_FATAL("An invalid norm type has been passed");
303 }
304
305 // Stores \f$diff = x-y\f$.
306 VectorType diff;
308 MathLib::LinAlg::axpy(diff, -1.0, y);
309
311
312 const double norm_x = MathLib::LinAlg::norm(x, norm_type);
313 if (norm_x > std::numeric_limits<double>::epsilon())
314 {
315 return norm_diff / norm_x;
316 }
317
318 // Both of norm_x and norm_diff are close to zero
319 if (norm_diff < std::numeric_limits<double>::epsilon())
320 {
321 return 1.0;
322 }
323
324 // Only norm_x is close to zero
325 return norm_diff / std::numeric_limits<double>::epsilon();
326}
327} // namespace MathLib::LinAlg
#define OGS_FATAL(...)
Definition Error.h:26
double norm(MatrixOrVector const &x, MathLib::VecNormType type)
Definition LinAlg.h:94
void linearSysNormalize(PETScMatrix const &, PETScMatrix &, PETScVector const &, PETScVector &)
Definition LinAlg.cpp:170
void finalizeAssembly(PETScMatrix &A)
Definition LinAlg.cpp:198
double norm1(PETScVector const &x)
Definition LinAlg.cpp:90
double normMax(PETScVector const &x)
Definition LinAlg.cpp:110
void copy(PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:37
void componentwiseDivide(PETScVector &w, PETScVector const &x, PETScVector const &y)
Definition LinAlg.cpp:81
void setLocalAccessibleVector(PETScVector const &x)
Definition LinAlg.cpp:27
void set(PETScVector &x, PetscScalar const a)
Definition LinAlg.cpp:32
void matMult(PETScMatrix const &A, PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:149
double computeRelativeNorm(VectorType const &x, VectorType const &y, MathLib::VecNormType norm_type)
Definition LinAlg.h:297
void matMultAdd(PETScMatrix const &A, PETScVector const &v1, PETScVector const &v2, PETScVector &v3)
Definition LinAlg.cpp:159
double norm2(PETScVector const &x)
Definition LinAlg.cpp:100
void scale(PETScVector &x, PetscScalar const a)
Definition LinAlg.cpp:44
void aypx(PETScVector &y, PetscScalar const a, PETScVector const &x)
Definition LinAlg.cpp:50
void axpy(PETScVector &y, PetscScalar const a, PETScVector const &x)
Definition LinAlg.cpp:57
void axpby(PETScVector &y, PetscScalar const a, PetscScalar const b, PETScVector const &x)
Definition LinAlg.cpp:64