Dune Core Modules (2.7.0)

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