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 
21 namespace MathLib
22 {
23 namespace LinAlg
24 {
25 // Vector
26 
28 {
30 }
31 
32 void set(PETScVector& x, double const a)
33 {
34  VecSet(x.getRawVector(), a);
35 }
36 
37 void copy(PETScVector const& x, PETScVector& y)
38 {
39  if (!y.getRawVector())
40  y.shallowCopy(x);
41  VecCopy(x.getRawVector(), y.getRawVector());
42 }
43 
44 void scale(PETScVector& x, double const a)
45 {
46  VecScale(x.getRawVector(), a);
47 }
48 
49 // y = a*y + X
50 void aypx(PETScVector& y, double const a, PETScVector const& x)
51 {
52  // TODO check sizes
53  VecAYPX(y.getRawVector(), a, x.getRawVector());
54 }
55 
56 // y = a*x + y
57 void axpy(PETScVector& y, double const a, PETScVector const& x)
58 {
59  // TODO check sizes
60  VecAXPY(y.getRawVector(), a, x.getRawVector());
61 }
62 
63 // y = a*x + b*y
64 void axpby(PETScVector& y, double const a, double const b, PETScVector const& x)
65 {
66  // TODO check sizes
67  VecAXPBY(y.getRawVector(), a, b, x.getRawVector());
68 }
69 
70 // Explicit specialization
71 // Computes w = x/y componentwise.
72 template <>
74  PETScVector const& y)
75 {
76  VecPointwiseDivide(w.getRawVector(), x.getRawVector(), y.getRawVector());
77 }
78 
79 // Explicit specialization
80 // Computes the Manhattan norm of x
81 template <>
82 double norm1(PETScVector const& x)
83 {
84  PetscScalar norm = 0.;
85  VecNorm(x.getRawVector(), NORM_1, &norm);
86  return norm;
87 }
88 
89 // Explicit specialization
90 // Computes the Euclidean norm of x
91 template <>
92 double norm2(PETScVector const& x)
93 {
94  PetscScalar norm = 0.;
95  VecNorm(x.getRawVector(), NORM_2, &norm);
96  return norm;
97 }
98 
99 // Explicit specialization
100 // Computes the Maximum norm of x
101 template <>
102 double normMax(PETScVector const& x)
103 {
104  PetscScalar norm = 0.;
105  VecNorm(x.getRawVector(), NORM_INFINITY, &norm);
106  return norm;
107 }
108 
109 // Matrix
110 
111 void copy(PETScMatrix const& A, PETScMatrix& B)
112 {
113  B = A;
114 }
115 
116 // A = a*A
117 void scale(PETScMatrix& A, double const a)
118 {
119  MatScale(A.getRawMatrix(), a);
120 }
121 
122 // Y = a*Y + X
123 void aypx(PETScMatrix& Y, double const a, PETScMatrix const& X)
124 {
125  // TODO check sizes
126  // TODO sparsity pattern, currently they are assumed to be different (slow)
127  MatAYPX(Y.getRawMatrix(), a, X.getRawMatrix(), DIFFERENT_NONZERO_PATTERN);
128 }
129 
130 // Y = a*X + Y
131 void axpy(PETScMatrix& Y, double const a, PETScMatrix const& X)
132 {
133  // TODO check sizes
134  // TODO sparsity pattern, currently they are assumed to be different (slow)
135  MatAXPY(Y.getRawMatrix(), a, X.getRawMatrix(), DIFFERENT_NONZERO_PATTERN);
136 }
137 
138 // Matrix and Vector
139 
140 // v3 = A*v1 + v2
141 void matMult(PETScMatrix const& A, PETScVector const& x, PETScVector& y)
142 {
143  // TODO check sizes
144  assert(&x != &y);
145  if (!y.getRawVector())
146  y.shallowCopy(x);
147  MatMult(A.getRawMatrix(), x.getRawVector(), y.getRawVector());
148 }
149 
150 // v3 = A*v1 + v2
151 void matMultAdd(PETScMatrix const& A, PETScVector const& v1,
152  PETScVector const& v2, PETScVector& v3)
153 {
154  // TODO check sizes
155  assert(&v1 != &v3);
156  if (!v3.getRawVector())
157  v3.shallowCopy(v1);
158  MatMultAdd(A.getRawMatrix(), v1.getRawVector(), v2.getRawVector(),
159  v3.getRawVector());
160 }
161 
163 {
164  A.finalizeAssembly(MAT_FINAL_ASSEMBLY);
165 }
166 
168 {
169  x.finalizeAssembly();
170 }
171 
172 } // namespace LinAlg
173 } // namespace MathLib
174 
175 // Sparse global EigenMatrix/EigenVector
176 // //////////////////////////////////////////
177 #else
178 
181 
182 namespace MathLib
183 {
184 namespace LinAlg
185 {
186 // Vector
187 
188 void setLocalAccessibleVector(EigenVector const& /*x*/) {}
189 
190 void set(EigenVector& x, double const a)
191 {
192  x.getRawVector().setConstant(a);
193 }
194 
195 void copy(EigenVector const& x, EigenVector& y)
196 {
197  y = x;
198 }
199 
200 void scale(EigenVector& x, double const a)
201 {
202  x.getRawVector() *= a;
203 }
204 
205 // y = a*y + X
206 void aypx(EigenVector& y, double const a, EigenVector const& x)
207 {
208  // TODO: does that break anything?
209  y.getRawVector() = a * y.getRawVector() + x.getRawVector();
210 }
211 
212 // y = a*x + y
213 void axpy(EigenVector& y, double const a, EigenVector const& x)
214 {
215  // TODO: does that break anything?
216  y.getRawVector() += a * x.getRawVector();
217 }
218 
219 // y = a*x + y
220 void axpby(EigenVector& y, double const a, double const b, EigenVector const& x)
221 {
222  // TODO: does that break anything?
223  y.getRawVector() = a * x.getRawVector() + b * y.getRawVector();
224 }
225 
226 // Explicit specialization
227 // Computes w = x/y componentwise.
228 template <>
229 void componentwiseDivide(EigenVector& w, EigenVector const& x,
230  EigenVector const& y)
231 {
232  w.getRawVector().noalias() =
233  x.getRawVector().cwiseQuotient(y.getRawVector());
234 }
235 
236 // Explicit specialization
237 // Computes the Manhattan norm of x
238 template <>
239 double norm1(EigenVector const& x)
240 {
241  return x.getRawVector().lpNorm<1>();
242 }
243 
244 // Explicit specialization
245 // Euclidean norm
246 template <>
247 double norm2(EigenVector const& x)
248 {
249  return x.getRawVector().norm();
250 }
251 
252 // Explicit specialization
253 // Computes the Maximum norm of x
254 template <>
255 double normMax(EigenVector const& x)
256 {
257  return x.getRawVector().lpNorm<Eigen::Infinity>();
258 }
259 
260 // Matrix
261 
262 void copy(EigenMatrix const& A, EigenMatrix& B)
263 {
264  B = A;
265 }
266 
267 // A = a*A
268 void scale(EigenMatrix& A, double const a)
269 {
270  // TODO: does that break anything?
271  A.getRawMatrix() *= a;
272 }
273 
274 // Y = a*Y + X
275 void aypx(EigenMatrix& Y, double const a, EigenMatrix const& X)
276 {
277  // TODO: does that break anything?
278  Y.getRawMatrix() = a * Y.getRawMatrix() + X.getRawMatrix();
279 }
280 
281 // Y = a*X + Y
282 void axpy(EigenMatrix& Y, double const a, EigenMatrix const& X)
283 {
284  // TODO: does that break anything?
285  Y.getRawMatrix() = a * X.getRawMatrix() + Y.getRawMatrix();
286 }
287 
288 // Matrix and Vector
289 
290 // v3 = A*v1 + v2
291 void matMult(EigenMatrix const& A, EigenVector const& x, EigenVector& y)
292 {
293  assert(&x != &y);
294  y.getRawVector() = A.getRawMatrix() * x.getRawVector();
295 }
296 
297 // v3 = A*v1 + v2
298 void matMultAdd(EigenMatrix const& A, EigenVector const& v1,
299  EigenVector const& v2, EigenVector& v3)
300 {
301  assert(&v1 != &v3);
302  // TODO: does that break anything?
303  v3.getRawVector() =
304  v2.getRawVector() + A.getRawMatrix() * v1.getRawVector();
305 }
306 
307 void finalizeAssembly(EigenMatrix& x)
308 {
309  x.getRawMatrix().makeCompressed();
310 }
311 
312 void finalizeAssembly(EigenVector& /*x*/) {}
313 
314 } // namespace LinAlg
315 
316 } // namespace MathLib
317 
318 #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.
Definition: PETScVector.h:171
double norm(MatrixOrVector const &x, MathLib::VecNormType type)
Definition: LinAlg.h:88
void axpy(PETScVector &y, double const a, PETScVector const &x)
Definition: LinAlg.cpp:57
void finalizeAssembly(PETScMatrix &A)
Definition: LinAlg.cpp:162
void scale(PETScVector &x, double const a)
Definition: LinAlg.cpp:44
double norm1(PETScVector const &x)
Definition: LinAlg.cpp:82
double normMax(PETScVector const &x)
Definition: LinAlg.cpp:102
void copy(PETScVector const &x, PETScVector &y)
Definition: LinAlg.cpp:37
void componentwiseDivide(PETScVector &w, PETScVector const &x, PETScVector const &y)
Definition: LinAlg.cpp:73
void setLocalAccessibleVector(PETScVector const &x)
Definition: LinAlg.cpp:27
void matMult(PETScMatrix const &A, PETScVector const &x, PETScVector &y)
Definition: LinAlg.cpp:141
void set(PETScVector &x, double const a)
Definition: LinAlg.cpp:32
void matMultAdd(PETScMatrix const &A, PETScVector const &v1, PETScVector const &v2, PETScVector &v3)
Definition: LinAlg.cpp:151
double norm2(PETScVector const &x)
Definition: LinAlg.cpp:92
void axpby(PETScVector &y, double const a, double const b, PETScVector const &x)
Definition: LinAlg.cpp:64
void aypx(PETScVector &y, double const a, PETScVector const &x)
Definition: LinAlg.cpp:50