OGS
LinAlg.cpp
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#include "LinAlg.h"
5
6// TODO reorder LinAlg function signatures?
7
8// Global PETScMatrix/PETScVector //////////////////////////////////////////
9#ifdef USE_PETSC
10
13
14namespace MathLib
15{
16namespace LinAlg
17{
18// Vector
19
24
25void set(PETScVector& x, PetscScalar const a)
26{
27 VecSet(x.getRawVector(), a);
28}
29
30void copy(PETScVector const& x, PETScVector& y)
31{
32 if (!y.getRawVector())
33 y.shallowCopy(x);
34 VecCopy(x.getRawVector(), y.getRawVector());
35}
36
37void scale(PETScVector& x, PetscScalar const a)
38{
39 VecScale(x.getRawVector(), a);
40}
41
42// y = a*y + X
43void aypx(PETScVector& y, PetscScalar const a, PETScVector const& x)
44{
45 // TODO check sizes
46 VecAYPX(y.getRawVector(), a, x.getRawVector());
47}
48
49// y = a*x + y
50void axpy(PETScVector& y, PetscScalar const a, PETScVector const& x)
51{
52 // TODO check sizes
53 VecAXPY(y.getRawVector(), a, x.getRawVector());
54}
55
56// y = a*x + b*y
57void axpby(PETScVector& y, PetscScalar const a, PetscScalar const b,
58 PETScVector const& x)
59{
60 // TODO check sizes
61 VecAXPBY(y.getRawVector(), a, b, x.getRawVector());
62}
63
64// Explicit specialization
65// Computes w = x/y componentwise.
66// \note that VecPointwiseDivide avoids to divide by values that are
67// identically zero such as
68// for (int i=0; i<n; i++)
69// {
70// w[i] = y[i] == 0.0 ? 0.0 : x[i] / y[i];
71// }
72//
73template <>
75 PETScVector const& y)
76{
77 VecPointwiseDivide(w.getRawVector(), x.getRawVector(), y.getRawVector());
78}
79
80// Explicit specialization
81// Computes the Manhattan norm of x
82template <>
83double norm1(PETScVector const& x)
84{
85 PetscScalar norm = 0.;
86 VecNorm(x.getRawVector(), NORM_1, &norm);
87 return norm;
88}
89
90// Explicit specialization
91// Computes the Euclidean norm of x
92template <>
93double norm2(PETScVector const& x)
94{
95 PetscScalar norm = 0.;
96 VecNorm(x.getRawVector(), NORM_2, &norm);
97 return norm;
98}
99
100// Explicit specialization
101// Computes the Maximum norm of x
102template <>
103double normMax(PETScVector const& x)
104{
105 PetscScalar norm = 0.;
106 VecNorm(x.getRawVector(), NORM_INFINITY, &norm);
107 return norm;
108}
109
110// Matrix
111
112void copy(PETScMatrix const& A, PETScMatrix& B)
113{
114 B = A;
115}
116
117// A = a*A
118void scale(PETScMatrix& A, PetscScalar const a)
119{
120 MatScale(A.getRawMatrix(), a);
121}
122
123// Y = a*Y + X
124void aypx(PETScMatrix& Y, PetscScalar const a, PETScMatrix const& X)
125{
126 // TODO check sizes
127 // TODO sparsity pattern, currently they are assumed to be different (slow)
128 MatAYPX(Y.getRawMatrix(), a, X.getRawMatrix(), DIFFERENT_NONZERO_PATTERN);
129}
130
131// Y = a*X + Y
132void axpy(PETScMatrix& Y, PetscScalar const a, PETScMatrix const& X)
133{
134 // TODO check sizes
135 // TODO sparsity pattern, currently they are assumed to be different (slow)
136 MatAXPY(Y.getRawMatrix(), a, X.getRawMatrix(), DIFFERENT_NONZERO_PATTERN);
137}
138
139// Matrix and Vector
140
141// v3 = A*v1 + v2
142void matMult(PETScMatrix const& A, PETScVector const& x, PETScVector& y)
143{
144 // TODO check sizes
145 assert(&x != &y);
146 if (!y.getRawVector())
147 y.shallowCopy(x);
148 MatMult(A.getRawMatrix(), x.getRawVector(), y.getRawVector());
149}
150
151// v3 = A*v1 + v2
152void matMultAdd(PETScMatrix const& A, PETScVector const& v1,
153 PETScVector const& v2, PETScVector& v3)
154{
155 // TODO check sizes
156 assert(&v1 != &v3);
157 if (!v3.getRawVector())
158 v3.shallowCopy(v1);
159 MatMultAdd(A.getRawMatrix(), v1.getRawVector(), v2.getRawVector(),
160 v3.getRawVector());
161}
162
163void linearSysNormalize(PETScMatrix const& /*A*/, PETScMatrix& /*new_A*/,
164 PETScVector const& /*b*/, PETScVector& /*new_b*/)
165{
166 // The following block is deactivated, because there is no tests yet for the
167 // normalization operation in PETSc. This will be a task for later.
168 /*
169 assert(&A != &new_A);
170 assert(&b != &new_b);
171
172 PetscInt n_rows(0);
173 PetscInt n_cols(0);
174 MatGetSize(A.getRawMatrix(), &n_rows, &n_cols);
175 // only when A matrix is not square
176 if (n_rows != n_cols)
177 {
178 // new_b = A^T * b
179 MatMultTranspose(A.getRawMatrix(), b.getRawVector(),
180 new_b.getRawVector());
181 // new_A = A^T * A
182 MatTranspose(A.getRawMatrix(), MAT_INITIAL_MATRIX,
183 &(new_A.getRawMatrix()));
184 }
185 */
186 OGS_FATAL(
187 "Normalization operation is not implemented yet for PETSc library! "
188 "Program terminated.");
189}
190
192{
193 A.finalizeAssembly(MAT_FINAL_ASSEMBLY);
194}
195
200
201} // namespace LinAlg
202} // namespace MathLib
203
204// Sparse global EigenMatrix/EigenVector
205// //////////////////////////////////////////
206#else
207
210
211namespace MathLib
212{
213namespace LinAlg
214{
215// Vector
216
217void setLocalAccessibleVector(EigenVector const& /*x*/) {}
218
219void set(EigenVector& x, double const a)
220{
221 x.getRawVector().setConstant(a);
222}
223
224void copy(EigenVector const& x, EigenVector& y)
225{
226 y = x;
227}
228
229void scale(EigenVector& x, double const a)
230{
231 x.getRawVector() *= a;
232}
233
234// y = a*y + X
235void aypx(EigenVector& y, double const a, EigenVector const& x)
236{
237 // TODO: does that break anything?
238 y.getRawVector() = a * y.getRawVector() + x.getRawVector();
239}
240
241// y = a*x + y
242void axpy(EigenVector& y, double const a, EigenVector const& x)
243{
244 // TODO: does that break anything?
245 y.getRawVector() += a * x.getRawVector();
246}
247
248// y = a*x + b*y
249void axpby(EigenVector& y, double const a, double const b, EigenVector const& x)
250{
251 // TODO: does that break anything?
252 y.getRawVector() = a * x.getRawVector() + b * y.getRawVector();
253}
254
255// Explicit specialization
256// Computes w = x/y componentwise.
257template <>
258void componentwiseDivide(EigenVector& w, EigenVector const& x,
259 EigenVector const& y)
260{
261 w.getRawVector().noalias() = x.getRawVector().binaryExpr(
262 y.getRawVector(),
263 [](auto const x, auto const y) { return y == 0 ? 0.0 : x / y; });
264}
265
266// Explicit specialization
267// Computes the Manhattan norm of x
268template <>
269double norm1(EigenVector const& x)
270{
271 return x.getRawVector().lpNorm<1>();
272}
273
274// Explicit specialization
275// Euclidean norm
276template <>
277double norm2(EigenVector const& x)
278{
279 return x.getRawVector().norm();
280}
281
282// Explicit specialization
283// Computes the Maximum norm of x
284template <>
285double normMax(EigenVector const& x)
286{
287 return x.getRawVector().lpNorm<Eigen::Infinity>();
288}
289
290// Matrix
291
292void copy(EigenMatrix const& A, EigenMatrix& B)
293{
294 B = A;
295}
296
297// A = a*A
298void scale(EigenMatrix& A, double const a)
299{
300 // TODO: does that break anything?
301 A.getRawMatrix() *= a;
302}
303
304// Y = a*Y + X
305void aypx(EigenMatrix& Y, double const a, EigenMatrix const& X)
306{
307 // TODO: does that break anything?
308 Y.getRawMatrix() = a * Y.getRawMatrix() + X.getRawMatrix();
309}
310
311// Y = a*X + Y
312void axpy(EigenMatrix& Y, double const a, EigenMatrix const& X)
313{
314 // TODO: does that break anything?
315 Y.getRawMatrix() = a * X.getRawMatrix() + Y.getRawMatrix();
316}
317
318// Matrix and Vector
319
320// v3 = A*v1 + v2
321void matMult(EigenMatrix const& A, EigenVector const& x, EigenVector& y)
322{
323 assert(&x != &y);
324 y.getRawVector() = A.getRawMatrix() * x.getRawVector();
325}
326
327// v3 = A*v1 + v2
328void matMultAdd(EigenMatrix const& A, EigenVector const& v1,
329 EigenVector const& v2, EigenVector& v3)
330{
331 assert(&v1 != &v3);
332 // TODO: does that break anything?
333 v3.getRawVector() =
334 v2.getRawVector() + A.getRawMatrix() * v1.getRawVector();
335}
336
337void linearSysNormalize(EigenMatrix const& A, EigenMatrix& new_A,
338 EigenVector const& b, EigenVector& new_b)
339{
340 // make sure that new_A and new_b are not the same memory
341 assert(&A != &new_A);
342 assert(&b != &new_b);
343
344 if (A.getRawMatrix().rows() == A.getRawMatrix().cols())
345 {
346 WARN(
347 "The number of rows and columns are the same for the LHS matrix."
348 "Are you sure you still need to normalize the LHS matrix and RHS "
349 "vector? ");
350 }
351
352 new_b.getRawVector() = A.getRawMatrix().transpose() * b.getRawVector();
353 new_A.getRawMatrix() = A.getRawMatrix().transpose() * A.getRawMatrix();
354}
355
356void finalizeAssembly(EigenMatrix& x)
357{
358 x.getRawMatrix().makeCompressed();
359}
360
361void finalizeAssembly(EigenVector& /*x*/) {}
362
363} // namespace LinAlg
364
365} // namespace MathLib
366
367#endif
#define OGS_FATAL(...)
Definition Error.h:19
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:34
Wrapper class for PETSc matrix routines for matrix.
Definition PETScMatrix.h:21
Mat & getRawMatrix()
Get matrix reference.
Definition PETScMatrix.h:75
void finalizeAssembly(const MatAssemblyType asm_type=MAT_FINAL_ASSEMBLY)
Perform MPI collection of assembled entries in buffer.
Definition PETScMatrix.h:56
Wrapper class for PETSc vector.
Definition PETScVector.h:28
void finalizeAssembly()
Perform MPI collection of assembled entries in buffer.
void setLocalAccessibleVector() const
void shallowCopy(const PETScVector &v)
PETSc_Vec & getRawVector()
Exposes the underlying PETSc vector.
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
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]