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> 62 static int refineStepsForHalf ()
67 static double refineWeight ()
84 static const int dimension = dim;
86 typedef MPIHelper::MPICommunicator MPICommunicatorType;
90 MPICommunicatorType comm = MPIHelper::getCommunicator() )
93 dgf_( rank( comm ), size( comm ) )
100 MPICommunicatorType comm = MPIHelper::getCommunicator() )
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 >
133 int numParameters ()
const 136 return dgf_.nofelparams;
137 else if( codim == dimension )
138 return dgf_.nofvtxparams;
144 template<
class Entity >
145 int numParameters (
const Entity & )
const 147 return numParameters< Entity::codimension >();
151 std::vector< double > ¶meter (
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 > ¶meter (
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 ) ];
173 bool haveBoundaryParameters ()
const 175 return dgf_.haveBndParameters;
179 template<
class GG,
class II >
186 const ReferenceElement< double, dimension > &refElem
187 = ReferenceElements< double, dimension >::general( entity.type() );
188 int corners = refElem.size( face, 1, dimension );
189 std::vector< unsigned int > bound( corners );
190 for(
int i = 0; i < corners; ++i )
192 const int k = refElem.subEntity( face, 1, i, dimension );
193 bound[ i ] = factory_.insertionIndex( entity.template subEntity< dimension >( k ) );
196 DuneGridFormatParser::facemap_t::key_type key( bound,
false );
197 const DuneGridFormatParser::facemap_t::const_iterator pos = dgf_.facemap.find( key );
198 if( pos != dgf_.facemap.end() )
199 return dgf_.facemap.find( key )->second.second;
206 void generate ( std::istream &input );
209 static int rank( MPICommunicatorType MPICOMM )
213 MPI_Comm_rank( MPICOMM, &rank );
219 static int size( MPICommunicatorType MPICOMM )
223 MPI_Comm_size( MPICOMM, &size );
232 #endif // #if ENABLE_UG 236 #endif // #ifndef DUNE_GRID_IO_FILE_DGFPARSER_DGFUG_HH bool noCopy_
Definition: dgfug.hh:50
bool noClosure() const
returns true if no closure should be used for UGGrid
Definition: dgfug.hh:42
Intersection of a mesh entities of codimension 0 ("elements") with a "neighboring" element or with th...
Definition: albertagrid/dgfparser.hh:26
Front-end for the grid manager of the finite element toolbox UG.
Definition: uggrid.hh:207
static const type & defaultValue()
default constructor
Definition: parser.hh:26
Wrapper class for entities.
Definition: common/entity.hh:61
size_t heapSize_
Definition: dgfug.hh:51
size_t boundarySegmentIndex() const
index of the boundary segment within the macro grid
Definition: common/intersection.hh:260
Provide a generic factory class for unstructured grids.
Definition: common/gridfactory.hh:263
std::string type
type of additional boundary parameters
Definition: parser.hh:23
Entity inside() const
return EntityPointer to the Entity on the inside of this intersection. That is the Entity where we st...
Definition: common/intersection.hh:286
bool noClosure_
Definition: dgfug.hh:49
UGGridParameterBlock(std::istream &input)
constructor taking istream
Definition: dgfug.cc:16
size_t heapSize() const
returns heap size used on construction of the grid
Definition: dgfug.hh:46
Grid abstract base classThis class is the base class for all grid implementations. Although no virtual functions are used we call it abstract since its methods do not contain an implementation but forward to the methods of the derived class via the Barton-Nackman trick.
Definition: common/grid.hh:388
Definition: common.hh:179
Include standard header files.
Definition: agrid.hh:59
The DuneGridFormatParser class: reads a DGF file and stores build information in vector structures us...
Definition: parser.hh:44
int indexInInside() const
Local index of codim 1 entity in the inside() entity where intersection is contained in...
Definition: common/intersection.hh:393
GridImp::template Codim< 0 >::Entity Entity
Type of entity that this Intersection belongs to.
Definition: common/intersection.hh:185
Common Grid parametersFor each grid implementation there is a set of parameters that can be passed vi...
Definition: gridparameter.hh:31
Some simple static information for a given GridType.
Definition: io/file/dgfparser/dgfparser.hh:54
bool noCopy() const
returns true if no copies are made for UGGrid elements
Definition: dgfug.hh:44
exception class for IO errors in the DGF parser
Definition: dgfexception.hh:12