Dune Core Modules (2.5.2)

gridptr.hh
1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 #ifndef DUNE_DGF_GRIDPTR_HH
4 #define DUNE_DGF_GRIDPTR_HH
5 
6 #include <cassert>
7 #include <iostream>
8 #include <string>
9 #include <vector>
10 #include <map>
11 #include <memory>
12 
13 //- Dune includes
16 
17 #include <dune/grid/common/gridenums.hh>
19 #include <dune/grid/common/intersection.hh>
20 
21 #include <dune/grid/io/file/dgfparser/dgfexception.hh>
22 #include <dune/grid/io/file/dgfparser/entitykey.hh>
23 #include <dune/grid/io/file/dgfparser/parser.hh>
24 
25 namespace Dune
26 {
27 
28  // External Forward Declarations
29  // -----------------------------
30 
31  template < class G >
32  struct DGFGridFactory;
33 
34  template< class GridImp, class IntersectionImp >
35  class Intersection;
36 
37 
38 
39  // GridPtr
40  // -------
41 
54  template< class GridType >
55  struct GridPtr
56  {
57  class mygrid_ptr : public std::shared_ptr< GridType >
58  {
59  typedef std::shared_ptr< GridType > base_t ;
60  // empty deleter to avoid deletion on release
61  typedef null_deleter< GridType > emptydeleter_t ;
62 
63  void removeObj()
64  {
65  // if use count is only 1 delete object
66  if( use_count() == 1 )
67  {
68  // delete point here, since we use the empty deleter
69  GridType* grd = release();
70  if( grd ) delete grd ;
71  }
72  }
73 
74  void assignObj( const mygrid_ptr& other )
75  {
76  removeObj();
77  base_t :: operator = ( other );
78  }
79  public:
80  using base_t :: get ;
81  using base_t :: swap ;
82  using base_t :: use_count ;
83 
84  // default constructor
85  mygrid_ptr() : base_t( ( GridType * ) 0, emptydeleter_t() ) {}
86  // copy constructor
87  mygrid_ptr( const mygrid_ptr& other ) { assignObj( other ); }
88  // constructor taking pointer
89  explicit mygrid_ptr( GridType* grd ) : base_t( grd, emptydeleter_t() ) {}
90 
91  // destructor
92  ~mygrid_ptr() { removeObj(); }
93 
94  // assigment operator
95  mygrid_ptr& operator = ( const mygrid_ptr& other )
96  {
97  assignObj( other );
98  return *this;
99  }
100 
101  // release pointer
102  GridType* release()
103  {
104  GridType* grd = this->get();
105  base_t ptr(( GridType * ) 0, emptydeleter_t() );
106  this->swap( ptr );
107  return grd ;
108  }
109  };
110 
111  typedef MPIHelper::MPICommunicator MPICommunicatorType;
112  static const int dimension = GridType::dimension;
113 
115  explicit GridPtr ( const std::string &filename,
116  MPICommunicatorType comm = MPIHelper::getCommunicator() )
117  : gridPtr_(),
118  elParam_(),
119  vtxParam_(),
120  bndParam_(),
121  bndId_(),
122  emptyParam_(),
123  nofElParam_( 0 ),
124  nofVtxParam_( 0 ),
125  haveBndParam_( false )
126  {
127  DGFGridFactory< GridType > dgfFactory( filename, comm );
128  initialize( dgfFactory );
129  }
130 
132  explicit GridPtr ( std::istream &input,
133  MPICommunicatorType comm = MPIHelper::getCommunicator() )
134  : gridPtr_(),
135  elParam_(),
136  vtxParam_(),
137  bndParam_(),
138  bndId_(),
139  emptyParam_(),
140  nofElParam_( 0 ),
141  nofVtxParam_( 0 ),
142  haveBndParam_( false )
143  {
144  DGFGridFactory< GridType > dgfFactory( input, comm );
145  initialize( dgfFactory );
146  }
147 
150  : gridPtr_(),
151  elParam_(),
152  vtxParam_(),
153  bndParam_(),
154  bndId_(),
155  emptyParam_(),
156  nofElParam_(0),
157  nofVtxParam_(0),
158  haveBndParam_( false )
159  {}
160 
162  explicit GridPtr( GridType *grd )
163  : gridPtr_(grd),
164  elParam_(),
165  vtxParam_(),
166  bndParam_(),
167  bndId_(),
168  emptyParam_(),
169  nofElParam_(0),
170  nofVtxParam_(0),
171  haveBndParam_( false )
172  {}
173 
175  GridPtr( const GridPtr &org )
176  : gridPtr_(org.gridPtr_),
177  elParam_(org.elParam_),
178  vtxParam_(org.vtxParam_),
179  bndParam_(org.bndParam_),
180  bndId_(org.bndId_),
181  emptyParam_( org.emptyParam_ ),
182  nofElParam_(org.nofElParam_),
183  nofVtxParam_(org.nofVtxParam_),
184  haveBndParam_(org.haveBndParam_)
185  {}
186 
188  GridPtr& operator= ( const GridPtr &org )
189  {
190  gridPtr_ = org.gridPtr_;
191  elParam_ = org.elParam_;
192  vtxParam_ = org.vtxParam_;
193  bndParam_ = org.bndParam_;
194  bndId_ = org.bndId_;
195  emptyParam_ = org.emptyParam_;
196 
197  nofElParam_ = org.nofElParam_;
198  nofVtxParam_ = org.nofVtxParam_;
199  haveBndParam_ = org.haveBndParam_;
200  return *this;
201  }
202 
204  GridPtr& operator = (GridType * grd)
205  {
206  gridPtr_ = mygrid_ptr( grd );
207  elParam_.resize(0);
208  vtxParam_.resize(0);
209  bndParam_.resize(0);
210  bndId_.resize(0);
211  emptyParam_.resize(0);
212 
213  nofVtxParam_ = 0;
214  nofElParam_ = 0;
215  haveBndParam_ = false;
216  return *this;
217  }
218 
220  GridType& operator*() {
221  return *gridPtr_;
222  }
223 
225  GridType* operator->() {
226  return gridPtr_.operator -> ();
227  }
228 
230  const GridType& operator*() const {
231  return *gridPtr_;
232  }
233 
235  const GridType* operator->() const {
236  return gridPtr_.operator -> ();
237  }
238 
240  GridType* release () { return gridPtr_.release(); }
241 
243  int nofParameters(int cdim) const {
244  switch (cdim) {
245  case 0 : return nofElParam_; break;
246  case GridType::dimension : return nofVtxParam_; break;
247  }
248  return 0;
249  }
250 
252  template <class Entity>
253  int nofParameters ( const Entity & ) const
254  {
255  return nofParameters( (int) Entity::codimension );
256  }
257 
259  template< class GridImp, class IntersectionImp >
260  int nofParameters ( const Intersection< GridImp, IntersectionImp > & intersection ) const
261  {
262  return parameters( intersection ).size();
263  }
264 
266  template <class Entity>
267  const std::vector< double > &parameters ( const Entity &entity ) const
268  {
269  typedef typename GridType::LevelGridView GridView;
270  GridView gridView = gridPtr_->levelGridView( 0 );
271  switch( (int)Entity::codimension )
272  {
273  case 0 :
274  if( nofElParam_ > 0 )
275  {
276  assert( (unsigned int)gridView.indexSet().index( entity ) < elParam_.size() );
277  return elParam_[ gridView.indexSet().index( entity ) ];
278  }
279  break;
280  case GridType::dimension :
281  if( nofVtxParam_ > 0 )
282  {
283  assert( (unsigned int)gridView.indexSet().index( entity ) < vtxParam_.size() );
284  return vtxParam_[ gridView.indexSet().index( entity ) ];
285  }
286  break;
287  }
288  return emptyParam_;
289  }
290 
292  template< class GridImp, class IntersectionImp >
294  {
295  // if no parameters given return empty vector
296  if ( !haveBndParam_ )
298 
299  return bndParam_[ intersection.boundarySegmentIndex() ];
300  }
301 
302  void loadBalance()
303  {
304  if ( gridPtr_->comm().size() == 1 )
305  return;
306  int nofParams = nofElParam_ + nofVtxParam_;
307  if ( haveBndParam_ )
308  nofParams++;
309  if ( gridPtr_->comm().max( nofParams ) > 0 )
310  {
311  DataHandle dh(*this);
312  gridPtr_->loadBalance( dh.interface() );
313  gridPtr_->communicate( dh.interface(), InteriorBorder_All_Interface,ForwardCommunication);
314  } else
315  {
316  gridPtr_->loadBalance();
317  }
318  }
319 
320  protected:
321  void initialize ( DGFGridFactory< GridType > &dgfFactory )
322  {
323  gridPtr_ = mygrid_ptr( dgfFactory.grid() );
324 
325  typedef typename GridType::LevelGridView GridView;
326  GridView gridView = gridPtr_->levelGridView( 0 );
327  const typename GridView::IndexSet &indexSet = gridView.indexSet();
328 
329  nofElParam_ = dgfFactory.template numParameters< 0 >();
330  nofVtxParam_ = dgfFactory.template numParameters< dimension >();
331  haveBndParam_ = dgfFactory.haveBoundaryParameters();
332 
333  if ( nofElParam_ > 0 )
334  elParam_.resize( indexSet.size(0) );
335  if ( nofVtxParam_ > 0 )
336  vtxParam_.resize( indexSet.size(dimension) );
337  bndId_.resize( indexSet.size(1) );
338  if ( haveBndParam_ )
339  bndParam_.resize( gridPtr_->numBoundarySegments() );
340 
341  const PartitionIteratorType partType = Interior_Partition;
342  typedef typename GridView::template Codim< 0 >::template Partition< partType >::Iterator Iterator;
343  const Iterator enditer = gridView.template end< 0, partType >();
344  for( Iterator iter = gridView.template begin< 0, partType >(); iter != enditer; ++iter )
345  {
346  const typename Iterator::Entity &el = *iter;
347  if ( nofElParam_ > 0 ) {
348  std::swap( elParam_[ indexSet.index(el) ], dgfFactory.parameter(el) );
349  assert( elParam_[ indexSet.index(el) ].size() == (size_t)nofElParam_ );
350  }
351  if ( nofVtxParam_ > 0 )
352  {
353  const unsigned int subEntities = el.subEntities(dimension);
354  for ( unsigned int v = 0; v < subEntities; ++v)
355  {
356  typename GridView::IndexSet::IndexType index = indexSet.subIndex(el,v,dimension);
357  if ( vtxParam_[ index ].empty() )
358  std::swap( vtxParam_[ index ], dgfFactory.parameter(el.template subEntity<dimension>(v) ) );
359  assert( vtxParam_[ index ].size() == (size_t)nofVtxParam_ );
360  }
361  }
362  if ( el.hasBoundaryIntersections() )
363  {
364  typedef typename GridView::IntersectionIterator IntersectionIterator;
365  const IntersectionIterator iend = gridView.iend(el);
366  for( IntersectionIterator iiter = gridView.ibegin(el); iiter != iend; ++iiter )
367  {
368  const typename IntersectionIterator::Intersection &inter = *iiter;
369  // dirty hack: check for "none" to make corner point grid work
370  if ( inter.boundary() && !inter.type().isNone() )
371  {
372  const int k = indexSet.subIndex(el,inter.indexInInside(),1);
373  bndId_[ k ] = dgfFactory.boundaryId( inter );
374  if ( haveBndParam_ )
375  bndParam_[ inter.boundarySegmentIndex() ] = dgfFactory.boundaryParameter( inter );
376  }
377  }
378  }
379  }
380  }
381 
382  template <class Entity>
383  std::vector< double > &params ( const Entity &entity )
384  {
385  typedef typename GridType::LevelGridView GridView;
386  GridView gridView = gridPtr_->levelGridView( 0 );
387  switch( (int)Entity::codimension )
388  {
389  case 0 :
390  if( nofElParam_ > 0 ) {
391  if ( gridView.indexSet().index( entity ) >= elParam_.size() )
392  elParam_.resize( gridView.indexSet().index( entity ) );
393  return elParam_[ gridView.indexSet().index( entity ) ];
394  }
395  break;
396  case GridType::dimension :
397  if( nofVtxParam_ > 0 ) {
398  if ( gridView.indexSet().index( entity ) >= vtxParam_.size() )
399  vtxParam_.resize( gridView.indexSet().index( entity ) );
400  return vtxParam_[ gridView.indexSet().index( entity ) ];
401  }
402  break;
403  }
404  return emptyParam_;
405  }
406 
407  void setNofParams( int cdim, int nofP )
408  {
409  switch (cdim) {
410  case 0 : nofElParam_ = nofP; break;
411  case GridType::dimension : nofVtxParam_ = nofP; break;
412  }
413  }
414 
415  struct DataHandle : public CommDataHandleIF<DataHandle,double>
416  {
417  DataHandle( GridPtr& gridPtr) :
418  gridPtr_(gridPtr),
419  idSet_(gridPtr->localIdSet())
420  {
421  typedef typename GridType::LevelGridView GridView;
422  GridView gridView = gridPtr_->levelGridView( 0 );
423  const typename GridView::IndexSet &indexSet = gridView.indexSet();
424 
425  const PartitionIteratorType partType = Interior_Partition;
426  typedef typename GridView::template Codim< 0 >::template Partition< partType >::Iterator Iterator;
427  const Iterator enditer = gridView.template end< 0, partType >();
428  for( Iterator iter = gridView.template begin< 0, partType >(); iter != enditer; ++iter )
429  {
430  const typename Iterator::Entity &el = *iter;
431  if ( gridPtr_.nofElParam_ > 0 )
432  std::swap( gridPtr_.elParam_[ indexSet.index(el) ], elData_[ idSet_.id(el) ] );
433  if ( gridPtr_.nofVtxParam_ > 0 )
434  {
435  for ( unsigned int v = 0; v < el.subEntities(dimension); ++v)
436  {
437  typename GridView::IndexSet::IndexType index = indexSet.subIndex(el,v,dimension);
438  if ( ! gridPtr_.vtxParam_[ index ].empty() )
439  std::swap( gridPtr_.vtxParam_[ index ], vtxData_[ idSet_.subId(el,v,dimension) ] );
440  }
441  }
442  }
443  }
444 
445  ~DataHandle()
446  {
447  typedef typename GridType::LevelGridView GridView;
448  GridView gridView = gridPtr_->levelGridView( 0 );
449  const typename GridView::IndexSet &indexSet = gridView.indexSet();
450 
451  if ( gridPtr_.nofElParam_ > 0 )
452  gridPtr_.elParam_.resize( indexSet.size(0) );
453  if ( gridPtr_.nofVtxParam_ > 0 )
454  gridPtr_.vtxParam_.resize( indexSet.size(dimension) );
455 
456  const PartitionIteratorType partType = All_Partition;
457  typedef typename GridView::template Codim< 0 >::template Partition< partType >::Iterator Iterator;
458  const Iterator enditer = gridView.template end< 0, partType >();
459  for( Iterator iter = gridView.template begin< 0, partType >(); iter != enditer; ++iter )
460  {
461  const typename Iterator::Entity &el = *iter;
462  if ( gridPtr_.nofElParam_ > 0 )
463  {
464  std::swap( gridPtr_.elParam_[ indexSet.index(el) ], elData_[ idSet_.id(el) ] );
465  assert( gridPtr_.elParam_[ indexSet.index(el) ].size() == (unsigned int)gridPtr_.nofElParam_ );
466  }
467  if ( gridPtr_.nofVtxParam_ > 0 )
468  {
469  for ( unsigned int v = 0; v < el.subEntities(dimension); ++v)
470  {
471  typename GridView::IndexSet::IndexType index = indexSet.subIndex(el,v,dimension);
472  if ( gridPtr_.vtxParam_[ index ].empty() )
473  std::swap( gridPtr_.vtxParam_[ index ], vtxData_[ idSet_.subId(el,v,dimension) ] );
474  assert( gridPtr_.vtxParam_[ index ].size() == (unsigned int)gridPtr_.nofVtxParam_ );
475  }
476  }
477  }
478  }
479 
480  CommDataHandleIF<DataHandle,double> &interface()
481  {
482  return *this;
483  }
484 
485  bool contains (int dim, int codim) const
486  {
487  return (codim==dim || codim==0);
488  }
489 
490  bool fixedSize (int dim, int codim) const
491  {
492  return false;
493  }
494 
495  template<class EntityType>
496  size_t size (const EntityType& e) const
497  {
498  return gridPtr_.nofParameters( (int) e.codimension);
499  }
500 
501  template<class MessageBufferImp, class EntityType>
502  void gather (MessageBufferImp& buff, const EntityType& e) const
503  {
504  const std::vector<double> &v = (e.codimension==0) ? elData_[idSet_.id(e)] : vtxData_[idSet_.id(e)];
505  const size_t s = v.size();
506  for (size_t i=0; i<s; ++i)
507  buff.write( v[i] );
508  assert( s == (size_t)gridPtr_.nofParameters(e.codimension) );
509  }
510 
511  template<class MessageBufferImp, class EntityType>
512  void scatter (MessageBufferImp& buff, const EntityType& e, size_t n)
513  {
514  std::vector<double> &v = (e.codimension==0) ? elData_[idSet_.id(e)] : vtxData_[idSet_.id(e)];
515  v.resize( n );
516  gridPtr_.setNofParams( e.codimension, n );
517  for (size_t i=0; i<n; ++i)
518  buff.read( v[i] );
519  }
520 
521  private:
522  typedef typename GridType::LocalIdSet IdSet;
523  GridPtr &gridPtr_;
524  const IdSet &idSet_;
525  mutable std::map< typename IdSet::IdType, std::vector<double> > elData_, vtxData_;
526  };
527 
528  // grid auto pointer
529  mutable mygrid_ptr gridPtr_;
530  // element and vertex parameters
531  std::vector< std::vector< double > > elParam_;
532  std::vector< std::vector< double > > vtxParam_;
533  std::vector< DGFBoundaryParameter::type > bndParam_;
534  std::vector< int > bndId_;
535  std::vector < double > emptyParam_;
536 
537  int nofElParam_;
538  int nofVtxParam_;
539  bool haveBndParam_;
540  }; // end of class GridPtr
541 
542 } // end namespace Dune
543 
544 #endif
Wrapper class for entities.
Definition: entity.hh:65
@ codimension
Know your own codimension.
Definition: entity.hh:107
Grid view abstract base class.
Definition: gridview.hh:60
Dune::Intersection< GridImp, IntersectionImp > Intersection
Type of Intersection this IntersectionIterator points to.
Definition: intersectioniterator.hh:105
Intersection of a mesh entity of codimension 0 ("element") with a "neighboring" element or with the d...
Definition: intersection.hh:162
size_t boundarySegmentIndex() const
index of the boundary segment within the macro grid
Definition: intersection.hh:264
MPI_Comm MPICommunicator
The type of the mpi communicator.
Definition: mpihelper.hh:173
static MPICommunicator getCommunicator()
get the default communicator
Definition: mpihelper.hh:181
Describes the parallel communication interface class for MessageBuffers and DataHandles.
Traits ::IntersectionIterator IntersectionIterator
type of the intersection iterator
Definition: gridview.hh:87
Traits ::IndexSet IndexSet
type of the index set
Definition: gridview.hh:81
const IndexSet & indexSet() const
obtain the index set
Definition: gridview.hh:173
PartitionIteratorType
Parameter to be used for the parallel level- and leaf iterators.
Definition: gridenums.hh:134
@ All_Partition
all entities
Definition: gridenums.hh:139
@ Interior_Partition
only interior entities
Definition: gridenums.hh:135
@ ForwardCommunication
communicate as given in InterfaceType
Definition: gridenums.hh:169
@ InteriorBorder_All_Interface
send interior and border, receive all entities
Definition: gridenums.hh:86
Helpers for dealing with MPI.
Dune namespace.
Definition: alignment.hh:11
This file implements the class shared_ptr (a reference counting pointer), for those systems that don'...
static const type & defaultValue()
default constructor
Definition: parser.hh:26
std::string type
type of additional boundary parameters
Definition: parser.hh:23
Class for constructing grids from DGF files.
Definition: gridptr.hh:56
const std::vector< double > & parameters(const Entity &entity) const
get parameters defined for each codim 0 und dim entity on the grid through the grid file
Definition: gridptr.hh:267
GridPtr & operator=(const GridPtr &org)
assignment of grid pointer
Definition: gridptr.hh:188
GridPtr(std::istream &input, MPICommunicatorType comm=MPIHelper::getCommunicator())
constructor given a std::istream
Definition: gridptr.hh:132
GridType * operator->()
return pointer to GridType instance
Definition: gridptr.hh:225
GridPtr(const GridPtr &org)
Copy constructor, copies internal auto pointer.
Definition: gridptr.hh:175
GridPtr()
Default constructor, creating empty GridPtr.
Definition: gridptr.hh:149
GridPtr(GridType *grd)
Constructor storing given pointer to internal auto pointer.
Definition: gridptr.hh:162
GridType * release()
release pointer from internal ownership
Definition: gridptr.hh:240
const GridType * operator->() const
return const pointer to GridType instance
Definition: gridptr.hh:235
int nofParameters(const Entity &) const
get parameters defined for given entity
Definition: gridptr.hh:253
GridPtr(const std::string &filename, MPICommunicatorType comm=MPIHelper::getCommunicator())
constructor given the name of a DGF file
Definition: gridptr.hh:115
const GridType & operator*() const
return const reference to GridType instance
Definition: gridptr.hh:230
int nofParameters(int cdim) const
get number of parameters defined for a given codimension
Definition: gridptr.hh:243
const DGFBoundaryParameter::type & parameters(const Intersection< GridImp, IntersectionImp > &intersection) const
get parameters for intersection
Definition: gridptr.hh:293
int nofParameters(const Intersection< GridImp, IntersectionImp > &intersection) const
get number of parameters defined for a given intersection
Definition: gridptr.hh:260
GridType & operator*()
return reference to GridType instance
Definition: gridptr.hh:220
implements the Deleter concept of shared_ptr without deleting anything
Definition: shared_ptr.hh:52
std::size_t fixedSize
The number of data items per index if it is fixed, 0 otherwise.
Definition: variablesizecommunicator.hh:238
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 8, 22:30, 2024)