5#ifndef DUNE_GRID_IO_FILE_DGFPARSER_DGFUG_HH
6#define DUNE_GRID_IO_FILE_DGFPARSER_DGFUG_HH
15#include <dune/common/exceptions.hh>
16#include <dune/common/fvector.hh>
17#include <dune/common/parallel/mpihelper.hh>
20#include <dune-grid-config.hh>
96 dgf_( rank( comm ), size( comm ) )
106 dgf_( rank( comm ), size( comm ) )
108 std::ifstream input( filename.c_str() );
110 DUNE_THROW(
DGFException,
"Error: Macrofile " << filename <<
" not found" );
121 template<
class GG,
class II >
124 return factory_.wasInserted( intersection );
128 template<
class GG,
class II >
135 template<
int codim >
139 return dgf_.nofelparams;
141 return dgf_.nofvtxparams;
147 template<
class Entity >
150 return numParameters< Entity::codimension >();
156 if( numParameters< 0 >() <= 0 )
158 DUNE_THROW( InvalidStateException,
159 "Calling DGFGridFactory::parameter is only allowed if there are parameters." );
161 return dgf_.elParams[ factory_.insertionIndex( element ) ];
167 if( numParameters< dimension >() <= 0 )
169 DUNE_THROW( InvalidStateException,
170 "Calling DGFGridFactory::parameter is only allowed if there are parameters." );
172 return dgf_.vtxParams[ factory_.insertionIndex( vertex ) ];
178 return dgf_.haveBndParameters;
182 template<
class GG,
class II >
189 auto refElem = referenceElement< double, dimension >( entity.type() );
190 int corners = refElem.size( face, 1,
dimension );
191 std::vector< unsigned int > bound( corners );
192 for(
int i = 0; i < corners; ++i )
194 const int k = refElem.subEntity( face, 1, i,
dimension );
195 bound[ i ] = factory_.insertionIndex( entity.template subEntity< dimension >( k ) );
198 DuneGridFormatParser::facemap_t::key_type key( bound,
false );
199 const DuneGridFormatParser::facemap_t::const_iterator pos = dgf_.facemap.find( key );
200 if( pos != dgf_.facemap.end() )
201 return dgf_.facemap.find( key )->second.second;
208 void generate ( std::istream &input );
215 MPI_Comm_rank( MPICOMM, &rank );
225 MPI_Comm_size( MPICOMM, &size );
231 GridFactory< UGGrid< dim > > factory_;
232 DuneGridFormatParser dgf_;
Include standard header files.
Definition agrid.hh:60
Definition dgfgridfactory.hh:38
MPIHelper::MPICommunicator MPICommunicatorType
Definition dgfgridfactory.hh:41
G Grid
Definition dgfgridfactory.hh:39
static const int dimension
Definition dgfgridfactory.hh:40
Intersection of a mesh entity of codimension 0 ("element") with a "neighboring" element or with the d...
Definition common/intersection.hh:164
int indexInInside() const
Local index of codim 1 entity in the inside() entity where intersection is contained in.
Definition common/intersection.hh:346
size_t boundarySegmentIndex() const
index of the boundary segment within the macro grid
Definition common/intersection.hh:236
Entity inside() const
return Entity on the inside of this intersection. That is the Entity where we started this.
Definition common/intersection.hh:250
GridImp::template Codim< 0 >::Entity Entity
Type of entity that this Intersection belongs to.
Definition common/intersection.hh:192
Wrapper class for entities.
Definition common/entity.hh:66
A Traits struct that collects all associated types of one implementation.
Definition common/grid.hh:411
Common Grid parameters.
Definition gridparameter.hh:35
exception class for IO errors in the DGF parser
Definition dgfexception.hh:16
Some simple static information for a given GridType.
Definition io/file/dgfparser/dgfparser.hh:56
size_t heapSize() const
returns heap size used on construction of the grid
Definition dgfug.hh:49
bool noCopy_
Definition dgfug.hh:53
size_t heapSize_
Definition dgfug.hh:54
bool noClosure_
Definition dgfug.hh:52
bool noCopy() const
returns true if no copies are made for UGGrid elements
Definition dgfug.hh:47
bool noClosure() const
returns true if no closure should be used for UGGrid
Definition dgfug.hh:45
static double refineWeight()
Definition dgfug.hh:70
static int refineStepsForHalf()
Definition dgfug.hh:65
const DGFBoundaryParameter::type & boundaryParameter(const Dune::Intersection< GG, II > &intersection) const
return invalid value
Definition dgfug.hh:183
std::vector< double > & parameter(const typename Grid::template Codim< dimension >::Entity &vertex)
return parameter for vertex
Definition dgfug.hh:165
DGFGridFactory(const std::string &filename, MPICommunicatorType comm=MPIHelper::getCommunicator())
constructor taking filename
Definition dgfug.hh:102
MPIHelper::MPICommunicator MPICommunicatorType
MPI communicator type.
Definition dgfug.hh:89
int numParameters(const Entity &) const
return number of parameters
Definition dgfug.hh:148
int boundaryId(const Dune::Intersection< GG, II > &intersection) const
will return boundary segment index
Definition dgfug.hh:129
UGGrid< dim > Grid
grid type
Definition dgfug.hh:85
bool wasInserted(const Dune::Intersection< GG, II > &intersection) const
please doc me
Definition dgfug.hh:122
DGFGridFactory(std::istream &input, MPICommunicatorType comm=MPIHelper::getCommunicator())
constructor taking istream
Definition dgfug.hh:92
int numParameters() const
return number of parameters
Definition dgfug.hh:136
bool haveBoundaryParameters() const
UGGrid does not support boundary parameters.
Definition dgfug.hh:176
std::vector< double > & parameter(const typename Grid::template Codim< 0 >::Entity &element)
return parameter for codim 0 entity
Definition dgfug.hh:154
Grid * grid()
return grid
Definition dgfug.hh:115
static const type & defaultValue()
default constructor
Definition parser.hh:28
std::string type
type of additional boundary parameters
Definition parser.hh:25
Front-end for the grid manager of the finite element toolbox UG3.
Definition uggrid.hh:208