dune-fem  2.4.1-rc
petscvector.hh
Go to the documentation of this file.
1 // vim: set expandtab ts=2 sw=2 sts=2:
2 #ifndef DUNE_FEM_PETSCVECTOR_HH
3 #define DUNE_FEM_PETSCVECTOR_HH
4 
6 
7 #if HAVE_PETSC
8 
14 
15 
16 namespace Dune
17 {
18 
19  namespace Fem
20  {
21 
22  // forward declarations
23  template< typename > class PetscDofBlock;
24  template< typename > class PetscDofProxy;
25  template< class DFSpace > class PetscVector;
26 
28  template < class DiscreteFunctionSpace, class Mapper >
29  class PetscManagedDofStorage
30  : public ManagedDofStorageImplementation< typename DiscreteFunctionSpace :: GridType,
31  Mapper,
32  PetscVector< DiscreteFunctionSpace > >
33  {
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;
39  protected:
40  DofArrayType myArray_;
41  public:
43  PetscManagedDofStorage( const DiscreteFunctionSpace& space,
44  const MapperType& mapper,
45  const std::string& name )
46  : BaseType( space.grid(), mapper, name, myArray_ ),
47  myArray_( space )
48  {
49  }
50  };
51 
52 
55  template< class DFS >
56  struct SpecialArrayFeatures< PetscVector< DFS > >
57  {
58  typedef PetscVector< DFS > ArrayType ;
60  typedef typename ArrayType :: value_type ValueType;
61 
63  static size_t used(const ArrayType & array)
64  {
65  return array.size() * sizeof(ValueType);
66  }
67 
69  static void setMemoryFactor(ArrayType & array, const double memFactor)
70  {
71  // do nothing here
72  }
73 
75  static void memMoveBackward(ArrayType& array, const int length,
76  const int oldStartIdx, const int newStartIdx)
77  {
78  DUNE_THROW(NotImplemented,"memMoveBackward is to be implemented");
79  }
80 
82  static void memMoveForward(ArrayType& array, const int length,
83  const int oldStartIdx, const int newStartIdx)
84  {
85  DUNE_THROW(NotImplemented,"memMoveForward is to be implemented");
86  }
87 
88  static void assign( ArrayType& array, const int newIndex, const int oldIndex )
89  {
90  /*
91  typedef typename ArrayType :: DofBlockPtrType DofBlockPtrType;
92  DofBlockPtrType newBlock = array.block( newIndex );
93  DofBlockPtrType oldBlock = array.block( oldIndex );
94 
95  const unsigned int blockSize = ArrayType :: blockSize;
96  for( unsigned int i = 0; i < blockSize; ++i )
97  (*newBlock)[ i ] = (*oldBlock)[ i ];
98  */
99  }
100  };
101 
102 
103  /* ========================================
104  * class PetscVector
105  * =======================================
106  */
107  /*
108  * This encapsules a PETSc Vec with ghosts.
109  * Some conceptual explanations:
110  * The PETSc vector, as modelled by this class, consists of three parts:
111  * 1) the whole PETSc vector, which might be distributed across several MPI processes.
112  * We call this the _global vector_
113  * 2) Each process holds a portion of this global vector, we call this part the
114  * _process vector_.
115  * 3) And there is a representation of the process vector, which also has 'ghost dof blocks'.
116  * We call this represenation the _ghosted vector_.
117  */
118  template< class DFSpace >
119  class PetscVector
120  {
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 >;
126  public:
127  typedef PetscSlaveDofProvider< DFSpace > PetscSlaveDofsType;
128  typedef PetscScalar value_type ;
129  typedef Vec DofContainerType;
130 
131  static const int blockSize = DFSpace :: localBlockSize;
132  typedef typename PetscSlaveDofsType :: PetscDofMappingType PetscDofMappingType;
133 
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;
141 
142  typedef DofIteratorType IteratorType;
143  typedef ConstDofIteratorType ConstIteratorType;
144  typedef typename DFSpace::RangeFieldType FieldType;
145  typedef int SizeType;
146 
147  typedef typename DFSpace::template CommDataHandle<void>::OperationType CommunicationOperationType;
148 
149  PetscVector ( const DFSpace& dfSpace )
150  : petscSlaveDofs_( dfSpace ),
151  memorySequence_( 0 ),
152  sequence_( 0 ),
153  communicateFlag_( false ),
154  localSize_( 0 ),
155  numGhosts_( 0 )
156  {
157  static_assert( CommunicationOperationType::value == DFCommunicationOperation::copy ||
158  CommunicationOperationType::value == DFCommunicationOperation::add,
159  "only copy/add are available communication operations for petsc");
160  // init vector
161  init();
162  }
163 
164  // TODO: think about sequence overflows...
165  PetscVector ( const ThisType &other )
166  : petscSlaveDofs_( other.petscSlaveDofs_.space() ),
167  memorySequence_( 0 ),
168  sequence_( 0 ),
169  communicateFlag_( false ),
170  localSize_( other.localSize_ ),
171  numGhosts_( other.numGhosts_ )
172  {
173  // assign vectors
174  assign( other );
175  }
176 
177  ~PetscVector ()
178  {
179  // destroy vectors
180  removeObj();
181  }
182 
183  size_t size () const { return dofMapping().size(); }
184  // static_cast< size_t >( localSize_ + numGhosts_ ) / blockSize; }
185 
186  void resize( const size_t newsize )
187  {
188  /*
189  std::vector< double > values( dofMapping().size() );
190  typedef typename std::vector< double > :: iterator iterator ;
191 
192  const DofIteratorType end = dend ();
193  iterator value = values.begin();
194  for( DofIteratorType it = dbegin(); it != end ; ++ it, ++value )
195  {
196  if( value == values.end() ) break ;
197  assert( value != values.end() );
198  *value = *it ;
199  }
200  */
201 
202  // TODO: keep old data stored in current vector
203  // remove old vectors
204  removeObj();
205 
206  // initialize new
207  init();
208 
209  /*
210  const size_t vsize = std::min( values.size(), size() );
211  DofIteratorType it = dbegin();
212  for( size_t i=0; i<vsize; ++ i, ++ it )
213  *it = values[ i ];
214  */
215 
216  hasBeenModified ();
217  }
218 
219  void reserve( const size_t capacity )
220  {
221  resize( capacity );
222  }
223 
224  void hasBeenModified () { ++sequence_; }
225 
226  void communicate ()
227  {
228  communicateFlag_ = true;
229  }
230 
231  // accessors for the underlying PETSc vectors
232  Vec* getVector ()
233  {
234  communicateIfNecessary();
235  return &vec_;
236  }
237 
238  const Vec* getVector () const
239  {
240  communicateIfNecessary();
241  return &vec_;
242  }
243 
244  // accessors for the underlying PETSc vectors
245  Vec& array ()
246  {
247  communicateIfNecessary();
248  return vec_;
249  }
250 
251  const Vec& array () const
252  {
253  communicateIfNecessary();
254  return vec_;
255  }
256 
257  Vec* getGhostedVector ()
258  {
259  communicateIfNecessary();
260  return &ghostedVec_;
261  }
262 
263  const Vec* getGhostedVector () const
264  {
265  communicateIfNecessary();
266  return &ghostedVec_;
267  }
268 
269  // force communication _now_
270  void communicateNow () const
271  {
272  communicateFlag_ = true;
273  ++sequence_;
274  communicateIfNecessary();
275  }
276 
277  DofBlockType operator[] ( const IndexType index ) { return DofBlockType( *this,index ); }
278  ConstDofBlockType operator[] ( const IndexType index ) const { return ConstDofBlockType( *this,index ); }
279 
280  ConstDofBlockPtrType block ( IndexType index ) const { return blockPtr( index ); }
281  DofBlockPtrType block ( IndexType index ) { return blockPtr( index ); }
282 
283  DofBlockPtrType blockPtr ( IndexType index )
284  {
285  assert( index < dofMapping().size() );
286  return DofBlockPtrType( typename DofBlockType::UnaryConstructorParamType( *this, index ) );
287  }
288 
289  ConstDofBlockPtrType blockPtr ( IndexType index ) const
290  {
291  assert( index < dofMapping().size() );
292  return ConstDofBlockPtrType( typename ConstDofBlockType::UnaryConstructorParamType( *this, index ) );
293  }
294 
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() ); }
299 
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() ); }
304 
305  void clear ()
306  {
307  ::Dune::Petsc::VecSet( *getVector(), 0. );
308  updateGhostRegions();
309  vectorIsUpToDateNow();
310  }
311 
312  PetscScalar operator* ( const ThisType &other ) const
313  {
314  PetscScalar ret;
315  ::Dune::Petsc::VecDot( *getVector(), *other.getVector(), &ret );
316  return ret;
317  }
318 
319  const ThisType& operator+= ( const ThisType &other )
320  {
321  ::Dune::Petsc::VecAXPY( *getVector(), 1., *other.getVector() );
322  updateGhostRegions();
323  vectorIsUpToDateNow();
324  return *this;
325  }
326 
327  const ThisType& operator-= ( const ThisType &other )
328  {
329  ::Dune::Petsc::VecAXPY( *getVector(), -1., *other.getVector() );
330  updateGhostRegions();
331  vectorIsUpToDateNow();
332  return *this;
333  }
334 
335  const ThisType& operator*= ( PetscScalar scalar )
336  {
337  ::Dune::Petsc::VecScale( *getVector(), scalar );
338  updateGhostRegions();
339  vectorIsUpToDateNow();
340  return *this;
341  }
342 
343  const ThisType& operator/= ( PetscScalar scalar )
344  {
345  assert( scalar != 0 );
346  return this->operator*=( 1./scalar );
347  }
348 
349  void axpy ( const PetscScalar &scalar, const ThisType &other )
350  {
351  ::Dune::Petsc::VecAXPY( *getVector(), scalar, *other.getVector() );
352  hasBeenModified();
353  }
354 
355  // debugging; comes in handy to call these 2 methods in gdb
356  // doit is only here to prevent the compiler from optimizing these calls away...
357  void printGlobal ( bool doit )
358  {
359  if( !doit )
360  return;
361  VecView( vec_, PETSC_VIEWER_STDOUT_WORLD );
362  }
363 
364  void printGhost ( bool doit)
365  {
366  if( !doit )
367  return;
368 
369  PetscScalar *array;
370  VecGetArray( ghostedVec_,&array );
371  for( int i=0; i < localSize_ + numGhosts_; i++ )
372  {
373  PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%D %G\n",i,PetscRealPart(array[i]));
374  }
375  VecRestoreArray( ghostedVec_, &array );
376 #if PETSC_VERSION_MAJOR <= 3 && PETSC_VERSION_MINOR < 5
377  PetscSynchronizedFlush( PETSC_COMM_WORLD );
378 #else
379  PetscSynchronizedFlush( PETSC_COMM_WORLD, PETSC_STDOUT );
380 #endif
381  }
382 
383  // assign from other given PetscVector
384  void assign( const ThisType& other )
385  {
386  // we want the 'other' to do all its communication right now before
387  // we start copying values from it
388  other.communicateIfNecessary();
389 
390  // Do the copying on the PETSc level
391  ::Dune::Petsc::VecDuplicate( other.vec_, &vec_ );
392  ::Dune::Petsc::VecCopy( other.vec_, vec_ );
393  ::Dune::Petsc::VecGhostGetLocalForm( vec_, &ghostedVec_ );
394 
395  updateGhostRegions();
396  }
397 
398  PetscVector& operator= ( const ThisType& other )
399  {
400  assign( other );
401  return *this;
402  }
403  protected:
404  // setup vector according to mapping sizes
405  void init()
406  {
407  // set up the DofMapping instance and all variables depending on it
408  localSize_ = dofMapping().numOwnedDofBlocks() * blockSize;
409  numGhosts_ = dofMapping().numSlaveBlocks() * blockSize;
410 
411  //std::cout << "PetscVector::init: "<< localSize_ << " " << numGhosts_ << std::endl;
412  assert( static_cast< size_t >( localSize_ + numGhosts_ ) == dofMapping().size() * blockSize );
413 
414  // set up the ghost array builder
415  typedef PetscGhostArrayBuilder< PetscSlaveDofsType, PetscDofMappingType > PetscGhostArrayBuilderType;
416  PetscGhostArrayBuilderType ghostArrayBuilder( petscSlaveDofs_, dofMapping() );
417  assert( int( ghostArrayBuilder.size() ) == dofMapping().numSlaveBlocks() );
418 
419  // finally, create the PETSc Vecs
420  ::Dune::Petsc::VecCreateGhost(
421  static_cast< PetscInt >( localSize_ ),
422  PETSC_DECIDE,
423  static_cast< PetscInt >( numGhosts_ ),
424  ghostArrayBuilder.array(),
425  &vec_
426  );
427  ::Dune::Petsc::VecGhostGetLocalForm( vec_, &ghostedVec_ );
428  }
429 
430  // delete vectors
431  void removeObj()
432  {
433  ::Dune::Petsc::VecGhostRestoreLocalForm( vec_, &ghostedVec_ );
434  ::Dune::Petsc::VecDestroy( &vec_ );
435  }
436 
437  PetscVector ();
438 
439  PetscDofMappingType& dofMapping () { return petscSlaveDofs_.dofMapping(); }
440  const PetscDofMappingType& dofMapping () const { return petscSlaveDofs_.dofMapping(); }
441 
442  void communicateIfNecessary () const
443  {
444  // communicate this process' values
445  if( communicateFlag_ && memorySequence_ < sequence_ )
446  {
447  if ( memorySequence_ < sequence_ )
448  {
449  if ( CommunicationOperationType::value == DFCommunicationOperation::add )
450  {
451  ::Dune::Petsc::VecGhostUpdateBegin( vec_, ADD_VALUES, SCATTER_REVERSE );
452  ::Dune::Petsc::VecGhostUpdateEnd( vec_, ADD_VALUES, SCATTER_REVERSE );
453  }
454  ::Dune::Petsc::VecGhostUpdateBegin( vec_, INSERT_VALUES, SCATTER_FORWARD );
455  ::Dune::Petsc::VecGhostUpdateEnd( vec_, INSERT_VALUES, SCATTER_FORWARD );
456 
457  memorySequence_ = sequence_;
458  }
459 
460  communicateFlag_ = false;
461  }
462  }
463 
464  // Updates the ghost dofs, obtains them from the owning process
465  void updateGhostRegions ()
466  {
467  ::Dune::Petsc::VecGhostUpdateBegin( vec_, INSERT_VALUES, SCATTER_FORWARD );
468  ::Dune::Petsc::VecGhostUpdateEnd( vec_, INSERT_VALUES, SCATTER_FORWARD );
469  }
470 
471  void vectorIsUpToDateNow () const
472  {
473  memorySequence_ = sequence_;
474  communicateFlag_ = false;
475  }
476 
477  /*
478  * data fields
479  */
480  PetscSlaveDofsType petscSlaveDofs_;
481  Vec vec_;
482  Vec ghostedVec_;
483 
484  mutable unsigned long memorySequence_; // represents the state of the PETSc vec in memory
485  mutable unsigned long sequence_; // represents the the modifications to the PETSc vec
486 
487  mutable bool communicateFlag_;
488  PetscInt localSize_;
489  PetscInt numGhosts_;
490 
491  };
492 
493  } // namespace Fem
494 
495 } // namespace Dune
496 
497 #endif // #if HAVE_PETSC
498 
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
discrete function space