DUNE-FEM (unstable)

spmatrix.hh
1#ifndef DUNE_FEM_SPMATRIX_HH
2#define DUNE_FEM_SPMATRIX_HH
3
4// system includes
5#include <algorithm>
6#include <array>
7#include <cstdlib>
8#include <iostream>
9#include <limits>
10#include <string>
11#include <utility>
12#include <vector>
13
14// local includes
15#include <dune/fem/function/adaptivefunction/adaptivefunction.hh>
16#include <dune/fem/misc/functor.hh>
17#include <dune/fem/operator/common/localmatrix.hh>
18#include <dune/fem/operator/common/localmatrixwrapper.hh>
19#include <dune/fem/io/file/asciiparser.hh>
20#include <dune/fem/io/parameter.hh>
21#include <dune/fem/operator/common/operator.hh>
22#include <dune/fem/operator/matrix/columnobject.hh>
23#include <dune/fem/space/mapper/nonblockmapper.hh>
24#include <dune/fem/storage/objectstack.hh>
25#include <dune/fem/operator/common/stencil.hh>
26#include <dune/fem/operator/matrix/functor.hh>
27#include <dune/fem/solver/parameter.hh>
28
29namespace Dune
30{
31
32 namespace Fem
33 {
34
36 template <class T, class IndexT = std::size_t,
37 class ValuesVector = std::vector< T >,
38 class IndicesVector = std::vector< IndexT > >
40 {
41 public:
43 typedef T field_type;
45 typedef IndexT size_type;
50
51 static const size_type defaultCol = std::numeric_limits<size_type>::max();
52 static const size_type zeroCol = std::numeric_limits<size_type>::max()-1;
53 static const int firstCol = 0;
54
55 SparseRowMatrix(const ThisType& ) = delete;
56
58 explicit SparseRowMatrix( const bool threading = true ) :
59 values_(0), columns_(0), rows_(0), dim_({{0,0}}), maxNzPerRow_(0), compressed_( false ), threading_( threading )
60 {}
61
64 SparseRowMatrix(const size_type rows, const size_type cols, const size_type nz, const bool threading = true ) :
65 values_(0), columns_(0), rows_(0), dim_({{0,0}}), maxNzPerRow_(0), compressed_( false ), threading_( threading )
66 {
67 reserve(rows,cols,nz);
68 }
69
71 void reserve(const size_type rows, const size_type cols, const size_type nz)
72 {
73 // if( (rows != dim_[0]) || (cols != dim_[1]) || (nz != maxNzPerRow_))
74 resize(rows,cols,nz);
75 clear();
76 }
77
79 template <class Stencil>
80 void fillPattern(const Stencil& stencil,
81 const size_type rowBlockSize,
82 const size_type colBlockSize )
83 {
84 const auto& sparsityPattern = stencil.globalStencil();
85 for( const auto& entry : sparsityPattern )
86 {
87 const auto& blockRow = entry.first;
88 const auto& blockColumnSet = entry.second;
89
90 // blocking of rows
91 const size_type nextRow = (blockRow + 1) * rowBlockSize;
92 for( size_type row = blockRow * rowBlockSize; row < nextRow; ++row )
93 {
94 size_type column = startRow( row );
95 for( const auto& blockColEntry : blockColumnSet )
96 {
97 size_type col = blockColEntry * colBlockSize;
98 for( size_type c = 0; c<colBlockSize; ++c, ++col, ++column )
99 {
100 assert( column < endRow( row ) );
101 columns_[ column ] = col ;
102 }
103 }
104 }
105 }
106 }
107
110 {
111 return dim_[0];
112 }
113
116 {
117 return dim_[1];
118 }
119
121 void set(const size_type row, const size_type col, const field_type val)
122 {
123 assert((col>=0) && (col < dim_[1]));
124 assert((row>=0) && (row < dim_[0]));
125
126 const size_type column = colIndex(row,col) ;
127 assert( column != defaultCol && column != zeroCol );
128
129 values_ [ column ] = val;
130 columns_[ column ] = col;
131 }
132
134 void add(const size_type row, const size_type col, const field_type val)
135 {
136 assert((col>=0) && (col < dim_[1]));
137 assert((row>=0) && (row < dim_[0]));
138
139 const size_type column = colIndex(row,col) ;
140 assert( column != defaultCol && column != zeroCol );
141
142 values_ [ column ] += val;
143 columns_[ column ] = col;
144 }
145
147 template<class ArgDFType, class DestDFType>
148 void apply(const ArgDFType& f, DestDFType& ret ) const
149 {
150 bool runParallel = threading_;
151
152 auto doApply = [this, &f, &ret, &runParallel]()
153 {
154 constexpr auto blockSize = ArgDFType::DiscreteFunctionSpaceType::localBlockSize;
155
156 // compute slice to be worked on in forward direction
157 const auto slice = runParallel ?
158 sliceBeginEnd( MPIManager::thread(), MPIManager::numThreads(), std::true_type() ) :
159 sliceBeginEnd( 0, 1, std::true_type() ) ;
160
161 // same as begin just with a row not necessarily zero
162 auto ret_it = ret.dofVector().find( slice.first );
163
164 const auto& fDofVec = f.dofVector();
165 for(size_type row = slice.first; row<slice.second; ++row)
166 {
167 const size_type endrow = endRow( row );
168 (*ret_it) = 0.0;
169 for(size_type col = startRow( row ); col<endrow; ++col)
170 {
171 const auto realCol = columns_[ col ];
172
173 if( ! compressed_ && ((realCol == defaultCol) || (realCol == zeroCol)) )
174 continue;
175
176 const auto blockNr = realCol / blockSize ;
177 const auto dofNr = realCol % blockSize ;
178 (*ret_it) += values_[ col ] * fDofVec[ blockNr ][ dofNr ];
179 }
180
181 ++ret_it;
182 }
183 };
184
185 if( runParallel )
186 {
187 try {
188 // execute in parallel
189 MPIManager :: run ( doApply );
190 }
191 catch ( const SingleThreadModeError& e )
192 {
193 runParallel = false;
194 }
195 }
196
197 // run serial if threading disabled or something else went wrong
198 if( ! runParallel )
199 {
200 doApply();
201 }
202 }
203
205 field_type get(const size_type row, const size_type col) const
206 {
207 assert((col>=0) && (col < dim_[1]));
208 assert((row>=0) && (row < dim_[0]));
209
210 const size_type endrow = endRow( row );
211 for( size_type i = startRow( row ); i<endrow; ++i )
212 {
213 if(columns_[ i ] == col)
214 return values_[ i ];
215 }
216 return 0;
217 }
218
220 void clear()
221 {
222 std::fill( values_.begin(), values_.end(), 0 );
223 for (auto &c : columns_) c = defaultCol;
224 }
225
227 void clearRow(const size_type row)
228 {
229 assert((row>=0) && (row < dim_[0]));
230
231 const size_type endrow = endRow( row );
232 for(size_type idx = startRow( row ); idx<endrow; ++idx )
233 {
234 values_[ idx ] = 0;
235 // if ( !compressed_ )
236 // columns_[idx] = zeroCol;
237 }
238 }
239
241 void scale(const size_type row, const size_type col, const field_type val)
242 {
243 assert((row>=0) && (row < rows()) );
244 assert((col>=0) && (col < cols()) );
245
246 const size_type column = colIndex(row,col) ;
247 assert( column != defaultCol && column != zeroCol );
248
249 // scale value
250 values_ [ column ] *= val;
251 }
252
256 {
257 return maxNzPerRow_;
258 }
259
263 {
264 return dim_[0] > 0 ? rows_[ dim_[0] ] : 0;
265 }
266
270 {
271 assert( row >= 0 && row < dim_[0] );
272 return endRow( row ) - startRow( row );
273 }
274
277 std::pair<const field_type, size_type> realValue(size_type index) const
278 {
279 return std::pair<const field_type, size_type>(values_[index], columns_[index]);
280 }
281
283 void print(std::ostream& s=std::cout, unsigned int offset=0) const
284 {
285 for(size_type row=0; row<dim_[0]; ++row)
286 {
287 const size_type endrow = endRow( row );
288 for( size_type pos = startRow( row ); pos<endrow; ++pos )
289 {
290 const auto rv(realValue(pos));
291 const auto column(rv.second);
292 const auto value(rv.first);
293 if((std::abs(value) > 1.e-15) && (column != defaultCol))
294 s << row+offset << " " << column+offset << " " << value << std::endl;
295 }
296 }
297 }
298
299 template <class SizeT, class NumericT >
300 void fillCSRStorage( std::vector< std::map<SizeT, NumericT> >& matrix ) const
301 {
302 matrix.resize( dim_[0] );
303
304 size_type thisCol = 0;
305 for(size_type row = 0; row<dim_[ 0 ]; ++row )
306 {
307 auto& matRow = matrix[ row ];
308 const size_type endrow = endRow( row );
309 for(size_type col = startRow( row ); col<endrow; ++col)
310 {
311 const size_type realCol = columns_[ col ];
312
313 if( ! compressed_ && (realCol == defaultCol || realCol == zeroCol) )
314 continue;
315
316 matRow[ realCol ] = values_[ thisCol ];
317 }
318 }
319 }
320
321 void compress ()
322 {
323 if( ! compressed_ && (dim_[0] != 0) && (dim_[1] != 0))
324 {
325 // determine first row nonZeros
326 size_type newpos = 0 ;
327 for( newpos = startRow( 0 ); newpos < endRow( 0 ); ++newpos )
328 {
329 if( columns_[ newpos ] == defaultCol )
330 break;
331 }
332
333 for( size_type row = 1; row<dim_[0]; ++row )
334 {
335 const size_type endrow = endRow( row );
336 size_type col = startRow( row );
337 // update new row start position
338 rows_[ row ] = newpos;
339 for(; col<endrow; ++col, ++newpos )
340 {
341 if( columns_[ col ] == defaultCol )
342 break ;
343
344 assert( newpos <= col );
345 values_ [ newpos ] = values_ [ col ];
346 columns_[ newpos ] = columns_[ col ];
347 }
348 }
349 rows_[ dim_[0] ] = newpos ;
350
351 // values_.resize( newpos );
352 // columns_.resize( newpos );
353 compressed_ = true ;
354 }
355 }
356
357 size_type startRow ( const size_type row ) const
358 {
359 return rows_[ row ];
360 }
361
362 size_type endRow ( const size_type row ) const
363 {
364 return rows_[ row+1 ];
365 }
366
367 std::tuple< ValuesVector&, IndicesVector&, IndicesVector& > exportCRS()
368 {
369 // only return internal data in compressed status
370 compress();
371 return std::tie(values_,columns_,rows_);
372 }
373
375 template<class DiagType, class ArgDFType, class DestDFType, class WType>
376 void forwardIterative(const DiagType& diagInv, const ArgDFType& b, const DestDFType& xold, DestDFType& xnew, const WType& w ) const
377 {
378 parallelIterative(diagInv, b, xold, xnew, w, std::true_type() );
379 }
380
382 template<class DiagType, class ArgDFType, class DestDFType, class WType>
383 void backwardIterative(const DiagType& diagInv, const ArgDFType& b, const DestDFType& xold, DestDFType& xnew, const WType& w ) const
384 {
385 parallelIterative(diagInv, b, xold, xnew, w, std::false_type() );
386 }
387
388 protected:
389 // return slice of rows for given thread in forward direction
390 std::pair< size_type, size_type > sliceBeginEnd(const size_type thread, const size_type numThreads, std::true_type ) const
391 {
392 const size_type sliceSize = this->rows() / numThreads;
393 const size_type sliceBegin = (thread * sliceSize) ;
394 const size_type sliceEnd = (thread == numThreads-1 ) ? this->rows(): (sliceBegin + sliceSize);
395 return std::make_pair( sliceBegin, sliceEnd );
396 }
397
398 // return slice of rows for given thread in backward direction
399 std::pair< size_type, size_type > sliceBeginEnd(const size_type thread, const size_type numThreads, std::false_type ) const
400 {
401 std::pair< size_type, size_type > beginEnd = sliceBeginEnd( thread, numThreads, std::true_type() );
402 return std::make_pair( beginEnd.second-1, beginEnd.first-1 );
403 }
404
405 template<class DiagType, class ArgDFType, class DestDFType, class WType, bool forward >
406 void parallelIterative(const DiagType& diagInv, const ArgDFType& b, const DestDFType& xold, DestDFType& xnew,
407 const WType& w, std::integral_constant<bool, forward> direction ) const
408 {
409 bool runParallel = threading_ && (&xold != &xnew) ; // only for Jacobi
410
411 auto doIterate = [this, &diagInv, &b, &xold, &xnew, &w, &direction, &runParallel] ()
412 {
413 // compute slice to be worked on depending on direction depending on parallel mode
414 const auto slice = runParallel ?
415 sliceBeginEnd( MPIManager::thread(), MPIManager::numThreads(), direction ) :
416 sliceBeginEnd( 0, 1, direction ) ;
417
418 doParallelIterative( diagInv.dofVector().find( slice.first ), // still O(1) operation just like for begin()
419 b.dofVector().find( slice.first ),
420 xold,
421 xnew.dofVector().find( slice.first ),
422 w,
423 slice.first, // row begin
424 slice.second, // row end
425 direction );
426 };
427
428 if( runParallel )
429 {
430 try {
431 // execute in parallel
432 MPIManager :: run ( doIterate );
433 }
434 catch ( const SingleThreadModeError& e )
435 {
436 runParallel = false;
437 }
438 }
439
440 // run serial if threading disabled or something else went wrong
441 if( ! runParallel )
442 {
443 doIterate();
444 }
445 }
446
448 template<class DiagIt, class ArgDFIt, class DestDFType, class DestDFIt,
449 class WType, bool forward>
450 void doParallelIterative(DiagIt diag, ArgDFIt bit, const DestDFType& xold, DestDFIt xit,
451 const WType& w,
452 size_type row, const size_type end,
453 std::integral_constant<bool, forward> ) const
454 {
455 constexpr auto blockSize = DestDFType::DiscreteFunctionSpaceType::localBlockSize;
456
457 const auto nextRow = [&diag, &xit, &bit](size_type &row)
458 {
459 if constexpr ( forward )
460 {
461 ++diag; ++xit; ++bit; ++row;
462 }
463 else
464 {
465 --diag; --xit; --bit; --row;
466 }
467 };
468
469 const auto& xOldVec = xold.dofVector();
470 for(; row != end; nextRow(row))
471 {
472 auto rhs = (*bit);
473
474 const size_type endrow = endRow( row );
475 for(size_type col = startRow( row ); col<endrow; ++col)
476 {
477 const auto realCol = columns_[ col ];
478
479 if( (realCol == row ) || (! compressed_ && ((realCol == defaultCol) || (realCol == zeroCol))) )
480 continue;
481
482 const auto blockNr = realCol / blockSize ;
483 const auto dofNr = realCol % blockSize ;
484
485 rhs -= values_[ col ] * xOldVec[ blockNr ][ dofNr ] ;
486 }
487
488 (*xit) = w * (rhs * (*diag));
489 }
490 }
493 {
494 constexpr auto colVal = defaultCol;
495 values_.resize( rows*nz , 0 );
496 columns_.resize( rows*nz , colVal );
497 rows_.resize( rows+1 , 0 );
498 rows_[ 0 ] = 0;
499 for( size_type i=1; i <= rows; ++i )
500 {
501 rows_[ i ] = rows_[ i-1 ] + nz ;
502 }
503 compressed_ = false;
504
505 dim_[0] = rows;
506 dim_[1] = cols;
507 maxNzPerRow_ = nz+firstCol;
508 }
509
512 {
513 assert((col>=0) && (col < dim_[1]));
514 assert((row>=0) && (row < dim_[0]));
515
516 const size_type endR = endRow( row );
517 size_type i = startRow( row );
518 // find local column or empty spot
519 for( ; i < endR; ++i )
520 {
521 if( columns_[ i ] == defaultCol || columns_[ i ] == zeroCol || columns_[ i ] == col )
522 {
523 return i;
524 }
525 }
526
527 assert(0);
528 DUNE_THROW( InvalidStateException, "Could not store entry in sparse matrix - no space available" );
529
530 // TODO: implement resize with 2*nz
531 std::abort();
532 return defaultCol;
533
534 /*
535 if(columns_[ i ] == col)
536 return i; // column already in matrix
537 else if( columns_[ i ] == defaultCol )
538 { // add this column at end of this row
539 ++nonZeros_[row];
540 return i;
541 }
542 else
543 {
544 std::abort();
545 // TODO re-implement
546 //
547 ++nonZeros_[row];
548 // must shift this row to add col at the position i
549 auto j = nz_-1; // last column
550 if (columns_[row*nz_+j] != defaultCol)
551 { // new space available - so resize
552 resize( rows(), cols(), (2 * nz_) );
553 j++;
554 }
555 for(;j>i;--j)
556 {
557 columns_[row*nz_+j] = columns_[row*nz_+j-1];
558 values_[row*nz_+j] = values_[row*nz_+j-1];
559 }
560 columns_[row*nz_+i] = col;
561 values_[row*nz_+i] = 0;
562 return i;
563 }
564 */
565 }
566
567 ValuesVector values_;
568 IndicesVector columns_;
569 IndicesVector rows_;
570
571 std::array<size_type,2> dim_;
572 size_type maxNzPerRow_;
573 bool compressed_;
574 const bool threading_;
575 };
576
577
578
580 template< class DomainSpace, class RangeSpace,
583 {
584 protected:
585 template< class MatrixObject >
586 struct LocalMatrixTraits;
587
588 template< class MatrixObject >
589 class LocalMatrix;
590
591 public:
592 typedef DomainSpace DomainSpaceType;
593 typedef RangeSpace RangeSpaceType;
594 typedef typename DomainSpaceType::EntityType DomainEntityType;
595 typedef typename RangeSpaceType::EntityType RangeEntityType;
596 typedef typename DomainSpaceType::EntityType ColumnEntityType;
597 typedef typename RangeSpaceType::EntityType RowEntityType;
598
599 typedef typename DomainSpaceType::BlockMapperType DomainBlockMapperType;
601 typedef typename RangeSpaceType::BlockMapperType RangeBlockMapperType;
603 typedef Matrix MatrixType;
604 typedef typename MatrixType::size_type size_type;
605 typedef typename MatrixType::field_type field_type;
606
608
609 static const size_type domainLocalBlockSize = DomainSpaceType::dimRange;
610 static const size_type rangeLocalBlockSize = RangeSpaceType::dimRange;
611
614
616
620 typedef LocalMatrixWrapper< LocalMatrixStackType > LocalMatrixType;
621 typedef ColumnObject< ThisType > LocalColumnObjectType;
622
624 SparseRowMatrixObject( const DomainSpaceType &domainSpace,
625 const RangeSpaceType &rangeSpace,
626 const SolverParameter& param = SolverParameter() )
627 : domainSpace_( domainSpace ),
628 rangeSpace_( rangeSpace ),
629 domainMapper_( domainSpace_.blockMapper() ),
630 rangeMapper_( rangeSpace_.blockMapper() ),
631 sequence_( -1 ),
632 matrix_( param.threading() ),
633 localMatrixStack_( *this )
634 {}
635
637 const DomainSpaceType& domainSpace() const
638 {
639 return domainSpace_;
640 }
641
643 const RangeSpaceType& rangeSpace() const
644 {
645 return rangeSpace_;
646 }
647
648 protected:
651 {
652 return matrix_;
653 }
654
655 void finalizeAssembly() const { const_cast< ThisType& > (*this).compress(); }
656
657 public:
660 {
661 finalizeAssembly();
662 return matrix_;
663 }
664
665
668 {
669 return new ObjectType( *this, domainSpace_, rangeSpace_, domainMapper_, rangeMapper_ );
670 }
671
676 [[deprecated("Use TemporaryLocalMatrix,LocalContribution and {get,add,set}LocalMatrix")]]
677 LocalMatrixType localMatrix( const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity ) const
678 {
679 return LocalMatrixType( localMatrixStack_, domainEntity, rangeEntity );
680 }
681
686 [[deprecated("Use TemporaryLocalMatrix,LocalContribution and {get,add,set}LocalMatrix")]]
687 LocalMatrixType localMatrix() const
688 {
689 return LocalMatrixType( localMatrixStack_ );
690 }
691
693 LocalColumnObjectType localColumn( const DomainEntityType &domainEntity ) const
694 {
695 return LocalColumnObjectType( *this, domainEntity );
696 }
697
698 void unitRow( const size_type row )
699 {
700 for( unsigned int i=0, r = row * domainLocalBlockSize; i<domainLocalBlockSize; ++i, ++r )
701 {
702 matrix_.clearRow( r );
703 matrix_.set( r, r, 1.0 );
704 }
705 }
706
707 template <class LocalBlock>
708 void addBlock( const size_type row, const size_type col, const LocalBlock& block )
709 {
710 std::array< size_type, rangeLocalBlockSize > rows;
711 std::array< size_type, domainLocalBlockSize > cols;
712 for( unsigned int i=0, r = row * domainLocalBlockSize, c = col * domainLocalBlockSize; i<domainLocalBlockSize; ++i, ++r, ++c )
713 {
714 rows[ i ] = r;
715 cols[ i ] = c;
716 }
717
718 for( unsigned int i=0; i<domainLocalBlockSize; ++i )
719 {
720 for( unsigned int j=0; j<domainLocalBlockSize; ++j )
721 {
722 matrix_.add( rows[ i ], cols[ j ], block[ i ][ j ]);
723 }
724 }
725 }
726
727 template <class LocalBlock>
728 void setBlock( const size_type row, const size_type col, const LocalBlock& block )
729 {
730 std::array< size_type, rangeLocalBlockSize > rows;
731 std::array< size_type, domainLocalBlockSize > cols;
732 for( unsigned int i=0, r = row * domainLocalBlockSize, c = col * domainLocalBlockSize; i<domainLocalBlockSize; ++i, ++r, ++c )
733 {
734 rows[ i ] = r;
735 cols[ i ] = c;
736 }
737
738 for( unsigned int i=0; i<domainLocalBlockSize; ++i )
739 {
740 for( unsigned int j=0; j<domainLocalBlockSize; ++j )
741 {
742 matrix_.set( rows[ i ], cols[ j ], block[ i ][ j ]);
743 }
744 }
745 }
746
747 template< class LocalMatrix >
748 void addLocalMatrix ( const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity, const LocalMatrix &localMat )
749 {
750 auto functor = [ &localMat, this ] ( std::pair< int, int > local, const std::pair< size_type, size_type >& global )
751 {
752 matrix_.add( global.first, global.second, localMat.get( local.first, local.second ) );
753 };
754
755 rangeMapper_.mapEach( rangeEntity, makePairFunctor( domainMapper_, domainEntity, functor ) );
756 }
757
758 template< class LocalMatrix, class Scalar >
759 void addScaledLocalMatrix ( const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity, const LocalMatrix &localMat, const Scalar &s )
760 {
761 auto functor = [ &localMat, &s, this ] ( std::pair< int, int > local, const std::pair< size_type, size_type >& global )
762 {
763 matrix_.add( global.first, global.second, s * localMat.get( local.first, local.second ) );
764 };
765
766 rangeMapper_.mapEach( rangeEntity, makePairFunctor( domainMapper_, domainEntity, functor ) );
767 }
768
769 template< class LocalMatrix >
770 void setLocalMatrix ( const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity, const LocalMatrix &localMat )
771 {
772 auto functor = [ &localMat, this ] ( std::pair< int, int > local, const std::pair< size_type, size_type >& global )
773 {
774 matrix_.set( global.first, global.second, localMat.get( local.first, local.second ) );
775 };
776
777 rangeMapper_.mapEach( rangeEntity, makePairFunctor( domainMapper_, domainEntity, functor ) );
778 }
779
780 template< class LocalMatrix >
781 void getLocalMatrix ( const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity, LocalMatrix &localMat ) const
782 {
783 auto functor = [ &localMat, this ] ( std::pair< int, int > local, const std::pair< size_type, size_type >& global )
784 {
785 localMat.set( local.first, local.second, matrix_.get( global.first, global.second ) );
786 };
787
788 rangeMapper_.mapEach( rangeEntity, makePairFunctor( domainMapper_, domainEntity, functor ) );
789 }
790
792 void clear()
793 {
794 matrix_.clear();
795 }
796
798 void compress() { matrix_.compress(); }
799
800 template <class Set>
801 void reserve (const std::vector< Set >& sparsityPattern )
802 {
803 reserve( StencilWrapper< DomainSpaceType,RangeSpaceType, Set >( sparsityPattern ) );
804 }
805
807 template <class Stencil>
808 void reserve(const Stencil &stencil, bool verbose = false )
809 {
810 if( sequence_ != domainSpace_.sequence() )
811 {
812 // if empty grid do nothing (can appear in parallel runs)
813 if( (domainSpace_.begin() != domainSpace_.end()) && (rangeSpace_.begin() != rangeSpace_.end()) )
814 {
815 // output some info
816 if( verbose )
817 {
818 std::cout << "Reserve Matrix with (" << rangeSpace_.size() << "," << domainSpace_.size()<< ")" << std::endl;
819 std::cout << "Max number of base functions = (" << rangeMapper_.maxNumDofs() << ","
820 << domainMapper_.maxNumDofs() << ")" << std::endl;
821 }
822
823 // reserve matrix
824 const auto nonZeros = std::max( static_cast<size_type>(stencil.maxNonZerosEstimate()*DomainSpaceType::localBlockSize),
825 matrix_.maxNzPerRow() );
826 matrix_.reserve( rangeSpace_.size(), domainSpace_.size(), nonZeros );
827 matrix_.fillPattern( stencil, RangeSpaceType::localBlockSize, DomainSpaceType::localBlockSize );
828 }
829 sequence_ = domainSpace_.sequence();
830 }
831 }
832
834 template< class DomainFunction, class RangeFunction >
835 void apply( const DomainFunction &arg, RangeFunction &dest ) const
836 {
837 // do matrix vector multiplication
838 matrix_.apply( arg, dest );
839 // communicate data
840 dest.communicate();
841 }
842
845 template < class DiscreteFunctionType >
847 {
848 assert( matrix_.rows() == matrix_.cols() );
849 const auto dofEnd = diag.dend();
850 size_type row = 0;
851 for( auto dofIt = diag.dbegin(); dofIt != dofEnd; ++dofIt, ++row )
852 {
853 assert( row < matrix_.rows() );
854 (*dofIt) = matrix_.get( row, row );
855 }
856 }
857
858 template <class Container>
859 void setUnitRows( const Container& unitRows, const Container& auxRows )
860 {
861 for (const auto& r : unitRows )
862 {
863 matrix_.clearRow(r);
864 matrix_.set(r, r, 1.0);
865 }
866
867 for (const auto& r : auxRows )
868 {
869 matrix_.clearRow(r);
870 // not sure if this is really needed,
871 // but for consistency with previous code
872 // matrix_.set(r, r, 0.0);
873 }
874 }
875
877 void resort()
878 {
879 DUNE_THROW(NotImplemented,"SpMatrixObject::resort is not implemented");
880 // this method does not even exist on SpMatrix!!!
881 // matrix_.resort();
882 }
883
884 protected:
885 const DomainSpaceType &domainSpace_;
886 const RangeSpaceType &rangeSpace_;
887 DomainMapperType domainMapper_ ;
888 RangeMapperType rangeMapper_ ;
889 int sequence_;
890 mutable MatrixType matrix_;
891 mutable LocalMatrixStackType localMatrixStack_;
892 };
893
894
895
897 template< class DomainSpace, class RangeSpace, class Matrix >
898 template< class MatrixObject >
899 struct SparseRowMatrixObject< DomainSpace, RangeSpace, Matrix >::LocalMatrixTraits
900 {
901 typedef DomainSpace DomainSpaceType;
902 typedef RangeSpace RangeSpaceType;
903
905
906 typedef typename SparseRowMatrixObjectType::template LocalMatrix< MatrixObject > LocalMatrixType;
907
908 typedef typename RangeSpaceType::RangeFieldType RangeFieldType;
909 typedef RangeFieldType LittleBlockType;
910
913 };
914
915
916
918 template< class DomainSpace, class RangeSpace, class Matrix >
919 template< class MatrixObject >
920 class SparseRowMatrixObject< DomainSpace, RangeSpace, Matrix >::LocalMatrix
921 : public LocalMatrixDefault< LocalMatrixTraits< MatrixObject > >
922 {
923 public:
925 typedef MatrixObject MatrixObjectType;
926
929
930 private:
932
933 public:
935 typedef typename MatrixObjectType::MatrixType MatrixType;
936
938 typedef typename Traits::RangeFieldType RangeFieldType;
939
942
944 typedef typename Traits::LittleBlockType LittleBlockType;
945
950
951 typedef std::vector< typename RangeMapperType::SizeType > RowIndicesType;
952 typedef std::vector< typename DomainMapperType::SizeType > ColumnIndicesType;
953
955 LocalMatrix( const MatrixObjectType &matrixObject,
956 const DomainSpaceType &domainSpace,
957 const RangeSpaceType &rangeSpace,
958 const DomainMapperType& domainMapper,
959 const RangeMapperType& rangeMapper )
961 matrix_( matrixObject.matrix() ),
962 domainMapper_( domainMapper ),
963 rangeMapper_( rangeMapper )
964 {}
965
966 LocalMatrix( const LocalMatrix & ) = delete;
967
968 void init( const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity )
969 {
970 // initialize base functions sets
971 BaseType::init( domainEntity, rangeEntity );
972 // rows are determined by the range space
973 rowIndices_.resize( rangeMapper_.numDofs( rangeEntity ) );
974 rangeMapper_.mapEach( rangeEntity, AssignFunctor< RowIndicesType >( rowIndices_ ) );
975 // columns are determined by the domain space
976 columnIndices_.resize( domainMapper_.numDofs( domainEntity ) );
977 domainMapper_.mapEach( domainEntity, AssignFunctor< ColumnIndicesType >( columnIndices_ ) );
978 }
979
981 size_type rows() const
982 {
983 return rowIndices_.size();
984 }
985
987 size_type columns() const
988 {
989 return columnIndices_.size();
990 }
991
993 void add(size_type localRow, size_type localCol, DofType value)
994 {
995 assert( value == value );
996 assert( (localRow >= 0) && (localRow < rows()) );
997 assert( (localCol >= 0) && (localCol < columns()) );
998 matrix_.add( rowIndices_[ localRow ], columnIndices_[ localCol ], value );
999 }
1000
1002 DofType get(size_type localRow, size_type localCol) const
1003 {
1004 assert( (localRow >= 0) && (localRow < rows()) );
1005 assert( (localCol >= 0) && (localCol < columns()) );
1006 return matrix_.get( rowIndices_[ localRow ], columnIndices_[ localCol ] );
1007 }
1008
1010 void set(size_type localRow, size_type localCol, DofType value)
1011 {
1012 assert( (localRow >= 0) && (localRow < rows()) );
1013 assert( (localCol >= 0) && (localCol < columns()) );
1014 matrix_.set( rowIndices_[ localRow ], columnIndices_[ localCol ], value );
1015 }
1016
1018 void unitRow(size_type localRow)
1019 {
1020 assert( (localRow >= 0) && (localRow < rows()) );
1021 matrix_.unitRow( rowIndices_[ localRow ] );
1022 }
1023
1025 void clearRow( size_type localRow )
1026 {
1027 assert( (localRow >= 0) && (localRow < rows()) );
1028 matrix_.clearRow( rowIndices_[localRow]);
1029 }
1030
1032 void clearCol( size_type localCol )
1033 {
1034 assert( (localCol >= 0) && (localCol < columns()) );
1035 matrix_.clearCol( columnIndices_[localCol] );
1036 }
1037
1039 void clear()
1040 {
1041 const size_type nrows = rows();
1042 for(size_type i=0; i < nrows; ++i )
1043 matrix_.clearRow( rowIndices_[ i ] );
1044 }
1045
1047 void resort()
1048 {
1049 DUNE_THROW(NotImplemented,"SpMatrixObject::LocalMatrix::resort is not implemented");
1050 //const size_type nrows = rows();
1051 //for(size_type i=0; i < nrows; ++i )
1052 //matrix_.resortRow( rowIndices_[ i ] );
1053 }
1054
1056 void scale( const DofType& value )
1057 {
1058 const size_type nrows = rows();
1059 const size_type ncols = columns();
1060 for(size_type i=0; i < nrows; ++i )
1061 {
1062 for( size_type j=0; j < ncols; ++j )
1063 {
1064 scale(i, j, value );
1065 }
1066 }
1067 }
1068
1069 protected:
1071 void scale(size_type localRow, size_type localCol, DofType value)
1072 {
1073 assert( (localRow >= 0) && (localRow < rows()) );
1074 assert( (localCol >= 0) && (localCol < columns()) );
1075 matrix_.scale( rowIndices_[ localRow ], columnIndices_[ localCol ], value );
1076 }
1077
1078 protected:
1079 MatrixType &matrix_;
1080 const DomainMapperType& domainMapper_;
1081 const RangeMapperType& rangeMapper_;
1082 RowIndicesType rowIndices_;
1083 ColumnIndicesType columnIndices_;
1084 };
1085
1086 } // namespace Fem
1087
1088} // namespace Dune
1089
1090#endif // #ifndef DUNE_FEM_SPMATRIX_HH
ConstDofIteratorType dbegin() const
Obtain the constant iterator pointing to the first dof.
Definition: discretefunction.hh:761
ConstDofIteratorType dend() const
Obtain the constant iterator pointing to the last dof.
Definition: discretefunction.hh:773
Default implementation for local matrix classes.
Definition: localmatrix.hh:287
Exception thrown when a code segment that is supposed to be only accessed in single thread mode is ac...
Definition: mpimanager.hh:43
LocalMatrix.
Definition: spmatrix.hh:922
void clearRow(size_type localRow)
set matrix row to zero
Definition: spmatrix.hh:1025
void add(size_type localRow, size_type localCol, DofType value)
add value to matrix entry
Definition: spmatrix.hh:993
MatrixObjectType::MatrixType MatrixType
type of matrix
Definition: spmatrix.hh:935
void clearCol(size_type localCol)
set matrix column to zero
Definition: spmatrix.hh:1032
void resort()
resort all global rows of matrix to have ascending numbering
Definition: spmatrix.hh:1047
size_type columns() const
return number of columns
Definition: spmatrix.hh:987
size_type rows() const
return number of rows
Definition: spmatrix.hh:981
Traits::RangeFieldType RangeFieldType
type of entries of little blocks
Definition: spmatrix.hh:938
DofType get(size_type localRow, size_type localCol) const
get matrix entry
Definition: spmatrix.hh:1002
MatrixObject MatrixObjectType
type of matrix object
Definition: spmatrix.hh:925
Traits::DomainMapperType DomainMapperType
type of nonblocked domain mapper
Definition: spmatrix.hh:947
Traits::RangeMapperType RangeMapperType
type of nonblocked domain mapper
Definition: spmatrix.hh:949
LocalMatrixTraits< MatrixObjectType > Traits
type of the traits
Definition: spmatrix.hh:928
LocalMatrix(const MatrixObjectType &matrixObject, const DomainSpaceType &domainSpace, const RangeSpaceType &rangeSpace, const DomainMapperType &domainMapper, const RangeMapperType &rangeMapper)
constructor
Definition: spmatrix.hh:955
void scale(const DofType &value)
scale local matrix with a certain value
Definition: spmatrix.hh:1056
RangeFieldType DofType
type of the DoFs
Definition: spmatrix.hh:941
void scale(size_type localRow, size_type localCol, DofType value)
scale matrix entry with value
Definition: spmatrix.hh:1071
void clear()
clear all entries belonging to local matrix
Definition: spmatrix.hh:1039
Traits::LittleBlockType LittleBlockType
type of little blocks
Definition: spmatrix.hh:944
void set(size_type localRow, size_type localCol, DofType value)
set matrix entry to value
Definition: spmatrix.hh:1010
void unitRow(size_type localRow)
set matrix row to zero except diagonla entry
Definition: spmatrix.hh:1018
SparseRowMatrixObject.
Definition: spmatrix.hh:583
const DomainSpaceType & domainSpace() const
get domain space (i.e. space that builds the rows)
Definition: spmatrix.hh:637
LocalMatrixType localMatrix() const
Definition: spmatrix.hh:687
MatrixType & exportMatrix() const
get reference to storage object
Definition: spmatrix.hh:659
void apply(const DomainFunction &arg, RangeFunction &dest) const
apply matrix to discrete function
Definition: spmatrix.hh:835
SparseRowMatrixObject(const DomainSpaceType &domainSpace, const RangeSpaceType &rangeSpace, const SolverParameter &param=SolverParameter())
construct matrix object
Definition: spmatrix.hh:624
LocalColumnObjectType localColumn(const DomainEntityType &domainEntity) const
get local column
Definition: spmatrix.hh:693
void extractDiagonal(DiscreteFunctionType &diag) const
Definition: spmatrix.hh:846
LocalMatrixType localMatrix(const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity) const
Definition: spmatrix.hh:677
ObjectType * newObject() const
interface method from LocalMatrixFactory
Definition: spmatrix.hh:667
void reserve(const Stencil &stencil, bool verbose=false)
reserve memory
Definition: spmatrix.hh:808
void compress()
compress matrix to a real CRS format
Definition: spmatrix.hh:798
void resort()
resort row numbering in matrix to have ascending numbering
Definition: spmatrix.hh:877
void clear()
clear matrix
Definition: spmatrix.hh:792
const RangeSpaceType & rangeSpace() const
get range space (i.e. space that builds the columns)
Definition: spmatrix.hh:643
MatrixType & matrix() const
get reference to storage object, for internal use
Definition: spmatrix.hh:650
SparseRowMatrix.
Definition: spmatrix.hh:40
void print(std::ostream &s=std::cout, unsigned int offset=0) const
print matrix
Definition: spmatrix.hh:283
size_type colIndex(size_type row, size_type col)
returns local col index for given global (row,col)
Definition: spmatrix.hh:511
std::pair< const field_type, size_type > realValue(size_type index) const
Definition: spmatrix.hh:277
SparseRowMatrix(const bool threading=true)
construct matrix of zero size
Definition: spmatrix.hh:58
void clearRow(const size_type row)
set all entries in row to zero
Definition: spmatrix.hh:227
IndexT size_type
matrix index type
Definition: spmatrix.hh:45
void add(const size_type row, const size_type col, const field_type val)
add value to row,col entry
Definition: spmatrix.hh:134
void apply(const ArgDFType &f, DestDFType &ret) const
ret = A*f
Definition: spmatrix.hh:148
size_type numNonZeros(size_type row) const
Definition: spmatrix.hh:269
void reserve(const size_type rows, const size_type cols, const size_type nz)
reserve memory for given rows, columns and number of non zeros
Definition: spmatrix.hh:71
void clear()
set all matrix entries to zero
Definition: spmatrix.hh:220
size_type numNonZeros() const
Definition: spmatrix.hh:262
void forwardIterative(const DiagType &diagInv, const ArgDFType &b, const DestDFType &xold, DestDFType &xnew, const WType &w) const
Apply Jacobi/SOR method.
Definition: spmatrix.hh:376
field_type get(const size_type row, const size_type col) const
return value of entry (row,col)
Definition: spmatrix.hh:205
void doParallelIterative(DiagIt diag, ArgDFIt bit, const DestDFType &xold, DestDFIt xit, const WType &w, size_type row, const size_type end, std::integral_constant< bool, forward >) const
Apply Jacobi/SOR method.
Definition: spmatrix.hh:450
void set(const size_type row, const size_type col, const field_type val)
set entry to value (also setting 0 will result in an entry)
Definition: spmatrix.hh:121
void resize(size_type rows, size_type cols, size_type nz)
resize matrix
Definition: spmatrix.hh:492
size_type rows() const
return number of rows
Definition: spmatrix.hh:109
size_type cols() const
return number of columns
Definition: spmatrix.hh:115
void scale(const size_type row, const size_type col, const field_type val)
scale all entries in row with a given value
Definition: spmatrix.hh:241
size_type maxNzPerRow() const
Definition: spmatrix.hh:255
void backwardIterative(const DiagType &diagInv, const ArgDFType &b, const DestDFType &xold, DestDFType &xnew, const WType &w) const
Apply Jacobi/SOR method.
Definition: spmatrix.hh:383
SparseRowMatrix(const size_type rows, const size_type cols, const size_type nz, const bool threading=true)
Definition: spmatrix.hh:64
ThisType MatrixBaseType
Definition: spmatrix.hh:49
void fillPattern(const Stencil &stencil, const size_type rowBlockSize, const size_type colBlockSize)
reserve memory for given rows, columns and number of non zeros
Definition: spmatrix.hh:80
T field_type
matrix field type
Definition: spmatrix.hh:43
a simple wrapper class for sparsity patterns provide as vector< set< size_t > >
Definition: stencil.hh:234
default implementation for a general operator stencil
Definition: stencil.hh:35
const GlobalStencilType & globalStencil() const
Return the full stencil.
Definition: stencil.hh:116
int maxNonZerosEstimate() const
Return an upper bound for the maximum number of non-zero entries in all rows.
Definition: stencil.hh:123
forward declaration
Definition: discretefunction.hh:51
A dense n x m matrix.
Definition: fmatrix.hh:117
Default exception if a function was called while the object is not in a valid state for that function...
Definition: exceptions.hh:375
A generic dynamic dense matrix.
Definition: matrix.hh:561
A::size_type size_type
Type for indices and sizes.
Definition: matrix.hh:577
typename Imp::BlockTraits< T >::field_type field_type
Export the type representing the underlying field.
Definition: matrix.hh:565
Default exception for dummy implementations.
Definition: exceptions.hh:357
#define DUNE_THROW(E,...)
Definition: exceptions.hh:314
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:485
typename Overloads::ScalarType< std::decay_t< V > >::type Scalar
Element type of some SIMD type.
Definition: interface.hh:235
Dune namespace.
Definition: alignedallocator.hh:13
LocalMatrixTraits.
Definition: spmatrix.hh:900
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Apr 19, 22:44, 2025)