OGS
LinAlg.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#pragma once
5
6#include <cassert>
7
8#include "BaseLib/Error.h"
9#include "LinAlgEnums.h"
10
11#ifdef USE_PETSC
12#include <petscsystypes.h>
13
14#include "PETSc/PETScVector.h"
15#endif
16
17namespace MathLib
18{
19
30namespace LinAlg
31{
32
33// Matrix or Vector
34
36template <typename MatrixOrVector>
37void copy(MatrixOrVector const& x, MatrixOrVector& y)
38{
39 y = x;
40}
41
43template <typename MatrixOrVector>
44void scale(MatrixOrVector& x, double const a)
45{
46 x *= a;
47}
48
50template <typename MatrixOrVector>
51void aypx(MatrixOrVector& y, double const a, MatrixOrVector const& x)
52{
53 y = a * y + x;
54}
55
57template <typename MatrixOrVector>
58void axpy(MatrixOrVector& y, double const a, MatrixOrVector const& x)
59{
60 y += a * x;
61}
62
64template <typename MatrixOrVector>
65void axpby(MatrixOrVector& y, double const a, double const b,
66 MatrixOrVector const& x)
67{
68 y = a * x + b * y;
69}
70
72template <typename MatrixOrVector>
73void componentwiseDivide(MatrixOrVector& w, MatrixOrVector const& x,
74 MatrixOrVector const& y);
75
77template <typename MatrixOrVector>
78double norm1(MatrixOrVector const& x);
79
81template <typename MatrixOrVector>
82double norm2(MatrixOrVector const& x);
83
85template <typename MatrixOrVector>
86double normMax(MatrixOrVector const& x);
87
88template <typename MatrixOrVector>
89double norm(MatrixOrVector const& x, MathLib::VecNormType type)
90{
91 switch (type)
92 {
94 return norm1(x);
96 return norm2(x);
98 return normMax(x);
99 default:
100 OGS_FATAL("Invalid norm type given.");
101 }
102}
103
104template <typename Matrix>
105void finalizeAssembly(Matrix& /*A*/);
106
107// Matrix and Vector
108
115template <typename Matrix, typename Vector>
116void matMult(Matrix const& A, Vector const& x, Vector& y)
117{
118 assert(&x != &y);
119 y = A * x;
120}
121
128template <typename Matrix, typename Vector>
129void matMultAdd(Matrix const& A, Vector const& v1, Vector const& v2, Vector& v3)
130{
131 assert(&v1 != &v3);
132 v3 = v2 + A * v1;
133}
134
135} // namespace LinAlg
136} // namespace MathLib
137
138// Global PETScMatrix/PETScVector //////////////////////////////////////////
139#ifdef USE_PETSC
140
141namespace MathLib
142{
143
144class PETScMatrix;
145class PETScVector;
146
147namespace LinAlg
148{
149
150// Vector
151
154void setLocalAccessibleVector(PETScVector const& x);
155
156void set(PETScVector& x, PetscScalar const a);
157
158void copy(PETScVector const& x, PETScVector& y);
159
160void scale(PETScVector& x, PetscScalar const a);
161
162// y = a*y + X
163void aypx(PETScVector& y, PetscScalar const a, PETScVector const& x);
164
165// y = a*x + y
166void axpy(PETScVector& y, PetscScalar const a, PETScVector const& x);
167
168// y = a*x + b*y
169void axpby(PETScVector& y, PetscScalar const a, PetscScalar const b,
170 PETScVector const& x);
171
172// Matrix
173
174void copy(PETScMatrix const& A, PETScMatrix& B);
175
176// A = a*A
177void scale(PETScMatrix& A, PetscScalar const a);
178
179// Y = a*Y + X
180void aypx(PETScMatrix& Y, PetscScalar const a, PETScMatrix const& X);
181
182// Y = a*X + Y
183void axpy(PETScMatrix& Y, PetscScalar const a, PETScMatrix const& X);
184
185// Matrix and Vector
186
187// v3 = A*v1 + v2
188void matMult(PETScMatrix const& A, PETScVector const& x, PETScVector& y);
189
190// y = A*x
191void matMultAdd(PETScMatrix const& A, PETScVector const& v1,
192 PETScVector const& v2, PETScVector& v3);
193
194// new_A = A^T * A
195// new_b = A^T * b
196void linearSysNormalize(PETScMatrix const& A, PETScMatrix& new_A,
197 PETScVector const& b, PETScVector& new_b);
198
199void finalizeAssembly(PETScMatrix& A);
200void finalizeAssembly(PETScVector& x);
201
202} // namespace LinAlg
203} // namespace MathLib
204
205// Sparse global EigenMatrix/EigenVector
206// //////////////////////////////////////////
207#else
208
209namespace MathLib
210{
211
212class EigenMatrix;
213class EigenVector;
214
215namespace LinAlg
216{
217
218// Vector
219
228void setLocalAccessibleVector(EigenVector const& x);
229
230void set(EigenVector& x, double const a);
231
232void copy(EigenVector const& x, EigenVector& y);
233
234void scale(EigenVector& x, double const a);
235
236// y = a*y + X
237void aypx(EigenVector& y, double const a, EigenVector const& x);
238
239// y = a*x + y
240void axpy(EigenVector& y, double const a, EigenVector const& x);
241
242// y = a*x + b*y
243void axpby(EigenVector& y, double const a, double const b,
244 EigenVector const& x);
245
246// Matrix
247
248void copy(EigenMatrix const& A, EigenMatrix& B);
249
250// A = a*A
251void scale(EigenMatrix& A, double const a);
252
253// Y = a*Y + X
254void aypx(EigenMatrix& Y, double const a, EigenMatrix const& X);
255
256// Y = a*X + Y
257void axpy(EigenMatrix& Y, double const a, EigenMatrix const& X);
258
259// Matrix and Vector
260
261// y = A*x
262void matMult(EigenMatrix const& A, EigenVector const& x, EigenVector& y);
263
264// v3 = A*v1 + v2
265void matMultAdd(EigenMatrix const& A, EigenVector const& v1,
266 EigenVector const& v2, EigenVector& v3);
267// new_A = A^T * A
268// new_b = A^T * b
269void linearSysNormalize(EigenMatrix const& A, EigenMatrix& new_A,
270 EigenVector const& b, EigenVector& new_b);
271void finalizeAssembly(EigenMatrix& x);
272void finalizeAssembly(EigenVector& A);
273
274} // namespace LinAlg
275
276} // namespace MathLib
277
278#endif
279
280namespace MathLib::LinAlg
281{
282
291template <typename VectorType>
292double computeRelativeNorm(VectorType const& x, VectorType const& y,
293 MathLib::VecNormType norm_type)
294{
295 if (norm_type == MathLib::VecNormType::INVALID)
296 {
297 OGS_FATAL("An invalid norm type has been passed");
298 }
299
300 // Stores \f$diff = x-y\f$.
301 VectorType diff;
302 MathLib::LinAlg::copy(x, diff);
303 MathLib::LinAlg::axpy(diff, -1.0, y);
304
305 const double norm_diff = MathLib::LinAlg::norm(diff, norm_type);
306
307 const double norm_x = MathLib::LinAlg::norm(x, norm_type);
308 if (norm_x > std::numeric_limits<double>::epsilon())
309 {
310 return norm_diff / norm_x;
311 }
312
313 // Both of norm_x and norm_diff are close to zero
314 if (norm_diff < std::numeric_limits<double>::epsilon())
315 {
316 return 1.0;
317 }
318
319 // Only norm_x is close to zero
320 return norm_diff / std::numeric_limits<double>::epsilon();
321}
322} // namespace MathLib::LinAlg
#define OGS_FATAL(...)
Definition Error.h:19
Global vector based on Eigen vector.
Definition EigenVector.h:19
Wrapper class for PETSc matrix routines for matrix.
Definition PETScMatrix.h:21
Wrapper class for PETSc vector.
Definition PETScVector.h:28
double norm(MatrixOrVector const &x, MathLib::VecNormType type)
Definition LinAlg.h:89
void linearSysNormalize(PETScMatrix const &, PETScMatrix &, PETScVector const &, PETScVector &)
Definition LinAlg.cpp:163
void finalizeAssembly(PETScMatrix &A)
Definition LinAlg.cpp:191
double norm1(PETScVector const &x)
Definition LinAlg.cpp:83
double normMax(PETScVector const &x)
Definition LinAlg.cpp:103
void copy(PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:30
void componentwiseDivide(PETScVector &w, PETScVector const &x, PETScVector const &y)
Definition LinAlg.cpp:74
void setLocalAccessibleVector(PETScVector const &x)
Definition LinAlg.cpp:20
void set(PETScVector &x, PetscScalar const a)
Definition LinAlg.cpp:25
void matMult(PETScMatrix const &A, PETScVector const &x, PETScVector &y)
Definition LinAlg.cpp:142
double computeRelativeNorm(VectorType const &x, VectorType const &y, MathLib::VecNormType norm_type)
Definition LinAlg.h:292
void matMultAdd(PETScMatrix const &A, PETScVector const &v1, PETScVector const &v2, PETScVector &v3)
Definition LinAlg.cpp:152
double norm2(PETScVector const &x)
Definition LinAlg.cpp:93
void scale(PETScVector &x, PetscScalar const a)
Definition LinAlg.cpp:37
void aypx(PETScVector &y, PetscScalar const a, PETScVector const &x)
Definition LinAlg.cpp:43
void axpy(PETScVector &y, PetscScalar const a, PETScVector const &x)
Definition LinAlg.cpp:50
void axpby(PETScVector &y, PetscScalar const a, PetscScalar const b, PETScVector const &x)
Definition LinAlg.cpp:57
double const GaussLegendre< 1 >::X[1]