DUNE-FEM (unstable)

boundaryidprovider.hh
1#ifndef DUNE_FEM_MISC_BOUNDARYIDPROVIDER_HH
2#define DUNE_FEM_MISC_BOUNDARYIDPROVIDER_HH
3
5// includes all grid declarations known to dune-fem
6#include <dune/fem/misc/griddeclaration.hh>
7#include <dune/fem/function/common/localcontribution.hh>
8
9namespace Dune
10{
11
12 // External Forward Declarations
13 // -----------------------------
14
15 template< class Grid >
16 struct HostGridAccess;
17
18
19
20 namespace Fem
21 {
22
23 // Internal Forward Declarations
24 // -----------------------------
25
26 template< class Grid >
27 struct BoundaryIdProvider;
28
29
30
31 // BoundaryIdProvider for AlbertaGrid
32 // ----------------------------------
33
34#if HAVE_ALBERTA
35 template< int dim, int dimW >
36 struct BoundaryIdProvider< AlbertaGrid< dim, dimW > >
37 {
38 typedef AlbertaGrid< dim, dimW > GridType;
39
40 template< class Intersection >
41 static int boundaryId ( const Intersection &intersection )
42 {
43 return intersection.impl().boundaryId();
44 }
45 };
46#endif // #if HAVE_ALBERTA
47
48
49
50 // BoundaryIdProvider for ALUGrid
51 // ------------------------------
52
53#if HAVE_DUNE_ALUGRID
54 template< int dim, int dimw, ALUGridElementType elType, ALUGridRefinementType refineType, class Comm >
55 struct BoundaryIdProvider< ALUGrid< dim, dimw, elType, refineType, Comm > >
56 {
57 typedef ALUGrid< dim, dimw, elType, refineType, Comm > GridType;
58
59 template< class Intersection >
60 static int boundaryId ( const Intersection &intersection )
61 {
62 return intersection.impl().boundaryId();
63 }
64 };
65#endif // #if HAVE_DUNE_ALUGRID
66
67
68
69 // BoundaryIdProvider for CacheItGrid
70 // ----------------------------------
71
72#if HAVE_DUNE_METAGRID
73 template< class HostGrid >
74 struct BoundaryIdProvider< CacheItGrid< HostGrid > >
75 {
76 typedef CacheItGrid< HostGrid > GridType;
77
78 template< class Intersection >
79 static int boundaryId ( const Intersection &intersection )
80 {
81 return BoundaryIdProvider< HostGrid >
82 ::boundaryId ( HostGridAccess< GridType >::hostIntersection( intersection ) );
83 }
84 };
85#endif // #if HAVE_DUNE_METAGRID
86
87
88
89 // BoundaryIdProvider for CartesianGrid
90 // ------------------------------------
91
92#if HAVE_DUNE_METAGRID
93 template< class HostGrid >
94 struct BoundaryIdProvider< CartesianGrid< HostGrid > >
95 {
96 typedef CartesianGrid< HostGrid > GridType;
97
98 template< class Intersection >
99 static int boundaryId ( const Intersection &intersection )
100 {
101 return BoundaryIdProvider< HostGrid >
102 ::boundaryId ( HostGridAccess< GridType >::hostIntersection( intersection ) );
103 }
104 };
105#endif // #if HAVE_DUNE_METAGRID
106
107
108
109 // BoundaryIdProvider for FilteredGrid
110 // -----------------------------------
111
112#if HAVE_DUNE_METAGRID
113 template< class HostGrid >
114 struct BoundaryIdProvider< FilteredGrid< HostGrid > >
115 {
116 typedef FilteredGrid< HostGrid > GridType;
117
118 // todo: FilteredGrid is a filtering grid and, hence, needs a specialized
119 // version of boundaryId.
120 template< class Intersection >
121 static int boundaryId ( const Intersection &intersection )
122 {
123 if( !HostGridAccess< GridType >::hostIntersection( intersection ).boundary() )
124 DUNE_THROW( NotImplemented, "BoundaryIdProvider for artificial boundaries of FilteredGrid not implemented." );
125 return BoundaryIdProvider< HostGrid >
126 ::boundaryId ( HostGridAccess< GridType >::hostIntersection( intersection ) );
127 }
128 };
129#endif // #if HAVE_DUNE_METAGRID
130
131
132
133 // BoundaryIdProvider for GeometryGrid
134 // -----------------------------------
135
136 template< class HostGrid, class CoordFunction, class Allocator >
137 struct BoundaryIdProvider< GeometryGrid< HostGrid, CoordFunction, Allocator > >
138 {
139 typedef GeometryGrid< HostGrid, CoordFunction, Allocator > GridType;
140
141 template< class Intersection >
142 static int boundaryId ( const Intersection &intersection )
143 {
144 return BoundaryIdProvider< HostGrid >
145 ::boundaryId ( HostGridAccess< GridType >::hostIntersection( intersection ) );
146 }
147 };
148
149
150
151 // BoundaryIdProvider for IdGrid
152 // -----------------------------
153
154#if HAVE_DUNE_METAGRID
155 template< class HostGrid >
156 struct BoundaryIdProvider< IdGrid< HostGrid > >
157 {
158 typedef IdGrid< HostGrid > GridType;
159
160 template< class Intersection >
161 static int boundaryId ( const Intersection &intersection )
162 {
163 return BoundaryIdProvider< HostGrid >
164 ::boundaryId ( HostGridAccess< GridType >::hostIntersection( intersection ) );
165 }
166 };
167#endif // #if HAVE_DUNE_METAGRID
168
169
170
171 // BoundaryIdProvider for OneDGrid
172 // -------------------------------
173
174 template<>
175 struct BoundaryIdProvider< OneDGrid >
176 {
177 typedef OneDGrid GridType;
178
179 template< class Intersection >
180 static int boundaryId ( const Intersection &intersection )
181 {
182 return intersection.boundarySegmentIndex();
183 }
184 };
185
186
187
188 // BoundaryIdProvider for ParallelGrid
189 // -----------------------------------
190
191#if HAVE_DUNE_METAGRID
192 template< class HostGrid >
193 struct BoundaryIdProvider< ParallelGrid< HostGrid > >
194 {
195 typedef ParallelGrid< HostGrid > GridType;
196
197 template< class Intersection >
198 static int boundaryId ( const Intersection &intersection )
199 {
200 return BoundaryIdProvider< HostGrid >
201 ::boundaryId ( HostGridAccess< GridType >::hostIntersection( intersection ) );
202 }
203 };
204#endif // #if HAVE_DUNE_METAGRID
205
206
207
208 // BoundaryIdProvider for SphereGrid
209 // ---------------------------------
210
211#if HAVE_DUNE_METAGRID
212 template< class HostGrid, class MapToSphere >
213 struct BoundaryIdProvider< SphereGrid< HostGrid, MapToSphere > >
214 {
215 typedef SphereGrid< HostGrid, MapToSphere > GridType;
216
217 template< class Intersection >
218 static int boundaryId ( const Intersection &intersection )
219 {
220 return BoundaryIdProvider< HostGrid >
221 ::boundaryId ( HostGridAccess< GridType >::hostIntersection( intersection ) );
222 }
223 };
224#endif // #if HAVE_DUNE_METAGRID
225
226
227
228 // BoundaryIdProvider for SPGrid
229 // -----------------------------
230
231#if HAVE_DUNE_SPGRID
232 template< class ct, int dim, template< int > class Strategy, class Comm >
233 struct BoundaryIdProvider< SPGrid< ct, dim, Strategy, Comm > >
234 {
235 typedef SPGrid< ct, dim, Strategy, Comm > GridType;
236
237 template< class Intersection >
238 static int boundaryId ( const Intersection &intersection )
239 {
240 return (intersection.boundary() ? (intersection.indexInInside()+1) : 0);
241 }
242 };
243#endif // #if HAVE_DUNE_SPGRID
244
245
246 // BoundaryIdProvider for PolygonGrid
247 // ----------------------------------
248
249#if HAVE_DUNE_POLYGONGRID
250 template< class ct >
251 struct BoundaryIdProvider< PolygonGrid< ct > >
252 {
253 typedef PolygonGrid< ct > GridType;
254
255 template< class Intersection >
256 static int boundaryId ( const Intersection &intersection )
257 {
258 return (intersection.boundary() ? (intersection.impl().boundaryId()) : 0);
259 }
260 };
261#endif // #if HAVE_OPM_GRID
262
263 // BoundaryIdProvider for PolygonGrid
264 // ----------------------------------
265
266#if HAVE_DUNE_P4ESTGRID
267 template< int dim, int dimworld, class ct >
268 struct BoundaryIdProvider< P4estGrid<dim, dimworld, ct > >
269 {
270 typedef P4estGrid<dim, dimworld, ct > GridType;
271
272 template< class Intersection >
273 static int boundaryId ( const Intersection &intersection )
274 {
275 return (intersection.boundary() ? (intersection.impl().boundaryId()) : 0);
276 }
277 };
278#endif // #if HAVE_OPM_GRID
279
280 // BoundaryIdProvider for PolyhedralGrid
281 // -------------------------------------
282
283#if HAVE_OPM_GRID
284 template< int dim, int dimworld, class ct >
285 struct BoundaryIdProvider< PolyhedralGrid< dim, dimworld, ct > >
286 {
287 typedef PolyhedralGrid< dim, dimworld, ct > GridType;
288
289 template< class Intersection >
290 static int boundaryId ( const Intersection &intersection )
291 {
292 return (intersection.boundary() ? (intersection.impl().boundaryId()) : 0);
293 }
294 };
295#endif // #if HAVE_OPM_GRID
296
297
298 // BoundaryIdProvider for UGGrid
299 // -----------------------------
300
301 template< int dim >
302 struct BoundaryIdProvider< UGGrid< dim > >
303 {
304 typedef UGGrid< dim > GridType;
305
306 template< class Intersection >
307 static int boundaryId ( const Intersection &intersection )
308 {
309 return intersection.boundarySegmentIndex();
310 }
311 };
312
313
314
315 // BoundaryIdProvider for YaspGrid
316 // -------------------------------
317
318 template< int dim, class CoordCont >
319 struct BoundaryIdProvider< YaspGrid< dim, CoordCont > >
320 {
321 typedef YaspGrid< dim, CoordCont > GridType;
322
323 template< class Intersection >
324 static int boundaryId ( const Intersection &intersection )
325 {
326 return (intersection.boundary() ? (intersection.indexInInside()+1) : 0);
327 }
328 };
329
330
331
332 // BoundaryIdProvider for general GridParts or GridViews
333 // -----------------------------------------------------
334
336 template< class GridView, class Intersection>
337 inline static int boundaryId ( const Intersection &intersection )
338 {
339 return Dune::Fem::BoundaryIdProvider< typename GridView::Grid > ::
340 boundaryId( intersection );
341 }
342
344 template< class GridView, class Intersection>
345 inline static int boundaryId ( const GridView&, const Intersection &intersection )
346 {
347 return Dune::Fem::BoundaryIdProvider< typename GridView::Grid > ::
348 boundaryId( intersection );
349 }
350
351
352 /* \brief BoundaryId inspection: project boundary ids to piecewise constant
353 function. Elements with more than one boundary segment will
354 contain the sum of all boundary ids.
355 */
356 template <class DiscreteFunction>
357 void projectBoundaryIds( DiscreteFunction& df )
358 {
359 if( df.space().order() != 0 ) // should be piecewise constant, and FV space
360 {
361 DUNE_THROW(InvalidStateException,"projectBoundaryIds: expect piecewise constant discrete function");
362 }
363
364 df.clear(); // reset all dofs
365
366 const auto& gridPart = df.space().gridPart();
367
368 typedef AddLocalContribution< DiscreteFunction > AddLocalContributionType;
369 AddLocalContributionType localDf( df );
370 for( const auto& element : df.space() )
371 {
372 auto guard = bindGuard( localDf, element );
373 for( const auto& intersection : intersections( gridPart, element ) )
374 {
375 if( intersection.boundary() )
376 {
377 localDf[ 0 ] += boundaryId( gridPart, intersection );
378 }
379 }
380 }
381 }
382
383 } // namespace Fem
384
385} // namespace Dune
386
387#endif // #ifndef DUNE_FEM_MISC_BOUNDARYIDPROVIDER_HH
A few common exception classes.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
Dune namespace.
Definition: alignedallocator.hh:13
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 21, 23:30, 2024)