2 #ifndef DUNE_FEM_PETSCLINEAROPERATOR_HH 3 #define DUNE_FEM_PETSCLINEAROPERATOR_HH 8 #include <dune/common/dynmatrix.hh> 25 #if defined HAVE_PETSC 35 struct PetscMatrixParameter
36 :
public MatrixParameter
38 typedef MatrixParameter BaseType;
40 PetscMatrixParameter(
const std::string keyPrefix =
"petscmatrix." )
41 : BaseType( keyPrefix )
49 template<
typename DomainFunction,
typename RangeFunction >
50 class PetscLinearOperator
51 :
public Fem::Operator< DomainFunction, RangeFunction >
53 typedef PetscLinearOperator< DomainFunction, RangeFunction > ThisType;
55 typedef Mat MatrixType;
56 typedef DomainFunction DomainFunctionType;
57 typedef RangeFunction RangeFunctionType;
58 typedef typename DomainFunctionType::RangeFieldType DomainFieldType;
59 typedef typename RangeFunctionType::RangeFieldType RangeFieldType;
60 typedef typename DomainFunctionType::DiscreteFunctionSpaceType DomainSpaceType;
61 typedef typename RangeFunctionType::DiscreteFunctionSpaceType RangeSpaceType;
63 typedef typename DomainSpaceType::GridPartType::template Codim< 0 >::EntityType RowEntityType;
64 typedef typename RangeSpaceType::GridPartType::template Codim< 0 >::EntityType ColumnEntityType;
66 const static size_t domainLocalBlockSize = DomainSpaceType::localBlockSize;
67 const static size_t rangeLocalBlockSize = RangeSpaceType::localBlockSize;
70 typedef PetscSlaveDofProvider< DomainSpaceType > RowPetscSlaveDofsType;
71 typedef PetscSlaveDofProvider< RangeSpaceType > ColPetscSlaveDofsType;
72 enum Status {statAssembled=0,statAdd=1,statInsert=2,statGet=3,statNothing=4};
75 typedef typename ColPetscSlaveDofsType :: PetscDofMappingType ColDofMappingType;
76 typedef typename RowPetscSlaveDofsType :: PetscDofMappingType RowDofMappingType;
81 struct LocalMatrixTraits
83 typedef typename DomainFunctionType::DiscreteFunctionSpaceType DomainSpaceType;
84 typedef typename RangeFunctionType::DiscreteFunctionSpaceType RangeSpaceType;
85 typedef LocalMatrix LocalMatrixType;
86 typedef typename RangeSpaceType::RangeFieldType RangeFieldType;
89 typedef RangeFieldType LittleBlockType;
93 typedef LocalMatrix ObjectType;
94 typedef ThisType LocalMatrixFactoryType;
98 typedef LocalMatrixWrapper< LocalMatrixStackType > LocalMatrixType;
99 typedef ColumnObject< ThisType > LocalColumnObjectType;
104 PetscLinearOperator (
const std::string &,
const DomainSpaceType &domainSpace,
const RangeSpaceType &rangeSpace,
105 const MatrixParameter& param = PetscMatrixParameter() )
106 : domainSpace_( domainSpace ),
107 rangeSpace_( rangeSpace ),
108 colSlaveDofs_( rangeSpace_ ),
109 rowSlaveDofs_( domainSpace_ ),
111 localMatrixStack_( *this ),
115 PetscLinearOperator (
const DomainSpaceType &domainSpace,
const RangeSpaceType &rangeSpace,
116 const MatrixParameter& param = PetscMatrixParameter() )
117 : domainSpace_( domainSpace ),
118 rangeSpace_( rangeSpace ),
119 colSlaveDofs_( rangeSpace_ ),
120 rowSlaveDofs_( domainSpace_ ),
122 localMatrixStack_( *this ),
128 ~PetscLinearOperator ()
130 if( status_ != statNothing )
131 ::Dune::Petsc::MatDestroy( &petscMatrix_ );
136 ::Dune::Petsc::MatAssemblyBegin( petscMatrix_, MAT_FINAL_ASSEMBLY );
137 ::Dune::Petsc::MatAssemblyEnd ( petscMatrix_, MAT_FINAL_ASSEMBLY );
138 status_ = statAssembled;
141 const DomainSpaceType& domainSpace ()
const {
return domainSpace_; }
142 const RangeSpaceType& rangeSpace ()
const {
return rangeSpace_; }
144 void apply (
const DomainFunctionType &arg, RangeFunctionType &dest )
const 146 ::Dune::Petsc::MatMult( petscMatrix_, *arg.petscVec() , *dest.petscVec() );
149 void operator() (
const DomainFunctionType &arg, RangeFunctionType &dest )
const 156 reserve( SimpleStencil<DomainSpaceType,RangeSpaceType>(0) );
160 template <
class StencilType>
161 void reserve (
const StencilType &stencil)
163 if(sequence_ != domainSpace().sequence())
168 const PetscInt localRows =
169 rowDofMapping().numOwnedDofBlocks() * domainLocalBlockSize;
170 const PetscInt localCols =
171 colDofMapping().numOwnedDofBlocks() * rangeLocalBlockSize;
173 assert( domainLocalBlockSize == rangeLocalBlockSize );
175 ::Dune::Petsc::MatCreate( &petscMatrix_ );
178 if( domainLocalBlockSize > 1 )
180 bs = domainLocalBlockSize ;
181 ::Dune::Petsc::MatSetType( petscMatrix_, MATBAIJ );
183 ::Dune::Petsc::MatSetBlockSize( petscMatrix_, bs );
187 ::Dune::Petsc::MatSetType( petscMatrix_, MATAIJ );
191 ::Dune::Petsc::MatSetSizes( petscMatrix_, localRows, localCols, PETSC_DETERMINE, PETSC_DETERMINE );
193 if (std::is_same< StencilType,SimpleStencil<DomainSpaceType,RangeSpaceType> >::value)
194 ::Dune::Petsc::MatSetUp( petscMatrix_, bs, stencil.maxNonZerosEstimate() );
197 std::vector<int> d_nnz(localRows/bs,0);
198 std::vector<int> o_nnz(localRows/bs,0);
199 typedef typename StencilType::GlobalStencilType GlobalStencilType;
200 typedef typename GlobalStencilType::const_iterator StencilIteratorType;
201 const GlobalStencilType &glStencil = stencil.globalStencil();
202 StencilIteratorType end = glStencil.end();
203 for ( StencilIteratorType it = glStencil.begin(); it != end; ++it)
205 int femIndex = it->first;
206 if ( rowDofMapping().isSlave( femIndex ) )
continue;
212 typedef typename StencilType::LocalStencilType LocalStencilType;
213 typedef typename LocalStencilType::const_iterator LocalStencilIteratorType;
214 LocalStencilIteratorType endLocal = it->second.end();
215 for ( LocalStencilIteratorType itLocal = it->second.begin(); itLocal != endLocal; ++itLocal)
217 if (!rowDofMapping().isSlave( *itLocal ))
231 int petscIndex = rowDofMapping().localSlaveMapping( femIndex );
232 assert( petscIndex >= 0 );
233 assert( petscIndex < (
int ) d_nnz.size() );
234 d_nnz[petscIndex] = nzDiag;
235 o_nnz[petscIndex] = nzOff;
237 ::Dune::Petsc::MatSetUp( petscMatrix_, bs, &d_nnz[0], &o_nnz[0] );
239 sequence_ = domainSpace().sequence();
242 status_ = statAssembled;
247 ::Dune::Petsc::MatZeroEntries( petscMatrix_ );
251 ObjectType* newObject()
const 253 return new ObjectType( *
this, domainSpace_, rangeSpace_ );
257 LocalMatrixType localMatrix (
const RowEntityType &rowEntity,
const ColumnEntityType &colEntity )
const 259 return LocalMatrixType(localMatrixStack_,rowEntity,colEntity);
261 LocalColumnObjectType localColumn(
const ColumnEntityType &colEntity )
const 263 return LocalColumnObjectType ( *
this, colEntity );
266 template<
class LocalMatrix >
267 void addLocalMatrix (
const RowEntityType &domainEntity,
const ColumnEntityType &rangeEntity,
const LocalMatrix &localMat )
269 assert( status_==statAssembled || status_==statAdd );
270 setStatus( statAdd );
272 std::vector< PetscInt > r, c;
273 setupIndices( rangeSpace_.blockMapper(), rowDofMapping(), rangeEntity, rangeLocalBlockSize, r);
274 setupIndices( domainSpace_.blockMapper(), colDofMapping(), domainEntity, domainLocalBlockSize, c);
276 std::vector< PetscScalar > v( r.size() * c.size() );
277 for( std::size_t i =0 ; i< r.size(); ++i )
278 for( std::size_t j =0; j< c.size(); ++j )
279 v[ i * c.size() +j ] = localMat.get( i, j );
281 ::Dune::Petsc::MatSetValues( petscMatrix_, r.size(), r.data(), c.size(), c.data(), v.data(), ADD_VALUES );
282 setStatus( statAssembled );
285 template<
class LocalMatrix,
class Scalar >
286 void addScaledLocalMatrix (
const RowEntityType &domainEntity,
const ColumnEntityType &rangeEntity,
const LocalMatrix &localMat,
const Scalar &s )
288 assert( status_==statAssembled || status_==statAdd );
289 setStatus( statAdd );
291 std::vector< PetscInt > r, c;
292 setupIndices( rangeSpace_.blockMapper(), rowDofMapping(), rangeEntity, rangeLocalBlockSize, r);
293 setupIndices( domainSpace_.blockMapper(), colDofMapping(), domainEntity, domainLocalBlockSize, c);
295 std::vector< PetscScalar > v( r.size() * c.size() );
296 for( std::size_t i =0 ; i< r.size(); ++i )
297 for( std::size_t j =0; j< c.size(); ++j )
298 v[ i * c.size() +j ] = s * localMat.get( i, j );
300 ::Dune::Petsc::MatSetValues( petscMatrix_, r.size(), r.data(), c.size(), c.data(), v.data(), ADD_VALUES );
301 setStatus( statAssembled );
304 template<
class LocalMatrix >
305 void setLocalMatrix (
const RowEntityType &domainEntity,
const ColumnEntityType &rangeEntity,
const LocalMatrix &localMat )
307 assert( status_==statAssembled || status_==statInsert );
308 setStatus( statInsert );
310 std::vector< PetscInt > r, c;
311 setupIndices( rangeSpace_.blockMapper(), rowDofMapping(), rangeEntity, rangeLocalBlockSize, r);
312 setupIndices( domainSpace_.blockMapper(), colDofMapping(), domainEntity, domainLocalBlockSize, c);
314 std::vector< PetscScalar > v( r.size() * c.size() );
315 for( std::size_t i =0 ; i< r.size(); ++i )
316 for( std::size_t j =0; j< c.size(); ++j )
317 v[ i * c.size() +j ] = localMat.get( i, j );
319 ::Dune::Petsc::MatSetValues( petscMatrix_, r.size(), r.data(), c.size(), c.data(), v.data(), INSERT_VALUES );
321 setStatus( statAssembled );
325 template<
class LocalMatrix >
326 void getLocalMatrix (
const RowEntityType &domainEntity,
const ColumnEntityType &rangeEntity, LocalMatrix &localMat )
const 328 assert( status_==statAssembled || status_==statGet );
329 setStatus( statGet );
331 std::vector< PetscInt > r, c;
332 setupIndices( rangeSpace_.blockMapper(), rowDofMapping(), rangeEntity, rangeLocalBlockSize, r);
333 setupIndices( domainSpace_.blockMapper(), colDofMapping(), domainEntity, domainLocalBlockSize, c);
335 std::vector< PetscScalar > v( r.size() * c.size() );
336 ::Dune::Petsc::MatGetValues( petscMatrix_, r.size(), r.data(), c.size(), c.data(), v.data() );
338 for( std::size_t i =0 ; i< r.size(); ++i )
339 for( std::size_t j =0; j< c.size(); ++j )
340 localMat.set( i, j, v[ i * c.size() +j ] );
342 setStatus( statAssembled );
348 ::Dune::Petsc::MatView( petscMatrix_, PETSC_VIEWER_STDOUT_WORLD );
358 Mat& petscMatrix ()
const {
return petscMatrix_; }
361 const RowDofMappingType& rowDofMapping()
const {
return rowSlaveDofs_.dofMapping(); }
363 const ColDofMappingType& colDofMapping()
const {
return colSlaveDofs_.dofMapping(); }
366 PetscLinearOperator ();
368 void setStatus(
const Status &newstatus)
const 371 if (status_ != statAssembled && status_ != newstatus)
373 ::Dune::Petsc::MatAssemblyBegin( petscMatrix_, MAT_FLUSH_ASSEMBLY );
374 ::Dune::Petsc::MatAssemblyEnd ( petscMatrix_, MAT_FLUSH_ASSEMBLY );
381 template<
typename DofMapper,
typename PetscMapping,
typename Entity >
382 static void setupIndices (
const DofMapper &dofMapper,
const PetscMapping &petscMapping,
const Entity &entity,
383 PetscInt blockSize, std::vector< PetscInt > &indices )
385 const int blockDofs = dofMapper.numDofs( entity );
386 const int numDofs = blockDofs * blockSize;
388 indices.resize( numDofs );
390 auto functor = [ &petscMapping, blockSize, &indices ] ( PetscInt local, PetscInt global )
392 const PetscInt block = petscMapping.globalMapping( global );
393 for( PetscInt b = 0; b < blockSize; ++b )
394 indices[ local *blockSize + b ] = block * blockSize + b;
399 dofMapper.mapEach( entity, functor );
406 const DomainSpaceType &domainSpace_;
407 const RangeSpaceType &rangeSpace_;
408 ColPetscSlaveDofsType colSlaveDofs_;
409 RowPetscSlaveDofsType rowSlaveDofs_;
412 mutable Mat petscMatrix_;
414 mutable LocalMatrixStackType localMatrixStack_;
415 mutable Status status_;
423 template<
typename DomainFunction,
typename RangeFunction >
424 class PetscLinearOperator< DomainFunction, RangeFunction >::LocalMatrix
425 :
public LocalMatrixDefault< typename PetscLinearOperator< DomainFunction, RangeFunction >::LocalMatrixTraits >
427 typedef LocalMatrix ThisType;
428 typedef LocalMatrixDefault< typename PetscLinearOperator< DomainFunction, RangeFunction >::LocalMatrixTraits > BaseType;
430 typedef PetscLinearOperator< DomainFunction, RangeFunction > PetscLinearOperatorType;
434 typedef PetscInt DofIndexType;
435 typedef std::vector< DofIndexType > IndexVectorType;
436 typedef typename DomainFunction::DiscreteFunctionSpaceType DomainSpaceType;
437 typedef typename RangeFunction::DiscreteFunctionSpaceType RangeSpaceType;
438 typedef typename DomainSpaceType::BasisFunctionSetType DomainBasisFunctionSetType;
439 typedef typename RangeSpaceType::BasisFunctionSetType RangeBasisFunctionSetType;
441 enum { rangeBlockSize = RangeSpaceType::localBlockSize };
442 enum { domainBlockSize = DomainSpaceType ::localBlockSize };
447 template<
typename PetscMapping >
448 struct PetscAssignFunctor
450 explicit PetscAssignFunctor (
const PetscMapping &petscMapping, IndexVectorType &indices )
451 : petscMapping_( petscMapping ),
455 template<
typename T >
456 void operator() (
const std::size_t localIndex, T globalIndex ) { indices_[ localIndex ] = petscMapping_.globalMapping( globalIndex ); }
459 const PetscMapping &petscMapping_;
460 IndexVectorType &indices_;
465 LocalMatrix (
const PetscLinearOperatorType &petscLinOp,
466 const DomainSpaceType &domainSpace,
467 const RangeSpaceType &rangeSpace )
468 : BaseType( domainSpace, rangeSpace ),
469 petscLinearOperator_( petscLinOp ),
477 if (status_ == statAdd)
479 PetscScalar v[rows()][columns()];
480 for (
int r=0;r<rows();++r)
482 for (
int c=0;c<columns();++c)
484 int idx = c + r*columns();
485 v[r][c] = values_[idx];
488 ::Dune::Petsc::MatSetValues( petscMatrix(), rows(), &(rowIndices_[0]), columns(), &(colIndices_[0]),
489 &v[0][0], ADD_VALUES );
493 void init (
const RowEntityType &domainEntity,
const ColumnEntityType &rangeEntity )
496 BaseType :: init( domainEntity, rangeEntity );
507 setupIndices( rangeSpace().blockMapper(), petscLinearOperator_.colDofMapping(), rangeEntity, rangeBlockSize, rowIndices_ );
510 setupIndices( domainSpace().blockMapper(), petscLinearOperator_.rowDofMapping(), domainEntity, domainBlockSize, colIndices_ );
513 values_.resize( columns()*rows(), 0. );
514 for (
int r=0;r<rows();++r)
516 for (
int c=0;c<columns();++c)
518 int idx = c + r*columns();
522 status_ = statAssembled;
523 petscLinearOperator_.setStatus(status_);
526 inline void add (
const int localRow,
const int localCol,
const RangeFieldType &value )
528 assert( status_==statAssembled || status_==statAdd );
530 petscLinearOperator_.setStatus(status_);
531 ::Dune::Petsc::MatSetValue( petscMatrix(), globalRowIndex( localRow ), globalColIndex( localCol ) , value, ADD_VALUES );
536 inline void set(
const int localRow,
const int localCol,
const RangeFieldType &value )
538 assert( status_==statAssembled || status_==statInsert );
539 status_ = statInsert;
540 petscLinearOperator_.setStatus(status_);
541 ::Dune::Petsc::MatSetValue( petscMatrix(), globalRowIndex( localRow ), globalColIndex( localCol ) , value, INSERT_VALUES );
547 void clearRow (
const int localRow )
549 assert( status_==statAssembled || status_==statInsert );
550 status_ = statInsert;
551 petscLinearOperator_.setStatus(status_);
552 const int col = this->columns();
553 for(
int localCol=0; localCol<col; ++localCol)
554 ::Dune::Petsc::MatSetValue( petscMatrix(), globalRowIndex( localRow ), globalColIndex( localCol ) ,
555 (localCol==localRow)?1:0., INSERT_VALUES );
563 inline const RangeFieldType
get (
const int localRow,
const int localCol )
const 565 assert( status_==statAssembled || status_==statGet );
567 petscLinearOperator_.setStatus(status_);
569 const int r[] = {globalRowIndex( localRow )};
570 const int c[] = {globalColIndex( localCol )};
571 ::Dune::Petsc::MatGetValues( petscMatrix(), 1, r, 1, c, v );
578 Mat& petscMatrix () {
return petscLinearOperator_.petscMatrix_; }
579 const Mat& petscMatrix ()
const {
return petscLinearOperator_.petscMatrix_; }
582 const int rows()
const {
return rowIndices_.size(); }
583 const int columns()
const {
return colIndices_.size(); }
587 DofIndexType globalRowIndex(
const int localRow )
const 589 assert( localRow < static_cast< int >( rowIndices_.size() ) );
590 return rowIndices_[ localRow ];
593 DofIndexType globalColIndex(
const int localCol )
const 595 assert( localCol < static_cast< int >( colIndices_.size() ) );
596 return colIndices_[ localCol ];
602 const PetscLinearOperatorType &petscLinearOperator_;
603 IndexVectorType rowIndices_;
604 IndexVectorType colIndices_;
605 std::vector<RangeFieldType> values_;
606 mutable Status status_;
614 #endif // #if defined HAVE_PETSC 616 #endif // #ifndef DUNE_FEM_PETSCLINEAROPERATOR_HH
Definition: coordinate.hh:4