Dune Core Modules (2.3.1)

dgfug.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_GRID_IO_FILE_DGFPARSER_DGFUG_HH
4#define DUNE_GRID_IO_FILE_DGFPARSER_DGFUG_HH
5
6//- C++ includes
7#include <fstream>
8#include <istream>
9#include <string>
10#include <vector>
11
12//- dune-common includes
16
17//- dune-grid includes
18#include <dune/grid/common/intersection.hh>
19#include <dune/grid/uggrid.hh>
20
21//- local includes
22#include "dgfparser.hh"
23#include "blocks/gridparameter.hh"
24
25
26namespace Dune
27{
28
29 namespace dgf
30 {
31
32 // UGGridParameterBlock
33 // --------------------
34
35 struct UGGridParameterBlock
36 : public GridParameterBlock
37 {
39 explicit UGGridParameterBlock ( std::istream &input );
40
42 bool noClosure () const { return noClosure_; }
44 bool noCopy () const { return noCopy_; }
46 size_t heapSize () const { return heapSize_; }
47
48 protected:
49 bool noClosure_; // no closure for UGGrid
50 bool noCopy_; // no copies for UGGrid
51 size_t heapSize_; // heap size for UGGrid
52 };
53
54 } // namespace dgf
55
56
57
58#if ENABLE_UG
59 template< int dim >
60 struct DGFGridInfo< UGGrid< dim > >
61 {
62 static int refineStepsForHalf ()
63 {
64 return 1;
65 }
66
67 static double refineWeight ()
68 {
69 return -1.;
70 }
71 };
72
73
74
75 // DGFGridFactory< UGGrid< dim > >
76 // -------------------------------
77
78 template< int dim >
79 struct DGFGridFactory< UGGrid< dim > >
80 {
82 typedef UGGrid< dim > Grid;
84 static const int dimension = dim;
86 typedef MPIHelper::MPICommunicator MPICommunicatorType;
87
89 explicit DGFGridFactory ( std::istream &input,
90 MPICommunicatorType comm = MPIHelper::getCommunicator() )
91 : grid_( 0 ),
92 factory_(),
93 dgf_( rank( comm ), size( comm ) )
94 {
95 generate( input );
96 }
97
99 explicit DGFGridFactory ( const std::string &filename,
100 MPICommunicatorType comm = MPIHelper::getCommunicator() )
101 : grid_( 0 ),
102 factory_(),
103 dgf_( rank( comm ), size( comm ) )
104 {
105 std::ifstream input( filename.c_str() );
106 if ( !input )
107 DUNE_THROW( DGFException, "Error: Macrofile " << filename << " not found" );
108 generate( input );
109 }
110
112 Grid *grid ()
113 {
114 return grid_;
115 }
116
118 template< class GG, class II >
119 bool wasInserted ( const Dune::Intersection< GG, II > &intersection ) const
120 {
121 return factory_.wasInserted( intersection );
122 }
123
125 template< class GG, class II >
126 int boundaryId ( const Dune::Intersection< GG, II > &intersection ) const
127 {
128 return intersection.boundarySegmentIndex();
129 }
130
132 template< int codim >
133 int numParameters () const
134 {
135 if( codim == 0 )
136 return dgf_.nofelparams;
137 else if( codim == dimension )
138 return dgf_.nofvtxparams;
139 else
140 return 0;
141 }
142
144 template< class Entity >
145 int numParameters ( const Entity & ) const
146 {
147 return numParameters< Entity::codimension >();
148 }
149
151 std::vector< double > &parameter ( const typename Grid::template Codim< 0 >::Entity &element )
152 {
153 if( numParameters< 0 >() <= 0 )
154 {
155 DUNE_THROW( InvalidStateException,
156 "Calling DGFGridFactory::parameter is only allowed if there are parameters." );
157 }
158 return dgf_.elParams[ factory_.insertionIndex( element ) ];
159 }
160
162 std::vector< double > &parameter ( const typename Grid::template Codim< dimension >::Entity &vertex )
163 {
164 if( numParameters< dimension >() <= 0 )
165 {
166 DUNE_THROW( InvalidStateException,
167 "Calling DGFGridFactory::parameter is only allowed if there are parameters." );
168 }
169 return dgf_.vtxParams[ factory_.insertionIndex( vertex ) ];
170 }
171
173 bool haveBoundaryParameters () const
174 {
175 return dgf_.haveBndParameters;
176 }
177
179 template< class GG, class II >
180 const DGFBoundaryParameter::type &boundaryParameter ( const Dune::Intersection< GG, II > &intersection ) const
181 {
182 typedef Dune::Intersection< GG, II > Intersection;
183 typename Intersection::EntityPointer inside = intersection.inside();
184 const typename Intersection::Entity &entity = *inside;
185 const int face = intersection.indexInInside();
186
187 const ReferenceElement< double, dimension > &refElem
189 int corners = refElem.size( face, 1, dimension );
190 std::vector< unsigned int > bound( corners );
191 for( int i = 0; i < corners; ++i )
192 {
193 const int k = refElem.subEntity( face, 1, i, dimension );
194 bound[ i ] = factory_.insertionIndex( *entity.template subEntity< dimension >( k ) );
195 }
196
197 DuneGridFormatParser::facemap_t::key_type key( bound, false );
198 const DuneGridFormatParser::facemap_t::const_iterator pos = dgf_.facemap.find( key );
199 if( pos != dgf_.facemap.end() )
200 return dgf_.facemap.find( key )->second.second;
201 else
203 }
204
205 private:
206 // create grid
207 void generate ( std::istream &input );
208
209 // return rank
210 static int rank( MPICommunicatorType MPICOMM )
211 {
212 int rank = 0;
213#if HAVE_MPI
214 MPI_Comm_rank( MPICOMM, &rank );
215#endif
216 return rank;
217 }
218
219 // return size
220 static int size( MPICommunicatorType MPICOMM )
221 {
222 int size = 1;
223#if HAVE_MPI
224 MPI_Comm_size( MPICOMM, &size );
225#endif
226 return size;
227 }
228
229 Grid *grid_;
230 GridFactory< UGGrid< dim > > factory_;
231 DuneGridFormatParser dgf_;
232 };
233#endif // #if ENABLE_UG
234
235} // namespace Dune
236
237#endif // #ifndef DUNE_GRID_IO_FILE_DGFPARSER_DGFUG_HH
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:363
GridImp::template Codim< 0 >::EntityPointer EntityPointer
Pointer to the type of entities that this Intersection belongs to.
Definition: intersection.hh:191
size_t boundarySegmentIndex() const
index of the boundary segment within the macro grid
Definition: intersection.hh:259
GridImp::template Codim< 0 >::Entity Entity
Type of entity that this Intersection belongs to.
Definition: intersection.hh:188
EntityPointer inside() const
return EntityPointer to the Entity on the inside of this intersection. That is the Entity where we st...
Definition: intersection.hh:273
MPI_Comm MPICommunicator
The type of the mpi communicator.
Definition: mpihelper.hh:174
static MPICommunicator getCommunicator()
get the default communicator
Definition: mpihelper.hh:182
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:244
Dune namespace.
Definition: alignment.hh:14
Helpers for dealing with MPI.
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
The UGGrid class.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 12, 23:30, 2024)