Go to the documentation of this file.
3 #ifndef DUNE_GRID_IO_FILE_DGFPARSER_DGFUG_HH
4 #define DUNE_GRID_IO_FILE_DGFPARSER_DGFUG_HH
13 #include <dune/common/exceptions.hh>
14 #include <dune/common/fvector.hh>
15 #include <dune/common/parallel/mpihelper.hh>
79 struct DGFGridFactory< UGGrid< dim > >
82 typedef UGGrid< dim >
Grid;
93 dgf_( rank( comm ), size( comm ) )
103 dgf_( rank( comm ), size( comm ) )
105 std::ifstream input( filename.c_str() );
107 DUNE_THROW( DGFException,
"Error: Macrofile " << filename <<
" not found" );
118 template<
class GG,
class II >
121 return factory_.wasInserted( intersection );
125 template<
class GG,
class II >
132 template<
int codim >
136 return dgf_.nofelparams;
138 return dgf_.nofvtxparams;
144 template<
class Entity >
147 return numParameters< Entity::codimension >();
151 std::vector< double > &
parameter (
const typename Grid::template Codim< 0 >::Entity &element )
153 if( numParameters< 0 >() <= 0 )
155 DUNE_THROW( InvalidStateException,
156 "Calling DGFGridFactory::parameter is only allowed if there are parameters." );
158 return dgf_.elParams[ factory_.insertionIndex( element ) ];
162 std::vector< double > &
parameter (
const typename Grid::template Codim< dimension >::Entity &
vertex )
164 if( numParameters< dimension >() <= 0 )
166 DUNE_THROW( InvalidStateException,
167 "Calling DGFGridFactory::parameter is only allowed if there are parameters." );
169 return dgf_.vtxParams[ factory_.insertionIndex(
vertex ) ];
175 return dgf_.haveBndParameters;
179 template<
class GG,
class II >
186 auto refElem = referenceElement< double, dimension >( entity.type() );
187 int corners = refElem.size( face, 1,
dimension );
188 std::vector< unsigned int > bound( corners );
189 for(
int i = 0; i < corners; ++i )
191 const int k = refElem.subEntity( face, 1, i,
dimension );
192 bound[ i ] = factory_.insertionIndex( entity.template subEntity< dimension >( k ) );
195 DuneGridFormatParser::facemap_t::key_type key( bound,
false );
196 const DuneGridFormatParser::facemap_t::const_iterator pos = dgf_.facemap.find( key );
197 if( pos != dgf_.facemap.end() )
198 return dgf_.facemap.find( key )->second.second;
205 void generate ( std::istream &input );
212 MPI_Comm_rank( MPICOMM, &rank );
222 MPI_Comm_size( MPICOMM, &size );
228 GridFactory< UGGrid< dim > > factory_;
229 DuneGridFormatParser dgf_;
231 #endif // #if ENABLE_UG
235 #endif // #ifndef DUNE_GRID_IO_FILE_DGFPARSER_DGFUG_HH
Entity inside() const
return Entity on the inside of this intersection. That is the Entity where we started this.
Definition: common/intersection.hh:260
int numParameters() const
Definition: dgfgridfactory.hh:113
Common Grid parameters.
Definition: gridparameter.hh:31
bool noClosure_
Definition: dgfug.hh:49
size_t boundarySegmentIndex() const
index of the boundary segment within the macro grid
Definition: common/intersection.hh:246
bool noClosure() const
returns true if no closure should be used for UGGrid
Definition: dgfug.hh:42
int indexInInside() const
Local index of codim 1 entity in the inside() entity where intersection is contained in.
Definition: common/intersection.hh:356
GridImp::template Codim< 0 >::Entity Entity
Type of entity that this Intersection belongs to.
Definition: common/intersection.hh:190
const DGFBoundaryParameter::type & boundaryParameter(const Intersection< GG, II > &intersection) const
Definition: dgfgridfactory.hh:163
bool noCopy() const
returns true if no copies are made for UGGrid elements
Definition: dgfug.hh:44
bool wasInserted(const Intersection &intersection) const
Definition: dgfgridfactory.hh:101
std::vector< double > & parameter(const Element &element)
Definition: dgfgridfactory.hh:129
static const type & defaultValue()
default constructor
Definition: parser.hh:26
static int refineStepsForHalf()
number of globalRefine steps needed to refuce h by 0.5
@ vertex
Definition: common.hh:179
G Grid
Definition: dgfgridfactory.hh:37
int boundaryId(const Intersection &intersection) const
Definition: dgfgridfactory.hh:107
MPIHelper::MPICommunicator MPICommunicatorType
Definition: dgfgridfactory.hh:39
bool noCopy_
Definition: dgfug.hh:50
std::string type
type of additional boundary parameters
Definition: parser.hh:23
UGGridParameterBlock(std::istream &input)
constructor taking istream
Definition: dgfug.cc:18
bool haveBoundaryParameters() const
Definition: dgfgridfactory.hh:156
Grid * grid()
Definition: dgfgridfactory.hh:95
static double refineWeight()
size_t heapSize_
Definition: dgfug.hh:51
Some simple static information for a given GridType.
Definition: io/file/dgfparser/dgfparser.hh:50
Front-end for the grid manager of the finite element toolbox UG.
Definition: uggrid.hh:230
Intersection of a mesh entity of codimension 0 ("element") with a "neighboring" element or with the d...
Definition: common/grid.hh:344
DGFGridFactory(std::istream &input, MPICommunicatorType comm=MPIHelper::getCommunicator())
Definition: dgfgridfactory.hh:47
size_t heapSize() const
returns heap size used on construction of the grid
Definition: dgfug.hh:46
Include standard header files.
Definition: agrid.hh:58
const static int dimension
Definition: dgfgridfactory.hh:38