Dune Core Modules (unstable)

matrixutils.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file LICENSE.md in module root
2// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
4// vi: set et ts=4 sw=2 sts=2:
5#ifndef DUNE_ISTL_MATRIXUTILS_HH
6#define DUNE_ISTL_MATRIXUTILS_HH
7
8#include <set>
9#include <vector>
10#include <limits>
17#include "istlexception.hh"
18
19namespace Dune
20{
21
22#ifndef DOYXGEN
23 template<typename B, typename A>
24 class BCRSMatrix;
25
26 template<typename K, int n, int m>
27 class FieldMatrix;
28
29 template<class T, class A>
30 class Matrix;
31#endif
32
46 template<class Matrix, std::size_t blocklevel, std::size_t l=blocklevel>
48 {
53 static void check([[maybe_unused]] const Matrix& mat)
54 {
55#ifdef DUNE_ISTL_WITH_CHECKING
56 typedef typename Matrix::ConstRowIterator Row;
57 typedef typename Matrix::ConstColIterator Entry;
58 for(Row row = mat.begin(); row!=mat.end(); ++row) {
59 Entry diagonal = row->find(row.index());
60 if(diagonal==row->end())
61 DUNE_THROW(ISTLError, "Missing diagonal value in row "<<row.index()
62 <<" at block recursion level "<<l-blocklevel);
63 else{
64 auto m = Impl::asMatrix(*diagonal);
65 CheckIfDiagonalPresent<decltype(m),blocklevel-1,l>::check(m);
66 }
67 }
68#endif
69 }
70 };
71
72 template<class Matrix, std::size_t l>
73 struct CheckIfDiagonalPresent<Matrix,0,l>
74 {
75 static void check(const Matrix& mat)
76 {
77 typedef typename Matrix::ConstRowIterator Row;
78 for(Row row = mat.begin(); row!=mat.end(); ++row) {
79 if(row->find(row.index())==row->end())
80 DUNE_THROW(ISTLError, "Missing diagonal value in row "<<row.index()
81 <<" at block recursion level "<<l);
82 }
83 }
84 };
85
86 template<typename FirstRow, typename... Args>
87 class MultiTypeBlockMatrix;
88
89 template<std::size_t blocklevel, std::size_t l, typename T1, typename... Args>
90 struct CheckIfDiagonalPresent<MultiTypeBlockMatrix<T1,Args...>,
91 blocklevel,l>
92 {
93 typedef MultiTypeBlockMatrix<T1,Args...> Matrix;
94
99 static void check(const Matrix& /* mat */)
100 {
101#ifdef DUNE_ISTL_WITH_CHECKING
102 // TODO Implement check
103#endif
104 }
105 };
106
118 template<class M>
119 inline auto countNonZeros(const M&,
120 [[maybe_unused]] typename std::enable_if_t<Dune::IsNumber<M>::value>* sfinae = nullptr)
121 {
122 return 1;
123 }
124
125 template<class M>
126 inline auto countNonZeros(const M& matrix,
127 [[maybe_unused]] typename std::enable_if_t<!Dune::IsNumber<M>::value>* sfinae = nullptr)
128 {
129 typename M::size_type nonZeros = 0;
130 for(auto&& row : matrix)
131 for(auto&& entry : row)
132 nonZeros += countNonZeros(entry);
133 return nonZeros;
134 }
135
136 /*
137 template<class M>
138 struct ProcessOnFieldsOfMatrix
139 */
140
142 namespace
143 {
144 struct CompPair {
145 template<class G,class M>
146 bool operator()(const std::pair<G,M>& p1, const std::pair<G,M>& p2) const
147 {
148 return p1.first<p2.first;
149 }
150 };
151
152 }
153 template<class M, class C>
154 void printGlobalSparseMatrix(const M& mat, C& ooc, std::ostream& os)
155 {
156 typedef typename C::ParallelIndexSet::const_iterator IIter;
157 typedef typename C::OwnerSet OwnerSet;
158 typedef typename C::ParallelIndexSet::GlobalIndex GlobalIndex;
159
160 GlobalIndex gmax=0;
161
162 for(IIter idx=ooc.indexSet().begin(), eidx=ooc.indexSet().end();
163 idx!=eidx; ++idx)
164 gmax=std::max(gmax,idx->global());
165
166 gmax=ooc.communicator().max(gmax);
167 ooc.buildGlobalLookup();
168
169 for(IIter idx=ooc.indexSet().begin(), eidx=ooc.indexSet().end();
170 idx!=eidx; ++idx) {
171 if(OwnerSet::contains(idx->local().attribute()))
172 {
173 typedef typename M::block_type Block;
174
175 std::set<std::pair<GlobalIndex,Block>,CompPair> entries;
176
177 // sort rows
178 typedef typename M::ConstColIterator CIter;
179 for(CIter c=mat[idx->local()].begin(), cend=mat[idx->local()].end();
180 c!=cend; ++c) {
181 const typename C::ParallelIndexSet::IndexPair* pair
182 =ooc.globalLookup().pair(c.index());
183 assert(pair);
184 entries.insert(std::make_pair(pair->global(), *c));
185 }
186
187 //wait until its the rows turn.
188 GlobalIndex rowidx = idx->global();
190 while(cur!=rowidx)
191 cur=ooc.communicator().min(rowidx);
192
193 // print rows
194 typedef typename std::set<std::pair<GlobalIndex,Block>,CompPair>::iterator SIter;
195 for(SIter s=entries.begin(), send=entries.end(); s!=send; ++s)
196 os<<idx->global()<<" "<<s->first<<" "<<s->second<<std::endl;
197
198
199 }
200 }
201
202 ooc.freeGlobalLookup();
203 // Wait until everybody is finished
205 while(cur!=ooc.communicator().min(cur)) ;
206 }
207
208 // Default implementation for scalar types
209 template<typename M>
210 struct MatrixDimension
211 {
212 static_assert(IsNumber<M>::value, "MatrixDimension is not implemented for this type!");
213
214 static auto rowdim(const M& A)
215 {
216 return 1;
217 }
218
219 static auto coldim(const M& A)
220 {
221 return 1;
222 }
223 };
224
225 // Default implementation for scalar types
226 template<typename B, typename TA>
227 struct MatrixDimension<Matrix<B,TA> >
228 {
229 using block_type = typename Matrix<B,TA>::block_type;
230 using size_type = typename Matrix<B,TA>::size_type;
231
232 static size_type rowdim (const Matrix<B,TA>& A, size_type i)
233 {
234 return MatrixDimension<block_type>::rowdim(A[i][0]);
235 }
236
237 static size_type coldim (const Matrix<B,TA>& A, size_type c)
238 {
239 return MatrixDimension<block_type>::coldim(A[0][c]);
240 }
241
242 static size_type rowdim (const Matrix<B,TA>& A)
243 {
244 size_type nn=0;
245 for (size_type i=0; i<A.N(); i++)
246 nn += rowdim(A,i);
247 return nn;
248 }
249
250 static size_type coldim (const Matrix<B,TA>& A)
251 {
252 size_type nn=0;
253 for (size_type i=0; i<A.M(); i++)
254 nn += coldim(A,i);
255 return nn;
256 }
257 };
258
259
260 template<typename B, typename TA>
261 struct MatrixDimension<BCRSMatrix<B,TA> >
262 {
263 typedef BCRSMatrix<B,TA> Matrix;
264 typedef typename Matrix::block_type block_type;
265 typedef typename Matrix::size_type size_type;
266
267 static size_type rowdim (const Matrix& A, size_type i)
268 {
269 const B* row = A.r[i].getptr();
270 if(row)
271 return MatrixDimension<block_type>::rowdim(*row);
272 else
273 return 0;
274 }
275
276 static size_type coldim (const Matrix& A, size_type c)
277 {
278 // find an entry in column c
279 if (A.nnz_ > 0)
280 {
281 for (size_type k=0; k<A.nnz_; k++) {
282 if (A.j_.get()[k] == c) {
283 return MatrixDimension<block_type>::coldim(A.a[k]);
284 }
285 }
286 }
287 else
288 {
289 for (size_type i=0; i<A.N(); i++)
290 {
291 size_type* j = A.r[i].getindexptr();
292 B* a = A.r[i].getptr();
293 for (size_type k=0; k<A.r[i].getsize(); k++)
294 if (j[k]==c) {
295 return MatrixDimension<block_type>::coldim(a[k]);
296 }
297 }
298 }
299
300 // not found
301 return 0;
302 }
303
304 static size_type rowdim (const Matrix& A){
305 size_type nn=0;
306 for (size_type i=0; i<A.N(); i++)
307 nn += rowdim(A,i);
308 return nn;
309 }
310
311 static size_type coldim (const Matrix& A){
312 typedef typename Matrix::ConstRowIterator ConstRowIterator;
313 typedef typename Matrix::ConstColIterator ConstColIterator;
314
315 // The following code has a complexity of nnz, and
316 // typically a very small constant.
317 //
318 std::vector<size_type> coldims(A.M(),
320
321 for (ConstRowIterator row=A.begin(); row!=A.end(); ++row)
322 for (ConstColIterator col=row->begin(); col!=row->end(); ++col)
323 // only compute blocksizes we don't already have
324 if (coldims[col.index()]==std::numeric_limits<size_type>::max())
325 coldims[col.index()] = MatrixDimension<block_type>::coldim(*col);
326
327 size_type sum = 0;
328 for (typename std::vector<size_type>::iterator it=coldims.begin();
329 it!=coldims.end(); ++it)
330 // skip rows for which no coldim could be determined
331 if ((*it)>=0)
332 sum += *it;
333
334 return sum;
335 }
336 };
337
338
339 template<typename B, int n, int m, typename TA>
340 struct MatrixDimension<BCRSMatrix<FieldMatrix<B,n,m> ,TA> >
341 {
342 typedef BCRSMatrix<FieldMatrix<B,n,m> ,TA> Matrix;
343 typedef typename Matrix::size_type size_type;
344
345 static size_type rowdim (const Matrix& /*A*/, size_type /*i*/)
346 {
347 return n;
348 }
349
350 static size_type coldim (const Matrix& /*A*/, size_type /*c*/)
351 {
352 return m;
353 }
354
355 static size_type rowdim (const Matrix& A) {
356 return A.N()*n;
357 }
358
359 static size_type coldim (const Matrix& A) {
360 return A.M()*m;
361 }
362 };
363
364 template<typename K, int n, int m>
365 struct MatrixDimension<FieldMatrix<K,n,m> >
366 {
367 typedef FieldMatrix<K,n,m> Matrix;
368 typedef typename Matrix::size_type size_type;
369
370 static size_type rowdim(const Matrix& /*A*/, size_type /*r*/)
371 {
372 return 1;
373 }
374
375 static size_type coldim(const Matrix& /*A*/, size_type /*r*/)
376 {
377 return 1;
378 }
379
380 static size_type rowdim(const Matrix& /*A*/)
381 {
382 return n;
383 }
384
385 static size_type coldim(const Matrix& /*A*/)
386 {
387 return m;
388 }
389 };
390
391 template <class T>
392 struct MatrixDimension<Dune::DynamicMatrix<T> >
393 {
394 typedef Dune::DynamicMatrix<T> MatrixType;
395 typedef typename MatrixType::size_type size_type;
396
397 static size_type rowdim(const MatrixType& /*A*/, size_type /*r*/)
398 {
399 return 1;
400 }
401
402 static size_type coldim(const MatrixType& /*A*/, size_type /*r*/)
403 {
404 return 1;
405 }
406
407 static size_type rowdim(const MatrixType& A)
408 {
409 return A.N();
410 }
411
412 static size_type coldim(const MatrixType& A)
413 {
414 return A.M();
415 }
416 };
417
418 template<typename K, int n, int m, typename TA>
419 struct MatrixDimension<Matrix<FieldMatrix<K,n,m>, TA> >
420 {
421 typedef Matrix<FieldMatrix<K,n,m>, TA> ThisMatrix;
422 typedef typename ThisMatrix::size_type size_type;
423
424 static size_type rowdim(const ThisMatrix& /*A*/, size_type /*r*/)
425 {
426 return n;
427 }
428
429 static size_type coldim(const ThisMatrix& /*A*/, size_type /*r*/)
430 {
431 return m;
432 }
433
434 static size_type rowdim(const ThisMatrix& A)
435 {
436 return A.N()*n;
437 }
438
439 static size_type coldim(const ThisMatrix& A)
440 {
441 return A.M()*m;
442 }
443 };
444
445 template<typename K, int n>
446 struct MatrixDimension<DiagonalMatrix<K,n> >
447 {
448 typedef DiagonalMatrix<K,n> Matrix;
449 typedef typename Matrix::size_type size_type;
450
451 static size_type rowdim(const Matrix& /*A*/, size_type /*r*/)
452 {
453 return 1;
454 }
455
456 static size_type coldim(const Matrix& /*A*/, size_type /*r*/)
457 {
458 return 1;
459 }
460
461 static size_type rowdim(const Matrix& /*A*/)
462 {
463 return n;
464 }
465
466 static size_type coldim(const Matrix& /*A*/)
467 {
468 return n;
469 }
470 };
471
472 template<typename K, int n>
473 struct MatrixDimension<ScaledIdentityMatrix<K,n> >
474 {
475 typedef ScaledIdentityMatrix<K,n> Matrix;
476 typedef typename Matrix::size_type size_type;
477
478 static size_type rowdim(const Matrix& /*A*/, size_type /*r*/)
479 {
480 return 1;
481 }
482
483 static size_type coldim(const Matrix& /*A*/, size_type /*r*/)
484 {
485 return 1;
486 }
487
488 static size_type rowdim(const Matrix& /*A*/)
489 {
490 return n;
491 }
492
493 static size_type coldim(const Matrix& /*A*/)
494 {
495 return n;
496 }
497 };
498
502 template<typename T>
503 struct IsMatrix
504 {
505 enum {
509 value = false
510 };
511 };
512
513 template<typename T>
514 struct IsMatrix<DenseMatrix<T> >
515 {
516 enum {
520 value = true
521 };
522 };
523
524
525 template<typename T, typename A>
526 struct IsMatrix<BCRSMatrix<T,A> >
527 {
528 enum {
532 value = true
533 };
534 };
535
536 template<typename T>
537 struct PointerCompare
538 {
539 bool operator()(const T* l, const T* r)
540 {
541 return *l < *r;
542 }
543 };
544
545}
546#endif
Construct a matrix with a dynamic size.
Definition: dynmatrix.hh:61
derive error class from the base class in common
Definition: istlexception.hh:19
ConstIterator class for sequential access.
Definition: matrix.hh:404
A generic dynamic dense matrix.
Definition: matrix.hh:561
A::size_type size_type
Type for indices and sizes.
Definition: matrix.hh:577
MatrixImp::DenseMatrixBase< T, A >::ConstIterator ConstRowIterator
Const iterator over the matrix rows.
Definition: matrix.hh:586
row_type::const_iterator ConstColIterator
Const iterator for the entries of each row.
Definition: matrix.hh:589
T block_type
Export the type representing the components.
Definition: matrix.hh:568
This file implements a quadratic diagonal matrix of fixed size.
This file implements a dense matrix with dynamic numbers of rows and columns.
Implements a matrix constructed from a given type representing a field and compile-time given number ...
#define DUNE_THROW(E,...)
Definition: exceptions.hh:312
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:484
auto countNonZeros(const M &, typename std::enable_if_t< Dune::IsNumber< M >::value > *sfinae=nullptr)
Get the number of nonzero fields in the matrix.
Definition: matrixutils.hh:119
Dune namespace.
Definition: alignedallocator.hh:13
Implements a scalar matrix view wrapper around an existing scalar.
This file implements a quadratic matrix of fixed size which is a multiple of the identity.
Check whether the a matrix has diagonal values on blocklevel recursion levels.
Definition: matrixutils.hh:48
static void check(const Matrix &mat)
Check whether the a matrix has diagonal values on blocklevel recursion levels.
Definition: matrixutils.hh:53
Test whether a type is an ISTL Matrix.
Definition: matrixutils.hh:504
@ value
True if T is an ISTL matrix.
Definition: matrixutils.hh:509
Whether this type acts as a scalar in the context of (hierarchically blocked) containers.
Definition: typetraits.hh:194
Traits for type conversions and type information.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 21, 23:30, 2024)