1
#ifndef DUNE_GRID_IO_FILE_DGFPARSER_DGFUG_HH
2
#define DUNE_GRID_IO_FILE_DGFPARSER_DGFUG_HH
10
//- dune-common includes
11
#include <dune/common/exceptions.hh>
12
#include <dune/common/fvector.hh>
13
#include <dune/common/mpihelper.hh>
15
//- dune-grid includes
16
#include <dune/grid/common/intersection.hh>
17
#include <dune/grid/uggrid.hh>
20
#include "dgfparser.hh"
21
#include "blocks/gridparameter.hh"
30
// UGGridParameterBlock
31
// --------------------
33
struct UGGridParameterBlock
34
: public GridParameterBlock
36
/** \brief constructor taking istream */
37
explicit UGGridParameterBlock ( std::istream &input );
39
/** \brief returns true if no closure should be used for UGGrid */
40
bool noClosure () const { return noClosure_; }
41
/** \brief returns true if no copies are made for UGGrid elements */
42
bool noCopy () const { return noCopy_; }
43
/** \brief returns heap size used on construction of the grid */
44
size_t heapSize () const { return heapSize_; }
47
bool noClosure_; // no closure for UGGrid
48
bool noCopy_; // no copies for UGGrid
49
size_t heapSize_; // heap size for UGGrid
58
struct DGFGridInfo< UGGrid< dim > >
60
static int refineStepsForHalf ()
65
static double refineWeight ()
73
// DGFGridFactory< UGGrid< dim > >
74
// -------------------------------
77
struct DGFGridFactory< UGGrid< dim > >
79
/** \brief grid type */
80
typedef UGGrid< dim > Grid;
81
/** \brief grid dimension */
82
static const int dimension = dim;
83
/** \brief MPI communicator type */
84
typedef MPIHelper::MPICommunicator MPICommunicatorType;
86
/** \brief constructor taking istream */
87
explicit DGFGridFactory ( std::istream &input,
88
MPICommunicatorType comm = MPIHelper::getCommunicator() )
91
dgf_( rank( comm ), size( comm ) )
96
/** \brief constructor taking filename */
97
explicit DGFGridFactory ( const std::string &filename,
98
MPICommunicatorType comm = MPIHelper::getCommunicator() )
101
dgf_( rank( comm ), size( comm ) )
103
std::ifstream input( filename.c_str() );
105
DUNE_THROW( DGFException, "Error: Macrofile " << filename << " not found" );
109
/** \brief return grid */
115
/** \brief please doc me */
116
template< class GG, template< class > class II >
117
bool wasInserted ( const Dune::Intersection< GG, II > &intersection ) const
119
return factory_.wasInserted( intersection );
122
/** \brief will return boundary segment index */
123
template < class GG, template< class > class II >
124
int boundaryId ( const Dune::Intersection< GG, II > &intersection ) const
126
return intersection.boundarySegmentIndex();
129
/** \brief return number of parameters */
130
template< int codim >
131
int numParameters () const
134
return dgf_.nofelparams;
135
else if( codim == dimension )
136
return dgf_.nofvtxparams;
141
/** \brief return number of parameters */
142
template< class Entity >
143
int numParameters ( const Entity & ) const
145
return numParameters< Entity::codimension >();
148
/** \brief return parameter for codim 0 entity */
149
std::vector< double > ¶meter ( const typename Grid::template Codim< 0 >::Entity &element )
151
if( numParameters< 0 >() <= 0 )
153
DUNE_THROW( InvalidStateException,
154
"Calling DGFGridFactory::parameter is only allowed if there are parameters." );
156
return dgf_.elParams[ factory_.insertionIndex( element ) ];
159
/** \brief return parameter for vertex */
160
std::vector< double > ¶meter ( const typename Grid::template Codim< dimension >::Entity &vertex )
162
if( numParameters< dimension >() <= 0 )
164
DUNE_THROW( InvalidStateException,
165
"Calling DGFGridFactory::parameter is only allowed if there are parameters." );
167
return dgf_.vtxParams[ factory_.insertionIndex( vertex ) ];
170
/** \brief UGGrid does not support boundary parameters */
171
bool haveBoundaryParameters () const
173
return dgf_.haveBndParameters;
176
/** \brief return invalid value */
177
template < class GG, template< class > class II >
178
const DGFBoundaryParameter::type &boundaryParameter ( const Dune::Intersection< GG, II > &intersection ) const
180
typedef Dune::Intersection< GG, II > Intersection;
181
typename Intersection::EntityPointer inside = intersection.inside();
182
const typename Intersection::Entity &entity = *inside;
183
const int face = intersection.indexInInside();
185
const GenericReferenceElement< double, dimension > &refElem
186
= GenericReferenceElements< double, dimension >::general( 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;
200
return DGFBoundaryParameter::defaultValue();
205
void generate ( std::istream &input );
208
static int rank( MPICommunicatorType MPICOMM )
212
MPI_Comm_rank( MPICOMM, &rank );
218
static int size( MPICommunicatorType MPICOMM )
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