grid.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 2; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=2 sw=2 sts=2:
3#ifndef DUNE_POLYHEDRALGRID_GRID_HH
4#define DUNE_POLYHEDRALGRID_GRID_HH
5
6#include <set>
7#include <vector>
8
9// Warning suppression for Dune includes.
11
12//- dune-common includes
13#include <dune/common/version.hh>
14#include <dune/common/parallel/mpihelper.hh>
15
16//- dune-grid includes
17#include <dune/grid/common/grid.hh>
18
19#include <dune/common/parallel/communication.hh>
20
21//- polyhedralgrid includes
29
30// Re-enable warnings.
32
33#include <opm/common/ErrorMacros.hpp>
34
36#include <opm/grid/cart_grid.h>
41
42#if HAVE_ECL_INPUT
43#include <opm/input/eclipse/EclipseState/Grid/EclipseGrid.hpp>
44#endif
45
46namespace Dune
47{
48
49
50 // PolyhedralGridFamily
51 // ------------
52
53 template< int dim, int dimworld, typename coord_t >
55 {
56 struct Traits
57 {
59
60 typedef coord_t ctype;
61
62 // type of data passed to entities, intersections, and iterators
63 // for PolyhedralGrid this is just an empty place holder
64 typedef const Grid* ExtraData;
65
66 typedef int Index ;
67
68 static const int dimension = dim;
69 static const int dimensionworld = dimworld;
70
71 typedef Dune::FieldVector< ctype, dimensionworld > GlobalCoordinate ;
72
77
78 typedef Dune::Intersection< const Grid, LeafIntersectionImpl > LeafIntersection;
79 typedef Dune::Intersection< const Grid, LevelIntersectionImpl > LevelIntersection;
80
81 typedef Dune::IntersectionIterator< const Grid, LeafIntersectionIteratorImpl, LeafIntersectionImpl > LeafIntersectionIterator;
82 typedef Dune::IntersectionIterator< const Grid, LevelIntersectionIteratorImpl, LevelIntersectionImpl > LevelIntersectionIterator;
83
85 typedef Dune::EntityIterator< 0, const Grid, HierarchicIteratorImpl > HierarchicIterator;
86
87 template< int codim >
88 struct Codim
89 {
91 typedef Dune::Geometry< dimension-codim, dimensionworld, const Grid, PolyhedralGridGeometry > Geometry;
92
94 typedef Dune::Geometry< dimension-codim, dimensionworld, const Grid, PolyhedralGridLocalGeometry > LocalGeometry;
95
97 typedef Dune::Entity< codim, dimension, const Grid, PolyhedralGridEntity > Entity;
98
101
102 //typedef Dune::EntitySeed< const Grid, PolyhedralGridEntitySeed< codim, const Grid > > EntitySeed;
104
105 template< PartitionIteratorType pitype >
107 {
109 typedef Dune::EntityIterator< codim, const Grid, LeafIteratorImpl > LeafIterator;
110
112 };
113
116 };
117
120
123
124 typedef Dune::MPIHelper::MPICommunicator MPICommunicator;
125 using Communication = Dune::Communication<MPICommunicator>;
126 using CollectiveCommunication = Dune::Communication<MPICommunicator>;
127
128 template< PartitionIteratorType pitype >
130 {
131 typedef Dune::GridView< PolyhedralGridViewTraits< dim, dimworld, ctype, pitype > > LeafGridView;
132 typedef Dune::GridView< PolyhedralGridViewTraits< dim, dimworld, ctype, pitype > > LevelGridView;
133 };
134
137 };
138 };
139
140
141
142 // PolyhedralGrid
143 // --------------
144
153 template < int dim, int dimworld, typename coord_t >
156 : public GridDefaultImplementation
157 < dim, dimworld, coord_t, PolyhedralGridFamily< dim, dimworld, coord_t > >
159 {
161
162 typedef GridDefaultImplementation
163 < dim, dimworld, coord_t, PolyhedralGridFamily< dim, dimworld, coord_t > > Base;
164
165 public:
167
168 protected:
170 {
171 inline void operator () ( UnstructuredGridType* grdPtr )
172 {
173 destroy_grid( grdPtr );
174 }
175 };
176
177 public:
178 typedef std::unique_ptr< UnstructuredGridType, UnstructuredGridDeleter > UnstructuredGridPtr;
179
181 allocateGrid ( std::size_t nCells, std::size_t nFaces, std::size_t nFaceNodes, std::size_t nCellFaces, std::size_t nNodes )
182 {
183 // Note that we here assign a grid of dimension dimworld in order to obtain global coordinates in the correct dimension
184 UnstructuredGridType *grid = allocate_grid( dimworld, nCells, nFaces, nFaceNodes, nCellFaces, nNodes );
185 if( !grid )
186 DUNE_THROW( GridError, "Unable to allocate grid" );
187 return UnstructuredGridPtr( grid );
188 }
189
190 static void
192 {
193 // get C pointer to UnstructuredGrid
194 UnstructuredGrid* ugPtr = ug.operator ->();
195
196 // compute geometric quantities like cell volume and face normals
197 compute_geometry( ugPtr );
198 }
199
208 typedef typename GridFamily::Traits Traits;
209
216 template< int codim >
217 struct Codim;
218
225 typedef typename Traits::HierarchicIterator HierarchicIterator;
227 typedef typename Traits::LeafIntersectionIterator LeafIntersectionIterator;
229 typedef typename Traits::LevelIntersectionIterator LevelIntersectionIterator;
230
237 template< PartitionIteratorType pitype >
239 {
240 typedef typename GridFamily::Traits::template Partition< pitype >::LevelGridView LevelGridView;
241 typedef typename GridFamily::Traits::template Partition< pitype >::LeafGridView LeafGridView;
242 };
243
247
261 typedef typename Traits::LeafIndexSet LeafIndexSet;
262
271 typedef typename Traits::LevelIndexSet LevelIndexSet;
272
283 typedef typename Traits::GlobalIdSet GlobalIdSet;
284
300 typedef typename Traits::LocalIdSet LocalIdSet;
301
308 typedef typename Traits::ctype ctype;
309
311 using Communication = typename Traits::Communication;
313
314 typedef typename Traits :: GlobalCoordinate GlobalCoordinate;
315
321#if HAVE_ECL_INPUT
327 explicit PolyhedralGrid ( const Opm::EclipseGrid& inputGrid,
328 const std::vector<double>& poreVolumes = std::vector<double> ())
329 : gridPtr_( createGrid( inputGrid, poreVolumes ) ),
330 grid_( *gridPtr_ ),
331 comm_( MPIHelper::getCommunicator() ),
332 leafIndexSet_( *this ),
333 globalIdSet_( *this ),
334 localIdSet_( *this ),
335 nBndSegments_( 0 )
336 {
337 init();
338 }
339#endif
340
346 explicit PolyhedralGrid ( const std::vector< int >& n,
347 const std::vector< double >& dx )
348 : gridPtr_( createGrid( n, dx ) ),
349 grid_( *gridPtr_ ),
350 comm_( MPIHelper::getCommunicator()),
351 leafIndexSet_( *this ),
352 globalIdSet_( *this ),
353 localIdSet_( *this ),
354 nBndSegments_( 0 )
355 {
356 init();
357 }
358
365 explicit PolyhedralGrid ( UnstructuredGridPtr &&gridPtr )
366 : gridPtr_( std::move( gridPtr ) ),
367 grid_( *gridPtr_ ),
368 comm_( MPIHelper::getCommunicator() ),
369 leafIndexSet_( *this ),
370 globalIdSet_( *this ),
371 localIdSet_( *this ),
372 nBndSegments_( 0 )
373 {
374 init();
375 }
376
384 explicit PolyhedralGrid ( const UnstructuredGridType& grid )
385 : gridPtr_(),
386 grid_( grid ),
387 comm_( MPIHelper::getCommunicator() ),
388 leafIndexSet_( *this ),
389 globalIdSet_( *this ),
390 localIdSet_( *this ),
391 nBndSegments_( 0 )
392 {
393 init();
394 }
395
400 operator const UnstructuredGridType& () const { return grid_; }
401
414 int maxLevel () const
415 {
416 return 0;
417 }
418
427 int size ( int /* level */, int codim ) const
428 {
429 return size( codim );
430 }
431
438 int size ( int codim ) const
439 {
440 if( codim == 0 )
441 {
442 return grid_.number_of_cells;
443 }
444 else if ( codim == 1 )
445 {
446 return grid_.number_of_faces;
447 }
448 else if ( codim == dim )
449 {
450 return grid_.number_of_nodes;
451 }
452 else
453 {
454 std::cerr << "Warning: codimension " << codim << " not available in PolyhedralGrid" << std::endl;
455 return 0;
456 }
457 }
458
467 int size ( int /* level */, GeometryType type ) const
468 {
469 return size( dim - type.dim() );
470 }
471
476 int size ( GeometryType type ) const
477 {
478 return size( dim - type.dim() );
479 }
480
487 size_t numBoundarySegments () const
488 {
489 return nBndSegments_;
490 }
493 template< int codim >
495 {
496 return leafbegin< codim, All_Partition >();
497 }
498
499 template< int codim >
501 {
502 return leafend< codim, All_Partition >();
503 }
504
505 template< int codim, PartitionIteratorType pitype >
506 typename Codim< codim >::template Partition< pitype >::LeafIterator
507 leafbegin () const
508 {
509 typedef typename Traits::template Codim< codim >::template Partition< pitype >::LeafIteratorImpl Impl;
510 return Impl( extraData(), true );
511 }
512
513 template< int codim, PartitionIteratorType pitype >
514 typename Codim< codim >::template Partition< pitype >::LeafIterator
515 leafend () const
516 {
517 typedef typename Traits::template Codim< codim >::template Partition< pitype >::LeafIteratorImpl Impl;
518 return Impl( extraData(), false );
519 }
520
521 template< int codim >
522 typename Codim< codim >::LevelIterator lbegin ( const int /* level */ ) const
523 {
524 return leafbegin< codim, All_Partition >();
525 }
526
527 template< int codim >
528 typename Codim< codim >::LevelIterator lend ( const int /* level */ ) const
529 {
530 return leafend< codim, All_Partition >();
531 }
532
533 template< int codim, PartitionIteratorType pitype >
534 typename Codim< codim >::template Partition< pitype >::LevelIterator
535 lbegin ( const int /* level */ ) const
536 {
537 return leafbegin< codim, pitype > ();
538 }
539
540 template< int codim, PartitionIteratorType pitype >
541 typename Codim< codim >::template Partition< pitype >::LevelIterator
542 lend ( const int /* level */ ) const
543 {
544 return leafend< codim, pitype > ();
545 }
546
547 const GlobalIdSet &globalIdSet () const
548 {
549 return globalIdSet_;
550 }
551
552 const LocalIdSet &localIdSet () const
553 {
554 return localIdSet_;
555 }
556
557 const LevelIndexSet &levelIndexSet ( int /* level */ ) const
558 {
559 return leafIndexSet();
560 }
561
563 {
564 return leafIndexSet_;
565 }
566
567 void globalRefine ( int /* refCount */ )
568 {
569 }
570
571 bool mark ( int /* refCount */, const typename Codim< 0 >::Entity& /* entity */ )
572 {
573 return false;
574 }
575
576 int getMark ( const typename Codim< 0 >::Entity& /* entity */) const
577 {
578 return false;
579 }
580
582 bool preAdapt ()
583 {
584 return false;
585 }
586
588 bool adapt ()
589 {
590 return false ;
591 }
592
597 template< class DataHandle >
598 bool adapt ( DataHandle & )
599 {
600 return false;
601 }
602
604 void postAdapt ()
605 {
606 }
607
615 int overlapSize ( int /* codim */) const
616 {
617 return 0;
618 }
619
624 int ghostSize( int codim ) const
625 {
626 return (codim == 0 ) ? 1 : 0;
627 }
628
634 int overlapSize ( int /* level */, int /* codim */ ) const
635 {
636 return 0;
637 }
638
644 int ghostSize ( int /* level */, int codim ) const
645 {
646 return ghostSize( codim );
647 }
648
662 template< class DataHandle>
663 void communicate ( DataHandle& /* dataHandle */,
664 InterfaceType /* interface */,
665 CommunicationDirection /* direction */,
666 int /* level */ ) const
667 {
668 OPM_THROW(std::runtime_error, "communicate not implemented for polyhedreal grid!");
669 }
670
683 template< class DataHandle>
684 void communicate ( DataHandle& /* dataHandle */,
685 InterfaceType /* interface */,
686 CommunicationDirection /* direction */ ) const
687 {
688 OPM_THROW(std::runtime_error, "communicate not implemented for polyhedreal grid!");
689 }
690
693 {
694 OPM_THROW(std::runtime_error, "switch to global view not implemented for polyhedreal grid!");
695 }
696
699 {
700 OPM_THROW(std::runtime_error, "switch to distributed view not implemented for polyhedreal grid!");
701 }
702
711 const CommunicationType &comm () const
712 {
713 return comm_;
714 }
715
716 // data handle interface different between geo and interface
717
728 {
729 return false ;
730 }
731
747 template< class DataHandle, class Data >
748 bool loadBalance ( CommDataHandleIF< DataHandle, Data >& /* datahandle */ )
749 {
750 return false;
751 }
752
767 template< class DofManager >
768 bool loadBalance ( DofManager& /* dofManager */ )
769 {
770 return false;
771 }
772
774 template< PartitionIteratorType pitype >
775 typename Partition< pitype >::LevelGridView levelGridView ( int /* level */ ) const
776 {
777 typedef typename Partition< pitype >::LevelGridView View;
778 typedef typename View::GridViewImp ViewImp;
779 return View( ViewImp( *this ) );
780 }
781
783 template< PartitionIteratorType pitype >
785 {
786 typedef typename Traits::template Partition< pitype >::LeafGridView View;
787 typedef typename View::GridViewImp ViewImp;
788 return View( ViewImp( *this ) );
789 }
790
792 LevelGridView levelGridView ( int /* level */ ) const
793 {
794 typedef typename LevelGridView::GridViewImp ViewImp;
795 return LevelGridView( ViewImp( *this ) );
796 }
797
800 {
801 typedef typename LeafGridView::GridViewImp ViewImp;
802 return LeafGridView( ViewImp( *this ) );
803 }
804
806 template< class EntitySeed >
808 entityPointer ( const EntitySeed &seed ) const
809 {
810 typedef typename Traits::template Codim< EntitySeed::codimension >::EntityPointer EntityPointer;
811 typedef typename Traits::template Codim< EntitySeed::codimension >::EntityPointerImpl EntityPointerImpl;
812 typedef typename Traits::template Codim< EntitySeed::codimension >::EntityImpl EntityImpl;
813 return EntityPointer( EntityPointerImpl( EntityImpl( extraData(), seed ) ) );
814 }
815
817 template< class EntitySeed >
818 typename Traits::template Codim< EntitySeed::codimension >::Entity
819 entity ( const EntitySeed &seed ) const
820 {
821 typedef typename Traits::template Codim< EntitySeed::codimension >::EntityImpl EntityImpl;
822 return EntityImpl( extraData(), seed );
823 }
824
838 void update ()
839 {
840 }
841
844 const std::array<int, 3>& logicalCartesianSize() const
845 {
846 return cartDims_;
847 }
848
849 const int* globalCell() const
850 {
851 assert( grid_.global_cell != 0 );
852 return grid_.global_cell;
853 }
854
855 const int* globalCellPtr() const
856 {
857 return grid_.global_cell;
858 }
859
860 void getIJK(const int c, std::array<int,3>& ijk) const
861 {
862 int gc = globalCell()[c];
863 ijk[0] = gc % logicalCartesianSize()[0]; gc /= logicalCartesianSize()[0];
864 ijk[1] = gc % logicalCartesianSize()[1];
865 ijk[2] = gc / logicalCartesianSize()[1];
866 }
867
872
873
882 template<class DataHandle>
883 void scatterData([[maybe_unused]] DataHandle& handle) const
884 {
885 OPM_THROW(std::runtime_error, "ScatterData not implemented for polyhedral grid!");
886 }
887
888 protected:
889#if HAVE_ECL_INPUT
890 UnstructuredGridType* createGrid( const Opm::EclipseGrid& inputGrid, const std::vector< double >& poreVolumes ) const
891 {
892 struct grdecl g;
893
894 g.dims[0] = inputGrid.getNX();
895 g.dims[1] = inputGrid.getNY();
896 g.dims[2] = inputGrid.getNZ();
897
898 std::vector<double> coord = inputGrid.getCOORD( );
899 std::vector<double> zcorn = inputGrid.getZCORN( );
900 std::vector<int> actnum = inputGrid.getACTNUM( );
901
902 g.coord = coord.data();
903 g.zcorn = zcorn.data();
904 g.actnum = actnum.data();
905
906 const double z_tolerance = inputGrid.isPinchActive() ? inputGrid.getPinchThresholdThickness() : 0.0;
907
908 if (!poreVolumes.empty() && (inputGrid.getMinpvMode() != Opm::MinpvMode::Inactive))
909 {
910 Opm::MinpvProcessor mp(g.dims[0], g.dims[1], g.dims[2]);
911 const std::vector<double>& minpvv = inputGrid.getMinpvVector();
912 // Currently the pinchProcessor is not used and only opmfil is supported
913 // The polyhedralgrid only only supports the opmfil option
914 //bool opmfil = inputGrid.getMinpvMode() == Opm::MinpvMode::OpmFIL;
915 bool opmfil = true;
916 const size_t cartGridSize = g.dims[0] * g.dims[1] * g.dims[2];
917 std::vector<double> thickness(cartGridSize);
918 for (size_t i = 0; i < cartGridSize; ++i) {
919 thickness[i] = inputGrid.getCellThickness(i);
920 }
921 mp.process(thickness, z_tolerance, inputGrid.getPinchMaxEmptyGap() , poreVolumes,
922 minpvv, actnum, opmfil, zcorn.data());
923 }
924
925 /*
926 if (!poreVolumes.empty() && (inputGrid.getMinpvMode() != Opm::MinpvMode::Inactive)) {
927 Opm::MinpvProcessor mp(g.dims[0], g.dims[1], g.dims[2]);
928 const double minpv_value = inputGrid.getMinpvValue();
929 // Currently the pinchProcessor is not used and only opmfil is supported
930 //bool opmfil = inputGrid.getMinpvMode() == Opm::MinpvMode::OpmFIL;
931 bool opmfil = true;
932 mp.process(poreVolumes, minpv_value, actnum, opmfil, zcorn.data());
933 }
934 */
935
936 UnstructuredGridType* cgrid = create_grid_cornerpoint(&g, z_tolerance);
937 if (!cgrid) {
938 OPM_THROW(std::runtime_error, "Failed to construct grid.");
939 }
940 return cgrid;
941 }
942#endif
943
944 UnstructuredGridType* createGrid( const std::vector< int >& n, const std::vector< double >& dx ) const
945 {
946 UnstructuredGridType* cgrid = nullptr ;
947 assert( int(n.size()) == dim );
948 if( dim == 2 )
949 {
950 cgrid = create_grid_cart2d( n[ 0 ], n[ 1 ], dx[ 0 ], dx[ 1 ] );
951 }
952 else if ( dim == 3 )
953 {
954 cgrid = create_grid_hexa3d( n[ 0 ], n[ 1 ], n[ 2 ], dx[ 0 ], dx[ 1 ], dx[ 2 ] );
955 }
956
957 //print_grid( cgrid );
958 if (!cgrid) {
959 OPM_THROW(std::runtime_error, "Failed to construct grid.");
960 }
961 return cgrid;
962 }
963
964 public:
965#if DUNE_VERSION_LT_REV(DUNE_GRID, 2, 7, 1)
966 using Base::getRealImplementation;
967#endif
968 typedef typename Traits :: ExtraData ExtraData;
969 ExtraData extraData () const { return this; }
970
971 template <class EntitySeed>
972 int corners( const EntitySeed& seed ) const
973 {
974 const int codim = EntitySeed :: codimension;
975 const int index = seed.index();
976 if (codim==0)
977 return cellVertices_[ index ].size();
978 if (codim==1)
979 return grid_.face_nodepos[ index+1 ] - grid_.face_nodepos[ index ];
980 if (codim==dim)
981 return 1;
982 return 0;
983 }
984
985 template <class EntitySeed>
987 corner ( const EntitySeed& seed, const int i ) const
988 {
989 const int codim = EntitySeed :: codimension;
990 if (codim==0)
991 {
992 const int coordIndex = GlobalCoordinate :: dimension * cellVertices_[ seed.index() ][ i ];
993 return copyToGlobalCoordinate( grid_.node_coordinates + coordIndex );
994 }
995 if (codim==1)
996 {
997 // for faces we need to swap vertices in 3d since in UnstructuredGrid
998 // those are ordered counter clockwise, for 2d this does not matter
999 // TODO: Improve this for performance reasons
1000 const int crners = corners( seed );
1001 const int crner = (crners == 4 && EntitySeed :: dimension == 3 && i > 1 ) ? 5 - i : i;
1002 const int faceVertex = grid_.face_nodes[ grid_.face_nodepos[seed.index() ] + crner ];
1003 return copyToGlobalCoordinate( grid_.node_coordinates + GlobalCoordinate :: dimension * faceVertex );
1004 }
1005 if (codim==dim)
1006 {
1007 const int coordIndex = GlobalCoordinate :: dimension * seed.index();
1008 return copyToGlobalCoordinate( grid_.node_coordinates + coordIndex );
1009 }
1010 return GlobalCoordinate( 0 );
1011 }
1012
1013 template <class EntitySeed>
1014 int subEntities( const EntitySeed& seed, const int codim ) const
1015 {
1016 const int index = seed.index();
1017 if( seed.codimension == 0 )
1018 {
1019 if (codim==0)
1020 return 1;
1021 if (codim==1)
1022 return grid_.cell_facepos[ index+1 ] - grid_.cell_facepos[ index ];
1023 if (codim==dim)
1024 return cellVertices_[ index ].size();
1025 }
1026 else if( seed.codimension == 1 )
1027 {
1028 if (codim==1)
1029 return 1;
1030 if (codim==dim)
1031 return grid_.face_nodepos[ index+1 ] - grid_.face_nodepos[ index ];
1032 }
1033 else if ( seed.codimension == dim )
1034 {
1035 return 1 ;
1036 }
1037
1038 return 0;
1039 }
1040
1041 template <int codim, class EntitySeedArg >
1042 typename Codim<codim>::EntitySeed
1043 subEntitySeed( const EntitySeedArg& baseSeed, const int i ) const
1044 {
1045 assert( codim >= EntitySeedArg::codimension );
1046 assert( i>= 0 && i<subEntities( baseSeed, codim ) );
1047 typedef typename Codim<codim>::EntitySeed EntitySeed;
1048
1049 // if codim equals entity seed codim just return same entity seed.
1050 if( codim == EntitySeedArg::codimension )
1051 {
1052 return EntitySeed( baseSeed.index() );
1053 }
1054
1055 if( EntitySeedArg::codimension == 0 )
1056 {
1057 if ( codim == 1 )
1058 {
1059 return EntitySeed( grid_.cell_faces[ grid_.cell_facepos[ baseSeed.index() ] + i ] );
1060 }
1061 else if ( codim == dim )
1062 {
1063 return EntitySeed( cellVertices_[ baseSeed.index() ][ i ] );
1064 }
1065 }
1066 else if ( EntitySeedArg::codimension == 1 && codim == dim )
1067 {
1068 return EntitySeed( grid_.face_nodes[ grid_.face_nodepos[ baseSeed.index() + i ] ]);
1069 }
1070
1071 DUNE_THROW(NotImplemented,"codimension not available");
1072 return EntitySeed();
1073 }
1074
1075 template <int codim>
1076 typename Codim<codim>::EntitySeed
1077 subEntitySeed( const typename Codim<1>::EntitySeed& faceSeed, const int i ) const
1078 {
1079 assert( i>= 0 && i<subEntities( faceSeed, codim ) );
1080 typedef typename Codim<codim>::EntitySeed EntitySeed;
1081 if ( codim == 1 )
1082 {
1083 return EntitySeed( faceSeed.index() );
1084 }
1085 else if ( codim == dim )
1086 {
1087 return EntitySeed( grid_.face_nodes[ grid_.face_nodepos[ faceSeed.index() ] + i ] );
1088 }
1089 else
1090 {
1091 DUNE_THROW(NotImplemented,"codimension not available");
1092 }
1093 }
1094
1095 bool hasBoundaryIntersections(const typename Codim<0>::EntitySeed& seed ) const
1096 {
1097 const int faces = subEntities( seed, 1 );
1098 for( int f=0; f<faces; ++f )
1099 {
1100 const auto faceSeed = this->template subEntitySeed<1>( seed, f );
1101 if( isBoundaryFace( faceSeed ) )
1102 return true;
1103 }
1104 return false;
1105 }
1106
1107 bool isBoundaryFace(const int face ) const
1108 {
1109 assert( face >= 0 && face < grid_.number_of_faces );
1110 const int facePos = 2 * face;
1111 return ((grid_.face_cells[ facePos ] < 0) || (grid_.face_cells[ facePos+1 ] < 0));
1112 }
1113
1114 bool isBoundaryFace(const typename Codim<1>::EntitySeed& faceSeed ) const
1115 {
1116 assert( faceSeed.isValid() );
1117 return isBoundaryFace( faceSeed.index() );
1118 }
1119
1120 int boundarySegmentIndex(const typename Codim<0>::EntitySeed& seed, const int face ) const
1121 {
1122 const auto faceSeed = this->template subEntitySeed<1>( seed, face );
1123 assert( faceSeed.isValid() );
1124 const int facePos = 2 * faceSeed.index();
1125 const int idx = std::min( grid_.face_cells[ facePos ], grid_.face_cells[ facePos+1 ]);
1126 // check that this is actually the boundary
1127 assert( idx < 0 );
1128 return -(idx+1); // +1 to include 0 boundary segment index
1129 }
1130
1131 const std::vector< GeometryType > &geomTypes ( const unsigned int codim ) const
1132 {
1133 static std::vector< GeometryType > emptyDummy;
1134 if (codim < geomTypes_.size())
1135 {
1136 return geomTypes_[codim];
1137 }
1138
1139 return emptyDummy;
1140 }
1141
1142 template < class Seed >
1143 GeometryType geometryType( const Seed& seed ) const
1144 {
1145 if( Seed::codimension == 0 )
1146 {
1147 assert(!geomTypes( Seed::codimension ).empty());
1148 return geomTypes( Seed::codimension )[ 0 ];
1149 }
1150 else
1151 {
1152 // 3d faces
1153 if( dim == 3 && Seed::codimension == 1 )
1154 {
1155 GeometryType face;
1156 const int nVx = corners( seed );
1157 if( nVx == 4 ) // quad face
1158 face = Dune::GeometryTypes::cube(2);
1159 else if( nVx == 3 ) // triangle face
1160 face = Dune::GeometryTypes::simplex(2);
1161 else // polygonal face
1162 face = Dune::GeometryTypes::none(2);
1163 return face;
1164 }
1165
1166 // for codim 2 and codim 3 there is only one option
1167 assert(!geomTypes( Seed::codimension ).empty());
1168 return geomTypes( Seed::codimension )[ 0 ];
1169 }
1170 }
1171
1172 int indexInInside( const typename Codim<0>::EntitySeed& seed, const int i ) const
1173 {
1174 return ( grid_.cell_facetag ) ? cartesianIndexInInside( seed, i ) : i;
1175 }
1176
1177 int cartesianIndexInInside( const typename Codim<0>::EntitySeed& seed, const int i ) const
1178 {
1179 assert( i>= 0 && i<subEntities( seed, 1 ) );
1180 return grid_.cell_facetag[ grid_.cell_facepos[ seed.index() ] + i ] ;
1181 }
1182
1183 typename Codim<0>::EntitySeed
1184 neighbor( const typename Codim<0>::EntitySeed& seed, const int i ) const
1185 {
1186 const int face = this->template subEntitySeed<1>( seed, i ).index();
1187 int nb = grid_.face_cells[ 2 * face ];
1188 if( nb == seed.index() )
1189 {
1190 nb = grid_.face_cells[ 2 * face + 1 ];
1191 }
1192
1193 typedef typename Codim<0>::EntitySeed EntitySeed;
1194 return EntitySeed( nb );
1195 }
1196
1197 int
1198 indexInOutside( const typename Codim<0>::EntitySeed& seed, const int i ) const
1199 {
1200 if( grid_.cell_facetag )
1201 {
1202 // if cell_facetag is present we assume pseudo Cartesian corner point case
1203 const int in_inside = cartesianIndexInInside( seed, i );
1204 return in_inside + ((in_inside % 2) ? -1 : 1);
1205 }
1206 else
1207 {
1208 typedef typename Codim<0>::EntitySeed EntitySeed;
1209 EntitySeed nb = neighbor( seed, i );
1210 const int faces = subEntities( seed, 1 );
1211 for( int face = 0; face<faces; ++ face )
1212 {
1213 if( neighbor( nb, face ).equals(seed) )
1214 {
1215 return indexInInside( nb, face );
1216 }
1217 }
1218 DUNE_THROW(InvalidStateException,"inverse intersection not found");
1219 return -1;
1220 }
1221 }
1222
1223 template <class EntitySeed>
1225 outerNormal( const EntitySeed& seed, const int i ) const
1226 {
1227 const int face = this->template subEntitySeed<1>( seed, i ).index();
1228 const int normalIdx = face * GlobalCoordinate :: dimension ;
1230 const int nb = grid_.face_cells[ 2*face ];
1231 if( nb != seed.index() )
1232 {
1233 normal *= -1.0;
1234 }
1235 return normal;
1236 }
1237
1238 template <class EntitySeed>
1240 unitOuterNormal( const EntitySeed& seed, const int i ) const
1241 {
1242 const int face = this->template subEntitySeed<1>( seed, i ).index();
1243 if( seed.index() == grid_.face_cells[ 2*face ] )
1244 {
1245 return unitOuterNormals_[ face ];
1246 }
1247 else
1248 {
1249 GlobalCoordinate normal = unitOuterNormals_[ face ];
1250 normal *= -1.0;
1251 return normal;
1252 }
1253 }
1254
1255 template <class EntitySeed>
1256 GlobalCoordinate centroids( const EntitySeed& seed ) const
1257 {
1258 if( ! seed.isValid() )
1259 return GlobalCoordinate( 0 );
1260
1261 const int index = GlobalCoordinate :: dimension * seed.index();
1262 const int codim = EntitySeed::codimension;
1263 assert( index >= 0 && index < size( codim ) * GlobalCoordinate :: dimension );
1264
1265 if( codim == 0 )
1266 {
1268 }
1269 else if ( codim == 1 )
1270 {
1272 }
1273 else if( codim == dim )
1274 {
1276 }
1277 else
1278 {
1279 DUNE_THROW(InvalidStateException,"codimension not implemented");
1280 return GlobalCoordinate( 0 );
1281 }
1282 }
1283
1284 GlobalCoordinate copyToGlobalCoordinate( const double* coords ) const
1285 {
1286 GlobalCoordinate coordinate;
1287 for( int i=0; i<GlobalCoordinate::dimension; ++i )
1288 {
1289 coordinate[ i ] = coords[ i ];
1290 }
1291 return coordinate;
1292 }
1293
1294 template <class EntitySeed>
1295 double volumes( const EntitySeed& seed ) const
1296 {
1297 static const int codim = EntitySeed::codimension;
1298 if( codim == dim || ! seed.isValid() )
1299 {
1300 return 1.0;
1301 }
1302 else
1303 {
1304 assert( seed.isValid() );
1305
1306 if( codim == 0 )
1307 {
1308 return grid_.cell_volumes[ seed.index() ];
1309 }
1310 else if ( codim == 1 )
1311 {
1312 return grid_.face_areas[ seed.index() ];
1313 }
1314 else
1315 {
1316 DUNE_THROW(InvalidStateException,"codimension not implemented");
1317 return 0.0;
1318 }
1319 }
1320 }
1321
1322 protected:
1323 void init ()
1324 {
1325 // copy Cartesian dimensions
1326 for( int i=0; i<3; ++i )
1327 {
1328 cartDims_[ i ] = grid_.cartdims[ i ];
1329 }
1330
1331 // setup list of cell vertices
1332 const int numCells = size( 0 );
1333
1334 cellVertices_.resize( numCells );
1335
1336 // sort vertices such that they comply with the dune reference cube
1337 if( grid_.cell_facetag )
1338 {
1339 typedef std::array<int, 3> KeyType;
1340 std::map< const KeyType, const int > vertexFaceTags;
1341 const int vertexFacePattern [8][3] = {
1342 { 0, 2, 4 }, // vertex 0
1343 { 1, 2, 4 }, // vertex 1
1344 { 0, 3, 4 }, // vertex 2
1345 { 1, 3, 4 }, // vertex 3
1346 { 0, 2, 5 }, // vertex 4
1347 { 1, 2, 5 }, // vertex 5
1348 { 0, 3, 5 }, // vertex 6
1349 { 1, 3, 5 } // vertex 7
1350 };
1351
1352 for( int i=0; i<8; ++i )
1353 {
1354 KeyType key; key.fill( 4 ); // default is 4 which is the first z coord (for the 2d case)
1355 for( int j=0; j<dim; ++j )
1356 {
1357 key[ j ] = vertexFacePattern[ i ][ j ];
1358 }
1359
1360 vertexFaceTags.insert( std::make_pair( key, i ) );
1361 }
1362
1363 for (int c = 0; c < numCells; ++c)
1364 {
1365 if( dim == 2 )
1366 {
1367 // for 2d Cartesian grids the face ordering is wrong
1368 int f = grid_.cell_facepos[ c ];
1369 std::swap( grid_.cell_faces[ f+1 ], grid_.cell_faces[ f+2 ] );
1370 std::swap( grid_.cell_facetag[ f+1 ], grid_.cell_facetag[ f+2 ] );
1371 }
1372
1373 typedef std::map<int,int> vertexmap_t;
1374 typedef typename vertexmap_t :: iterator iterator;
1375
1376 std::vector< vertexmap_t > cell_pts( dim*2 );
1377
1378 for (int hf=grid_.cell_facepos[ c ]; hf < grid_.cell_facepos[c+1]; ++hf)
1379 {
1380 const int f = grid_.cell_faces[ hf ];
1381 const int faceTag = grid_.cell_facetag[ hf ];
1382
1383 for( int nodepos=grid_.face_nodepos[f]; nodepos<grid_.face_nodepos[f+1]; ++nodepos )
1384 {
1385 const int node = grid_.face_nodes[ nodepos ];
1386 iterator it = cell_pts[ faceTag ].find( node );
1387 if( it == cell_pts[ faceTag ].end() )
1388 {
1389 cell_pts[ faceTag ].insert( std::make_pair( node, 1 ) );
1390 }
1391 else
1392 {
1393 // increase vertex reference counter
1394 (*it).second++;
1395 }
1396 }
1397 }
1398
1399 typedef std::map< int, std::set<int> > vertexlist_t;
1400 vertexlist_t vertexList;
1401
1402 for( int faceTag = 0; faceTag<dim*2; ++faceTag )
1403 {
1404 for( iterator it = cell_pts[ faceTag ].begin(),
1405 end = cell_pts[ faceTag ].end(); it != end; ++it )
1406 {
1407
1408 // only consider vertices with one appearance
1409 if( (*it).second == 1 )
1410 {
1411 vertexList[ (*it).first ].insert( faceTag );
1412 }
1413 }
1414 }
1415
1416 assert( int(vertexList.size()) == ( dim == 2 ? 4 : 8) );
1417
1418 cellVertices_[ c ].resize( vertexList.size() );
1419 for( auto it = vertexList.begin(), end = vertexList.end(); it != end; ++it )
1420 {
1421 assert( (*it).second.size() == dim );
1422 KeyType key; key.fill( 4 ); // fill with 4 which is the first z coord
1423
1424 std::copy( (*it).second.begin(), (*it).second.end(), key.begin() );
1425 auto vx = vertexFaceTags.find( key );
1426 assert( vx != vertexFaceTags.end() );
1427 if( vx != vertexFaceTags.end() )
1428 {
1429 if( (*vx).second >= int(cellVertices_[ c ].size()) )
1430 cellVertices_[ c ].resize( (*vx).second+1 );
1431 // store node number on correct local position
1432 cellVertices_[ c ][ (*vx).second ] = (*it).first ;
1433 }
1434 }
1435 }
1436
1437 // if face_tag is available we assume that the elements follow a cube-like structure
1438 geomTypes_.resize(dim + 1);
1439 GeometryType tmp;
1440 for (int codim = 0; codim <= dim; ++codim)
1441 {
1442 tmp = Dune::GeometryTypes::cube(dim - codim);
1443 geomTypes_[codim].push_back(tmp);
1444 }
1445 }
1446 else // if ( grid_.cell_facetag )
1447 {
1448 int maxVx = 0 ;
1449 int minVx = std::numeric_limits<int>::max();
1450
1451 for (int c = 0; c < numCells; ++c)
1452 {
1453 std::set<int> cell_pts;
1454 for (int hf=grid_.cell_facepos[ c ]; hf < grid_.cell_facepos[c+1]; ++hf)
1455 {
1456 int f = grid_.cell_faces[ hf ];
1457 const int* fnbeg = grid_.face_nodes + grid_.face_nodepos[f];
1458 const int* fnend = grid_.face_nodes + grid_.face_nodepos[f+1];
1459 cell_pts.insert(fnbeg, fnend);
1460 }
1461
1462 cellVertices_[ c ].resize( cell_pts.size() );
1463 std::copy(cell_pts.begin(), cell_pts.end(), cellVertices_[ c ].begin() );
1464 maxVx = std::max( maxVx, int( cell_pts.size() ) );
1465 minVx = std::min( minVx, int( cell_pts.size() ) );
1466 }
1467
1468 if( minVx == maxVx && maxVx == 4 )
1469 {
1470 for (int c = 0; c < numCells; ++c)
1471 {
1472 assert( cellVertices_[ c ].size() == 4 );
1473 GlobalCoordinate center( 0 );
1474 GlobalCoordinate p[ dim+1 ];
1475 for( int i=0; i<dim+1; ++i )
1476 {
1477 const int vertex = cellVertices_[ c ][ i ];
1478
1479 for( int d=0; d<dim; ++d )
1480 {
1481 center[ d ] += grid_.node_coordinates[ vertex*dim + d ];
1482 p[ i ][ d ] = grid_.node_coordinates[ vertex*dim + d ];
1483 }
1484 }
1485 center *= 0.25;
1486 for( int d=0; d<dim; ++d )
1487 {
1488 grid_.cell_centroids[ c*dim + d ] = center[ d ];
1489 }
1490
1491 Dune::GeometryType simplex;
1492 simplex = Dune::GeometryTypes::simplex(dim);
1493
1494 typedef Dune::AffineGeometry< ctype, dim, dimworld> AffineGeometryType;
1495 AffineGeometryType geometry( simplex, p );
1496 grid_.cell_volumes[ c ] = geometry.volume();
1497 }
1498 }
1499
1500 // check face normals
1501 {
1502 const int faces = grid_.number_of_faces;
1503 for( int face = 0 ; face < faces; ++face )
1504 {
1505 const int a = grid_.face_cells[ 2*face ];
1506 const int b = grid_.face_cells[ 2*face + 1 ];
1507
1508 assert( a >=0 || b >=0 );
1509
1510 if( grid_.face_areas[ face ] < 0 )
1511 std::abort();
1512
1513 GlobalCoordinate centerDiff( 0 );
1514 if( b >= 0 )
1515 {
1516 for( int d=0; d<dimworld; ++d )
1517 {
1518 centerDiff[ d ] = grid_.cell_centroids[ b*dimworld + d ];
1519 }
1520 }
1521 else
1522 {
1523 for( int d=0; d<dimworld; ++d )
1524 {
1525 centerDiff[ d ] = grid_.face_centroids[ face*dimworld + d ];
1526 }
1527 }
1528
1529 if( a >= 0 )
1530 {
1531 for( int d=0; d<dimworld; ++d )
1532 {
1533 centerDiff[ d ] -= grid_.cell_centroids[ a*dimworld + d ];
1534 }
1535 }
1536 else
1537 {
1538 for( int d=0; d<dimworld; ++d )
1539 {
1540 centerDiff[ d ] -= grid_.face_centroids[ face*dimworld + d ];
1541 }
1542 }
1543
1544 GlobalCoordinate normal( 0 );
1545 for( int d=0; d<dimworld; ++d )
1546 {
1547 normal[ d ] = grid_.face_normals[ face*dimworld + d ];
1548 }
1549
1550 if( centerDiff.two_norm() < 1e-10 )
1551 std::abort();
1552
1553 // if diff and normal point in different direction, flip faces
1554 if( centerDiff * normal < 0 )
1555 {
1556 grid_.face_cells[ 2*face ] = b;
1557 grid_.face_cells[ 2*face + 1 ] = a;
1558 }
1559 }
1560 }
1561
1562 bool allSimplex = true ;
1563 bool allCube = true ;
1564
1565 for (int c = 0; c < numCells; ++c)
1566 {
1567 const int nVx = cellVertices_[ c ].size();
1568 if( nVx != 4 )
1569 {
1570 allSimplex = false;
1571 }
1572
1573 if( nVx != 8 )
1574 {
1575 allCube = false;
1576 }
1577 }
1578 // Propogate the cell geometry type to all codimensions
1579 geomTypes_.resize(dim + 1);
1580 GeometryType tmp;
1581 for (int codim = 0; codim <= dim; ++codim)
1582 {
1583 if( allSimplex )
1584 {
1585 tmp = Dune::GeometryTypes::simplex(dim - codim);
1586 geomTypes_[ codim ].push_back( tmp );
1587 }
1588 else if ( allCube )
1589 {
1590 tmp = Dune::GeometryTypes::cube(dim - codim);
1591 geomTypes_[ codim ].push_back( tmp );
1592 }
1593 else
1594 {
1595 tmp = Dune::GeometryTypes::none(dim - codim);
1596 geomTypes_[ codim ].push_back( tmp );
1597 }
1598 }
1599
1600 } // end else of ( grid_.cell_facetag )
1601
1602 nBndSegments_ = 0;
1604 for( int face = 0; face < grid_.number_of_faces; ++face )
1605 {
1606 const int normalIdx = face * GlobalCoordinate :: dimension ;
1608 normal /= normal.two_norm();
1609 unitOuterNormals_[ face ] = normal;
1610
1611 if( isBoundaryFace( face ) )
1612 {
1613 // increase number if boundary segments
1614 ++nBndSegments_;
1615 const int facePos = 2 * face ;
1616 // store negative number to indicate boundary
1617 // the abstract value is the segment index
1618 if( grid_.face_cells[ facePos ] < 0 )
1619 {
1620 grid_.face_cells[ facePos ] = -nBndSegments_;
1621 }
1622 else if ( grid_.face_cells[ facePos+1 ] < 0 )
1623 {
1624 grid_.face_cells[ facePos+1 ] = -nBndSegments_;
1625 }
1626 }
1627 }
1628 }
1629
1630 void print( std::ostream& out, const UnstructuredGridType& grid ) const
1631 {
1632 const int numCells = grid.number_of_cells;
1633 for( int c=0; c<numCells; ++c )
1634 {
1635 out << "cell " << c << " : faces = " << std::endl;
1636 for (int hf=grid.cell_facepos[ c ]; hf < grid.cell_facepos[c+1]; ++hf)
1637 {
1638 int f = grid_.cell_faces[ hf ];
1639 const int* fnbeg = grid_.face_nodes + grid_.face_nodepos[f];
1640 const int* fnend = grid_.face_nodes + grid_.face_nodepos[f+1];
1641 out << f << " vx = " ;
1642 while( fnbeg != fnend )
1643 {
1644 out << *fnbeg << " ";
1645 ++fnbeg;
1646 }
1647 out << std::endl;
1648 }
1649 out << std::endl;
1650
1651 const auto& vx = cellVertices_[ c ];
1652 out << "cell " << c << " : vertices = ";
1653 for( size_t i=0; i<vx.size(); ++i )
1654 out << vx[ i ] << " ";
1655 out << std::endl;
1656 }
1657
1658 }
1659
1660 protected:
1663
1665 std::array< int, 3 > cartDims_;
1666 std::vector< std::vector< GeometryType > > geomTypes_;
1667 std::vector< std::vector< int > > cellVertices_;
1668
1669 std::vector< GlobalCoordinate > unitOuterNormals_;
1670
1674
1676
1677 private:
1678 // no copying
1679 PolyhedralGrid ( const PolyhedralGrid& );
1680 };
1681
1682
1683
1684 // PolyhedralGrid::Codim
1685 // -------------
1686
1687 template< int dim, int dimworld, typename coord_t >
1688 template< int codim >
1690 : public Base::template Codim< codim >
1691 {
1700 typedef typename Traits::template Codim< codim >::Entity Entity;
1701
1706 typedef typename Traits::template Codim< codim >::EntityPointer EntityPointer;
1707
1721 typedef typename Traits::template Codim< codim >::Geometry Geometry;
1722
1731 typedef typename Traits::template Codim< codim >::LocalGeometry LocalGeometry;
1732
1738 template< PartitionIteratorType pitype >
1740 {
1741 typedef typename Traits::template Codim< codim >
1744 typedef typename Traits::template Codim< codim >
1747 };
1748
1757
1766
1768 };
1769
1770} // namespace Dune
1771
1775
1776#endif // #ifndef DUNE_POLYHEDRALGRID_GRID_HH
void destroy_grid(struct UnstructuredGrid *g)
struct UnstructuredGrid * allocate_grid(size_t ndims, size_t ncells, size_t nfaces, size_t nfacenodes, size_t ncellfaces, size_t nnodes)
struct UnstructuredGrid * create_grid_cart2d(int nx, int ny, double dx, double dy)
struct UnstructuredGrid * create_grid_hexa3d(int nx, int ny, int nz, double dx, double dy, double dz)
Definition: entity.hh:152
Definition: entityseed.hh:16
Definition: geometry.hh:245
Definition: idset.hh:18
Definition: indexset.hh:24
Definition: intersectioniterator.hh:16
Definition: intersection.hh:22
Definition: iterator.hh:21
Definition: geometry.hh:264
identical grid wrapper
Definition: grid.hh:159
void postAdapt()
Definition: grid.hh:604
Traits::ctype ctype
type of vector coordinates (e.g., double)
Definition: grid.hh:308
int ghostSize(int, int codim) const
obtain size of ghost region for a grid level
Definition: grid.hh:644
bool loadBalance(CommDataHandleIF< DataHandle, Data > &)
rebalance the load each process has to handle
Definition: grid.hh:748
Partition< All_Partition >::LevelGridView LevelGridView
View types for All_Partition.
Definition: grid.hh:245
Traits::template Codim< EntitySeed::codimension >::EntityPointer entityPointer(const EntitySeed &seed) const
obtain EntityPointer from EntitySeed.
Definition: grid.hh:808
Traits::template Codim< EntitySeed::codimension >::Entity entity(const EntitySeed &seed) const
obtain EntityPointer from EntitySeed.
Definition: grid.hh:819
static void computeGeometry(UnstructuredGridPtr &ug)
Definition: grid.hh:191
bool preAdapt()
Definition: grid.hh:582
std::vector< std::vector< int > > cellVertices_
Definition: grid.hh:1667
GlobalCoordinate outerNormal(const EntitySeed &seed, const int i) const
Definition: grid.hh:1225
int ghostSize(int codim) const
obtain size of ghost region for the leaf grid
Definition: grid.hh:624
bool loadBalance(DofManager &)
rebalance the load each process has to handle
Definition: grid.hh:768
Traits::GlobalIdSet GlobalIdSet
type of global id set
Definition: grid.hh:283
int maxLevel() const
obtain maximal grid level
Definition: grid.hh:414
Codim< codim >::template Partition< pitype >::LeafIterator leafend() const
Definition: grid.hh:515
Codim< codim >::LevelIterator lbegin(const int) const
Definition: grid.hh:522
std::array< int, 3 > cartDims_
Definition: grid.hh:1665
const LeafIndexSet & leafIndexSet() const
Definition: grid.hh:562
const int * globalCellPtr() const
Definition: grid.hh:855
Traits::LevelIndexSet LevelIndexSet
type of level index set
Definition: grid.hh:271
LevelGridView levelGridView(int) const
View for a grid level for All_Partition.
Definition: grid.hh:792
std::unique_ptr< UnstructuredGridType, UnstructuredGridDeleter > UnstructuredGridPtr
Definition: grid.hh:178
PolyhedralGrid(const std::vector< int > &n, const std::vector< double > &dx)
constructor
Definition: grid.hh:346
bool isBoundaryFace(const typename Codim< 1 >::EntitySeed &faceSeed) const
Definition: grid.hh:1114
static UnstructuredGridPtr allocateGrid(std::size_t nCells, std::size_t nFaces, std::size_t nFaceNodes, std::size_t nCellFaces, std::size_t nNodes)
Definition: grid.hh:181
size_t nBndSegments_
Definition: grid.hh:1675
Codim< codim >::LeafIterator leafbegin() const
Definition: grid.hh:494
Codim< codim >::LevelIterator lend(const int) const
Definition: grid.hh:528
ExtraData extraData() const
Definition: grid.hh:969
Partition< pitype >::LeafGridView leafGridView() const
View for the leaf grid.
Definition: grid.hh:784
void print(std::ostream &out, const UnstructuredGridType &grid) const
Definition: grid.hh:1630
int getMark(const typename Codim< 0 >::Entity &) const
Definition: grid.hh:576
Traits::HierarchicIterator HierarchicIterator
iterator over the grid hierarchy
Definition: grid.hh:225
Codim< 0 >::EntitySeed neighbor(const typename Codim< 0 >::EntitySeed &seed, const int i) const
Definition: grid.hh:1184
CommunicationType comm_
Definition: grid.hh:1664
Codim< codim >::LeafIterator leafend() const
Definition: grid.hh:500
int size(int, int codim) const
obtain number of entites on a level
Definition: grid.hh:427
Partition< pitype >::LevelGridView levelGridView(int) const
View for a grid level.
Definition: grid.hh:775
GlobalCoordinate copyToGlobalCoordinate(const double *coords) const
Definition: grid.hh:1284
void update()
update grid caches
Definition: grid.hh:838
void getIJK(const int c, std::array< int, 3 > &ijk) const
Definition: grid.hh:860
Partition< All_Partition >::LeafGridView LeafGridView
Definition: grid.hh:246
void switchToDistributedView()
Switch to the distributed view.
Definition: grid.hh:698
int overlapSize(int, int) const
obtain size of overlap region for a grid level
Definition: grid.hh:634
int indexInOutside(const typename Codim< 0 >::EntitySeed &seed, const int i) const
Definition: grid.hh:1198
std::vector< std::vector< GeometryType > > geomTypes_
Definition: grid.hh:1666
int size(int codim) const
obtain number of leaf entities
Definition: grid.hh:438
void scatterData(DataHandle &handle) const
Moves data from the global (all data on process) view to the distributed view.
Definition: grid.hh:883
LeafIndexSet leafIndexSet_
Definition: grid.hh:1671
int indexInInside(const typename Codim< 0 >::EntitySeed &seed, const int i) const
Definition: grid.hh:1172
Codim< codim >::EntitySeed subEntitySeed(const EntitySeedArg &baseSeed, const int i) const
Definition: grid.hh:1043
Traits::LocalIdSet LocalIdSet
type of local id set
Definition: grid.hh:300
int corners(const EntitySeed &seed) const
Definition: grid.hh:972
const GlobalIdSet & globalIdSet() const
Definition: grid.hh:547
typename Traits::Communication Communication
communicator with all other processes having some part of the grid
Definition: grid.hh:311
std::vector< GlobalCoordinate > unitOuterNormals_
Definition: grid.hh:1669
Codim< codim >::template Partition< pitype >::LeafIterator leafbegin() const
Definition: grid.hh:507
int subEntities(const EntitySeed &seed, const int codim) const
Definition: grid.hh:1014
UnstructuredGridType * createGrid(const std::vector< int > &n, const std::vector< double > &dx) const
Definition: grid.hh:944
const CommunicationType & comm() const
obtain CollectiveCommunication object
Definition: grid.hh:711
int overlapSize(int) const
obtain size of overlap region for the leaf grid
Definition: grid.hh:615
const std::vector< GeometryType > & geomTypes(const unsigned int codim) const
Definition: grid.hh:1131
Traits::LevelIntersectionIterator LevelIntersectionIterator
iterator over intersections with other entities on the same level
Definition: grid.hh:229
Communication CommunicationType
Definition: grid.hh:312
GridFamily::Traits Traits
type of the grid traits
Definition: grid.hh:208
PolyhedralGrid(UnstructuredGridPtr &&gridPtr)
constructor
Definition: grid.hh:365
Codim< codim >::EntitySeed subEntitySeed(const typename Codim< 1 >::EntitySeed &faceSeed, const int i) const
Definition: grid.hh:1077
Traits::LeafIntersectionIterator LeafIntersectionIterator
iterator over intersections with other entities on the leaf level
Definition: grid.hh:227
bool adapt(DataHandle &)
Definition: grid.hh:598
Codim< codim >::template Partition< pitype >::LevelIterator lend(const int) const
Definition: grid.hh:542
UnstructuredGrid UnstructuredGridType
Definition: grid.hh:166
bool hasBoundaryIntersections(const typename Codim< 0 >::EntitySeed &seed) const
Definition: grid.hh:1095
double volumes(const EntitySeed &seed) const
Definition: grid.hh:1295
int size(int, GeometryType type) const
obtain number of entites on a level
Definition: grid.hh:467
void communicate(DataHandle &, InterfaceType, CommunicationDirection) const
communicate information on leaf entities
Definition: grid.hh:684
Codim< codim >::template Partition< pitype >::LevelIterator lbegin(const int) const
Definition: grid.hh:535
UnstructuredGridPtr gridPtr_
Definition: grid.hh:1661
const LocalIdSet & localIdSet() const
Definition: grid.hh:552
GlobalIdSet globalIdSet_
Definition: grid.hh:1672
size_t numBoundarySegments() const
obtain number of leaf entities
Definition: grid.hh:487
void switchToGlobalView()
Switch to the global view.
Definition: grid.hh:692
GeometryType geometryType(const Seed &seed) const
Definition: grid.hh:1143
const std::array< int, 3 > & logicalCartesianSize() const
Definition: grid.hh:844
bool loadBalance()
rebalance the load each process has to handle
Definition: grid.hh:727
PolyhedralGrid(const UnstructuredGridType &grid)
constructor
Definition: grid.hh:384
LeafGridView leafGridView() const
View for the leaf grid for All_Partition.
Definition: grid.hh:799
int boundarySegmentIndex(const typename Codim< 0 >::EntitySeed &seed, const int face) const
Definition: grid.hh:1120
bool isBoundaryFace(const int face) const
Definition: grid.hh:1107
int size(GeometryType type) const
returns the number of boundary segments within the macro grid
Definition: grid.hh:476
int cartesianIndexInInside(const typename Codim< 0 >::EntitySeed &seed, const int i) const
Definition: grid.hh:1177
bool mark(int, const typename Codim< 0 >::Entity &)
Definition: grid.hh:571
Traits::GlobalCoordinate GlobalCoordinate
Definition: grid.hh:314
const LevelIndexSet & levelIndexSet(int) const
Definition: grid.hh:557
void communicate(DataHandle &, InterfaceType, CommunicationDirection, int) const
communicate information on a grid level
Definition: grid.hh:663
bool adapt()
Definition: grid.hh:588
Traits::ExtraData ExtraData
Definition: grid.hh:968
const int * globalCell() const
Definition: grid.hh:849
const UnstructuredGridType & grid_
Definition: grid.hh:1662
GlobalCoordinate unitOuterNormal(const EntitySeed &seed, const int i) const
Definition: grid.hh:1240
GlobalCoordinate centroids(const EntitySeed &seed) const
Definition: grid.hh:1256
LocalIdSet localIdSet_
Definition: grid.hh:1673
Traits::LeafIndexSet LeafIndexSet
type of leaf index set
Definition: grid.hh:261
void globalRefine(int)
Definition: grid.hh:567
void init()
Definition: grid.hh:1323
GlobalCoordinate corner(const EntitySeed &seed, const int i) const
Definition: grid.hh:987
Transform a corner-point grid ZCORN field to account for MINPV processing.
Definition: MinpvProcessor.hpp:35
struct UnstructuredGrid * create_grid_cornerpoint(const struct grdecl *in, double tol)
void compute_geometry(struct UnstructuredGrid *g)
The namespace Dune is the main namespace for all Dune code.
Definition: common/CartesianIndexMapper.hpp:10
int faceTag(const Dune::CpGrid &grid, const Dune::cpgrid::Cell2FacesRow::iterator &cell_face)
Get Eclipse Cartesian tag of a face.
int numCells(const Dune::CpGrid &grid)
Get the number of cells of a grid.
STL namespace.
Traits::template Codim< codim >::template Partition< pitype >::LeafIterator LeafIterator
Definition: grid.hh:1743
Traits::template Codim< codim >::template Partition< pitype >::LevelIterator LevelIterator
Definition: grid.hh:1746
traits structure containing types for a codimension
Definition: grid.hh:1691
Partition< All_Partition >::LeafIterator LeafIterator
type of level iterator
Definition: grid.hh:1756
Traits::template Codim< codim >::Entity Entity
type of entity
Definition: grid.hh:1700
Traits::template Codim< codim >::EntityPointer EntityPointer
type of entity pointer
Definition: grid.hh:1706
Traits::template Codim< codim >::LocalGeometry LocalGeometry
type of local geometry
Definition: grid.hh:1731
Partition< All_Partition >::LevelIterator LevelIterator
type of leaf iterator
Definition: grid.hh:1765
Traits::template Codim< codim >::Geometry Geometry
type of world geometry
Definition: grid.hh:1721
Types for GridView.
Definition: grid.hh:239
GridFamily::Traits::template Partition< pitype >::LeafGridView LeafGridView
Definition: grid.hh:241
GridFamily::Traits::template Partition< pitype >::LevelGridView LevelGridView
Definition: grid.hh:240
void operator()(UnstructuredGridType *grdPtr)
Definition: grid.hh:171
PolyhedralGridIterator< codim, const Grid, pitype > LeafIteratorImpl
Definition: grid.hh:108
Dune::EntityIterator< codim, const Grid, LeafIteratorImpl > LeafIterator
Definition: grid.hh:109
LeafIterator LevelIterator
Definition: grid.hh:111
Dune::Geometry< dimension-codim, dimensionworld, const Grid, PolyhedralGridLocalGeometry > LocalGeometry
Definition: grid.hh:94
Entity EntityPointer
Definition: grid.hh:100
PolyhedralGridGeometry< dimension-codim, dimensionworld, const Grid > GeometryImpl
Definition: grid.hh:90
Dune::Geometry< dimension-codim, dimensionworld, const Grid, PolyhedralGridGeometry > Geometry
Definition: grid.hh:91
Partition< All_Partition >::LeafIterator LeafIterator
Definition: grid.hh:114
EntityImpl EntityPointerImpl
Definition: grid.hh:99
PolyhedralGridEntitySeed< codim, const Grid > EntitySeed
Definition: grid.hh:103
Dune::Entity< codim, dimension, const Grid, PolyhedralGridEntity > Entity
Definition: grid.hh:97
PolyhedralGridLocalGeometry< dimension-codim, dimensionworld, const Grid > LocalGeometryImpl
Definition: grid.hh:93
Partition< All_Partition >::LevelIterator LevelIterator
Definition: grid.hh:115
PolyhedralGridEntity< codim, dimension, const Grid > EntityImpl
Definition: grid.hh:96
Dune::GridView< PolyhedralGridViewTraits< dim, dimworld, ctype, pitype > > LeafGridView
Definition: grid.hh:131
Dune::GridView< PolyhedralGridViewTraits< dim, dimworld, ctype, pitype > > LevelGridView
Definition: grid.hh:132
PolyhedralGridIterator< 0, const Grid, All_Partition > HierarchicIteratorImpl
Definition: grid.hh:84
Dune::EntityIterator< 0, const Grid, HierarchicIteratorImpl > HierarchicIterator
Definition: grid.hh:85
Dune::IntersectionIterator< const Grid, LevelIntersectionIteratorImpl, LevelIntersectionImpl > LevelIntersectionIterator
Definition: grid.hh:82
GlobalIdSet LocalIdSet
Definition: grid.hh:122
PolyhedralGridIndexSet< dim, dimworld, ctype > LeafIndexSet
Definition: grid.hh:118
Dune::Communication< MPICommunicator > Communication
Definition: grid.hh:125
Dune::IntersectionIterator< const Grid, LeafIntersectionIteratorImpl, LeafIntersectionImpl > LeafIntersectionIterator
Definition: grid.hh:81
PolyhedralGridIntersection< const Grid > LevelIntersectionImpl
Definition: grid.hh:74
PolyhedralGridIntersectionIterator< const Grid > LeafIntersectionIteratorImpl
Definition: grid.hh:75
Dune::Intersection< const Grid, LeafIntersectionImpl > LeafIntersection
Definition: grid.hh:78
PolyhedralGridIntersectionIterator< const Grid > LevelIntersectionIteratorImpl
Definition: grid.hh:76
static const int dimension
Definition: grid.hh:68
const Grid * ExtraData
Definition: grid.hh:64
int Index
Definition: grid.hh:66
Partition< All_Partition >::LevelGridView LevelGridView
Definition: grid.hh:136
Partition< All_Partition >::LeafGridView LeafGridView
Definition: grid.hh:135
PolyhedralGrid< dim, dimworld, coord_t > Grid
Definition: grid.hh:58
PolyhedralGridIdSet< dim, dimworld, ctype > GlobalIdSet
Definition: grid.hh:121
Dune::FieldVector< ctype, dimensionworld > GlobalCoordinate
Definition: grid.hh:71
Dune::Communication< MPICommunicator > CollectiveCommunication
Definition: grid.hh:126
static const int dimensionworld
Definition: grid.hh:69
PolyhedralGridIntersection< const Grid > LeafIntersectionImpl
Definition: grid.hh:73
Dune::MPIHelper::MPICommunicator MPICommunicator
Definition: grid.hh:124
coord_t ctype
Definition: grid.hh:60
Dune::Intersection< const Grid, LevelIntersectionImpl > LevelIntersection
Definition: grid.hh:79
PolyhedralGridIndexSet< dim, dimworld, ctype > LevelIndexSet
Definition: grid.hh:119
Definition: grid.hh:55
Definition: UnstructuredGrid.h:99
int * face_nodes
Definition: UnstructuredGrid.h:121
int * face_nodepos
Definition: UnstructuredGrid.h:127
int number_of_cells
Definition: UnstructuredGrid.h:109
double * face_areas
Definition: UnstructuredGrid.h:173
int * cell_faces
Definition: UnstructuredGrid.h:146
int number_of_faces
Definition: UnstructuredGrid.h:111
double * cell_centroids
Definition: UnstructuredGrid.h:192
int * cell_facepos
Definition: UnstructuredGrid.h:152
int * cell_facetag
Definition: UnstructuredGrid.h:244
double * cell_volumes
Definition: UnstructuredGrid.h:197
int * face_cells
Definition: UnstructuredGrid.h:138
double * face_centroids
Definition: UnstructuredGrid.h:168
int number_of_nodes
Definition: UnstructuredGrid.h:113
int cartdims[3]
Definition: UnstructuredGrid.h:227
int * global_cell
Definition: UnstructuredGrid.h:214
double * face_normals
Definition: UnstructuredGrid.h:184
double * node_coordinates
Definition: UnstructuredGrid.h:160
Definition: preprocess.h:56
const double * coord
Definition: preprocess.h:58
int dims[3]
Definition: preprocess.h:57
const double * zcorn
Definition: preprocess.h:59
const int * actnum
Definition: preprocess.h:60
double b(const Dune::FieldVector< ct, dimworld > &x, double t)
Definition: transportproblem2.hh:25