dune-fem  2.4.1-rc
petscdofblock.hh
Go to the documentation of this file.
1 // vim: set expandtab ts=2 sw=2 sts=2:
2 #ifndef DUNE_FEM_PETSCDOFBLOCK_HH
3 #define DUNE_FEM_PETSCDOFBLOCK_HH
4 
5 #include <algorithm>
6 #include <memory>
7 
9 
10 #if HAVE_PETSC
11 
14 
16 
17 namespace Dune
18 {
19 
20  namespace Fem
21  {
22  template < class PVector >
23  class PetscDofBlock;
24 
25  /* ========================================
26  * class PetscDofBlock::DofProxy
27  * =======================================
28  */
29  template< class PVector >
30  class PetscDofProxy
31  {
32  public:
33  typedef PVector PetscVectorType;
34  typedef typename PetscDofBlock< PetscVectorType >::DofProxy ThisType;
35  typedef typename PetscDofBlock< PetscVectorType >::IndexType IndexType;
36 
37  static const unsigned int blockSize = PetscVectorType::blockSize;
38 
39  // This is needed to put DofProxies in STL (or STL-like) containers...
40  PetscDofProxy ( PetscScalar s = 0 )
41  : petscVector_( 0 ),
42  blockIndex_( 0 ),
43  indexInBlock_( 0 )
44  {}
45 
46  // this is called by a friend
47  PetscDofProxy ( PetscVectorType &petscVector, IndexType blockIndex, PetscInt indexInBlock )
48  : petscVector_( &petscVector ),
49  blockIndex_( blockIndex ),
50  indexInBlock_( indexInBlock )
51  {}
52 
53  // Does not change the value of the DoF which this proxy, but it assigns this proxy instance...
54  void assign ( const ThisType &other )
55  {
56  petscVector_ = other.petscVector_;
57  blockIndex_ = other.blockIndex_;
58  indexInBlock_ = other.indexInBlock_;
59  }
60 
61  const ThisType& operator= ( ThisType &other )
62  {
63  setValue( other.getValue() );
64  return *this;
65  }
66 
67  const ThisType& operator= ( PetscScalar val )
68  {
69  setValue( val );
70  return *this;
71  }
72 
73  const ThisType& operator*= ( const ThisType& other )
74  {
75  PetscScalar value = getValue() * other.getValue();
76  setValue( value );
77  return *this;
78  }
79 
80  const ThisType& operator*= ( const PetscScalar scalar )
81  {
82  PetscScalar value = getValue() * scalar ;
83  setValue( value );
84  return *this;
85  }
86 
87  const ThisType& operator+= ( const ThisType &other )
88  {
89  setValue( other.getValue(), ADD_VALUES );
90  return *this;
91  }
92 
93  const ThisType& operator+= ( const PetscScalar scalar )
94  {
95  setValue( scalar, ADD_VALUES );
96  return *this;
97  }
98 
99  const ThisType& operator-= ( const ThisType &other )
100  {
101  setValue( -other.getValue(), ADD_VALUES );
102  return *this;
103  }
104 
105  const ThisType& operator-= ( const PetscScalar scalar )
106  {
107  setValue( -scalar, ADD_VALUES );
108  return *this;
109  }
110 
111  // conversion operators
112  operator PetscScalar () const
113  {
114  return valid() ? getValue() : PetscScalar( 0 );
115  }
116 
117  bool valid() const { return bool( petscVector_ ); }
118 
119  protected:
120  PetscScalar getValue () const
121  {
122  assert( valid() );
123  PetscInt index = blockSize*petscVector().dofMapping().localSlaveMapping( blockIndex_ ) + indexInBlock_;
124  PetscScalar ret;
125  ::Dune::Petsc::VecGetValues( *petscVector().getGhostedVector(), 1, &index, &ret );
126  return ret;
127  }
128 
129  void setValue ( const PetscScalar &val, InsertMode mode = INSERT_VALUES )
130  {
131  PetscInt index = blockSize*petscVector().dofMapping().localSlaveMapping( blockIndex_ ) + indexInBlock_;
132  ::Dune::Petsc::VecSetValue( *( petscVector().getGhostedVector() ), index, val, mode );
133  petscVector().hasBeenModified();
134  }
135 
136  PetscVectorType& petscVector ()
137  {
138  assert( petscVector_ );
139  return *petscVector_;
140  }
141 
142  const PetscVectorType& petscVector () const
143  {
144  assert( petscVector_ );
145  return *petscVector_;
146  }
147 
148  // data fields
149  PetscVectorType *petscVector_;
150  IndexType blockIndex_;
151  PetscInt indexInBlock_;
152 
153  };
154 
155 
156  // TODO: give this the same interface as Dune::Fem::ReferenceBlockVectorBlock?
157  /* ========================================
158  * class PetscDofBlock
159  * =======================================
160  */
161  template< class PVector >
162  class PetscDofBlock
163  {
164  typedef PetscDofBlock< PVector > ThisType;
165 
166  public:
167  typedef PVector PetscVectorType;
168  typedef typename PetscVectorType::PetscDofMappingType PetscDofMappingType;
169  typedef typename PetscDofMappingType::IndexType IndexType;
170 
171  static const unsigned int blockSize = PetscVectorType::blockSize;
172 
173  typedef PetscDofProxy< PVector > DofProxy;
174  class DofIterator;
175 
176  // Needed so that we can put this class in a Dune::Envelope
177  typedef std::pair< PetscVectorType&, IndexType > UnaryConstructorParamType;
178 
179  // this is the ctor to be used
180  PetscDofBlock ( PetscVectorType &petscVector, IndexType blockIndex )
181  : petscVector_( petscVector ),
182  blockIndex_( blockIndex )
183  {}
184 
185  // this is the ctor to be used
186  PetscDofBlock ( const PetscDofBlock& other )
187  : petscVector_( other.petscVector_ ),
188  blockIndex_( other.blockIndex_ )
189  {}
190 
191  // ..or this one, which is semantically equivalent to the above ctor
192  explicit PetscDofBlock ( UnaryConstructorParamType arg )
193  : petscVector_( arg.first ),
194  blockIndex_( arg.second )
195  {}
196 
197  const PetscDofBlock& operator*= ( const PetscScalar value )
198  {
199  for( unsigned int i=0; i<blockSize; ++i )
200  (*this)[ i ] *= value ;
201  return *this;
202  }
203 
204  private:
205  PetscDofBlock ();
206  // The standard copy ctor is okay... and needed.
207 
208  public:
209  IndexType size() const { return blockSize; }
210 
211  DofProxy operator [] ( unsigned int index )
212  {
213  assert( index < blockSize );
214  return DofProxy( petscVector_, blockIndex_, index );
215  }
216 
217  const DofProxy operator [] ( unsigned int index ) const
218  {
219  assert( index < blockSize );
220  return DofProxy( petscVector_, blockIndex_, index );
221  }
222 
223  private:
224  /*
225  * data fields
226  */
227  PetscVectorType &petscVector_;
228  IndexType blockIndex_;
229  };
230 
232  //template< class Traits, class PVector >
233  template< class Traits, class PVector >
234  inline OutStreamInterface< Traits >&
235  operator<< ( OutStreamInterface< Traits > &out,
236  const PetscDofProxy< PVector >& value )
237  {
238  // write to stream
239  out << PetscScalar( value );
240  return out;
241  }
242 
244  //template< class Traits, class PVector >
245  template< class Traits, class PVector >
246  inline InStreamInterface< Traits >&
247  operator>> ( InStreamInterface< Traits > &in,
248  PetscDofProxy< PVector > value )
249  {
250  PetscScalar val;
251  // get value from stream
252  in >> val;
253  // write value to discrete function
254  value = val;
255  return in;
256  }
257 
258 
259  /*
260  * This is almost a bidirectional iterator but does not completely satisfy the required
261  * interface (see the C++2003 standard, 24.1.4) [no default ctor, no operator->].
262  */
263 
264  /* ========================================
265  * class PetscDofBlock::DofIterator
266  * =======================================
267  */
268  template< class PVector >
269  class PetscDofBlock< PVector >::DofIterator
270  : public std::iterator< std::input_iterator_tag, PetscScalar >
271  {
272  typedef typename PetscDofBlock< PVector >::DofIterator ThisType;
273  typedef PetscDofBlock< PVector > DofBlockType;
274 
275  typedef std::iterator< std::input_iterator_tag, PetscScalar > BaseType;
276 
277  // TODO: get rid of this! we don't like shared pointers. Own a real instance instead!
278  typedef std::shared_ptr< DofBlockType > DofBlockSharedPointer;
279 
280  public:
281  typedef PVector PetscVectorType;
282  typedef typename DofBlockType::DofProxy value_type; // (this overrides the 2nd template parameter of BaseType...)
283 
284  // standard ctor
285  DofIterator ( PetscVectorType &petscVector, unsigned int blockIndex, PetscInt indexInBlock = 0 )
286  : petscVector_( petscVector ),
287  blockIndex_( blockIndex ),
288  indexInBlock_( indexInBlock )
289  {
290  // blockIndex == size denotes the end iterator
291  assert( blockIndex <= petscVector_.dofMapping().size() );
292 
293  // Is this not the end iterator?
294  if( blockIndex < petscVector_.dofMapping().size() )
295  {
296  resetBlockPtr();
297  }
298  }
299 
300  bool operator== ( const ThisType &other ) const
301  {
302  return blockIndex_ == other.blockIndex_ &&
303  indexInBlock_ == other.indexInBlock_;
304  }
305 
306  bool operator!= ( const ThisType &other ) const { return !this->operator==( other ); };
307 
308  value_type operator* () { return block()[ indexInBlock_]; }
309  const value_type operator* () const { return block()[ indexInBlock_ ]; }
310 
311  // prefix increment
312  ThisType& operator++ () { increment(); return *this; }
313 
314  // prefix decrement
315  ThisType& operator-- () { decrement(); return *this; }
316 
317  private:
318  // forbidden
319  DofIterator ();
320 
321  PetscVectorType& petscVector () { return petscVector_; }
322 
323  void increment ()
324  {
325  ++indexInBlock_;
326  if( static_cast< unsigned int >( indexInBlock_ ) >= DofBlockType::blockSize )
327  {
328  ++blockIndex_;
329  indexInBlock_ = 0;
330  if( static_cast< unsigned int >( blockIndex_ ) < petscVector().dofMapping().size() )
331  resetBlockPtr();
332  }
333  }
334 
335  void decrement ()
336  {
337  assert( blockIndex_ > 0 );
338  if( indexInBlock_ == 0 )
339  {
340  indexInBlock_ = DofBlockType::blockSize - 1;
341  --blockIndex_;
342  resetBlockPtr();
343  }
344  else
345  {
346  --blockIndex_;
347  }
348  }
349 
350  // TODO: iterator should own the block - _no_ new...
351  void resetBlockPtr () { blockPtr_.reset( new DofBlockType( *petscVector().block( blockIndex_ ) ) ); }
352 
353  DofBlockType& block () const { return *blockPtr_.get(); }
354 
355  /*
356  * data fields
357  */
358  PetscVectorType &petscVector_;
359  unsigned int blockIndex_;
360  PetscInt indexInBlock_;
361  DofBlockSharedPointer blockPtr_;
362  };
363 
364  } // namespace Fem
365 
366 } // namespace Dune
367 
368 #endif // #if HAVE_PETSC
369 
370 #endif // DUNE_FEM_PETSCDOFBLOCK_HH
InStreamInterface< StreamTraits > & operator>>(InStreamInterface< StreamTraits > &in, DiscreteFunctionInterface< Impl > &df)
read a discrete function from an input stream
Definition: discretefunction_inline.hh:395
bool operator!=(const Double &a, const Double &b)
Definition: double.hh:629
bool operator==(const Double &a, const Double &b)
Definition: double.hh:589
Definition: coordinate.hh:4
Double operator*(const Double &a, const Double &b)
Definition: double.hh:495