3#ifndef DUNE_GRID_IO_FILE_DGFPARSER_DGFUG_HH 
    4#define DUNE_GRID_IO_FILE_DGFPARSER_DGFUG_HH 
   18#include <dune/grid/common/intersection.hh> 
   22#include "dgfparser.hh" 
   23#include "blocks/gridparameter.hh" 
   35    struct UGGridParameterBlock
 
   36      : 
public GridParameterBlock
 
   39      explicit UGGridParameterBlock ( std::istream &input );
 
   42      bool noClosure ()
 const { 
return noClosure_; }
 
   44      bool noCopy ()
 const { 
return noCopy_; }
 
   46      size_t heapSize ()
 const { 
return heapSize_; }
 
   60  struct DGFGridInfo< UGGrid< dim > >
 
   79  struct DGFGridFactory< UGGrid< dim > >
 
   82    typedef UGGrid< dim > Grid;
 
   84    static const int dimension = dim;
 
   89    explicit DGFGridFactory ( std::istream &input,
 
   93        dgf_( rank( comm ), size( comm ) )
 
   99    explicit DGFGridFactory ( 
const std::string &filename,
 
  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 )
 
  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 )
 
  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
 
  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 );
 
  229    GridFactory< UGGrid< dim > > factory_;
 
  230    DuneGridFormatParser dgf_;
 
Intersection of a mesh entities of codimension 0 ("elements") with a "neighboring" element or with th...
Definition: intersection.hh:161
 
int indexInInside() const
Local index of codim 1 entity in the inside() entity where intersection is contained in.
Definition: intersection.hh:393
 
size_t boundarySegmentIndex() const
index of the boundary segment within the macro grid
Definition: intersection.hh:260
 
Entity inside() const
return EntityPointer to the Entity on the inside of this intersection. That is the Entity where we st...
Definition: intersection.hh:286
 
GridImp::template Codim< 0 >::Entity Entity
Type of entity that this Intersection belongs to.
Definition: intersection.hh:185
 
MPI_Comm MPICommunicator
The type of the mpi communicator.
Definition: mpihelper.hh:173
 
static MPICommunicator getCommunicator()
get the default communicator
Definition: mpihelper.hh:181
 
A few common exception classes.
 
Implements a vector constructed from a given type representing a field and a compile-time given size.
 
#define DUNE_THROW(E, m)
Definition: exceptions.hh:243
 
Helpers for dealing with MPI.
 
Dune namespace.
Definition: alignment.hh:10
 
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:484