OGS
LinAlg.cpp
Go to the documentation of this file.
1
11#include "LinAlg.h"
12
13// TODO reorder LinAlg function signatures?
14
15// Global PETScMatrix/PETScVector //////////////////////////////////////////
16#ifdef USE_PETSC
17
20
21namespace MathLib
22{
23namespace LinAlg
24{
25// Vector
26
31
32void set(PETScVector& x, PetscScalar const a)
33{
34 VecSet(x.getRawVector(), a);
35}
36
37void copy(PETScVector const& x, PETScVector& y)
38{
39 if (!y.getRawVector())
40 y.shallowCopy(x);
42}
43
44void scale(PETScVector& x, PetscScalar const a)
45{
46 VecScale(x.getRawVector(), a);
47}
48
49// y = a*y + X
50void aypx(PETScVector& y, PetscScalar const a, PETScVector const& x)
51{
52 // TODO check sizes
54}
55
56// y = a*x + y
57void axpy(PETScVector& y, PetscScalar const a, PETScVector const& x)
58{
59 // TODO check sizes
61}
62
63// y = a*x + b*y
64void axpby(PETScVector& y, PetscScalar const a, PetscScalar const b,
65 PETScVector const& x)
66{
67 // TODO check sizes
68 VecAXPBY(y.getRawVector(), a, b, x.getRawVector());
69}
70
71// Explicit specialization
72// Computes w = x/y componentwise.
73// \note that VecPointwiseDivide avoids to divide by values that are
74// identically zero such as
75// for (int i=0; i<n; i++)
76// {
77// w[i] = y[i] == 0.0 ? 0.0 : x[i] / y[i];
78// }
79//
80template <>
86
87// Explicit specialization
88// Computes the Manhattan norm of x
89template <>
90double norm1(PETScVector const& x)
91{
92 PetscScalar norm = 0.;
94 return norm;
95}
96
97// Explicit specialization
98// Computes the Euclidean norm of x
99template <>
100double norm2(PETScVector const& x)
101{
102 PetscScalar norm = 0.;
104 return norm;
105}
106
107// Explicit specialization
108// Computes the Maximum norm of x
109template <>
110double normMax(PETScVector const& x)
111{
112 PetscScalar norm = 0.;
114 return norm;
115}
116
117// Matrix
118
119void copy(PETScMatrix const& A, PETScMatrix& B)
120{
121 B = A;
122}
123
124// A = a*A
126{
127 MatScale(A.getRawMatrix(), a);
128}
129
130// Y = a*Y + X
131void aypx(PETScMatrix& Y, PetscScalar const a, PETScMatrix const& X)
132{
133 // TODO check sizes
134 // TODO sparsity pattern, currently they are assumed to be different (slow)
136}
137
138// Y = a*X + Y
139void axpy(PETScMatrix& Y, PetscScalar const a, PETScMatrix const& X)
140{
141 // TODO check sizes
142 // TODO sparsity pattern, currently they are assumed to be different (slow)
144}
145
146// Matrix and Vector
147
148// v3 = A*v1 + v2
149void matMult(PETScMatrix const& A, PETScVector const& x, PETScVector& y)
150{
151 // TODO check sizes
152 assert(&x != &y);
153 if (!y.getRawVector())
154 y.shallowCopy(x);
156}
157
158// v3 = A*v1 + v2
159void matMultAdd(PETScMatrix const& A, PETScVector const& v1,
160 PETScVector const& v2, PETScVector& v3)
161{
162 // TODO check sizes
163 assert(&v1 != &v3);
164 if (!v3.getRawVector())
165 v3.shallowCopy(v1);
166 MatMultAdd(A.getRawMatrix(), v1.getRawVector(), v2.getRawVector(),
167 v3.getRawVector());
168}
169
174
179
180} // namespace LinAlg
181} // namespace MathLib
182
183// Sparse global EigenMatrix/EigenVector
184// //////////////////////////////////////////
185#else
186
189
190namespace MathLib
191{
192namespace LinAlg
193{
194// Vector
195
196void setLocalAccessibleVector(EigenVector const& /*x*/) {}
197
198void set(EigenVector& x, double const a)
199{
200 x.getRawVector().setConstant(a);
201}
202
203void copy(EigenVector const& x, EigenVector& y)
204{
205 y = x;
206}
207
208void scale(EigenVector& x, double const a)
209{
210 x.getRawVector() *= a;
211}
212
213// y = a*y + X
214void aypx(EigenVector& y, double const a, EigenVector const& x)
215{
216 // TODO: does that break anything?
217 y.getRawVector() = a * y.getRawVector() + x.getRawVector();
218}
219
220// y = a*x + y
221void axpy(EigenVector& y, double const a, EigenVector const& x)
222{
223 // TODO: does that break anything?
224 y.getRawVector() += a * x.getRawVector();
225}
226
227// y = a*x + b*y
228void axpby(EigenVector& y, double const a, double const b, EigenVector const& x)
229{
230 // TODO: does that break anything?
231 y.getRawVector() = a * x.getRawVector() + b * y.getRawVector();
232}
233
234// Explicit specialization
235// Computes w = x/y componentwise.
236template <>
237void componentwiseDivide(EigenVector& w, EigenVector const& x,
238 EigenVector const& y)
239{
240 w.getRawVector().noalias() = x.getRawVector().binaryExpr(
241 y.getRawVector(),
242 [](auto const x, auto const y) { return y == 0 ? 0.0 : x / y; });
243}
244
245// Explicit specialization
246// Computes the Manhattan norm of x
247template <>
248double norm1(EigenVector const& x)
249{
250 return x.getRawVector().lpNorm<1>();
251}
252
253// Explicit specialization
254// Euclidean norm
255template <>
256double norm2(EigenVector const& x)
257{
258 return x.getRawVector().norm();
259}
260
261// Explicit specialization
262// Computes the Maximum norm of x
263template <>
264double normMax(EigenVector const& x)
265{
266 return x.getRawVector().lpNorm<Eigen::Infinity>();
267}
268
269// Matrix
270
271void copy(EigenMatrix const& A, EigenMatrix& B)
272{
273 B = A;
274}
275
276// A = a*A
277void scale(EigenMatrix& A, double const a)
278{
279 // TODO: does that break anything?
280 A.getRawMatrix() *= a;
281}
282
283// Y = a*Y + X
284void aypx(EigenMatrix& Y, double const a, EigenMatrix const& X)
285{
286 // TODO: does that break anything?
287 Y.getRawMatrix() = a * Y.getRawMatrix() + X.getRawMatrix();
288}
289
290// Y = a*X + Y
291void axpy(EigenMatrix& Y, double const a, EigenMatrix const& X)
292{
293 // TODO: does that break anything?
294 Y.getRawMatrix() = a * X.getRawMatrix() + Y.getRawMatrix();
295}
296
297// Matrix and Vector
298
299// v3 = A*v1 + v2
300void matMult(EigenMatrix const& A, EigenVector const& x, EigenVector& y)
301{
302 assert(&x != &y);
303 y.getRawVector() = A.getRawMatrix() * x.getRawVector();
304}
305
306// v3 = A*v1 + v2
307void matMultAdd(EigenMatrix const& A, EigenVector const& v1,
308 EigenVector const& v2, EigenVector& v3)
309{
310 assert(&v1 != &v3);
311 // TODO: does that break anything?
312 v3.getRawVector() =
313 v2.getRawVector() + A.getRawMatrix() * v1.getRawVector();
314}
315
316void finalizeAssembly(EigenMatrix& x)
317{
318 x.getRawMatrix().makeCompressed();
319}
320
321void finalizeAssembly(EigenVector& /*x*/) {}
322
323} // namespace LinAlg
324
325} // namespace MathLib
326
327#endif
Declaration of class PETScMatrix, which provides an interface to PETSc matrix routines.
Declaration of class PETScVector, which provides an interface to PETSc vector routines.
Wrapper class for PETSc matrix routines for matrix.
Definition PETScMatrix.h:32
Mat & getRawMatrix()
Get matrix reference.
Definition PETScMatrix.h:86
void finalizeAssembly(const MatAssemblyType asm_type=MAT_FINAL_ASSEMBLY)
Perform MPI collection of assembled entries in buffer.
Definition PETScMatrix.h:67
Wrapper class for PETSc vector.
Definition PETScVector.h:33
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:94
void finalizeAssembly(PETScMatrix &A)
Definition LinAlg.cpp:170
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
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