Dune Core Modules (unstable)

gridptr.hh
1 // SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file LICENSE.md in module root
2 // SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
4 // vi: set et ts=4 sw=2 sts=2:
5 #ifndef DUNE_DGF_GRIDPTR_HH
6 #define DUNE_DGF_GRIDPTR_HH
7 
8 #include <cassert>
9 #include <cctype>
10 
11 #include <array>
12 #include <iostream>
13 #include <map>
14 #include <memory>
15 #include <string>
16 #include <type_traits>
17 #include <vector>
18 
19 //- Dune includes
22 
23 #include <dune/grid/common/gridenums.hh>
25 #include <dune/grid/common/intersection.hh>
26 #include <dune/grid/common/partitionset.hh>
27 #include <dune/grid/common/rangegenerators.hh>
28 
29 #include <dune/grid/io/file/dgfparser/dgfexception.hh>
30 #include <dune/grid/io/file/dgfparser/entitykey.hh>
31 #include <dune/grid/io/file/dgfparser/parser.hh>
32 
33 #include <dune/grid/io/file/gmshreader.hh>
34 
35 namespace Dune
36 {
37 
38  // External Forward Declarations
39  // -----------------------------
40 
41  template < class G >
42  struct DGFGridFactory;
43 
44  template< class GridImp, class IntersectionImp >
45  class Intersection;
46 
47 
48 
49  // GridPtr
50  // -------
51 
64  template< class GridType >
65  struct GridPtr
66  {
67  class mygrid_ptr : public std::shared_ptr< GridType >
68  {
69  typedef std::shared_ptr< GridType > base_t ;
70  // empty deleter to avoid deletion on release
71  typedef null_deleter< GridType > emptydeleter_t ;
72 
73  void removeObj()
74  {
75  // if use count is only 1 delete object
76  if( use_count() == 1 )
77  {
78  // delete point here, since we use the empty deleter
79  GridType* grd = release();
80  if( grd ) delete grd ;
81  }
82  }
83 
84  void assignObj( const mygrid_ptr& other )
85  {
86  removeObj();
87  base_t :: operator = ( other );
88  }
89  public:
90  using base_t :: get ;
91  using base_t :: swap ;
92  using base_t :: use_count ;
93 
94  // default constructor
95  mygrid_ptr() : base_t( ( GridType * ) 0, emptydeleter_t() ) {}
96  // copy constructor
97  mygrid_ptr( const mygrid_ptr& other ) : base_t(nullptr) { assignObj( other ); }
98  // constructor taking pointer
99  explicit mygrid_ptr( GridType* grd ) : base_t( grd, emptydeleter_t() ) {}
100 
101  // destructor
102  ~mygrid_ptr() { removeObj(); }
103 
104  // assignment operator
105  mygrid_ptr& operator = ( const mygrid_ptr& other )
106  {
107  assignObj( other );
108  return *this;
109  }
110 
111  // release pointer
112  GridType* release()
113  {
114  GridType* grd = this->get();
115  base_t ptr(( GridType * ) 0, emptydeleter_t() );
116  this->swap( ptr );
117  return grd ;
118  }
119  };
120 
121  protected:
122  std::string getFileExtension( const std::string& filename ) const
123  {
124  // extract file extension
125  auto extpos = filename.find_last_of(".");
126  std::string ext;
127  if( extpos != std::string::npos)
128  ext = filename.substr( extpos + 1 );
129 
130  // convert all letters to lower case
131  for( auto& item : ext )
132  item = std::tolower( item );
133  return ext;
134  }
135 
136  // read gmsh file if dimension world <= 3
137  void readGmsh( const std::string& filename, std::integral_constant< bool, true > )
138  {
139  GridFactory<GridType> gridFactory;
140  std::vector<int> boundaryIDs;
141  std::vector<int> elementsIDs;
142  GmshReader<GridType>::read(gridFactory,filename,boundaryIDs,elementsIDs);
143  initialize( gridFactory, boundaryIDs,elementsIDs);
144  }
145 
146  // if dimension world > 3 throw GridError
147  void readGmsh( const std::string& filename, std::integral_constant< bool, false > )
148  {
149  DUNE_THROW(GridError, "GmshReader requires dimWorld <= 3." );
150  }
151 
152  public:
153 
154  typedef MPIHelper::MPICommunicator MPICommunicatorType;
155  static const int dimension = GridType::dimension;
156 
158  explicit GridPtr ( const std::string &filename,
159  MPICommunicatorType comm = MPIHelper::getCommunicator() )
160  : gridPtr_(),
161  elParam_(),
162  vtxParam_(),
163  bndParam_(),
164  bndId_(),
165  emptyParam_(),
166  nofElParam_( 0 ),
167  nofVtxParam_( 0 ),
168  haveBndParam_( false )
169  {
170  std::string fileExt = getFileExtension( filename );
171 
172  if( fileExt == "dgf" )
173  {
174  DGFGridFactory< GridType > dgfFactory( filename, comm );
175  initialize( dgfFactory );
176  }
177  else if( fileExt == "msh" )
178  {
179  // Gmsh reader only compiles for dimworld <= 3
180  readGmsh( filename, std::integral_constant< bool, GridType::dimensionworld <= 3 > () );
181  }
182  else if( fileExt == "amc" || fileExt == "2d" || fileExt == "3d" )
183  {
184  // TODO: AlbertaReader
185  DUNE_THROW( NotImplemented, "GridPtr: file format '" << fileExt << "' not supported yet!" );
186  }
187  else if( fileExt == "vtu" )
188  {
189  // TODO: vtu/vtk reader
190  DUNE_THROW( NotImplemented, "GridPtr: file format '" << fileExt << "' not supported yet!" );
191  }
192  else
193  {
194  DUNE_THROW( NotImplemented, "GridPtr: file format '" << fileExt << "' not supported yet!" );
195  }
196  }
197 
199  explicit GridPtr ( std::istream &input,
200  MPICommunicatorType comm = MPIHelper::getCommunicator() )
201  : gridPtr_(),
202  elParam_(),
203  vtxParam_(),
204  bndParam_(),
205  bndId_(),
206  emptyParam_(),
207  nofElParam_( 0 ),
208  nofVtxParam_( 0 ),
209  haveBndParam_( false )
210  {
211  // input stream only works for DGF format right now
212  DGFGridFactory< GridType > dgfFactory( input, comm );
213  initialize( dgfFactory );
214  }
215 
218  : gridPtr_(),
219  elParam_(),
220  vtxParam_(),
221  bndParam_(),
222  bndId_(),
223  emptyParam_(),
224  nofElParam_(0),
225  nofVtxParam_(0),
226  haveBndParam_( false )
227  {}
228 
230  explicit GridPtr( GridType *grd )
231  : gridPtr_(grd),
232  elParam_(),
233  vtxParam_(),
234  bndParam_(),
235  bndId_(),
236  emptyParam_(),
237  nofElParam_(0),
238  nofVtxParam_(0),
239  haveBndParam_( false )
240  {}
241 
243  GridPtr( const GridPtr &org ) = default;
244 
246  GridPtr& operator= ( const GridPtr &org )
247  {
248  gridPtr_ = org.gridPtr_;
249  elParam_ = org.elParam_;
250  vtxParam_ = org.vtxParam_;
251  bndParam_ = org.bndParam_;
252  bndId_ = org.bndId_;
253  emptyParam_ = org.emptyParam_;
254 
255  nofElParam_ = org.nofElParam_;
256  nofVtxParam_ = org.nofVtxParam_;
257  haveBndParam_ = org.haveBndParam_;
258  return *this;
259  }
260 
262  GridPtr& operator = (GridType * grd)
263  {
264  gridPtr_ = mygrid_ptr( grd );
265  elParam_.resize(0);
266  vtxParam_.resize(0);
267  bndParam_.resize(0);
268  bndId_.resize(0);
269  emptyParam_.resize(0);
270 
271  nofVtxParam_ = 0;
272  nofElParam_ = 0;
273  haveBndParam_ = false;
274  return *this;
275  }
276 
278  GridType& operator*() {
279  return *gridPtr_;
280  }
281 
283  GridType* operator->() {
284  return gridPtr_.operator -> ();
285  }
286 
288  const GridType& operator*() const {
289  return *gridPtr_;
290  }
291 
293  const GridType* operator->() const {
294  return gridPtr_.operator -> ();
295  }
296 
298  GridType* release () { return gridPtr_.release(); }
299 
301  int nofParameters(int cdim) const {
302  switch (cdim) {
303  case 0 : return nofElParam_; break;
304  case GridType::dimension : return nofVtxParam_; break;
305  }
306  return 0;
307  }
308 
310  template <class Entity>
311  int nofParameters ( const Entity & ) const
312  {
313  return nofParameters( (int) Entity::codimension );
314  }
315 
317  template< class GridImp, class IntersectionImp >
318  int nofParameters ( const Intersection< GridImp, IntersectionImp > & intersection ) const
319  {
320  return parameters( intersection ).size();
321  }
322 
324  template <class Entity>
325  const std::vector< double > &parameters ( const Entity &entity ) const
326  {
327  typedef typename GridType::LevelGridView GridView;
328  GridView gridView = gridPtr_->levelGridView( 0 );
329  switch( (int)Entity::codimension )
330  {
331  case 0 :
332  if( nofElParam_ > 0 )
333  {
334  assert( (unsigned int)gridView.indexSet().index( entity ) < elParam_.size() );
335  return elParam_[ gridView.indexSet().index( entity ) ];
336  }
337  break;
338  case GridType::dimension :
339  if( nofVtxParam_ > 0 )
340  {
341  assert( (unsigned int)gridView.indexSet().index( entity ) < vtxParam_.size() );
342  return vtxParam_[ gridView.indexSet().index( entity ) ];
343  }
344  break;
345  }
346  return emptyParam_;
347  }
348 
350  template< class GridImp, class IntersectionImp >
352  {
353  // if no parameters given return empty vector
354  if ( !haveBndParam_ )
356 
357  return bndParam_[ intersection.boundarySegmentIndex() ];
358  }
359 
360  void communicate ()
361  {
362  if( gridPtr_->comm().size() > 1 )
363  {
364  DataHandle dh(*this);
365  gridPtr_->levelGridView( 0 ).communicate( dh.interface(), InteriorBorder_All_Interface,ForwardCommunication );
366  }
367  }
368 
369  void loadBalance ()
370  {
371  if( gridPtr_->comm().size() > 1 )
372  {
373  DataHandle dh(*this);
374  gridPtr_->loadBalance( dh.interface() );
375  gridPtr_->levelGridView( 0 ).communicate( dh.interface(), InteriorBorder_All_Interface,ForwardCommunication );
376  }
377  }
378 
379  protected:
380  template< class Range >
381  static bool isEmpty ( Range &&range )
382  {
383  return range.begin() == range.end();
384  }
385 
386  void initialize ( DGFGridFactory< GridType > &dgfFactory )
387  {
388  gridPtr_ = mygrid_ptr( dgfFactory.grid() );
389 
390  const auto gridView = gridPtr_->levelGridView( 0 );
391  const auto &indexSet = gridView.indexSet();
392 
393  nofElParam_ = dgfFactory.template numParameters< 0 >();
394  nofVtxParam_ = dgfFactory.template numParameters< dimension >();
395  haveBndParam_ = dgfFactory.haveBoundaryParameters();
396 
397  std::array< int, 3 > nofParams = {{ nofElParam_, nofVtxParam_, static_cast< int >( haveBndParam_ ) }};
398  gridView.comm().max( nofParams.data(), nofParams.size() );
399 
400  // empty grids have no parameters associated
401  if( isEmpty( elements( gridView, Partitions::interiorBorder ) ) )
402  {
403  nofElParam_ = nofParams[ 0 ];
404  nofVtxParam_ = nofParams[ 1 ];
405  }
406 
407  // boundary parameters may be empty
408  haveBndParam_ = static_cast< bool >( nofParams[ 2 ] );
409 
410  if( (nofElParam_ != nofParams[ 0 ]) || (nofVtxParam_ != nofParams[ 1 ]) )
411  DUNE_THROW( DGFException, "Number of parameters differs between processes" );
412 
413  elParam_.resize( nofElParam_ > 0 ? indexSet.size( 0 ) : 0 );
414  vtxParam_.resize( nofVtxParam_ > 0 ? indexSet.size( dimension ) : 0 );
415 
416  bndId_.resize( indexSet.size( 1 ) );
417  if( haveBndParam_ )
418  bndParam_.resize( gridPtr_->numBoundarySegments() );
419 
420  for( const auto &element : elements( gridView, Partitions::interiorBorder ) )
421  {
422  if( nofElParam_ > 0 )
423  {
424  std::swap( elParam_[ indexSet.index( element ) ], dgfFactory.parameter( element ) );
425  assert( elParam_[ indexSet.index( element ) ].size() == static_cast< std::size_t >( nofElParam_ ) );
426  }
427 
428  if( nofVtxParam_ > 0 )
429  {
430  for( unsigned int v = 0, n = element.subEntities( dimension ); v < n; ++v )
431  {
432  const auto index = indexSet.subIndex( element, v, dimension );
433  if( vtxParam_[ index ].empty() )
434  std::swap( vtxParam_[ index ], dgfFactory.parameter( element.template subEntity< dimension >( v ) ) );
435  assert( vtxParam_[ index ].size() == static_cast< std::size_t >( nofVtxParam_ ) );
436  }
437  }
438 
439  if( element.hasBoundaryIntersections() )
440  {
441  for( const auto &intersection : intersections( gridView, element ) )
442  {
443  // dirty hack: check for "none" to make corner point grid work
444  if( !intersection.boundary() || intersection.type().isNone() )
445  continue;
446 
447  const auto k = indexSet.subIndex( element, intersection.indexInInside(), 1 );
448  bndId_[ k ] = dgfFactory.boundaryId( intersection );
449  if( haveBndParam_ )
450  bndParam_[ intersection.boundarySegmentIndex() ] = dgfFactory.boundaryParameter( intersection );
451  }
452  }
453  }
454  }
455 
456  void initialize ( GridFactory< GridType > &factory,
457  std::vector<int>& boundaryIds,
458  std::vector<int>& elementIds )
459  {
460  gridPtr_ = mygrid_ptr( factory.createGrid().release() );
461 
462  const auto& gridView = gridPtr_->leafGridView();
463  const auto& indexSet = gridView.indexSet();
464 
465  nofElParam_ = elementIds.empty() ? 0 : 1 ;
466  nofVtxParam_ = 0;
467  haveBndParam_ = boundaryIds.empty() ? 0 : 1 ;
468 
469  std::array< int, 3 > nofParams = {{ nofElParam_, nofVtxParam_, static_cast< int >( haveBndParam_ ) }};
470  gridView.comm().max( nofParams.data(), nofParams.size() );
471 
472  // empty grids have no parameters associated
473  if( isEmpty( elements( gridView, Partitions::interiorBorder ) ) )
474  {
475  nofElParam_ = nofParams[ 0 ];
476  }
477 
478  // boundary parameters may be empty
479  haveBndParam_ = static_cast< bool >( nofParams[ 2 ] );
480 
481  // Reorder boundary IDs according to the insertion index
482  if(!boundaryIds.empty() || !elementIds.empty() )
483  {
484  bndParam_.resize( boundaryIds.size() );
485  elParam_.resize( elementIds.size(), std::vector<double>(1) );
486  for(const auto& entity : elements( gridView ))
487  {
488  elParam_[ indexSet.index( entity ) ][ 0 ] = elementIds[ factory.insertionIndex( entity ) ];
489  if( haveBndParam_ )
490  {
491  for(const auto& intersection : intersections( gridView,entity) )
492  {
493  if(intersection.boundary())
494  {
495  // DGFBoundaryParameter::type is of type string.
496  bndParam_[intersection.boundarySegmentIndex()] = std::to_string(boundaryIds[factory.insertionIndex(intersection)]);
497  }
498  }
499  }
500  }
501  }
502  }
503 
504  template <class Entity>
505  std::vector< double > &params ( const Entity &entity )
506  {
507  const auto gridView = gridPtr_->levelGridView( 0 );
508  switch( (int)Entity::codimension )
509  {
510  case 0 :
511  if( nofElParam_ > 0 ) {
512  if ( gridView.indexSet().index( entity ) >= elParam_.size() )
513  elParam_.resize( gridView.indexSet().index( entity ) );
514  return elParam_[ gridView.indexSet().index( entity ) ];
515  }
516  break;
517  case GridType::dimension :
518  if( nofVtxParam_ > 0 ) {
519  if ( gridView.indexSet().index( entity ) >= vtxParam_.size() )
520  vtxParam_.resize( gridView.indexSet().index( entity ) );
521  return vtxParam_[ gridView.indexSet().index( entity ) ];
522  }
523  break;
524  }
525  return emptyParam_;
526  }
527 
528  void setNofParams( int cdim, int nofP )
529  {
530  switch (cdim) {
531  case 0 : nofElParam_ = nofP; break;
532  case GridType::dimension : nofVtxParam_ = nofP; break;
533  }
534  }
535 
536  struct DataHandle
537  : public CommDataHandleIF< DataHandle, char >
538  {
539  explicit DataHandle ( GridPtr &gridPtr )
540  : gridPtr_( gridPtr ), idSet_( gridPtr->localIdSet() )
541  {
542  const auto gridView = gridPtr_->levelGridView( 0 );
543  const auto &indexSet = gridView.indexSet();
544 
545  for( const auto &element : elements( gridView, Partitions::interiorBorder ) )
546  {
547  if( gridPtr_.nofElParam_ > 0 )
548  std::swap( gridPtr_.elParam_[ indexSet.index( element ) ], elData_[ idSet_.id( element ) ] );
549 
550  if( gridPtr_.nofVtxParam_ > 0 )
551  {
552  for( unsigned int v = 0, n = element.subEntities( dimension ); v < n; ++v )
553  {
554  const auto index = indexSet.subIndex( element, v, dimension );
555  if ( !gridPtr_.vtxParam_[ index ].empty() )
556  std::swap( gridPtr_.vtxParam_[ index ], vtxData_[ idSet_.subId( element, v, dimension ) ] );
557  }
558  }
559 
560  if( element.hasBoundaryIntersections() )
561  {
562  for( const auto &intersection : intersections( gridView, element ) )
563  {
564  // dirty hack: check for "none" to make corner point grid work
565  if( !intersection.boundary() || intersection.type().isNone() )
566  continue;
567 
568  const int i = intersection.indexInInside();
569  auto &bndData = bndData_[ idSet_.subId( element, i, 1 ) ];
570  bndData.first = gridPtr_.bndId_[ indexSet.subIndex( element, i, 1 ) ];
571  if( gridPtr_.haveBndParam_ )
572  std::swap( bndData.second, gridPtr_.bndParam_[ intersection.boundarySegmentIndex() ] );
573  }
574  }
575  }
576  }
577 
578  DataHandle ( const DataHandle & ) = delete;
579  DataHandle ( DataHandle && ) = delete;
580 
581  ~DataHandle ()
582  {
583  const auto gridView = gridPtr_->levelGridView( 0 );
584  const auto &indexSet = gridView.indexSet();
585 
586  if( gridPtr_.nofElParam_ > 0 )
587  gridPtr_.elParam_.resize( indexSet.size( 0 ) );
588  if( gridPtr_.nofVtxParam_ > 0 )
589  gridPtr_.vtxParam_.resize( indexSet.size( dimension ) );
590  gridPtr_.bndId_.resize( indexSet.size( 1 ) );
591  if( gridPtr_.haveBndParam_ )
592  gridPtr_.bndParam_.resize( gridPtr_->numBoundarySegments() );
593 
594  for( const auto &element : elements( gridView, Partitions::all ) )
595  {
596  if( gridPtr_.nofElParam_ > 0 )
597  {
598  std::swap( gridPtr_.elParam_[ indexSet.index( element ) ], elData_[ idSet_.id( element ) ] );
599  assert( gridPtr_.elParam_[ indexSet.index( element ) ].size() == static_cast< std::size_t >( gridPtr_.nofElParam_ ) );
600  }
601 
602  if( gridPtr_.nofVtxParam_ > 0 )
603  {
604  for( unsigned int v = 0; v < element.subEntities( dimension ); ++v )
605  {
606  const auto index = indexSet.subIndex( element, v, dimension );
607  if( gridPtr_.vtxParam_[ index ].empty() )
608  std::swap( gridPtr_.vtxParam_[ index ], vtxData_[ idSet_.subId( element, v, dimension ) ] );
609  assert( gridPtr_.vtxParam_[ index ].size() == static_cast< std::size_t >( gridPtr_.nofVtxParam_ ) );
610  }
611  }
612 
613  if( element.hasBoundaryIntersections() )
614  {
615  for( const auto &intersection : intersections( gridView, element ) )
616  {
617  // dirty hack: check for "none" to make corner point grid work
618  if( !intersection.boundary() || intersection.type().isNone() )
619  continue;
620 
621  const int i = intersection.indexInInside();
622  auto &bndData = bndData_[ idSet_.subId( element, i, 1 ) ];
623  gridPtr_.bndId_[ indexSet.subIndex( element, i, 1 ) ] = bndData.first;
624  if( gridPtr_.haveBndParam_ )
625  std::swap( bndData.second, gridPtr_.bndParam_[ intersection.boundarySegmentIndex() ] );
626  }
627  }
628  }
629  }
630 
631  CommDataHandleIF< DataHandle, char > &interface () { return *this; }
632 
633  bool contains ( int dim, int codim ) const
634  {
635  assert( dim == dimension );
636  // do not use a switch statement, because dimension == 1 is possible
637  return (codim == 1) || ((codim == dimension) && (gridPtr_.nofVtxParam_ > 0)) || ((codim == 0) && (gridPtr_.nofElParam_ > 0));
638  }
639 
640  bool fixedSize (int /* dim */, int /* codim */) const { return false; }
641 
642  template< class Entity >
643  std::size_t size ( const Entity &entity ) const
644  {
645  std::size_t totalSize = 0;
646 
647  // do not use a switch statement, because dimension == 1 is possible
648  if( (Entity::codimension == 0) && (gridPtr_.nofElParam_ > 0) )
649  {
650  assert( elData_[ idSet_.id( entity ) ].size() == static_cast< std::size_t >( gridPtr_.nofElParam_ ) );
651  for( double &v : elData_[ idSet_.id( entity ) ] )
652  totalSize += dataSize( v );
653  }
654 
655  if( (Entity::codimension == dimension) && (gridPtr_.nofVtxParam_ > 0) )
656  {
657  assert( vtxData_[ idSet_.id( entity ) ].size() == static_cast< std::size_t >( gridPtr_.nofVtxParam_ ) );
658  for( double &v : vtxData_[ idSet_.id( entity ) ] )
659  totalSize += dataSize( v );
660  }
661 
662  if( Entity::codimension == 1 )
663  {
664  const auto bndData = bndData_.find( idSet_.id( entity ) );
665  if( bndData != bndData_.end() )
666  totalSize += dataSize( bndData->second.first ) + dataSize( bndData->second.second );
667  }
668 
669  return totalSize;
670  }
671 
672  template< class Buffer, class Entity >
673  void gather ( Buffer &buffer, const Entity &entity ) const
674  {
675  // do not use a switch statement, because dimension == 1 is possible
676  if( (Entity::codimension == 0) && (gridPtr_.nofElParam_ > 0) )
677  {
678  assert( elData_[ idSet_.id( entity ) ].size() == static_cast< std::size_t >( gridPtr_.nofElParam_ ) );
679  for( double &v : elData_[ idSet_.id( entity ) ] )
680  write( buffer, v );
681  }
682 
683  if( (Entity::codimension == dimension) && (gridPtr_.nofVtxParam_ > 0) )
684  {
685  assert( vtxData_[ idSet_.id( entity ) ].size() == static_cast< std::size_t >( gridPtr_.nofVtxParam_ ) );
686  for( double &v : vtxData_[ idSet_.id( entity ) ] )
687  write( buffer, v );
688  }
689 
690  if( Entity::codimension == 1 )
691  {
692  const auto bndData = bndData_.find( idSet_.id( entity ) );
693  if( bndData != bndData_.end() )
694  {
695  write( buffer, bndData->second.first );
696  write( buffer, bndData->second.second );
697  }
698  }
699  }
700 
701  template< class Buffer, class Entity >
702  void scatter ( Buffer &buffer, const Entity &entity, std::size_t n )
703  {
704  // do not use a switch statement, because dimension == 1 is possible
705  if( (Entity::codimension == 0) && (gridPtr_.nofElParam_ > 0) )
706  {
707  auto &p = elData_[ idSet_.id( entity ) ];
708  p.resize( gridPtr_.nofElParam_ );
709  for( double &v : p )
710  read( buffer, v, n );
711  }
712 
713  if( (Entity::codimension == dimension) && (gridPtr_.nofVtxParam_ > 0) )
714  {
715  auto &p = vtxData_[ idSet_.id( entity ) ];
716  p.resize( gridPtr_.nofVtxParam_ );
717  for( double &v : p )
718  read( buffer, v, n );
719  }
720 
721  if( (Entity::codimension == 1) && (n > 0) )
722  {
723  auto &bndData = bndData_[ idSet_.id( entity ) ];
724  read( buffer, bndData.first, n );
725  read( buffer, bndData.second, n );
726  }
727 
728  assert( n == 0 );
729  }
730 
731  private:
732  template< class T >
733  static std::enable_if_t< std::is_trivially_copyable< T >::value, std::size_t > dataSize ( const T & /* value */ )
734  {
735  return sizeof( T );
736  }
737 
738  static std::size_t dataSize ( const std::string &s )
739  {
740  return dataSize( s.size() ) + s.size();
741  }
742 
743  template< class Buffer, class T >
744  static std::enable_if_t< std::is_trivially_copyable< T >::value > write ( Buffer &buffer, const T &value )
745  {
746  std::array< char, sizeof( T ) > bytes;
747  std::memcpy( bytes.data(), &value, sizeof( T ) );
748  for( char &b : bytes )
749  buffer.write( b );
750  }
751 
752  template< class Buffer >
753  static void write ( Buffer &buffer, const std::string &s )
754  {
755  write( buffer, s.size() );
756  for( const char &c : s )
757  buffer.write( c );
758  }
759 
760  template< class Buffer, class T >
761  static std::enable_if_t< std::is_trivially_copyable< T >::value > read ( Buffer &buffer, T &value, std::size_t &n )
762  {
763  assert( n >= sizeof( T ) );
764  n -= sizeof( T );
765 
766  std::array< char, sizeof( T ) > bytes;
767  for( char &b : bytes )
768  buffer.read( b );
769  std::memcpy( &value, bytes.data(), sizeof( T ) );
770  }
771 
772  template< class Buffer >
773  static void read ( Buffer &buffer, std::string &s, std::size_t &n )
774  {
775  std::size_t size;
776  read( buffer, size, n );
777  s.resize( size );
778 
779  assert( n >= size );
780  n -= size;
781 
782  for( char &c : s )
783  buffer.read( c );
784  }
785 
786  GridPtr &gridPtr_;
787  const typename GridType::LocalIdSet &idSet_;
788  mutable std::map< typename GridType::LocalIdSet::IdType, std::vector< double > > elData_, vtxData_;
789  mutable std::map< typename GridType::LocalIdSet::IdType, std::pair< int, DGFBoundaryParameter::type > > bndData_;
790  };
791 
792  // grid auto pointer
793  mutable mygrid_ptr gridPtr_;
794  // element and vertex parameters
795  std::vector< std::vector< double > > elParam_;
796  std::vector< std::vector< double > > vtxParam_;
797  std::vector< DGFBoundaryParameter::type > bndParam_;
798  std::vector< int > bndId_;
799  std::vector< double > emptyParam_;
800 
801  int nofElParam_;
802  int nofVtxParam_;
803  bool haveBndParam_;
804  }; // end of class GridPtr
805 
806 } // end namespace Dune
807 
808 #endif
void scatter(MessageBufferImp &buff, const EntityType &e, size_t n)
unpack data from message buffer to user.
Definition: datahandleif.hh:143
void gather(MessageBufferImp &buff, const EntityType &e) const
pack data from user to message buffer
Definition: datahandleif.hh:129
Wrapper class for entities.
Definition: entity.hh:66
constexpr static int codimension
Know your own codimension.
Definition: entity.hh:106
static std::unique_ptr< Grid > read(const std::string &fileName, bool verbose=true, bool insertBoundarySegments=true)
Definition: gmshreader.hh:935
Base class for exceptions in Dune grid modules.
Definition: exceptions.hh:20
Provide a generic factory class for unstructured grids.
Definition: gridfactory.hh:275
Grid view abstract base class.
Definition: gridview.hh:66
Intersection of a mesh entity of codimension 0 ("element") with a "neighboring" element or with the d...
Definition: intersection.hh:164
size_t boundarySegmentIndex() const
index of the boundary segment within the macro grid
Definition: intersection.hh:236
MPI_Comm MPICommunicator
The type of the mpi communicator.
Definition: mpihelper.hh:192
static MPICommunicator getCommunicator()
get the default communicator
Definition: mpihelper.hh:200
Default exception for dummy implementations.
Definition: exceptions.hh:263
Describes the parallel communication interface class for MessageBuffers and DataHandles.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
const IndexSet & indexSet() const
obtain the index set
Definition: gridview.hh:177
@ ForwardCommunication
communicate as given in InterfaceType
Definition: gridenums.hh:171
@ InteriorBorder_All_Interface
send interior and border, receive all entities
Definition: gridenums.hh:88
concept Intersection
Model of an intersection.
Definition: intersection.hh:23
concept GridView
Model of a grid view.
Definition: gridview.hh:81
concept Entity
Model of a grid entity.
Definition: entity.hh:107
Helpers for dealing with MPI.
constexpr All all
PartitionSet for all partitions.
Definition: partitionset.hh:295
constexpr InteriorBorder interiorBorder
PartitionSet for the interior and border partitions.
Definition: partitionset.hh:286
Dune namespace.
Definition: alignedallocator.hh:13
constexpr std::bool_constant<(sizeof...(II)==0)> empty(std::integer_sequence< T, II... >)
Checks whether the sequence is empty.
Definition: integersequence.hh:80
constexpr std::bool_constant<((II==value)||...)> contains(std::integer_sequence< T, II... >, std::integral_constant< T, value >)
Checks whether or not a given sequence contains a value.
Definition: integersequence.hh:137
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
constexpr auto get(std::integer_sequence< T, II... >, std::integral_constant< std::size_t, pos >={})
Return the entry at position pos of the given sequence.
Definition: integersequence.hh:22
This file implements several utilities related to std::shared_ptr.
static const type & defaultValue()
default constructor
Definition: parser.hh:28
std::string type
type of additional boundary parameters
Definition: parser.hh:25
Class for constructing grids from DGF files.
Definition: gridptr.hh:66
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:325
GridPtr & operator=(const GridPtr &org)
assignment of grid pointer
Definition: gridptr.hh:246
GridPtr(std::istream &input, MPICommunicatorType comm=MPIHelper::getCommunicator())
constructor given a std::istream
Definition: gridptr.hh:199
GridType * operator->()
return pointer to GridType instance
Definition: gridptr.hh:283
GridPtr()
Default constructor, creating empty GridPtr.
Definition: gridptr.hh:217
GridPtr(GridType *grd)
Constructor storing given pointer to internal auto pointer.
Definition: gridptr.hh:230
GridType * release()
release pointer from internal ownership
Definition: gridptr.hh:298
const GridType * operator->() const
return const pointer to GridType instance
Definition: gridptr.hh:293
int nofParameters(const Entity &) const
get parameters defined for given entity
Definition: gridptr.hh:311
GridPtr(const std::string &filename, MPICommunicatorType comm=MPIHelper::getCommunicator())
constructor given the name of a DGF file
Definition: gridptr.hh:158
const GridType & operator*() const
return const reference to GridType instance
Definition: gridptr.hh:288
int nofParameters(int cdim) const
get number of parameters defined for a given codimension
Definition: gridptr.hh:301
const DGFBoundaryParameter::type & parameters(const Intersection< GridImp, IntersectionImp > &intersection) const
get parameters for intersection
Definition: gridptr.hh:351
int nofParameters(const Intersection< GridImp, IntersectionImp > &intersection) const
get number of parameters defined for a given intersection
Definition: gridptr.hh:318
GridPtr(const GridPtr &org)=default
Copy constructor, copies internal auto pointer.
GridType & operator*()
return reference to GridType instance
Definition: gridptr.hh:278
implements the Deleter concept of shared_ptr without deleting anything
Definition: shared_ptr.hh:49
std::size_t fixedSize
The number of data items per index if it is fixed, 0 otherwise.
Definition: variablesizecommunicator.hh:264
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (Apr 24, 22:30, 2024)