3#ifndef DUNE_DGFPARSERALU_HH
4#define DUNE_DGFPARSERALU_HH
10#include <dune/grid/io/file/dgfparser/dgfparser.hh>
11#include <dune/grid/io/file/dgfparser/parser.hh>
12#include <dune/grid/io/file/dgfparser/blocks/projection.hh>
14#include <dune/grid/common/intersection.hh>
22 template<
class Gr
idImp,
class IntersectionImp >
32 struct DGFGridInfo< ALUCubeGrid< 3, dimw > >
39 struct DGFGridInfo< ALUSimplexGrid< 3, dimw > >
46 struct DGFGridInfo< ALUSimplexGrid< 2, dimw > >
53 struct DGFGridInfo< ALUCubeGrid< 2, dimw > >
60 struct DGFGridInfo<
Dune::ALUConformGrid< 2, dimw > >
66 template<
int dimg,
int dimw, ALUGr
idElementType eltype, ALUGr
idRefinementType refinementtype,
class Comm >
67 struct DGFGridInfo<
Dune::ALUGrid< dimg, dimw, eltype, refinementtype, Comm > >
69 static int refineStepsForHalf () {
return ( refinementtype == conforming ) ? dimg : 1; }
70 static double refineWeight () {
return ( refinementtype == conforming ) ? 0.5 : 1.0/(std::pow( 2.0,
double(dimg))); }
84 typedef typename Grid::template Codim<0>::Entity Element;
85 typedef typename Grid::template Codim<dimension>::Entity Vertex;
93 explicit DGFBaseFactory ( MPICommunicatorType comm )
95 dgf_( rank(comm), size(comm) )
103 template<
class Intersection >
104 bool wasInserted (
const Intersection &intersection )
const
106 return factory_.wasInserted( intersection );
109 template<
class GG,
class II >
110 int boundaryId (
const Intersection< GG, II > & intersection )
const
115 const int face = intersection.indexInInside();
117 const ReferenceElement< double, dimension > & refElem =
119 int corners = refElem.size( face, 1, dimension );
120 std :: vector< unsigned int > bound( corners );
121 for(
int i=0; i < corners; ++i )
123 const int k = refElem.subEntity( face, 1, i, dimension );
124 bound[ i ] = factory_.insertionIndex( *entity.template subEntity< dimension >( k ) );
127 DuneGridFormatParser::facemap_t::key_type key( bound,
false );
128 const DuneGridFormatParser::facemap_t::const_iterator pos = dgf_.facemap.find( key );
129 if( pos != dgf_.facemap.end() )
130 return dgf_.facemap.find( key )->second.first;
132 return (intersection.boundary() ? 1 : 0);
135 template<
class GG,
class II >
137 boundaryParameter (
const Intersection< GG, II > & intersection )
const
142 const int face = intersection.indexInInside();
144 const ReferenceElement< double, dimension > & refElem =
146 int corners = refElem.size( face, 1, dimension );
147 std :: vector< unsigned int > bound( corners );
148 for(
int i=0; i < corners; ++i )
150 const int k = refElem.subEntity( face, 1, i, dimension );
151 bound[ i ] = factory_.insertionIndex( *entity.template subEntity< dimension >( k ) );
154 DuneGridFormatParser::facemap_t::key_type key( bound,
false );
155 const DuneGridFormatParser::facemap_t::const_iterator pos = dgf_.facemap.find( key );
156 if( pos != dgf_.facemap.end() )
157 return dgf_.facemap.find( key )->second.second;
162 template<
int codim >
163 int numParameters ()
const
166 return dgf_.nofelparams;
167 else if( codim == dimension )
168 return dgf_.nofvtxparams;
174 bool haveBoundaryParameters ()
const
176 return dgf_.haveBndParameters;
179 std::vector< double > ¶meter (
const Element &element )
181 if( numParameters< 0 >() <= 0 )
184 "Calling DGFGridFactory::parameter is only allowed if there are parameters." );
186 return dgf_.elParams[ factory_.insertionIndex( element ) ];
189 std::vector< double > ¶meter (
const Vertex &vertex )
191 if( numParameters< dimension >() <= 0 )
194 "Calling DGFGridFactory::parameter is only allowed if there are parameters." );
196 return dgf_.vtxParams[ factory_.insertionIndex( vertex ) ];
202 MPICommunicatorType communicator,
203 const std::string &filename );
207 MPICommunicatorType communicator,
208 const std::string &filename );
210 static Grid* callDirectly(
const char* gridname,
212 const char *filename,
213 MPICommunicatorType communicator )
215 #if ALU3DGRID_PARALLEL
217 std :: stringstream tmps;
218 tmps << filename <<
"." << rank;
219 const std :: string &tmp = tmps.str();
222 if( fileExists( tmp.c_str() ) )
223 return new Grid( tmp.c_str(), communicator );
228 if( fileExists( filename ) )
229 return new Grid( filename , communicator );
233 DUNE_THROW( GridError,
"Unable to create " << gridname <<
" from '"
234 << filename <<
"'." );
237 dwarn <<
"WARNING: P[" << rank <<
"]: Creating empty grid!" << std::endl;
240 return new Grid( communicator );
242 static bool fileExists (
const char *fileName )
244 std :: ifstream testfile( fileName );
250 static int rank( MPICommunicatorType MPICOMM )
254 MPI_Comm_rank( MPICOMM, &rank );
258 static int size( MPICommunicatorType MPICOMM )
262 MPI_Comm_size( MPICOMM, &size );
267 GridFactory factory_;
268 DuneGridFormatParser dgf_;
272 template <
int dimw >
273 struct DGFGridFactory< ALUSimplexGrid<3,dimw> > :
274 public DGFBaseFactory< ALUSimplexGrid<3,dimw> >
276 typedef ALUSimplexGrid<3,dimw> DGFGridType;
277 typedef DGFBaseFactory< DGFGridType > BaseType;
278 typedef typename BaseType :: MPICommunicatorType MPICommunicatorType;
280 using BaseType :: grid_;
281 using BaseType :: callDirectly;
283 explicit DGFGridFactory ( std::istream &input,
285 : DGFBaseFactory< DGFGridType >( comm )
290 DUNE_THROW( DGFException,
"Error resetting input stream." );
291 generate( input, comm );
294 explicit DGFGridFactory (
const std::string &filename,
296 : DGFBaseFactory< DGFGridType >( comm )
298 std::ifstream input( filename.c_str() );
300 bool fileFound = input.is_open() ;
302 fileFound = generate( input, comm, filename );
305 grid_ = callDirectly(
"ALUSimplexGrid< 3 , 3 >", this->rank( comm ), filename.c_str(), comm );
309 bool generate( std::istream &file, MPICommunicatorType comm,
const std::string &filename =
"" );
313 template <
int dimw >
314 struct DGFGridFactory< ALUCubeGrid<3, dimw> > :
315 public DGFBaseFactory< ALUCubeGrid<3, dimw> >
317 typedef ALUCubeGrid<3,dimw> DGFGridType;
318 typedef DGFBaseFactory< DGFGridType > BaseType;
319 typedef typename BaseType :: MPICommunicatorType MPICommunicatorType;
321 using BaseType :: grid_;
322 using BaseType :: callDirectly;
324 explicit DGFGridFactory ( std::istream &input,
326 : DGFBaseFactory< DGFGridType >( comm )
331 DUNE_THROW( DGFException,
"Error resetting input stream." );
332 generate( input, comm );
335 explicit DGFGridFactory (
const std::string &filename,
337 : DGFBaseFactory< DGFGridType >( comm )
339 std::ifstream input( filename.c_str() );
340 bool fileFound = input.is_open() ;
342 fileFound = generate( input, comm, filename );
345 grid_ = callDirectly(
"ALUCubeGrid< 3 , 3 >", this->rank( comm ), filename.c_str(), comm );
349 bool generate( std::istream &file, MPICommunicatorType comm,
const std::string &filename =
"" );
352 template < ALUGr
idElementType eltype, ALUGr
idRefinementType refinementtype,
class Comm >
353 struct DGFGridFactory< ALUGrid<3,3, eltype, refinementtype, Comm > > :
354 public DGFBaseFactory< ALUGrid<3,3, eltype, refinementtype, Comm > >
356 typedef ALUGrid<3,3, eltype, refinementtype, Comm > DGFGridType;
357 typedef DGFBaseFactory< DGFGridType > BaseType;
358 typedef typename BaseType :: MPICommunicatorType MPICommunicatorType;
360 using BaseType :: grid_;
361 using BaseType :: callDirectly;
363 explicit DGFGridFactory ( std::istream &input,
370 DUNE_THROW( DGFException,
"Error resetting input stream." );
371 generate( input, comm );
374 explicit DGFGridFactory (
const std::string &filename,
378 std::ifstream input( filename.c_str() );
379 bool fileFound = input.is_open() ;
381 fileFound = generate( input, comm, filename );
384 grid_ = callDirectly(
"ALUGrid< 3, 3, eltype, ref, comm >", this->rank( comm ), filename.c_str(), comm );
388 bool generate( std::istream &file, MPICommunicatorType comm,
const std::string &filename =
"" );
392 struct DGFGridFactory< ALUConformGrid<2,dimw> > :
393 public DGFBaseFactory< ALUConformGrid<2,dimw> >
395 typedef DGFBaseFactory< ALUConformGrid<2,dimw> > Base;
396 typedef typename Base:: MPICommunicatorType MPICommunicatorType;
397 typedef typename Base::Grid Grid;
399 typedef typename Base::GridFactory GridFactory;
401 explicit DGFGridFactory ( std::istream &input,
408 DUNE_THROW( DGFException,
"Error resetting input stream." );
409 generate( input, comm );
412 explicit DGFGridFactory (
const std::string &filename,
414 : DGFBaseFactory< ALUConformGrid<2,dimw> >()
416 std::ifstream input( filename.c_str() );
418 DUNE_THROW( DGFException,
"Macrofile '" << filename <<
"' not found." );
419 if( !generate( input, comm, filename ) )
421 if( Base::fileExists( filename.c_str() ) )
422 Base::grid_ =
new Grid( filename );
424 DUNE_THROW( GridError,
"Unable to create a 2d ALUGrid from '" << filename <<
"'." );
428 bool generate( std::istream &file, MPICommunicatorType comm,
const std::string &filename =
"" );
430 using Base::factory_;
435 struct DGFGridFactory< ALUSimplexGrid<2,dimw> > :
436 public DGFBaseFactory< ALUSimplexGrid<2,dimw> >
438 typedef DGFBaseFactory< ALUSimplexGrid<2,dimw> > Base;
439 typedef typename Base::MPICommunicatorType MPICommunicatorType;
440 typedef typename Base::Grid Grid;
442 typedef typename Base::GridFactory GridFactory;
444 explicit DGFGridFactory ( std::istream &input,
451 DUNE_THROW(DGFException,
"Error resetting input stream." );
452 generate( input, comm );
455 explicit DGFGridFactory (
const std::string &filename,
457 : DGFBaseFactory< ALUSimplexGrid<2,dimw> >()
459 std::ifstream input( filename.c_str() );
461 DUNE_THROW( DGFException,
"Macrofile '" << filename <<
"' not found." );
462 if( !generate( input, comm, filename ) )
464 if( Base::fileExists( filename.c_str() ) )
465 Base::grid_ =
new Grid( filename );
467 DUNE_THROW( GridError,
"Unable to create a 2d ALUGrid from '" << filename <<
"'." );
472 bool generate( std::istream &file, MPICommunicatorType comm,
const std::string &filename =
"" );
474 using Base::factory_;
479 struct DGFGridFactory< ALUCubeGrid<2,dimw> > :
480 public DGFBaseFactory< ALUCubeGrid<2,dimw> >
482 typedef DGFBaseFactory< ALUCubeGrid<2,dimw> > Base;
483 typedef typename Base::MPICommunicatorType MPICommunicatorType;
484 typedef typename Base::Grid Grid;
486 typedef typename Base::GridFactory GridFactory;
488 explicit DGFGridFactory ( std::istream &input,
495 DUNE_THROW(DGFException,
"Error resetting input stream." );
496 generate( input, comm );
499 explicit DGFGridFactory (
const std::string &filename,
501 : DGFBaseFactory< ALUCubeGrid<2,dimw> >()
503 std::ifstream input( filename.c_str() );
505 DUNE_THROW( DGFException,
"Macrofile '" << filename <<
"' not found." );
506 if( !generate( input, comm, filename ) )
508 if( Base::fileExists( filename.c_str() ) )
509 Base::grid_ =
new Grid( filename );
511 DUNE_THROW( GridError,
"Unable to create a 2d ALUGrid from '" << filename <<
"'." );
516 bool generate( std::istream &file, MPICommunicatorType comm,
const std::string &filename =
"" );
518 using Base::factory_;
522 template <
int dimw, ALUGr
idElementType eltype, ALUGr
idRefinementType refinementtype,
class Comm >
523 struct DGFGridFactory< ALUGrid<2, dimw, eltype, refinementtype, Comm > > :
524 public DGFBaseFactory< ALUGrid< 2, dimw, eltype, refinementtype, Comm > >
526 typedef ALUGrid< 2, dimw, eltype, refinementtype, Comm > DGFGridType;
527 typedef DGFBaseFactory< DGFGridType > BaseType;
528 typedef typename BaseType :: MPICommunicatorType MPICommunicatorType;
529 typedef typename BaseType::Grid Grid;
531 typedef typename BaseType::GridFactory GridFactory;
533 explicit DGFGridFactory ( std::istream &input,
540 DUNE_THROW(DGFException,
"Error resetting input stream." );
541 generate( input, comm );
544 explicit DGFGridFactory (
const std::string &filename,
548 std::ifstream input( filename.c_str() );
550 DUNE_THROW( DGFException,
"Macrofile '" << filename <<
"' not found." );
551 if( !generate( input, comm, filename ) )
553 if( BaseType::fileExists( filename.c_str() ) )
554 grid_ =
new Grid( filename );
556 DUNE_THROW( GridError,
"Unable to create a 2d ALUGrid from '" << filename <<
"'." );
561 bool generate( std::istream &file, MPICommunicatorType comm,
const std::string &filename =
"" );
562 using BaseType::grid_;
563 using BaseType::factory_;
564 using BaseType::dgf_;
Provides base classes for ALUGrid.
@ dimension
The dimension of the grid.
Definition: grid.hh:400
Intersection of a mesh entities of codimension 0 ("elements") with a "neighboring" element or with th...
Definition: intersection.hh:161
GridImp::template Codim< 0 >::EntityPointer EntityPointer
Pointer to the type of entities that this Intersection belongs to.
Definition: intersection.hh:191
GridImp::template Codim< 0 >::Entity Entity
Type of entity that this Intersection belongs to.
Definition: intersection.hh:188
MPI_Comm MPICommunicator
The type of the mpi communicator.
Definition: mpihelper.hh:174
static MPICommunicator getCommunicator()
get the default communicator
Definition: mpihelper.hh:182
#define DUNE_THROW(E, m)
Definition: exceptions.hh:244
DWarnType dwarn(std::cerr)
Stream for warnings indicating problems.
Definition: stdstreams.hh:160
Dune namespace.
Definition: alignment.hh:14
ALUGridElementType
basic element types for ALUGrid
Definition: declaration.hh:18
static const type & defaultValue()
default constructor
Definition: parser.hh:26
std::string type
type of additional boundary parameters
Definition: parser.hh:23
static double refineWeight()
static int refineStepsForHalf()
number of globalRefine steps needed to refuce h by 0.5
static const ReferenceElement< ctype, dim > & general(const GeometryType &type)
get general reference elements
Definition: referenceelements.hh:568