dune-fem  2.4.1-rc
petscghostarraybuilder.hh
Go to the documentation of this file.
1 // vim: set expandtab ts=2 sw=2 sts=2:
2 #ifndef DUNE_FEM_PETSCGHOSTARRAYBUILDER_HH
3 #define DUNE_FEM_PETSCGHOSTARRAYBUILDER_HH
4 
5 #include <dune/grid/common/datahandleif.hh>
6 
8 
9 #if HAVE_PETSC
10 
14 
15 
16 namespace Dune
17 {
18 
19  namespace Fem
20  {
21 
22  /*
23  * forward declarations
24  */
25  template< typename, typename > class PetscDofCommunicator;
26 
27 
28 
29  /* =================================================
30  * class PetscGhostArrayBuilder
31  *
32  * This class is responsable to build the array of indices that indicates the position
33  * of ghost dofs in the PETSc vec. Once an instance is initialized properly, one can
34  * access the array using array()
35  * =================================================
36  */
37  template< typename SlaveProvider, typename PDofMapping >
38  class PetscGhostArrayBuilder
39  {
40  typedef PetscGhostArrayBuilder< SlaveProvider, PDofMapping > ThisType;
41 
42  typedef std::vector< std::pair< int, int > > StorageType;
43  struct StorageSort
44  {
45  template< typename IT > bool operator () ( IT a, IT b ) { return a.first < b.first; }
46  };
47 
48  public:
49  typedef SlaveProvider SlaveProviderType;
50  typedef typename SlaveProviderType::DiscreteFunctionSpaceType DiscreteFunctionSpaceType;
51  typedef PDofMapping PetscDofMappingType;
52  typedef PetscInt value_type;
53  static const int blockSize = DiscreteFunctionSpaceType::localBlockSize;
54 
55 
56  explicit PetscGhostArrayBuilder ( SlaveProviderType& slaveDofs,
57  PetscDofMappingType& petscDofMapping )
58  : array_( 0 ),
59  arrayNeedsUpdate_( true )
60  {
61  if( slaveDofs.space().gridPart().comm().size() > 1 ) // YaspGrid and lagrange polorder = 2 fails unless we skip this in serial code
62  {
63  try
64  {
65  PetscDofCommunicator< DiscreteFunctionSpaceType, ThisType > handle( slaveDofs.space(), petscDofMapping, *this );
66  slaveDofs.space().gridPart().communicate( handle, DiscreteFunctionSpaceType::GridPartType::indexSetInterfaceType, ForwardCommunication );
67  }
68  catch( const Dune::Exception &e )
69  {
70  std::cerr << e << std::endl;
71  std::cerr << "Exception thrown in: " << __FILE__ << " line: " << __LINE__ << std::endl;
72  abort();
73  }
74  }
75 
76  // TODO: This is a corner case, but it might be possible that
77  // everything works fine but this assertion fails
78  //
79  // the test is being removed because it does not work correctly with
80  // dg spaces - needs to be modified first
81  // assert( !hasDuplicatedGhosts());
82  }
83 
84  ~PetscGhostArrayBuilder ()
85  {
86  delete[] array_;
87  }
88 
89  void clear ()
90  {
91  vector_.clear();
92  arrayNeedsUpdate_ = true;
93  }
94 
95  bool hasDuplicatedGhosts () const
96  {
97  std::map< int, bool > map;
98  for( StorageType::const_iterator it = vector_.begin(); it != vector_.end(); ++it )
99  {
100  if( map.find( it->second ) != map.end() )
101  std::cerr << it->second << " found more than once on rk " << MPIManager::rank() << "\n";
102 
103  map[ it->second ] = true;
104  }
105  return ( map.size() < vector_.size() );
106  }
107 
108  value_type* array ()
109  {
110  updateArray();
111  return array_;
112  }
113 
114  value_type* array () const
115  {
116  updateArray();
117  return array_;
118  }
119 
120  const value_type& operator[] ( size_t i ) const { return vector_[ i ].second; }
121 
122  value_type& operator[] ( size_t i )
123  {
124  arrayNeedsUpdate_ = true;
125  return vector_[ i ].second;
126  }
127 
128  void push_back ( int index, int value )
129  {
130  vector_.push_back( std::make_pair( index, value ) );
131  arrayNeedsUpdate_ = true;
132  }
133 
134  size_t size () const { return vector_.size(); }
135 
136  private:
137  PetscGhostArrayBuilder ( const ThisType& );
138  ThisType& operator= ( const ThisType & );
139  PetscGhostArrayBuilder ();
140 
141  void updateArray () const
142  {
143  if( arrayNeedsUpdate_ )
144  {
145  delete[] array_;
146  array_ = new value_type[ blockSize * size() ];
147  std::sort( vector_.begin(), vector_.end(), StorageSort() );
148  for( size_t block = 0; block < size(); ++block )
149  {
150  for( size_t inBlock = 0; inBlock < static_cast< size_t >( blockSize ); ++inBlock )
151  {
152  array_[ block*blockSize + inBlock ] = ( vector_[ block ].second ) * blockSize + inBlock;
153  }
154  }
155  arrayNeedsUpdate_ = false;
156  }
157  }
158 
159 
160  mutable StorageType vector_;
161  mutable value_type *array_;
162  mutable bool arrayNeedsUpdate_;
163  };
164 
165  /* =================================================
166  * class PetscDofCommunicator
167  *
168  * This is a helper class for PetscGhostArrayBuilder which is responsable for the
169  * communication necessary to build the ghost array.
170  * =================================================
171  */
172  template< typename DFSpace, typename PGhostArrayBuilder >
173  class PetscDofCommunicator
174  : public CommDataHandleIF< PetscDofCommunicator< DFSpace, PGhostArrayBuilder >, int >
175  {
176  typedef PetscDofCommunicator< DFSpace, PGhostArrayBuilder > ThisType;
177  typedef CommDataHandleIF< PetscDofCommunicator< DFSpace, PGhostArrayBuilder >, int > BaseType;
178 
179  public:
180 
181  typedef PGhostArrayBuilder PetscGhostArrayBuilderType;
182  typedef DFSpace DiscreteFunctionSpaceType;
183  typedef typename BaseType::DataType DataType;
184  typedef PetscSlaveDofProvider< DiscreteFunctionSpaceType > PetscSlaveDofProviderType ;
185  typedef typename PetscSlaveDofProviderType :: PetscDofMappingType PetscDofMappingType;
186  typedef typename DiscreteFunctionSpaceType::BlockMapperType BlockMapperType;
187 
188  /*
189  * methods
190  */
191  PetscDofCommunicator ( const DiscreteFunctionSpaceType &space,
192  PetscDofMappingType &petscDofMapping,
193  PetscGhostArrayBuilderType& ghostArrayBuilder )
194  : space_( space ),
195  mapper_( space_.blockMapper() ),
196  petscDofMapping_( petscDofMapping ),
197  ghostArrayBuilder_( ghostArrayBuilder ),
198  duplicates_()
199  {}
200 
201  bool contains ( int dim, int codim ) const
202  {
203  return mapper_.contains( codim );
204  }
205 
206  bool fixedsize ( int dim, int codim ) const
207  {
208  return false;
209  }
210 
211  template< typename MessageBuffer, typename Entity >
212  void gather ( MessageBuffer &buffer, const Entity &entity ) const
213  {
214  const int numDofs = mapper_.numEntityDofs( entity );
215  int counter = 0;
216  std::vector< int > slaveDofsWithIndex( 2*numDofs, -1 );
217 
218  std::vector< size_t > globalDofs;
219  globalDofs.resize( numDofs );
220  AssignFunctor< std::vector< size_t> > functor( globalDofs );
221  mapper_.mapEachEntityDof( entity, functor );
222 
223  for( int i = 0; i < numDofs; ++i )
224  {
225  size_t &dof = globalDofs[i];
226  assert( size_t( dof ) < petscDofMapping_.size() );
227  const int petscDof = petscDofMapping_.localSlaveMapping( dof );
228 
229  // does this process own the dof?
230  if( petscDof < petscDofMapping_.numOwnedDofBlocks() )
231  {
232  const int petscGlobalDof = petscDofMapping_.processStartIndex() + petscDof;
233 
234  // write the entity dof index and the PETSc dof
235  slaveDofsWithIndex[ counter++ ] = i;
236  slaveDofsWithIndex[ counter++ ] = petscGlobalDof;
237  }
238  }
239  assert( counter % 2 == 0 );
240 
241  // write the number of communicated slave dofs first
242  buffer.write( counter );
243 
244  assert( size_t( counter ) <= slaveDofsWithIndex.size() );
245  for( int i = 0; i < counter; ++i )
246  {
247  assert( slaveDofsWithIndex[ i ] != -1 );
248  buffer.write( slaveDofsWithIndex[ i ] );
249  }
250  }
251 
252  template< typename MessageBuffer, typename EntityType >
253  void scatter ( MessageBuffer &buffer, const EntityType &entity, size_t n )
254  {
255  bool isDuplicate = false;
256  int idx = space_.indexSet().index(entity);
257  if (duplicates_.find(idx) != duplicates_.end())
258  isDuplicate = true;
259  else
260  duplicates_.insert(idx);
261 
262  int twiceNumSlaveDofs;
263  buffer.read( twiceNumSlaveDofs );
264 
265  assert( twiceNumSlaveDofs % 2 == 0 );
266  assert( static_cast< size_t >( twiceNumSlaveDofs ) + 1 == n );
267  assert( twiceNumSlaveDofs + 1 == int( n ) );
268 
269 
270  const int numDofs = mapper_.numEntityDofs( entity );
271  std::vector< size_t > globalDofs;
272  globalDofs.resize( numDofs );
273  AssignFunctor< std::vector< size_t> > functor( globalDofs );
274  mapper_.mapEachEntityDof( entity, functor );
275 
276  for( size_t i = 0; i < size_t( twiceNumSlaveDofs/2 ); ++i )
277  {
278  int petscDof, index;
279  buffer.read( index );
280 
281  size_t &globalDof = globalDofs[ index ];
282  assert( petscDofMapping_.isSlave( globalDof ) );
283 
284  buffer.read( petscDof );
285  if ( !isDuplicate )
286  ghostArrayBuilder_.push_back( globalDof, petscDof );
287  }
288  }
289 
291  template< typename Entity >
292  size_t size ( const Entity &entity ) const
293  {
294  size_t size = 0;
295 
296  const int numDofs = mapper_.numEntityDofs( entity );
297 
298  std::vector< size_t > globalDofs;
299  globalDofs.resize( numDofs );
300  AssignFunctor< std::vector< size_t> > functor( globalDofs );
301  mapper_.mapEachEntityDof( entity, functor );
302 
303  for( int i = 0; i < numDofs; ++i )
304  {
305  if( !petscDofMapping_.isSlave( globalDofs[i] ) )
306  ++size;
307  }
308  assert( size == 0 || size == size_t( numDofs ) );
309 
310  return 2*size + 1;
311  }
312 
313  private:
314  PetscDofCommunicator ();
315  PetscDofCommunicator ( const ThisType& );
316  ThisType operator= ( const ThisType& );
317 
318 
319  const DiscreteFunctionSpaceType &space_;
320  const BlockMapperType &mapper_;
321  PetscDofMappingType &petscDofMapping_;
322  PetscGhostArrayBuilderType& ghostArrayBuilder_;
323  std::set<int> duplicates_;
324 
325  };
326 
327 
328  } // namespace Fem
329 
330 } // namespace Dune
331 
332 #endif // #if HAVE_PETSC
333 
334 #endif // DUNE_FEM_PETSCGHOSTARRAYBUILDER_HH
static int rank()
Definition: mpimanager.hh:116
Definition: coordinate.hh:4