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 const bool edge_conformal = false)
330 : gridPtr_ { createGrid(inputGrid, poreVolumes, static_cast<int>(edge_conformal)) }
331 , grid_ { *gridPtr_ }
332 , comm_ { MPIHelper::getCommunicator() }
333 , leafIndexSet_ { *this }
334 , globalIdSet_ { *this }
335 , localIdSet_ { *this }
336 , nBndSegments_ { 0 }
337 {
338 init();
339 }
340#endif
341
347 explicit PolyhedralGrid ( const std::vector< int >& n,
348 const std::vector< double >& dx )
349 : gridPtr_( createGrid( n, dx ) ),
350 grid_( *gridPtr_ ),
351 comm_( MPIHelper::getCommunicator()),
352 leafIndexSet_( *this ),
353 globalIdSet_( *this ),
354 localIdSet_( *this ),
355 nBndSegments_( 0 )
356 {
357 init();
358 }
359
366 explicit PolyhedralGrid ( UnstructuredGridPtr &&gridPtr )
367 : gridPtr_( std::move( gridPtr ) ),
368 grid_( *gridPtr_ ),
369 comm_( MPIHelper::getCommunicator() ),
370 leafIndexSet_( *this ),
371 globalIdSet_( *this ),
372 localIdSet_( *this ),
373 nBndSegments_( 0 )
374 {
375 init();
376 }
377
385 explicit PolyhedralGrid ( const UnstructuredGridType& grid )
386 : gridPtr_(),
387 grid_( grid ),
388 comm_( MPIHelper::getCommunicator() ),
389 leafIndexSet_( *this ),
390 globalIdSet_( *this ),
391 localIdSet_( *this ),
392 nBndSegments_( 0 )
393 {
394 init();
395 }
396
401 operator const UnstructuredGridType& () const { return grid_; }
402
415 int maxLevel () const
416 {
417 return 0;
418 }
419
428 int size ( int /* level */, int codim ) const
429 {
430 return size( codim );
431 }
432
439 int size ( int codim ) const
440 {
441 if( codim == 0 )
442 {
443 return grid_.number_of_cells;
444 }
445 else if ( codim == 1 )
446 {
447 return grid_.number_of_faces;
448 }
449 else if ( codim == dim )
450 {
451 return grid_.number_of_nodes;
452 }
453 else
454 {
455 std::cerr << "Warning: codimension " << codim << " not available in PolyhedralGrid" << std::endl;
456 return 0;
457 }
458 }
459
468 int size ( int /* level */, GeometryType type ) const
469 {
470 return size( dim - type.dim() );
471 }
472
477 int size ( GeometryType type ) const
478 {
479 return size( dim - type.dim() );
480 }
481
488 size_t numBoundarySegments () const
489 {
490 return nBndSegments_;
491 }
494 template< int codim >
496 {
497 return leafbegin< codim, All_Partition >();
498 }
499
500 template< int codim >
502 {
503 return leafend< codim, All_Partition >();
504 }
505
506 template< int codim, PartitionIteratorType pitype >
507 typename Codim< codim >::template Partition< pitype >::LeafIterator
508 leafbegin () const
509 {
510 typedef typename Traits::template Codim< codim >::template Partition< pitype >::LeafIteratorImpl Impl;
511 return Impl( extraData(), true );
512 }
513
514 template< int codim, PartitionIteratorType pitype >
515 typename Codim< codim >::template Partition< pitype >::LeafIterator
516 leafend () const
517 {
518 typedef typename Traits::template Codim< codim >::template Partition< pitype >::LeafIteratorImpl Impl;
519 return Impl( extraData(), false );
520 }
521
522 template< int codim >
523 typename Codim< codim >::LevelIterator lbegin ( const int /* level */ ) const
524 {
525 return leafbegin< codim, All_Partition >();
526 }
527
528 template< int codim >
529 typename Codim< codim >::LevelIterator lend ( const int /* level */ ) const
530 {
531 return leafend< codim, All_Partition >();
532 }
533
534 template< int codim, PartitionIteratorType pitype >
535 typename Codim< codim >::template Partition< pitype >::LevelIterator
536 lbegin ( const int /* level */ ) const
537 {
538 return leafbegin< codim, pitype > ();
539 }
540
541 template< int codim, PartitionIteratorType pitype >
542 typename Codim< codim >::template Partition< pitype >::LevelIterator
543 lend ( const int /* level */ ) const
544 {
545 return leafend< codim, pitype > ();
546 }
547
548 const GlobalIdSet &globalIdSet () const
549 {
550 return globalIdSet_;
551 }
552
553 const LocalIdSet &localIdSet () const
554 {
555 return localIdSet_;
556 }
557
558 const LevelIndexSet &levelIndexSet ( int /* level */ ) const
559 {
560 return leafIndexSet();
561 }
562
564 {
565 return leafIndexSet_;
566 }
567
568 void globalRefine ( int /* refCount */ )
569 {
570 }
571
572 bool mark ( int /* refCount */, const typename Codim< 0 >::Entity& /* entity */ )
573 {
574 return false;
575 }
576
577 int getMark ( const typename Codim< 0 >::Entity& /* entity */) const
578 {
579 return false;
580 }
581
583 bool preAdapt ()
584 {
585 return false;
586 }
587
589 bool adapt ()
590 {
591 return false ;
592 }
593
598 template< class DataHandle >
599 bool adapt ( DataHandle & )
600 {
601 return false;
602 }
603
605 void postAdapt ()
606 {
607 }
608
616 int overlapSize ( int /* codim */) const
617 {
618 return 0;
619 }
620
625 int ghostSize( int codim ) const
626 {
627 return (codim == 0 ) ? 1 : 0;
628 }
629
635 int overlapSize ( int /* level */, int /* codim */ ) const
636 {
637 return 0;
638 }
639
645 int ghostSize ( int /* level */, int codim ) const
646 {
647 return ghostSize( codim );
648 }
649
663 template< class DataHandle>
664 void communicate ( DataHandle& /* dataHandle */,
665 InterfaceType /* interface */,
666 CommunicationDirection /* direction */,
667 int /* level */ ) const
668 {
669 OPM_THROW(std::runtime_error, "communicate not implemented for polyhedreal grid!");
670 }
671
684 template< class DataHandle>
685 void communicate ( DataHandle& /* dataHandle */,
686 InterfaceType /* interface */,
687 CommunicationDirection /* direction */ ) const
688 {
689 OPM_THROW(std::runtime_error, "communicate not implemented for polyhedreal grid!");
690 }
691
694 {
695 OPM_THROW(std::runtime_error, "switch to global view not implemented for polyhedreal grid!");
696 }
697
700 {
701 OPM_THROW(std::runtime_error, "switch to distributed view not implemented for polyhedreal grid!");
702 }
703
712 const CommunicationType &comm () const
713 {
714 return comm_;
715 }
716
717 // data handle interface different between geo and interface
718
729 {
730 return false ;
731 }
732
748 template< class DataHandle, class Data >
749 bool loadBalance ( CommDataHandleIF< DataHandle, Data >& /* datahandle */ )
750 {
751 return false;
752 }
753
768 template< class DofManager >
769 bool loadBalance ( DofManager& /* dofManager */ )
770 {
771 return false;
772 }
773
775 template< PartitionIteratorType pitype >
776 typename Partition< pitype >::LevelGridView levelGridView ( int /* level */ ) const
777 {
778 typedef typename Partition< pitype >::LevelGridView View;
779 typedef typename View::GridViewImp ViewImp;
780 return View( ViewImp( *this ) );
781 }
782
784 template< PartitionIteratorType pitype >
786 {
787 typedef typename Traits::template Partition< pitype >::LeafGridView View;
788 typedef typename View::GridViewImp ViewImp;
789 return View( ViewImp( *this ) );
790 }
791
793 LevelGridView levelGridView ( int /* level */ ) const
794 {
795 typedef typename LevelGridView::GridViewImp ViewImp;
796 return LevelGridView( ViewImp( *this ) );
797 }
798
801 {
802 typedef typename LeafGridView::GridViewImp ViewImp;
803 return LeafGridView( ViewImp( *this ) );
804 }
805
807 template< class EntitySeed >
809 entityPointer ( const EntitySeed &seed ) const
810 {
811 typedef typename Traits::template Codim< EntitySeed::codimension >::EntityPointer EntityPointer;
812 typedef typename Traits::template Codim< EntitySeed::codimension >::EntityPointerImpl EntityPointerImpl;
813 typedef typename Traits::template Codim< EntitySeed::codimension >::EntityImpl EntityImpl;
814 return EntityPointer( EntityPointerImpl( EntityImpl( extraData(), seed ) ) );
815 }
816
818 template< class EntitySeed >
819 typename Traits::template Codim< EntitySeed::codimension >::Entity
820 entity ( const EntitySeed &seed ) const
821 {
822 typedef typename Traits::template Codim< EntitySeed::codimension >::EntityImpl EntityImpl;
823 return EntityImpl( extraData(), seed );
824 }
825
839 void update ()
840 {
841 }
842
845 const std::array<int, 3>& logicalCartesianSize() const
846 {
847 return cartDims_;
848 }
849
850 const int* globalCell() const
851 {
852 assert( grid_.global_cell != 0 );
853 return grid_.global_cell;
854 }
855
856 const int* globalCellPtr() const
857 {
858 return grid_.global_cell;
859 }
860
861 void getIJK(const int c, std::array<int,3>& ijk) const
862 {
863 int gc = globalCell()[c];
864 ijk[0] = gc % logicalCartesianSize()[0]; gc /= logicalCartesianSize()[0];
865 ijk[1] = gc % logicalCartesianSize()[1];
866 ijk[2] = gc / logicalCartesianSize()[1];
867 }
868
873
874
883 template<class DataHandle>
884 void scatterData([[maybe_unused]] DataHandle& handle) const
885 {
886 OPM_THROW(std::runtime_error, "ScatterData not implemented for polyhedral grid!");
887 }
888
889 protected:
890#if HAVE_ECL_INPUT
892 createGrid(const Opm::EclipseGrid& inputGrid,
893 const std::vector<double>& poreVolumes,
894 const bool edge_conformal) const
895 {
896 grdecl g{};
897
898 g.dims[0] = inputGrid.getNX();
899 g.dims[1] = inputGrid.getNY();
900 g.dims[2] = inputGrid.getNZ();
901
902 std::vector<double> coord = inputGrid.getCOORD( );
903 std::vector<double> zcorn = inputGrid.getZCORN( );
904 std::vector<int> actnum = inputGrid.getACTNUM( );
905
906 g.coord = coord.data();
907 g.zcorn = zcorn.data();
908 g.actnum = actnum.data();
909
910 const double z_tolerance = inputGrid.isPinchActive() ? inputGrid.getPinchThresholdThickness() : 0.0;
911
912 if (!poreVolumes.empty() && (inputGrid.getMinpvMode() != Opm::MinpvMode::Inactive))
913 {
914 Opm::MinpvProcessor mp(g.dims[0], g.dims[1], g.dims[2]);
915 const std::vector<double>& minpvv = inputGrid.getMinpvVector();
916 // Currently the pinchProcessor is not used and only opmfil is supported
917 // The polyhedralgrid only only supports the opmfil option
918 //bool opmfil = inputGrid.getMinpvMode() == Opm::MinpvMode::OpmFIL;
919 bool opmfil = true;
920 const size_t cartGridSize = g.dims[0] * g.dims[1] * g.dims[2];
921 std::vector<double> thickness(cartGridSize);
922 for (size_t i = 0; i < cartGridSize; ++i) {
923 thickness[i] = inputGrid.getCellThickness(i);
924 }
925 mp.process(thickness, z_tolerance, inputGrid.getPinchMaxEmptyGap() , poreVolumes,
926 minpvv, actnum, opmfil, zcorn.data());
927 }
928
929 /*
930 if (!poreVolumes.empty() && (inputGrid.getMinpvMode() != Opm::MinpvMode::Inactive)) {
931 Opm::MinpvProcessor mp(g.dims[0], g.dims[1], g.dims[2]);
932 const double minpv_value = inputGrid.getMinpvValue();
933 // Currently the pinchProcessor is not used and only opmfil is supported
934 //bool opmfil = inputGrid.getMinpvMode() == Opm::MinpvMode::OpmFIL;
935 bool opmfil = true;
936 mp.process(poreVolumes, minpv_value, actnum, opmfil, zcorn.data());
937 }
938 */
939
941 (&g, z_tolerance, static_cast<int>(edge_conformal));
942
943 if (cgrid == nullptr) {
944 OPM_THROW(std::runtime_error, "Failed to construct grid.");
945 }
946
947 return cgrid;
948 }
949#endif
950
951 UnstructuredGridType* createGrid( const std::vector< int >& n, const std::vector< double >& dx ) const
952 {
953 UnstructuredGridType* cgrid = nullptr ;
954 assert( int(n.size()) == dim );
955 if( dim == 2 )
956 {
957 cgrid = create_grid_cart2d( n[ 0 ], n[ 1 ], dx[ 0 ], dx[ 1 ] );
958 }
959 else if ( dim == 3 )
960 {
961 cgrid = create_grid_hexa3d( n[ 0 ], n[ 1 ], n[ 2 ], dx[ 0 ], dx[ 1 ], dx[ 2 ] );
962 }
963
964 //print_grid( cgrid );
965 if (!cgrid) {
966 OPM_THROW(std::runtime_error, "Failed to construct grid.");
967 }
968 return cgrid;
969 }
970
971 public:
972 typedef typename Traits :: ExtraData ExtraData;
973 ExtraData extraData () const { return this; }
974
975 template <class EntitySeed>
976 int corners( const EntitySeed& seed ) const
977 {
978 const int codim = EntitySeed :: codimension;
979 const int index = seed.index();
980 if (codim==0)
981 return cellVertices_[ index ].size();
982 if (codim==1)
983 return grid_.face_nodepos[ index+1 ] - grid_.face_nodepos[ index ];
984 if (codim==dim)
985 return 1;
986 return 0;
987 }
988
989 template <class EntitySeed>
991 corner ( const EntitySeed& seed, const int i ) const
992 {
993 const int codim = EntitySeed :: codimension;
994 if (codim==0)
995 {
996 const int coordIndex = GlobalCoordinate :: dimension * cellVertices_[ seed.index() ][ i ];
997 return copyToGlobalCoordinate( grid_.node_coordinates + coordIndex );
998 }
999 if (codim==1)
1000 {
1001 // for faces we need to swap vertices in 3d since in UnstructuredGrid
1002 // those are ordered counter clockwise, for 2d this does not matter
1003 // TODO: Improve this for performance reasons
1004 const int crners = corners( seed );
1005 const int crner = (crners == 4 && EntitySeed :: dimension == 3 && i > 1 ) ? 5 - i : i;
1006 const int faceVertex = grid_.face_nodes[ grid_.face_nodepos[seed.index() ] + crner ];
1007 return copyToGlobalCoordinate( grid_.node_coordinates + GlobalCoordinate :: dimension * faceVertex );
1008 }
1009 if (codim==dim)
1010 {
1011 const int coordIndex = GlobalCoordinate :: dimension * seed.index();
1012 return copyToGlobalCoordinate( grid_.node_coordinates + coordIndex );
1013 }
1014 return GlobalCoordinate( 0 );
1015 }
1016
1017 template <class EntitySeed>
1018 int subEntities( const EntitySeed& seed, const int codim ) const
1019 {
1020 const int index = seed.index();
1021 if( seed.codimension == 0 )
1022 {
1023 if (codim==0)
1024 return 1;
1025 if (codim==1)
1026 return grid_.cell_facepos[ index+1 ] - grid_.cell_facepos[ index ];
1027 if (codim==dim)
1028 return cellVertices_[ index ].size();
1029 }
1030 else if( seed.codimension == 1 )
1031 {
1032 if (codim==1)
1033 return 1;
1034 if (codim==dim)
1035 return grid_.face_nodepos[ index+1 ] - grid_.face_nodepos[ index ];
1036 }
1037 else if ( seed.codimension == dim )
1038 {
1039 return 1 ;
1040 }
1041
1042 return 0;
1043 }
1044
1045 template <int codim, class EntitySeedArg >
1046 typename Codim<codim>::EntitySeed
1047 subEntitySeed( const EntitySeedArg& baseSeed, const int i ) const
1048 {
1049 assert( codim >= EntitySeedArg::codimension );
1050 assert( i>= 0 && i<subEntities( baseSeed, codim ) );
1051 typedef typename Codim<codim>::EntitySeed EntitySeed;
1052
1053 // if codim equals entity seed codim just return same entity seed.
1054 if( codim == EntitySeedArg::codimension )
1055 {
1056 return EntitySeed( baseSeed.index() );
1057 }
1058
1059 if( EntitySeedArg::codimension == 0 )
1060 {
1061 if ( codim == 1 )
1062 {
1063 return EntitySeed( grid_.cell_faces[ grid_.cell_facepos[ baseSeed.index() ] + i ] );
1064 }
1065 else if ( codim == dim )
1066 {
1067 return EntitySeed( cellVertices_[ baseSeed.index() ][ i ] );
1068 }
1069 }
1070 else if ( EntitySeedArg::codimension == 1 && codim == dim )
1071 {
1072 return EntitySeed( grid_.face_nodes[ grid_.face_nodepos[ baseSeed.index() + i ] ]);
1073 }
1074
1075 DUNE_THROW(NotImplemented,"codimension not available");
1076 return EntitySeed();
1077 }
1078
1079 template <int codim>
1080 typename Codim<codim>::EntitySeed
1081 subEntitySeed( const typename Codim<1>::EntitySeed& faceSeed, const int i ) const
1082 {
1083 assert( i>= 0 && i<subEntities( faceSeed, codim ) );
1084 typedef typename Codim<codim>::EntitySeed EntitySeed;
1085 if ( codim == 1 )
1086 {
1087 return EntitySeed( faceSeed.index() );
1088 }
1089 else if ( codim == dim )
1090 {
1091 return EntitySeed( grid_.face_nodes[ grid_.face_nodepos[ faceSeed.index() ] + i ] );
1092 }
1093 else
1094 {
1095 DUNE_THROW(NotImplemented,"codimension not available");
1096 }
1097 }
1098
1099 bool hasBoundaryIntersections(const typename Codim<0>::EntitySeed& seed ) const
1100 {
1101 const int faces = subEntities( seed, 1 );
1102 for( int f=0; f<faces; ++f )
1103 {
1104 const auto faceSeed = this->template subEntitySeed<1>( seed, f );
1105 if( isBoundaryFace( faceSeed ) )
1106 return true;
1107 }
1108 return false;
1109 }
1110
1111 bool isBoundaryFace(const int face ) const
1112 {
1113 assert( face >= 0 && face < grid_.number_of_faces );
1114 const int facePos = 2 * face;
1115 return ((grid_.face_cells[ facePos ] < 0) || (grid_.face_cells[ facePos+1 ] < 0));
1116 }
1117
1118 bool isBoundaryFace(const typename Codim<1>::EntitySeed& faceSeed ) const
1119 {
1120 assert( faceSeed.isValid() );
1121 return isBoundaryFace( faceSeed.index() );
1122 }
1123
1124 int boundarySegmentIndex(const typename Codim<0>::EntitySeed& seed, const int face ) const
1125 {
1126 const auto faceSeed = this->template subEntitySeed<1>( seed, face );
1127 assert( faceSeed.isValid() );
1128 const int facePos = 2 * faceSeed.index();
1129 const int idx = std::min( grid_.face_cells[ facePos ], grid_.face_cells[ facePos+1 ]);
1130 // check that this is actually the boundary
1131 assert( idx < 0 );
1132 return -(idx+1); // +1 to include 0 boundary segment index
1133 }
1134
1135 const std::vector< GeometryType > &geomTypes ( const unsigned int codim ) const
1136 {
1137 static std::vector< GeometryType > emptyDummy;
1138 if (codim < geomTypes_.size())
1139 {
1140 return geomTypes_[codim];
1141 }
1142
1143 return emptyDummy;
1144 }
1145
1146 template < class Seed >
1147 GeometryType geometryType( const Seed& seed ) const
1148 {
1149 if( Seed::codimension == 0 )
1150 {
1151 assert(!geomTypes( Seed::codimension ).empty());
1152 return geomTypes( Seed::codimension )[ 0 ];
1153 }
1154 else
1155 {
1156 // 3d faces
1157 if( dim == 3 && Seed::codimension == 1 )
1158 {
1159 GeometryType face;
1160 const int nVx = corners( seed );
1161 if( nVx == 4 ) // quad face
1162 face = Dune::GeometryTypes::cube(2);
1163 else if( nVx == 3 ) // triangle face
1164 face = Dune::GeometryTypes::simplex(2);
1165 else // polygonal face
1166 face = Dune::GeometryTypes::none(2);
1167 return face;
1168 }
1169
1170 // for codim 2 and codim 3 there is only one option
1171 assert(!geomTypes( Seed::codimension ).empty());
1172 return geomTypes( Seed::codimension )[ 0 ];
1173 }
1174 }
1175
1176 int indexInInside( const typename Codim<0>::EntitySeed& seed, const int i ) const
1177 {
1178 return ( grid_.cell_facetag ) ? cartesianIndexInInside( seed, i ) : i;
1179 }
1180
1181 int cartesianIndexInInside( const typename Codim<0>::EntitySeed& seed, const int i ) const
1182 {
1183 assert( i>= 0 && i<subEntities( seed, 1 ) );
1184 return grid_.cell_facetag[ grid_.cell_facepos[ seed.index() ] + i ] ;
1185 }
1186
1187 typename Codim<0>::EntitySeed
1188 neighbor( const typename Codim<0>::EntitySeed& seed, const int i ) const
1189 {
1190 const int face = this->template subEntitySeed<1>( seed, i ).index();
1191 int nb = grid_.face_cells[ 2 * face ];
1192 if( nb == seed.index() )
1193 {
1194 nb = grid_.face_cells[ 2 * face + 1 ];
1195 }
1196
1197 typedef typename Codim<0>::EntitySeed EntitySeed;
1198 return EntitySeed( nb );
1199 }
1200
1201 int
1202 indexInOutside( const typename Codim<0>::EntitySeed& seed, const int i ) const
1203 {
1204 if( grid_.cell_facetag )
1205 {
1206 // if cell_facetag is present we assume pseudo Cartesian corner point case
1207 const int in_inside = cartesianIndexInInside( seed, i );
1208 return in_inside + ((in_inside % 2) ? -1 : 1);
1209 }
1210 else
1211 {
1212 typedef typename Codim<0>::EntitySeed EntitySeed;
1213 EntitySeed nb = neighbor( seed, i );
1214 const int faces = subEntities( seed, 1 );
1215 for( int face = 0; face<faces; ++ face )
1216 {
1217 if( neighbor( nb, face ).equals(seed) )
1218 {
1219 return indexInInside( nb, face );
1220 }
1221 }
1222 DUNE_THROW(InvalidStateException,"inverse intersection not found");
1223 return -1;
1224 }
1225 }
1226
1227 template <class EntitySeed>
1229 outerNormal( const EntitySeed& seed, const int i ) const
1230 {
1231 const int face = this->template subEntitySeed<1>( seed, i ).index();
1232 const int normalIdx = face * GlobalCoordinate :: dimension ;
1234 const int nb = grid_.face_cells[ 2*face ];
1235 if( nb != seed.index() )
1236 {
1237 normal *= -1.0;
1238 }
1239 return normal;
1240 }
1241
1242 template <class EntitySeed>
1244 unitOuterNormal( const EntitySeed& seed, const int i ) const
1245 {
1246 const int face = this->template subEntitySeed<1>( seed, i ).index();
1247 if( seed.index() == grid_.face_cells[ 2*face ] )
1248 {
1249 return unitOuterNormals_[ face ];
1250 }
1251 else
1252 {
1253 GlobalCoordinate normal = unitOuterNormals_[ face ];
1254 normal *= -1.0;
1255 return normal;
1256 }
1257 }
1258
1259 template <class EntitySeed>
1260 GlobalCoordinate centroids( const EntitySeed& seed ) const
1261 {
1262 if( ! seed.isValid() )
1263 return GlobalCoordinate( 0 );
1264
1265 const int index = GlobalCoordinate :: dimension * seed.index();
1266 const int codim = EntitySeed::codimension;
1267 assert( index >= 0 && index < size( codim ) * GlobalCoordinate :: dimension );
1268
1269 if( codim == 0 )
1270 {
1272 }
1273 else if ( codim == 1 )
1274 {
1276 }
1277 else if( codim == dim )
1278 {
1280 }
1281 else
1282 {
1283 DUNE_THROW(InvalidStateException,"codimension not implemented");
1284 return GlobalCoordinate( 0 );
1285 }
1286 }
1287
1288 GlobalCoordinate copyToGlobalCoordinate( const double* coords ) const
1289 {
1290 GlobalCoordinate coordinate;
1291 for( int i=0; i<GlobalCoordinate::dimension; ++i )
1292 {
1293 coordinate[ i ] = coords[ i ];
1294 }
1295 return coordinate;
1296 }
1297
1298 template <class EntitySeed>
1299 double volumes( const EntitySeed& seed ) const
1300 {
1301 static const int codim = EntitySeed::codimension;
1302 if( codim == dim || ! seed.isValid() )
1303 {
1304 return 1.0;
1305 }
1306 else
1307 {
1308 assert( seed.isValid() );
1309
1310 if( codim == 0 )
1311 {
1312 return grid_.cell_volumes[ seed.index() ];
1313 }
1314 else if ( codim == 1 )
1315 {
1316 return grid_.face_areas[ seed.index() ];
1317 }
1318 else
1319 {
1320 DUNE_THROW(InvalidStateException,"codimension not implemented");
1321 return 0.0;
1322 }
1323 }
1324 }
1325
1326 protected:
1327 void init ()
1328 {
1329 // copy Cartesian dimensions
1330 for( int i=0; i<3; ++i )
1331 {
1332 cartDims_[ i ] = grid_.cartdims[ i ];
1333 }
1334
1335 // setup list of cell vertices
1336 const int numCells = size( 0 );
1337
1338 cellVertices_.resize( numCells );
1339
1340 // sort vertices such that they comply with the dune reference cube
1341 if( grid_.cell_facetag )
1342 {
1343 typedef std::array<int, 3> KeyType;
1344 std::map< const KeyType, const int > vertexFaceTags;
1345 const int vertexFacePattern [8][3] = {
1346 { 0, 2, 4 }, // vertex 0
1347 { 1, 2, 4 }, // vertex 1
1348 { 0, 3, 4 }, // vertex 2
1349 { 1, 3, 4 }, // vertex 3
1350 { 0, 2, 5 }, // vertex 4
1351 { 1, 2, 5 }, // vertex 5
1352 { 0, 3, 5 }, // vertex 6
1353 { 1, 3, 5 } // vertex 7
1354 };
1355
1356 for( int i=0; i<8; ++i )
1357 {
1358 KeyType key; key.fill( 4 ); // default is 4 which is the first z coord (for the 2d case)
1359 for( int j=0; j<dim; ++j )
1360 {
1361 key[ j ] = vertexFacePattern[ i ][ j ];
1362 }
1363
1364 vertexFaceTags.insert( std::make_pair( key, i ) );
1365 }
1366
1367 for (int c = 0; c < numCells; ++c)
1368 {
1369 if( dim == 2 )
1370 {
1371 // for 2d Cartesian grids the face ordering is wrong
1372 int f = grid_.cell_facepos[ c ];
1373 std::swap( grid_.cell_faces[ f+1 ], grid_.cell_faces[ f+2 ] );
1374 std::swap( grid_.cell_facetag[ f+1 ], grid_.cell_facetag[ f+2 ] );
1375 }
1376
1377 typedef std::map<int,int> vertexmap_t;
1378 typedef typename vertexmap_t :: iterator iterator;
1379
1380 std::vector< vertexmap_t > cell_pts( dim*2 );
1381
1382 for (unsigned hf = grid_.cell_facepos[ c ]; hf < grid_.cell_facepos[c+1]; ++hf)
1383 {
1384 const int f = grid_.cell_faces[ hf ];
1385 const int faceTag = grid_.cell_facetag[ hf ];
1386
1387 for (unsigned nodepos = grid_.face_nodepos[f]; nodepos < grid_.face_nodepos[f+1]; ++nodepos )
1388 {
1389 const int node = grid_.face_nodes[ nodepos ];
1390 iterator it = cell_pts[ faceTag ].find( node );
1391 if( it == cell_pts[ faceTag ].end() )
1392 {
1393 cell_pts[ faceTag ].insert( std::make_pair( node, 1 ) );
1394 }
1395 else
1396 {
1397 // increase vertex reference counter
1398 (*it).second++;
1399 }
1400 }
1401 }
1402
1403 typedef std::map< int, std::set<int> > vertexlist_t;
1404 vertexlist_t vertexList;
1405
1406 for( int faceTag = 0; faceTag<dim*2; ++faceTag )
1407 {
1408 for( iterator it = cell_pts[ faceTag ].begin(),
1409 end = cell_pts[ faceTag ].end(); it != end; ++it )
1410 {
1411
1412 // only consider vertices with one appearance
1413 if( (*it).second == 1 )
1414 {
1415 vertexList[ (*it).first ].insert( faceTag );
1416 }
1417 }
1418 }
1419
1420 assert( int(vertexList.size()) == ( dim == 2 ? 4 : 8) );
1421
1422 cellVertices_[ c ].resize( vertexList.size() );
1423 for( auto it = vertexList.begin(), end = vertexList.end(); it != end; ++it )
1424 {
1425 assert( (*it).second.size() == dim );
1426 KeyType key; key.fill( 4 ); // fill with 4 which is the first z coord
1427
1428 std::copy( (*it).second.begin(), (*it).second.end(), key.begin() );
1429 auto vx = vertexFaceTags.find( key );
1430 assert( vx != vertexFaceTags.end() );
1431 if( vx != vertexFaceTags.end() )
1432 {
1433 if( (*vx).second >= int(cellVertices_[ c ].size()) )
1434 cellVertices_[ c ].resize( (*vx).second+1 );
1435 // store node number on correct local position
1436 cellVertices_[ c ][ (*vx).second ] = (*it).first ;
1437 }
1438 }
1439 }
1440
1441 // if face_tag is available we assume that the elements follow a cube-like structure
1442 geomTypes_.resize(dim + 1);
1443 GeometryType tmp;
1444 for (int codim = 0; codim <= dim; ++codim)
1445 {
1446 tmp = Dune::GeometryTypes::cube(dim - codim);
1447 geomTypes_[codim].push_back(tmp);
1448 }
1449 }
1450 else // if ( grid_.cell_facetag )
1451 {
1452 int maxVx = 0 ;
1453 int minVx = std::numeric_limits<int>::max();
1454
1455 for (int c = 0; c < numCells; ++c)
1456 {
1457 std::set<int> cell_pts;
1458 for (unsigned hf = grid_.cell_facepos[ c ]; hf < grid_.cell_facepos[c+1]; ++hf)
1459 {
1460 int f = grid_.cell_faces[ hf ];
1461 const int* fnbeg = grid_.face_nodes + grid_.face_nodepos[f];
1462 const int* fnend = grid_.face_nodes + grid_.face_nodepos[f+1];
1463 cell_pts.insert(fnbeg, fnend);
1464 }
1465
1466 cellVertices_[ c ].resize( cell_pts.size() );
1467 std::copy(cell_pts.begin(), cell_pts.end(), cellVertices_[ c ].begin() );
1468 maxVx = std::max( maxVx, int( cell_pts.size() ) );
1469 minVx = std::min( minVx, int( cell_pts.size() ) );
1470 }
1471
1472 if( minVx == maxVx && maxVx == 4 )
1473 {
1474 for (int c = 0; c < numCells; ++c)
1475 {
1476 assert( cellVertices_[ c ].size() == 4 );
1477 GlobalCoordinate center( 0 );
1478 GlobalCoordinate p[ dim+1 ];
1479 for( int i=0; i<dim+1; ++i )
1480 {
1481 const int vertex = cellVertices_[ c ][ i ];
1482
1483 for( int d=0; d<dim; ++d )
1484 {
1485 center[ d ] += grid_.node_coordinates[ vertex*dim + d ];
1486 p[ i ][ d ] = grid_.node_coordinates[ vertex*dim + d ];
1487 }
1488 }
1489 center *= 0.25;
1490 for( int d=0; d<dim; ++d )
1491 {
1492 grid_.cell_centroids[ c*dim + d ] = center[ d ];
1493 }
1494
1495 Dune::GeometryType simplex;
1496 simplex = Dune::GeometryTypes::simplex(dim);
1497
1498 typedef Dune::AffineGeometry< ctype, dim, dimworld> AffineGeometryType;
1499 AffineGeometryType geometry( simplex, p );
1500 grid_.cell_volumes[ c ] = geometry.volume();
1501 }
1502 }
1503
1504 // check face normals
1505 {
1506 const int faces = grid_.number_of_faces;
1507 for( int face = 0 ; face < faces; ++face )
1508 {
1509 const int a = grid_.face_cells[ 2*face ];
1510 const int b = grid_.face_cells[ 2*face + 1 ];
1511
1512 assert( a >=0 || b >=0 );
1513
1514 if( grid_.face_areas[ face ] < 0 )
1515 std::abort();
1516
1517 GlobalCoordinate centerDiff( 0 );
1518 if( b >= 0 )
1519 {
1520 for( int d=0; d<dimworld; ++d )
1521 {
1522 centerDiff[ d ] = grid_.cell_centroids[ b*dimworld + d ];
1523 }
1524 }
1525 else
1526 {
1527 for( int d=0; d<dimworld; ++d )
1528 {
1529 centerDiff[ d ] = grid_.face_centroids[ face*dimworld + d ];
1530 }
1531 }
1532
1533 if( a >= 0 )
1534 {
1535 for( int d=0; d<dimworld; ++d )
1536 {
1537 centerDiff[ d ] -= grid_.cell_centroids[ a*dimworld + d ];
1538 }
1539 }
1540 else
1541 {
1542 for( int d=0; d<dimworld; ++d )
1543 {
1544 centerDiff[ d ] -= grid_.face_centroids[ face*dimworld + d ];
1545 }
1546 }
1547
1548 GlobalCoordinate normal( 0 );
1549 for( int d=0; d<dimworld; ++d )
1550 {
1551 normal[ d ] = grid_.face_normals[ face*dimworld + d ];
1552 }
1553
1554 if( centerDiff.two_norm() < 1e-10 )
1555 std::abort();
1556
1557 // if diff and normal point in different direction, flip faces
1558 if( centerDiff * normal < 0 )
1559 {
1560 grid_.face_cells[ 2*face ] = b;
1561 grid_.face_cells[ 2*face + 1 ] = a;
1562 }
1563 }
1564 }
1565
1566 bool allSimplex = true ;
1567 bool allCube = true ;
1568
1569 for (int c = 0; c < numCells; ++c)
1570 {
1571 const int nVx = cellVertices_[ c ].size();
1572 if( nVx != 4 )
1573 {
1574 allSimplex = false;
1575 }
1576
1577 if( nVx != 8 )
1578 {
1579 allCube = false;
1580 }
1581 }
1582 // Propogate the cell geometry type to all codimensions
1583 geomTypes_.resize(dim + 1);
1584 GeometryType tmp;
1585 for (int codim = 0; codim <= dim; ++codim)
1586 {
1587 if( allSimplex )
1588 {
1589 tmp = Dune::GeometryTypes::simplex(dim - codim);
1590 geomTypes_[ codim ].push_back( tmp );
1591 }
1592 else if ( allCube )
1593 {
1594 tmp = Dune::GeometryTypes::cube(dim - codim);
1595 geomTypes_[ codim ].push_back( tmp );
1596 }
1597 else
1598 {
1599 tmp = Dune::GeometryTypes::none(dim - codim);
1600 geomTypes_[ codim ].push_back( tmp );
1601 }
1602 }
1603
1604 } // end else of ( grid_.cell_facetag )
1605
1606 nBndSegments_ = 0;
1608 for( int face = 0; face < grid_.number_of_faces; ++face )
1609 {
1610 const int normalIdx = face * GlobalCoordinate :: dimension ;
1612 normal /= normal.two_norm();
1613 unitOuterNormals_[ face ] = normal;
1614
1615 if( isBoundaryFace( face ) )
1616 {
1617 // increase number if boundary segments
1618 ++nBndSegments_;
1619 const int facePos = 2 * face ;
1620 // store negative number to indicate boundary
1621 // the abstract value is the segment index
1622 if( grid_.face_cells[ facePos ] < 0 )
1623 {
1624 grid_.face_cells[ facePos ] = -nBndSegments_;
1625 }
1626 else if ( grid_.face_cells[ facePos+1 ] < 0 )
1627 {
1628 grid_.face_cells[ facePos+1 ] = -nBndSegments_;
1629 }
1630 }
1631 }
1632 }
1633
1634 void print( std::ostream& out, const UnstructuredGridType& grid ) const
1635 {
1636 const int numCells = grid.number_of_cells;
1637 for( int c=0; c<numCells; ++c )
1638 {
1639 out << "cell " << c << " : faces = " << std::endl;
1640 for (int hf=grid.cell_facepos[ c ]; hf < grid.cell_facepos[c+1]; ++hf)
1641 {
1642 int f = grid_.cell_faces[ hf ];
1643 const int* fnbeg = grid_.face_nodes + grid_.face_nodepos[f];
1644 const int* fnend = grid_.face_nodes + grid_.face_nodepos[f+1];
1645 out << f << " vx = " ;
1646 while( fnbeg != fnend )
1647 {
1648 out << *fnbeg << " ";
1649 ++fnbeg;
1650 }
1651 out << std::endl;
1652 }
1653 out << std::endl;
1654
1655 const auto& vx = cellVertices_[ c ];
1656 out << "cell " << c << " : vertices = ";
1657 for( size_t i=0; i<vx.size(); ++i )
1658 out << vx[ i ] << " ";
1659 out << std::endl;
1660 }
1661
1662 }
1663
1664 protected:
1667
1669 std::array< int, 3 > cartDims_;
1670 std::vector< std::vector< GeometryType > > geomTypes_;
1671 std::vector< std::vector< int > > cellVertices_;
1672
1673 std::vector< GlobalCoordinate > unitOuterNormals_;
1674
1678
1680
1681 private:
1682 // no copying
1683 PolyhedralGrid ( const PolyhedralGrid& );
1684 };
1685
1686
1687
1688 // PolyhedralGrid::Codim
1689 // -------------
1690
1691 template< int dim, int dimworld, typename coord_t >
1692 template< int codim >
1694 : public Base::template Codim< codim >
1695 {
1704 typedef typename Traits::template Codim< codim >::Entity Entity;
1705
1710 typedef typename Traits::template Codim< codim >::EntityPointer EntityPointer;
1711
1725 typedef typename Traits::template Codim< codim >::Geometry Geometry;
1726
1735 typedef typename Traits::template Codim< codim >::LocalGeometry LocalGeometry;
1736
1742 template< PartitionIteratorType pitype >
1744 {
1745 typedef typename Traits::template Codim< codim >
1748 typedef typename Traits::template Codim< codim >
1751 };
1752
1761
1770
1772 };
1773
1774} // namespace Dune
1775
1779
1780#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:605
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:645
bool loadBalance(CommDataHandleIF< DataHandle, Data > &)
rebalance the load each process has to handle
Definition: grid.hh:749
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:809
Traits::template Codim< EntitySeed::codimension >::Entity entity(const EntitySeed &seed) const
obtain EntityPointer from EntitySeed.
Definition: grid.hh:820
static void computeGeometry(UnstructuredGridPtr &ug)
Definition: grid.hh:191
bool preAdapt()
Definition: grid.hh:583
std::vector< std::vector< int > > cellVertices_
Definition: grid.hh:1671
GlobalCoordinate outerNormal(const EntitySeed &seed, const int i) const
Definition: grid.hh:1229
int ghostSize(int codim) const
obtain size of ghost region for the leaf grid
Definition: grid.hh:625
bool loadBalance(DofManager &)
rebalance the load each process has to handle
Definition: grid.hh:769
Traits::GlobalIdSet GlobalIdSet
type of global id set
Definition: grid.hh:283
int maxLevel() const
obtain maximal grid level
Definition: grid.hh:415
Codim< codim >::template Partition< pitype >::LeafIterator leafend() const
Definition: grid.hh:516
Codim< codim >::LevelIterator lbegin(const int) const
Definition: grid.hh:523
std::array< int, 3 > cartDims_
Definition: grid.hh:1669
const LeafIndexSet & leafIndexSet() const
Definition: grid.hh:563
const int * globalCellPtr() const
Definition: grid.hh:856
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:793
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:347
bool isBoundaryFace(const typename Codim< 1 >::EntitySeed &faceSeed) const
Definition: grid.hh:1118
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:1679
Codim< codim >::LeafIterator leafbegin() const
Definition: grid.hh:495
Codim< codim >::LevelIterator lend(const int) const
Definition: grid.hh:529
ExtraData extraData() const
Definition: grid.hh:973
Partition< pitype >::LeafGridView leafGridView() const
View for the leaf grid.
Definition: grid.hh:785
void print(std::ostream &out, const UnstructuredGridType &grid) const
Definition: grid.hh:1634
int getMark(const typename Codim< 0 >::Entity &) const
Definition: grid.hh:577
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:1188
CommunicationType comm_
Definition: grid.hh:1668
Codim< codim >::LeafIterator leafend() const
Definition: grid.hh:501
int size(int, int codim) const
obtain number of entites on a level
Definition: grid.hh:428
Partition< pitype >::LevelGridView levelGridView(int) const
View for a grid level.
Definition: grid.hh:776
GlobalCoordinate copyToGlobalCoordinate(const double *coords) const
Definition: grid.hh:1288
void update()
update grid caches
Definition: grid.hh:839
void getIJK(const int c, std::array< int, 3 > &ijk) const
Definition: grid.hh:861
Partition< All_Partition >::LeafGridView LeafGridView
Definition: grid.hh:246
void switchToDistributedView()
Switch to the distributed view.
Definition: grid.hh:699
int overlapSize(int, int) const
obtain size of overlap region for a grid level
Definition: grid.hh:635
int indexInOutside(const typename Codim< 0 >::EntitySeed &seed, const int i) const
Definition: grid.hh:1202
std::vector< std::vector< GeometryType > > geomTypes_
Definition: grid.hh:1670
int size(int codim) const
obtain number of leaf entities
Definition: grid.hh:439
void scatterData(DataHandle &handle) const
Moves data from the global (all data on process) view to the distributed view.
Definition: grid.hh:884
LeafIndexSet leafIndexSet_
Definition: grid.hh:1675
int indexInInside(const typename Codim< 0 >::EntitySeed &seed, const int i) const
Definition: grid.hh:1176
Codim< codim >::EntitySeed subEntitySeed(const EntitySeedArg &baseSeed, const int i) const
Definition: grid.hh:1047
Traits::LocalIdSet LocalIdSet
type of local id set
Definition: grid.hh:300
int corners(const EntitySeed &seed) const
Definition: grid.hh:976
const GlobalIdSet & globalIdSet() const
Definition: grid.hh:548
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:1673
Codim< codim >::template Partition< pitype >::LeafIterator leafbegin() const
Definition: grid.hh:508
int subEntities(const EntitySeed &seed, const int codim) const
Definition: grid.hh:1018
UnstructuredGridType * createGrid(const std::vector< int > &n, const std::vector< double > &dx) const
Definition: grid.hh:951
const CommunicationType & comm() const
obtain CollectiveCommunication object
Definition: grid.hh:712
int overlapSize(int) const
obtain size of overlap region for the leaf grid
Definition: grid.hh:616
const std::vector< GeometryType > & geomTypes(const unsigned int codim) const
Definition: grid.hh:1135
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:366
Codim< codim >::EntitySeed subEntitySeed(const typename Codim< 1 >::EntitySeed &faceSeed, const int i) const
Definition: grid.hh:1081
Traits::LeafIntersectionIterator LeafIntersectionIterator
iterator over intersections with other entities on the leaf level
Definition: grid.hh:227
bool adapt(DataHandle &)
Definition: grid.hh:599
Codim< codim >::template Partition< pitype >::LevelIterator lend(const int) const
Definition: grid.hh:543
UnstructuredGrid UnstructuredGridType
Definition: grid.hh:166
bool hasBoundaryIntersections(const typename Codim< 0 >::EntitySeed &seed) const
Definition: grid.hh:1099
double volumes(const EntitySeed &seed) const
Definition: grid.hh:1299
int size(int, GeometryType type) const
obtain number of entites on a level
Definition: grid.hh:468
void communicate(DataHandle &, InterfaceType, CommunicationDirection) const
communicate information on leaf entities
Definition: grid.hh:685
Codim< codim >::template Partition< pitype >::LevelIterator lbegin(const int) const
Definition: grid.hh:536
UnstructuredGridPtr gridPtr_
Definition: grid.hh:1665
const LocalIdSet & localIdSet() const
Definition: grid.hh:553
GlobalIdSet globalIdSet_
Definition: grid.hh:1676
size_t numBoundarySegments() const
obtain number of leaf entities
Definition: grid.hh:488
void switchToGlobalView()
Switch to the global view.
Definition: grid.hh:693
GeometryType geometryType(const Seed &seed) const
Definition: grid.hh:1147
const std::array< int, 3 > & logicalCartesianSize() const
Definition: grid.hh:845
bool loadBalance()
rebalance the load each process has to handle
Definition: grid.hh:728
PolyhedralGrid(const UnstructuredGridType &grid)
constructor
Definition: grid.hh:385
LeafGridView leafGridView() const
View for the leaf grid for All_Partition.
Definition: grid.hh:800
int boundarySegmentIndex(const typename Codim< 0 >::EntitySeed &seed, const int face) const
Definition: grid.hh:1124
bool isBoundaryFace(const int face) const
Definition: grid.hh:1111
int size(GeometryType type) const
returns the number of boundary segments within the macro grid
Definition: grid.hh:477
int cartesianIndexInInside(const typename Codim< 0 >::EntitySeed &seed, const int i) const
Definition: grid.hh:1181
bool mark(int, const typename Codim< 0 >::Entity &)
Definition: grid.hh:572
Traits::GlobalCoordinate GlobalCoordinate
Definition: grid.hh:314
const LevelIndexSet & levelIndexSet(int) const
Definition: grid.hh:558
void communicate(DataHandle &, InterfaceType, CommunicationDirection, int) const
communicate information on a grid level
Definition: grid.hh:664
bool adapt()
Definition: grid.hh:589
Traits::ExtraData ExtraData
Definition: grid.hh:972
const int * globalCell() const
Definition: grid.hh:850
const UnstructuredGridType & grid_
Definition: grid.hh:1666
GlobalCoordinate unitOuterNormal(const EntitySeed &seed, const int i) const
Definition: grid.hh:1244
GlobalCoordinate centroids(const EntitySeed &seed) const
Definition: grid.hh:1260
LocalIdSet localIdSet_
Definition: grid.hh:1677
Traits::LeafIndexSet LeafIndexSet
type of leaf index set
Definition: grid.hh:261
void globalRefine(int)
Definition: grid.hh:568
void init()
Definition: grid.hh:1327
GlobalCoordinate corner(const EntitySeed &seed, const int i) const
Definition: grid.hh:991
Transform a corner-point grid ZCORN field to account for MINPV processing.
Definition: MinpvProcessor.hpp:35
void compute_geometry(struct UnstructuredGrid *g)
struct UnstructuredGrid * create_grid_cornerpoint(const struct grdecl *in, double tol, int edge_conformal)
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:1747
Traits::template Codim< codim >::template Partition< pitype >::LevelIterator LevelIterator
Definition: grid.hh:1750
traits structure containing types for a codimension
Definition: grid.hh:1695
Partition< All_Partition >::LeafIterator LeafIterator
type of level iterator
Definition: grid.hh:1760
Traits::template Codim< codim >::Entity Entity
type of entity
Definition: grid.hh:1704
Traits::template Codim< codim >::EntityPointer EntityPointer
type of entity pointer
Definition: grid.hh:1710
Traits::template Codim< codim >::LocalGeometry LocalGeometry
type of local geometry
Definition: grid.hh:1735
Partition< All_Partition >::LevelIterator LevelIterator
type of leaf iterator
Definition: grid.hh:1769
Traits::template Codim< codim >::Geometry Geometry
type of world geometry
Definition: grid.hh:1725
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:101
int * face_nodes
Definition: UnstructuredGrid.h:123
grid_size_t * face_nodepos
Definition: UnstructuredGrid.h:129
int number_of_cells
Definition: UnstructuredGrid.h:111
double * face_areas
Definition: UnstructuredGrid.h:175
int * cell_faces
Definition: UnstructuredGrid.h:148
int number_of_faces
Definition: UnstructuredGrid.h:113
double * cell_centroids
Definition: UnstructuredGrid.h:194
int * cell_facetag
Definition: UnstructuredGrid.h:246
double * cell_volumes
Definition: UnstructuredGrid.h:199
int * face_cells
Definition: UnstructuredGrid.h:140
grid_size_t * cell_facepos
Definition: UnstructuredGrid.h:154
double * face_centroids
Definition: UnstructuredGrid.h:170
int number_of_nodes
Definition: UnstructuredGrid.h:115
int cartdims[3]
Definition: UnstructuredGrid.h:229
int * global_cell
Definition: UnstructuredGrid.h:216
double * face_normals
Definition: UnstructuredGrid.h:186
double * node_coordinates
Definition: UnstructuredGrid.h:162
Definition: preprocess.h:56
int dims[3]
Definition: preprocess.h:57
double b(const Dune::FieldVector< ct, dimworld > &, double)
Definition: transportproblem2.hh:25