dune-grid 2.10
Loading...
Searching...
No Matches
dgfug.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file LICENSE.md in module root
2// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
4// vi: set et ts=4 sw=2 sts=2:
5#ifndef DUNE_GRID_IO_FILE_DGFPARSER_DGFUG_HH
6#define DUNE_GRID_IO_FILE_DGFPARSER_DGFUG_HH
7
8//- C++ includes
9#include <fstream>
10#include <istream>
11#include <string>
12#include <vector>
13
14//- dune-common includes
15#include <dune/common/exceptions.hh>
16#include <dune/common/fvector.hh>
17#include <dune/common/parallel/mpihelper.hh>
18
19//- dune-grid includes
20#include <dune-grid-config.hh> // HAVE_DUNE_UGGRID
22#include <dune/grid/uggrid.hh>
23
24//- local includes
25#include "dgfparser.hh"
27
28
29namespace Dune
30{
31
32 namespace dgf
33 {
34
35 // UGGridParameterBlock
36 // --------------------
37
39 : public GridParameterBlock
40 {
42 explicit UGGridParameterBlock ( std::istream &input );
43
45 bool noClosure () const { return noClosure_; }
47 bool noCopy () const { return noCopy_; }
49 size_t heapSize () const { return heapSize_; }
50
51 protected:
52 bool noClosure_; // no closure for UGGrid
53 bool noCopy_; // no copies for UGGrid
54 size_t heapSize_; // heap size for UGGrid
55 };
56
57 } // namespace dgf
58
59
60
61#if HAVE_DUNE_UGGRID
62 template< int dim >
63 struct DGFGridInfo< UGGrid< dim > >
64 {
65 static int refineStepsForHalf ()
66 {
67 return 1;
68 }
69
70 static double refineWeight ()
71 {
72 return -1.;
73 }
74 };
75
76
77
78 // DGFGridFactory< UGGrid< dim > >
79 // -------------------------------
80
81 template< int dim >
82 struct DGFGridFactory< UGGrid< dim > >
83 {
87 static const int dimension = dim;
89 typedef MPIHelper::MPICommunicator MPICommunicatorType;
90
92 explicit DGFGridFactory ( std::istream &input,
93 MPICommunicatorType comm = MPIHelper::getCommunicator() )
94 : grid_( 0 ),
95 factory_(),
96 dgf_( rank( comm ), size( comm ) )
97 {
98 generate( input );
99 }
100
102 explicit DGFGridFactory ( const std::string &filename,
103 MPICommunicatorType comm = MPIHelper::getCommunicator() )
104 : grid_( 0 ),
105 factory_(),
106 dgf_( rank( comm ), size( comm ) )
107 {
108 std::ifstream input( filename.c_str() );
109 if ( !input )
110 DUNE_THROW( DGFException, "Error: Macrofile " << filename << " not found" );
111 generate( input );
112 }
113
116 {
117 return grid_;
118 }
119
121 template< class GG, class II >
122 bool wasInserted ( const Dune::Intersection< GG, II > &intersection ) const
123 {
124 return factory_.wasInserted( intersection );
125 }
126
128 template< class GG, class II >
129 int boundaryId ( const Dune::Intersection< GG, II > &intersection ) const
130 {
131 return intersection.boundarySegmentIndex();
132 }
133
135 template< int codim >
136 int numParameters () const
137 {
138 if( codim == 0 )
139 return dgf_.nofelparams;
140 else if( codim == dimension )
141 return dgf_.nofvtxparams;
142 else
143 return 0;
144 }
145
147 template< class Entity >
148 int numParameters ( const Entity & ) const
149 {
150 return numParameters< Entity::codimension >();
151 }
152
154 std::vector< double > &parameter ( const typename Grid::template Codim< 0 >::Entity &element )
155 {
156 if( numParameters< 0 >() <= 0 )
157 {
158 DUNE_THROW( InvalidStateException,
159 "Calling DGFGridFactory::parameter is only allowed if there are parameters." );
160 }
161 return dgf_.elParams[ factory_.insertionIndex( element ) ];
162 }
163
165 std::vector< double > &parameter ( const typename Grid::template Codim< dimension >::Entity &vertex )
166 {
167 if( numParameters< dimension >() <= 0 )
168 {
169 DUNE_THROW( InvalidStateException,
170 "Calling DGFGridFactory::parameter is only allowed if there are parameters." );
171 }
172 return dgf_.vtxParams[ factory_.insertionIndex( vertex ) ];
173 }
174
177 {
178 return dgf_.haveBndParameters;
179 }
180
182 template< class GG, class II >
184 {
186 typename Intersection::Entity entity = intersection.inside();
187 const int face = intersection.indexInInside();
188
189 auto refElem = referenceElement< double, dimension >( entity.type() );
190 int corners = refElem.size( face, 1, dimension );
191 std::vector< unsigned int > bound( corners );
192 for( int i = 0; i < corners; ++i )
193 {
194 const int k = refElem.subEntity( face, 1, i, dimension );
195 bound[ i ] = factory_.insertionIndex( entity.template subEntity< dimension >( k ) );
196 }
197
198 DuneGridFormatParser::facemap_t::key_type key( bound, false );
199 const DuneGridFormatParser::facemap_t::const_iterator pos = dgf_.facemap.find( key );
200 if( pos != dgf_.facemap.end() )
201 return dgf_.facemap.find( key )->second.second;
202 else
204 }
205
206 private:
207 // create grid
208 void generate ( std::istream &input );
209
210 // return rank
211 static int rank( MPICommunicatorType MPICOMM )
212 {
213 int rank = 0;
214#if HAVE_MPI
215 MPI_Comm_rank( MPICOMM, &rank );
216#endif
217 return rank;
218 }
219
220 // return size
221 static int size( MPICommunicatorType MPICOMM )
222 {
223 int size = 1;
224#if HAVE_MPI
225 MPI_Comm_size( MPICOMM, &size );
226#endif
227 return size;
228 }
229
230 Grid *grid_;
231 GridFactory< UGGrid< dim > > factory_;
232 DuneGridFormatParser dgf_;
233 };
234#endif // #if HAVE_DUNE_UGGRID
235
236} // namespace Dune
237
238#endif // #ifndef DUNE_GRID_IO_FILE_DGFPARSER_DGFUG_HH
The UGGrid class.
Include standard header files.
Definition agrid.hh:60
Definition dgfgridfactory.hh:38
MPIHelper::MPICommunicator MPICommunicatorType
Definition dgfgridfactory.hh:41
G Grid
Definition dgfgridfactory.hh:39
static const int dimension
Definition dgfgridfactory.hh:40
Intersection of a mesh entity of codimension 0 ("element") with a "neighboring" element or with the d...
Definition common/intersection.hh:164
int indexInInside() const
Local index of codim 1 entity in the inside() entity where intersection is contained in.
Definition common/intersection.hh:346
size_t boundarySegmentIndex() const
index of the boundary segment within the macro grid
Definition common/intersection.hh:236
Entity inside() const
return Entity on the inside of this intersection. That is the Entity where we started this.
Definition common/intersection.hh:250
GridImp::template Codim< 0 >::Entity Entity
Type of entity that this Intersection belongs to.
Definition common/intersection.hh:192
Wrapper class for entities.
Definition common/entity.hh:66
A Traits struct that collects all associated types of one implementation.
Definition common/grid.hh:411
Common Grid parameters.
Definition gridparameter.hh:35
exception class for IO errors in the DGF parser
Definition dgfexception.hh:16
Some simple static information for a given GridType.
Definition io/file/dgfparser/dgfparser.hh:56
Definition dgfug.hh:40
size_t heapSize() const
returns heap size used on construction of the grid
Definition dgfug.hh:49
bool noCopy_
Definition dgfug.hh:53
size_t heapSize_
Definition dgfug.hh:54
bool noClosure_
Definition dgfug.hh:52
bool noCopy() const
returns true if no copies are made for UGGrid elements
Definition dgfug.hh:47
bool noClosure() const
returns true if no closure should be used for UGGrid
Definition dgfug.hh:45
static double refineWeight()
Definition dgfug.hh:70
static int refineStepsForHalf()
Definition dgfug.hh:65
const DGFBoundaryParameter::type & boundaryParameter(const Dune::Intersection< GG, II > &intersection) const
return invalid value
Definition dgfug.hh:183
std::vector< double > & parameter(const typename Grid::template Codim< dimension >::Entity &vertex)
return parameter for vertex
Definition dgfug.hh:165
DGFGridFactory(const std::string &filename, MPICommunicatorType comm=MPIHelper::getCommunicator())
constructor taking filename
Definition dgfug.hh:102
MPIHelper::MPICommunicator MPICommunicatorType
MPI communicator type.
Definition dgfug.hh:89
int numParameters(const Entity &) const
return number of parameters
Definition dgfug.hh:148
int boundaryId(const Dune::Intersection< GG, II > &intersection) const
will return boundary segment index
Definition dgfug.hh:129
UGGrid< dim > Grid
grid type
Definition dgfug.hh:85
bool wasInserted(const Dune::Intersection< GG, II > &intersection) const
please doc me
Definition dgfug.hh:122
DGFGridFactory(std::istream &input, MPICommunicatorType comm=MPIHelper::getCommunicator())
constructor taking istream
Definition dgfug.hh:92
int numParameters() const
return number of parameters
Definition dgfug.hh:136
bool haveBoundaryParameters() const
UGGrid does not support boundary parameters.
Definition dgfug.hh:176
std::vector< double > & parameter(const typename Grid::template Codim< 0 >::Entity &element)
return parameter for codim 0 entity
Definition dgfug.hh:154
Grid * grid()
return grid
Definition dgfug.hh:115
static const type & defaultValue()
default constructor
Definition parser.hh:28
std::string type
type of additional boundary parameters
Definition parser.hh:25
Front-end for the grid manager of the finite element toolbox UG3.
Definition uggrid.hh:208