2 #ifndef DUNE_FEM_PETSCVECTOR_HH 3 #define DUNE_FEM_PETSCVECTOR_HH 23 template<
typename >
class PetscDofBlock;
24 template<
typename >
class PetscDofProxy;
25 template<
class DFSpace >
class PetscVector;
28 template <
class DiscreteFunctionSpace,
class Mapper >
29 class PetscManagedDofStorage
30 :
public ManagedDofStorageImplementation< typename DiscreteFunctionSpace :: GridType,
32 PetscVector< DiscreteFunctionSpace > >
34 typedef typename DiscreteFunctionSpace :: GridType GridType;
35 typedef PetscSlaveDofProvider< DiscreteFunctionSpace > PetscSlaveDofProviderType;
36 typedef Mapper MapperType ;
37 typedef PetscVector< DiscreteFunctionSpace > DofArrayType ;
38 typedef ManagedDofStorageImplementation< GridType, MapperType, DofArrayType > BaseType;
40 DofArrayType myArray_;
44 const MapperType& mapper,
45 const std::string& name )
46 : BaseType( space.grid(), mapper, name, myArray_ ),
56 struct SpecialArrayFeatures< PetscVector< DFS > >
58 typedef PetscVector< DFS > ArrayType ;
60 typedef typename ArrayType :: value_type ValueType;
63 static size_t used(
const ArrayType & array)
65 return array.size() *
sizeof(ValueType);
69 static void setMemoryFactor(ArrayType & array,
const double memFactor)
75 static void memMoveBackward(ArrayType& array,
const int length,
76 const int oldStartIdx,
const int newStartIdx)
78 DUNE_THROW(NotImplemented,
"memMoveBackward is to be implemented");
82 static void memMoveForward(ArrayType& array,
const int length,
83 const int oldStartIdx,
const int newStartIdx)
85 DUNE_THROW(NotImplemented,
"memMoveForward is to be implemented");
88 static void assign( ArrayType& array,
const int newIndex,
const int oldIndex )
118 template<
class DFSpace >
121 typedef PetscVector< DFSpace > ThisType;
122 friend class PetscDofBlock< ThisType >;
123 friend class PetscDofBlock< const ThisType >;
124 friend class PetscDofProxy< ThisType >;
125 friend class PetscDofProxy< const ThisType >;
127 typedef PetscSlaveDofProvider< DFSpace > PetscSlaveDofsType;
128 typedef PetscScalar value_type ;
129 typedef Vec DofContainerType;
131 static const int blockSize = DFSpace :: localBlockSize;
132 typedef typename PetscSlaveDofsType :: PetscDofMappingType PetscDofMappingType;
134 typedef PetscDofBlock< ThisType > DofBlockType;
135 typedef PetscDofBlock< const ThisType > ConstDofBlockType;
136 typedef typename DofBlockType::DofIterator DofIteratorType;
137 typedef typename ConstDofBlockType::DofIterator ConstDofIteratorType;
138 typedef Envelope< DofBlockType > DofBlockPtrType;
139 typedef Envelope< ConstDofBlockType > ConstDofBlockPtrType;
140 typedef typename DofBlockType::IndexType IndexType;
142 typedef DofIteratorType IteratorType;
143 typedef ConstDofIteratorType ConstIteratorType;
144 typedef typename DFSpace::RangeFieldType FieldType;
145 typedef int SizeType;
147 typedef typename DFSpace::template CommDataHandle<void>::OperationType CommunicationOperationType;
149 PetscVector (
const DFSpace& dfSpace )
150 : petscSlaveDofs_( dfSpace ),
151 memorySequence_( 0 ),
153 communicateFlag_( false ),
159 "only copy/add are available communication operations for petsc");
165 PetscVector (
const ThisType &other )
166 : petscSlaveDofs_( other.petscSlaveDofs_.space() ),
167 memorySequence_( 0 ),
169 communicateFlag_( false ),
170 localSize_( other.localSize_ ),
171 numGhosts_( other.numGhosts_ )
183 size_t size ()
const {
return dofMapping().size(); }
186 void resize(
const size_t newsize )
219 void reserve(
const size_t capacity )
224 void hasBeenModified () { ++sequence_; }
228 communicateFlag_ =
true;
234 communicateIfNecessary();
238 const Vec* getVector ()
const 240 communicateIfNecessary();
247 communicateIfNecessary();
251 const Vec& array ()
const 253 communicateIfNecessary();
257 Vec* getGhostedVector ()
259 communicateIfNecessary();
263 const Vec* getGhostedVector ()
const 265 communicateIfNecessary();
270 void communicateNow ()
const 272 communicateFlag_ =
true;
274 communicateIfNecessary();
277 DofBlockType operator[] (
const IndexType index ) {
return DofBlockType( *
this,index ); }
278 ConstDofBlockType operator[] (
const IndexType index )
const {
return ConstDofBlockType( *
this,index ); }
280 ConstDofBlockPtrType block ( IndexType index )
const {
return blockPtr( index ); }
281 DofBlockPtrType block ( IndexType index ) {
return blockPtr( index ); }
283 DofBlockPtrType blockPtr ( IndexType index )
285 assert( index < dofMapping().size() );
286 return DofBlockPtrType(
typename DofBlockType::UnaryConstructorParamType( *
this, index ) );
289 ConstDofBlockPtrType blockPtr ( IndexType index )
const 291 assert( index < dofMapping().size() );
292 return ConstDofBlockPtrType(
typename ConstDofBlockType::UnaryConstructorParamType( *
this, index ) );
295 DofIteratorType dbegin () {
return DofIteratorType( *
this, 0, 0 ); }
296 ConstDofIteratorType dbegin ()
const {
return ConstDofIteratorType( *
this, 0, 0 ); }
297 DofIteratorType dend() {
return DofIteratorType( *
this, dofMapping().size() ); }
298 ConstDofIteratorType dend()
const {
return ConstDofIteratorType( *
this, dofMapping().size() ); }
300 DofIteratorType begin () {
return DofIteratorType( *
this, 0, 0 ); }
301 ConstDofIteratorType begin ()
const {
return ConstDofIteratorType( *
this, 0, 0 ); }
302 DofIteratorType end() {
return DofIteratorType( *
this, dofMapping().size() ); }
303 ConstDofIteratorType end()
const {
return ConstDofIteratorType( *
this, dofMapping().size() ); }
307 ::Dune::Petsc::VecSet( *getVector(), 0. );
308 updateGhostRegions();
309 vectorIsUpToDateNow();
312 PetscScalar
operator* (
const ThisType &other )
const 315 ::Dune::Petsc::VecDot( *getVector(), *other.getVector(), &ret );
319 const ThisType& operator+= (
const ThisType &other )
321 ::Dune::Petsc::VecAXPY( *getVector(), 1., *other.getVector() );
322 updateGhostRegions();
323 vectorIsUpToDateNow();
327 const ThisType& operator-= (
const ThisType &other )
329 ::Dune::Petsc::VecAXPY( *getVector(), -1., *other.getVector() );
330 updateGhostRegions();
331 vectorIsUpToDateNow();
335 const ThisType& operator*= ( PetscScalar scalar )
337 ::Dune::Petsc::VecScale( *getVector(), scalar );
338 updateGhostRegions();
339 vectorIsUpToDateNow();
343 const ThisType& operator/= ( PetscScalar scalar )
345 assert( scalar != 0 );
346 return this->operator*=( 1./scalar );
349 void axpy (
const PetscScalar &scalar,
const ThisType &other )
351 ::Dune::Petsc::VecAXPY( *getVector(), scalar, *other.getVector() );
357 void printGlobal (
bool doit )
361 VecView( vec_, PETSC_VIEWER_STDOUT_WORLD );
364 void printGhost (
bool doit)
370 VecGetArray( ghostedVec_,&array );
371 for(
int i=0; i < localSize_ + numGhosts_; i++ )
373 PetscSynchronizedPrintf(PETSC_COMM_WORLD,
"%D %G\n",i,PetscRealPart(array[i]));
375 VecRestoreArray( ghostedVec_, &array );
376 #if PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 5 377 PetscSynchronizedFlush( PETSC_COMM_WORLD );
379 PetscSynchronizedFlush( PETSC_COMM_WORLD, PETSC_STDOUT );
384 void assign(
const ThisType& other )
388 other.communicateIfNecessary();
391 ::Dune::Petsc::VecDuplicate( other.vec_, &vec_ );
392 ::Dune::Petsc::VecCopy( other.vec_, vec_ );
393 ::Dune::Petsc::VecGhostGetLocalForm( vec_, &ghostedVec_ );
395 updateGhostRegions();
398 PetscVector& operator= (
const ThisType& other )
408 localSize_ = dofMapping().numOwnedDofBlocks() * blockSize;
409 numGhosts_ = dofMapping().numSlaveBlocks() * blockSize;
412 assert( static_cast< size_t >( localSize_ + numGhosts_ ) == dofMapping().size() * blockSize );
415 typedef PetscGhostArrayBuilder< PetscSlaveDofsType, PetscDofMappingType > PetscGhostArrayBuilderType;
416 PetscGhostArrayBuilderType ghostArrayBuilder( petscSlaveDofs_, dofMapping() );
417 assert(
int( ghostArrayBuilder.size() ) == dofMapping().numSlaveBlocks() );
420 ::Dune::Petsc::VecCreateGhost(
421 static_cast< PetscInt >( localSize_ ),
423 static_cast< PetscInt >( numGhosts_ ),
424 ghostArrayBuilder.array(),
427 ::Dune::Petsc::VecGhostGetLocalForm( vec_, &ghostedVec_ );
433 ::Dune::Petsc::VecGhostRestoreLocalForm( vec_, &ghostedVec_ );
434 ::Dune::Petsc::VecDestroy( &vec_ );
439 PetscDofMappingType& dofMapping () {
return petscSlaveDofs_.dofMapping(); }
440 const PetscDofMappingType& dofMapping ()
const {
return petscSlaveDofs_.dofMapping(); }
442 void communicateIfNecessary ()
const 445 if( communicateFlag_ && memorySequence_ < sequence_ )
447 if ( memorySequence_ < sequence_ )
451 ::Dune::Petsc::VecGhostUpdateBegin( vec_, ADD_VALUES, SCATTER_REVERSE );
452 ::Dune::Petsc::VecGhostUpdateEnd( vec_, ADD_VALUES, SCATTER_REVERSE );
454 ::Dune::Petsc::VecGhostUpdateBegin( vec_, INSERT_VALUES, SCATTER_FORWARD );
455 ::Dune::Petsc::VecGhostUpdateEnd( vec_, INSERT_VALUES, SCATTER_FORWARD );
457 memorySequence_ = sequence_;
460 communicateFlag_ =
false;
465 void updateGhostRegions ()
467 ::Dune::Petsc::VecGhostUpdateBegin( vec_, INSERT_VALUES, SCATTER_FORWARD );
468 ::Dune::Petsc::VecGhostUpdateEnd( vec_, INSERT_VALUES, SCATTER_FORWARD );
471 void vectorIsUpToDateNow ()
const 473 memorySequence_ = sequence_;
474 communicateFlag_ =
false;
480 PetscSlaveDofsType petscSlaveDofs_;
484 mutable unsigned long memorySequence_;
485 mutable unsigned long sequence_;
487 mutable bool communicateFlag_;
497 #endif // #if HAVE_PETSC 499 #endif // DUNE_FEM_PETSCVECTOR_HH Definition: commoperations.hh:261
void axpy(const T &a, const T &x, T &y)
Definition: space/basisfunctionset/functor.hh:37
Definition: coordinate.hh:4
Double operator*(const Double &a, const Double &b)
Definition: double.hh:495
Definition: commoperations.hh:261