Dune Core Modules (unstable)

dgfug.hh
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
18
19//- dune-grid includes
20#include <dune-grid-config.hh> // HAVE_DUNE_UGGRID
21#include <dune/grid/common/intersection.hh>
22#include <dune/grid/uggrid.hh>
23
24//- local includes
25#include "dgfparser.hh"
26#include "blocks/gridparameter.hh"
27
28
29namespace Dune
30{
31
32 namespace dgf
33 {
34
35 // UGGridParameterBlock
36 // --------------------
37
38 struct UGGridParameterBlock
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 {
85 typedef UGGrid< dim > Grid;
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
115 Grid *grid ()
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
176 bool haveBoundaryParameters () const
177 {
178 return dgf_.haveBndParameters;
179 }
180
182 template< class GG, class II >
183 const DGFBoundaryParameter::type &boundaryParameter ( const Dune::Intersection< GG, II > &intersection ) const
184 {
185 typedef Dune::Intersection< GG, II > Intersection;
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
Intersection of a mesh entity of codimension 0 ("element") with a "neighboring" element or with the d...
Definition: intersection.hh:164
int indexInInside() const
Local index of codim 1 entity in the inside() entity where intersection is contained in.
Definition: intersection.hh:346
size_t boundarySegmentIndex() const
index of the boundary segment within the macro grid
Definition: intersection.hh:236
Entity inside() const
return Entity on the inside of this intersection. That is the Entity where we started this.
Definition: intersection.hh:250
GridImp::template Codim< 0 >::Entity Entity
Type of entity that this Intersection belongs to.
Definition: intersection.hh:192
MPI_Comm MPICommunicator
The type of the mpi communicator.
Definition: mpihelper.hh:192
static MPICommunicator getCommunicator()
get the default communicator
Definition: mpihelper.hh:200
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,...)
Definition: exceptions.hh:312
constexpr GeometryType vertex
GeometryType representing a vertex.
Definition: type.hh:492
Helpers for dealing with MPI.
Dune namespace.
Definition: alignedallocator.hh:13
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
static const type & defaultValue()
default constructor
Definition: parser.hh:28
std::string type
type of additional boundary parameters
Definition: parser.hh:25
static double refineWeight()
static int refineStepsForHalf()
number of globalRefine steps needed to refuce h by 0.5
The UGGrid class.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 22, 23:30, 2024)