DUNE-FEM (unstable)

dirichletconstraints.hh
1/**************************************************************************
2
3 The dune-fem module is a module of DUNE (see www.dune-project.org).
4 It is based on the dune-grid interface library
5 extending the grid interface by a number of discretization algorithms
6 for solving non-linear systems of partial differential equations.
7
8 Copyright (C) 2003 - 2015 Robert Kloefkorn
9 Copyright (C) 2003 - 2010 Mario Ohlberger
10 Copyright (C) 2004 - 2015 Andreas Dedner
11 Copyright (C) 2005 Adrian Burri
12 Copyright (C) 2005 - 2015 Mirko Kraenkel
13 Copyright (C) 2006 - 2015 Christoph Gersbacher
14 Copyright (C) 2006 - 2015 Martin Nolte
15 Copyright (C) 2011 - 2015 Tobias Malkmus
16 Copyright (C) 2012 - 2015 Stefan Girke
17 Copyright (C) 2013 - 2015 Claus-Justus Heine
18 Copyright (C) 2013 - 2014 Janick Gerstenberger
19 Copyright (C) 2013 Sven Kaulman
20 Copyright (C) 2013 Tom Ranner
21 Copyright (C) 2015 Marco Agnese
22 Copyright (C) 2015 Martin Alkaemper
23
24
25 The dune-fem module is free software; you can redistribute it and/or
26 modify it under the terms of the GNU General Public License as
27 published by the Free Software Foundation; either version 2 of
28 the License, or (at your option) any later version.
29
30 The dune-fem module is distributed in the hope that it will be useful,
31 but WITHOUT ANY WARRANTY; without even the implied warranty of
32 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
33 GNU General Public License for more details.
34
35 You should have received a copy of the GNU General Public License along
36 with this program; if not, write to the Free Software Foundation, Inc.,
37 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
38
39**************************************************************************/
40#ifndef DUNE_DIRICHLETCONSTRAINTS_HH
41#define DUNE_DIRICHLETCONSTRAINTS_HH
42
43#include <dune/fem/function/common/scalarproducts.hh>
44#include <dune/fem/function/localfunction/const.hh>
45#include <dune/fem/operator/common/temporarylocalmatrix.hh>
46#include <dune/fem/common/bindguard.hh>
47#include <dune/fem/common/coordinate.hh>
48#include <dune/fem/common/localcontribution.hh>
49#include <dune/fem/space/common/localinterpolation.hh>
50#include <dune/fem/space/mapper/nonblockmapper.hh>
51#include <dune/fem/misc/hasboundaryintersection.hh>
52
53namespace Dune {
54
55namespace Fem {
56class HasLocalFunction;
57}
58
59template < class Model, class DiscreteFunctionSpace, bool useIdentity >
60class DirichletConstraints
61{
62public:
63 enum Operation { set = 0, sub = 1, add = 2 };
64 typedef Model ModelType;
65 typedef DiscreteFunctionSpace DiscreteFunctionSpaceType;
66 typedef typename DiscreteFunctionSpaceType::DomainType DomainType;
67 typedef typename DiscreteFunctionSpaceType::RangeType RangeType;
68 typedef typename DiscreteFunctionSpaceType::JacobianRangeType JacobianRangeType;
69 typedef typename DiscreteFunctionSpaceType::HessianRangeType HessianRangeType;
70
72 typedef typename DiscreteFunctionSpaceType :: GridPartType GridPartType;
74 typedef typename DiscreteFunctionSpaceType :: GridType GridType;
75
76 // types for boundary treatment
77 // ----------------------------
78 static const int localBlockSize = DiscreteFunctionSpaceType :: localBlockSize ;
79 static const int dimRange = DiscreteFunctionSpaceType::FunctionSpaceType::dimRange;
80 static_assert( localBlockSize <= dimRange,
81 "local block size of the space must be less than or equal to the dimension of the range of the function space.");
82 typedef std::array<int,localBlockSize> DirichletBlock;
83 typedef std::vector< DirichletBlock > DirichletBlockVector;
84
85 class BoundaryWrapper
86 {
87 public:
88 typedef typename DiscreteFunctionSpaceType::EntityType EntityType;
90 typedef typename DiscreteFunctionSpaceType::DomainType DomainType;
91 typedef typename DiscreteFunctionSpace::RangeType RangeType;
92 typedef typename DiscreteFunctionSpace::JacobianRangeType JacobianRangeType;
93 typedef typename DiscreteFunctionSpace::HessianRangeType HessianRangeType;
94
95 private:
96 const ModelType& impl_;
97 const EntityType& entity_;
98 const int order_;
99 int bndId_;
100
101 public:
102 static const int dimRange = RangeType::dimension;
103 BoundaryWrapper( const ModelType& impl, const EntityType& entity, const int order, int bndId )
104 : impl_( impl ), entity_(entity), order_(order), bndId_(bndId) {}
105 bool valid() const { return true; }
106 const EntityType& entity() const { return entity_; }
107 const int order () const { return order_; }
108 template <class Point>
109 void evaluate( const Point& x, RangeType& ret ) const
110 { impl_.dirichlet(bndId_,Dune::Fem::coordinate(x),ret); }
111 template <class Point>
112 void jacobian( const Point& x, JacobianRangeType& ret ) const
113 { DUNE_THROW(Dune::NotImplemented,"rhs jacobian not implemented"); }
114 };
115
116 DirichletConstraints( ModelType &model, const DiscreteFunctionSpaceType& space )
117 : model_(model),
118 space_( space ),
119 dirichletBlocks_(),
120 // mark DoFs on the Dirichlet boundary
121 hasDirichletDofs_( false ),
122 sequence_( -1 )
123 {}
124
132 template < class DiscreteFunctionType >
133 void operator ()( const DiscreteFunctionType& u, DiscreteFunctionType& w ) const
134 {
135 updateDirichletDofs();
136
137 // if Dirichlet Dofs have been found, treat them
138 if( hasDirichletDofs_ )
139 {
140 typedef typename DiscreteFunctionType :: DofIteratorType DofIteratorType ;
141 typedef typename DiscreteFunctionType :: ConstDofIteratorType ConstDofIteratorType ;
142
143 ConstDofIteratorType uIt = u.dbegin();
144 DofIteratorType wIt = w.dbegin();
145
146 // loop over all blocks
147 const unsigned int blocks = space_.blockMapper().size();
148 for( unsigned int blockDof = 0; blockDof < blocks ; ++ blockDof )
149 {
150 for( int l = 0; l < localBlockSize ; ++ l, ++ wIt, ++ uIt )
151 {
152 if( dirichletBlocks_[ blockDof ][l] )
153 {
154 // copy dofs of the block
155 assert( uIt != u.dend() );
156 assert( wIt != w.dend() );
157 (*wIt) = (*uIt);
158 }
159 }
160 }
161 }
162 }
170 template < class DiscreteFunctionType >
171 void operator ()( const typename DiscreteFunctionType::RangeType& value, DiscreteFunctionType& w ) const
172 {
173 updateDirichletDofs();
174
175 // if Dirichlet Dofs have been found, treat them
176 if( hasDirichletDofs_ )
177 {
178 typedef typename DiscreteFunctionType :: DofIteratorType DofIteratorType ;
179
180 DofIteratorType wIt = w.dbegin();
181
182 // loop over all blocks
183 const unsigned int blocks = space_.blockMapper().size();
184 for( unsigned int blockDof = 0; blockDof < blocks ; ++ blockDof )
185 {
186 for( int l = 0; l < localBlockSize ; ++ l, ++ wIt )
187 {
188 if( dirichletBlocks_[ blockDof ][l] )
189 {
190 // copy dofs of the block
191 assert( wIt != w.dend() );
192 (*wIt) = value[l];
193 }
194 }
195 }
196 }
197 }
198
205 template < class DiscreteFunctionType >
206 void operator ()( DiscreteFunctionType& w, Operation op=Operation::set ) const
207 {
208 if (op != Operation::set)
209 {
211 return (*this)(u,w,op);
212 }
213
214 updateDirichletDofs();
215
216 if( hasDirichletDofs_ )
217 {
218 typedef typename DiscreteFunctionSpaceType :: IteratorType IteratorType;
219 typedef typename IteratorType :: Entity EntityType;
220
221 Dune::Fem::SetSelectedLocalContribution< DiscreteFunctionType > wLocal( w );
222 Dune::Fem::LocalInterpolation< DiscreteFunctionSpaceType > interpolation( space_ );
223 for( const EntityType &entity : space_ )
224 {
225 model_.init(entity);
226
227 // bind local contribution to entity
228 auto wGuard = Dune::Fem::bindGuard( wLocal, entity );
229 // bind interpolation to entity
230 auto iGuard = bindGuard( interpolation, entity );
231
232 // interpolate dirichlet dofs
233 dirichletDofTreatment( interpolation, wLocal );
234 }
235 }
236 }
237 template < class GridFunctionType, class DiscreteFunctionType,
238 typename = std::enable_if_t< std::is_base_of<Dune::Fem::HasLocalFunction, GridFunctionType>::value > >
239 void operator ()( const GridFunctionType &u,
240 DiscreteFunctionType& w, Operation op=Operation::set ) const
241 {
242 updateDirichletDofs();
243
244 if( hasDirichletDofs_ )
245 {
246 Dune::Fem::ConstLocalFunction< GridFunctionType > uLocal( u );
247 Dune::Fem::SetSelectedLocalContribution< DiscreteFunctionType > wLocal( w );
248 Dune::Fem::LocalInterpolation< DiscreteFunctionSpaceType > interpolation( space_ );
249
250 for( const auto &entity : space_ )
251 {
252 auto guard = Dune::Fem::bindGuard( std::tie(uLocal,wLocal), entity );
253 // bind interpolation to entity
254 auto iGuard = bindGuard( interpolation, entity );
255
256 // interpolate dirichlet dofs
257 if (op == Operation::sub || op == Operation::add)
258 model_.init(entity);
259 dirichletDofTreatment( interpolation, uLocal, wLocal, op );
260 }
261 }
262 }
263
272 template <class LinearOperator>
273 void applyToOperator( LinearOperator& linearOperator ) const
274 {
275 updateDirichletDofs();
276
277 typedef typename DiscreteFunctionSpaceType :: IteratorType IteratorType;
278 typedef typename IteratorType :: Entity EntityType;
279
280 // if Dirichlet Dofs have been found, treat them
281 if( hasDirichletDofs_ )
282 {
283 typedef typename DiscreteFunctionSpaceType :: BlockMapperType BlockMapperType;
285
286 std::vector<std::size_t> globalBlockDofs;
287 std::vector<std::size_t> globalDofs;
288
289 // storage for unit rows and auxiliary rows, should be unique and sorted
290 std::set<std::size_t> unitRows;
291 std::set<std::size_t> auxRows;
292
293 const auto& space = space_; // linearOperator.rangeSpace();
294 // get auxiliary dof structure (for parallel runs) /*@LST0S@*/
295 //const auto &auxiliaryDofs = space.auxiliaryDofs();
296
297 const IteratorType end = space_.end();
298 for( IteratorType it = space_.begin(); it != end; ++it )
299 {
300 const EntityType &entity = *it;
301
302 // get number of basis functions
303 const int localBlocks = space.blockMapper().numDofs( entity );
304
305 // map local to global dofs
306 globalBlockDofs.resize(localBlocks);
307 // obtain all DofBlocks for this element
308 space.blockMapper().map( entity, globalBlockDofs );
309
310 // obtain all non-blocked dofs
311 globalDofs.resize(localBlocks * localBlockSize);
312 mapper.map( entity, globalDofs );
313
314 // counter for all local dofs (i.e. localBlockDof * localBlockSize + ... )
315 int localDof = 0;
316 // iterate over face dofs and set unit row
317 for( int localBlockDof = 0 ; localBlockDof < localBlocks; ++ localBlockDof )
318 {
319 int global = globalBlockDofs[localBlockDof];
320 // all diagonals were set to one anyway independent of auxDof or not
321 // so don't need this distinction anymore
322 // std::set<std::size_t>& rows = auxiliaryDofs.contains( global ) ? auxRows : unitRows;
323 std::set<std::size_t>& rows = useIdentity ? unitRows : auxRows;
324 for( int l = 0; l < localBlockSize; ++ l, ++ localDof )
325 {
326 if( dirichletBlocks_[global][l] )
327 {
328 // push non-blocked dof
329 rows.insert( globalDofs[ localDof ] );
330 }
331 }
332 }
333 } // end for elements
334
335 // set unit and auxiliary rows at once
336 linearOperator.setUnitRows( unitRows, auxRows );
337 }
338 }
339
340 const DirichletBlockVector &dirichletBlocks() const
341 {
342 updateDirichletDofs();
343 return dirichletBlocks_;
344 }
345
346protected:
348 template< class LocalInterpolationType, class LocalFunctionType >
349 void dirichletDofTreatment( const LocalInterpolationType& interpolation, LocalFunctionType &wLocal) const
350 {
351 // get entity
352 const typename LocalFunctionType::EntityType &entity = wLocal.entity();
353
354 // get number of Lagrange Points
355 const int localBlocks = space_.blockMapper().numDofs( entity );
356
357 // map local to global BlockDofs
358 std::vector<std::size_t> globalBlockDofs(localBlocks);
359 space_.blockMapper().map( entity, globalBlockDofs );
360 std::vector<typename LocalFunctionType::RangeFieldType> values( localBlocks*localBlockSize );
361
362 int localDof = 0;
363
364 // iterate over face dofs and set unit row
365 for( int localBlock = 0 ; localBlock < localBlocks; ++ localBlock )
366 {
367 // store result to dof vector
368 int global = globalBlockDofs[localBlock];
369 for( int l = 0; l < localBlockSize ; ++ l, ++localDof )
370 {
371 if( dirichletBlocks_[ global ][l] )
372 {
373 interpolation( BoundaryWrapper(model_, entity, wLocal.order(), dirichletBlocks_[global][l]), values );
374 // store result
375 assert( (unsigned int)localDof < wLocal.size() );
376 wLocal[ localDof ] = values[ localDof ];
377 }
378 }
379 }
380 }
381
382 template< class LocalInterpolationType, class GridLocalFunctionType, class LocalFunctionType >
383 void dirichletDofTreatment( const LocalInterpolationType& interpolation,
384 const GridLocalFunctionType &uLocal,
385 LocalFunctionType &wLocal, Operation op ) const
386 {
387 // get entity
388 const typename LocalFunctionType::EntityType &entity = wLocal.entity();
389
390 // get number of Lagrange Points
391 const int localBlocks = space_.blockMapper().numDofs( entity );
392
393 typedef typename DiscreteFunctionSpaceType::BlockMapperType::GlobalKeyType GlobalKeyType;
394
395 // map local to global BlockDofs
396 std::vector< GlobalKeyType > globalBlockDofs(localBlocks);
397 space_.blockMapper().map(entity,globalBlockDofs);
398 std::vector<double> values( localBlocks*localBlockSize,0 );
399 std::vector<double> valuesModel( localBlocks*localBlockSize );
400
401 if constexpr ( useIdentity )
402 interpolation( uLocal, values );
403
404 int localDof = 0;
405
406 // iterate over face dofs and set unit row
407 for( int localBlock = 0 ; localBlock < localBlocks; ++ localBlock )
408 {
409 // store result to dof vector
410 int global = globalBlockDofs[localBlock];
411 for( int l = 0; l < localBlockSize ; ++ l, ++localDof )
412 {
413 if( dirichletBlocks_[ global ][l] )
414 {
415 if (op == Operation::sub)
416 {
417 interpolation(BoundaryWrapper(model_, entity, wLocal.order(), dirichletBlocks_[global][l]), valuesModel);
418 values[ localDof ] -= valuesModel[ localDof ];
419 }
420 else if (op == Operation::add)
421 {
422 interpolation(BoundaryWrapper(model_, entity, wLocal.order(), dirichletBlocks_[global][l]), valuesModel);
423 values[ localDof ] += valuesModel[ localDof ];
424 }
425 // store result
426 assert( (unsigned int)localDof < wLocal.size() );
427 wLocal[ localDof ] = values[ localDof ];
428 }
429 }
430 }
431 }
432
433protected:
434 // detect all DoFs on the Dirichlet boundary
435 void updateDirichletDofs() const
436 {
437 if( sequence_ != space_.sequence() )
438 {
439 // only start search if Dirichlet boundary is present
440 if( ! model_.hasDirichletBoundary() )
441 {
442 hasDirichletDofs_ = false ;
443 return ;
444 }
445
446 // resize flag vector with number of blocks and reset flags
447 const int blocks = space_.blockMapper().size() ;
448 dirichletBlocks_.resize( blocks );
449 for( int i=0; i<blocks; ++i )
450 dirichletBlocks_[ i ].fill(0);
451
452 typedef typename DiscreteFunctionSpaceType :: IteratorType IteratorType;
453 typedef typename IteratorType :: Entity EntityType;
454
455 bool hasDirichletBoundary = false;
456 const IteratorType end = space_.end();
457 for( IteratorType it = space_.begin(); it != end; ++it )
458 {
459 const EntityType &entity = *it;
460 // if entity has boundary intersections
461 if( Dune::Fem::HasBoundaryIntersection<GridPartType>::apply(entity) )
462 {
463 hasDirichletBoundary |= searchEntityDirichletDofs( entity, model_ );
464 }
465 }
466
467 // update sequence number
468 sequence_ = space_.sequence();
469 if( space_.gridPart().comm().size() > 1 )
470 {
471 try
472 {
473 DirichletBuilder handle( *this, space_ , space_.blockMapper() );
474 space_.gridPart().communicate
475 ( handle, GridPartType::indexSetInterfaceType, ForwardCommunication );
476 }
477 // catch possible exceptions here to have a clue where it happend
478 catch( const Exception &e )
479 {
480 std::cerr << e << std::endl;
481 std::cerr << "Exception thrown in: " << __FILE__ << " line:" << __LINE__ << std::endl;
482 abort();
483 }
484 hasDirichletDofs_ = space_.gridPart().grid().comm().max( hasDirichletBoundary );
485 }
486 else
487 {
488 hasDirichletDofs_ = hasDirichletBoundary;
489 }
490 }
491 }
492
493 // detect all DoFs on the Dirichlet boundary of the given entity
494 template< class EntityType >
495 bool searchEntityDirichletDofs( const EntityType &entity, ModelType& model ) const
496 {
497 typedef typename DiscreteFunctionSpaceType :: GridPartType GridPartType;
498
499 typedef typename GridPartType :: IntersectionIteratorType
500 IntersectionIteratorType;
501
502 const GridPartType &gridPart = space_.gridPart();
503
504 // default is false
505 bool hasDirichletBoundary = false;
506
507 //map local to global BlockDofs
508 std::vector<size_t> globalBlockDofs(space_.blockMapper().numDofs(entity));
509 space_.blockMapper().map(entity,globalBlockDofs);
510
511 std::vector<char> subEntityFilter(space_.blockMapper().numDofs(entity));
512
513 IntersectionIteratorType it = gridPart.ibegin( entity );
514 const IntersectionIteratorType endit = gridPart.iend( entity );
515 for( ; it != endit; ++it )
516 {
517 typedef typename IntersectionIteratorType :: Intersection IntersectionType;
518 const IntersectionType& intersection = *it;
519
520 // if intersection is with boundary, adjust data
521 if( intersection.boundary() )
522 {
523 // get dirichlet information from model
524 DirichletBlock block;
525 block.fill(0);
526 std::array<int,dimRange> dblock;
527
528 model.init(intersection);
529 const bool isDirichletIntersection = model.isDirichletIntersection( intersection, dblock );
530 if (isDirichletIntersection)
531 {
532 // get face number of boundary intersection
533 const int face = intersection.indexInInside();
534 space_.blockMapper().onSubEntity(entity,face,1,subEntityFilter);
535 for( unsigned int i=0;i<globalBlockDofs.size();++i)
536 {
537 assert( i < subEntityFilter.size());
538 if ( subEntityFilter[i]==0) continue;
539 // mark global DoF number
540 for(int k = 0; k < localBlockSize; ++k)
541 {
542 int comp = subEntityFilter[i]-1+k;
543 assert( comp < dimRange );
544 dirichletBlocks_[ globalBlockDofs[ i ] ][k] = dblock[comp];
545 }
546 // we have Dirichlet values
547 hasDirichletBoundary = true ;
548 }
549 }
550 }
551 }
552
553 return hasDirichletBoundary;
554 }
555
556 ModelType &model_;
557 const DiscreteFunctionSpaceType& space_;
558 mutable DirichletBlockVector dirichletBlocks_;
559 mutable bool hasDirichletDofs_ ;
560 mutable int sequence_ ;
561
562 class DirichletBuilder;
563};
564
565template < class Model, class Space, bool useIdentity >
566class DirichletConstraints< Model,Space, useIdentity > :: DirichletBuilder
567 : public CommDataHandleIF< DirichletBuilder, int >
568{
569public:
570 typedef Space SpaceType;
571 typedef typename SpaceType::BlockMapperType MapperType;
572
573 enum { nCodim = SpaceType :: GridType :: dimension + 1 };
574
575public:
576 typedef int DataType;
577
578 const int myRank_;
579 const int mySize_;
580
581 typedef DirichletConstraints< Model,Space, useIdentity > DirichletType;
582 const DirichletType &dirichlet_;
583
584 const SpaceType &space_;
585 const MapperType &mapper_;
586
587 static const int blockSize = SpaceType::localBlockSize;
588
589public:
590 DirichletBuilder( const DirichletType &dirichlet,
591 const SpaceType &space,
592 const MapperType &mapper )
593 : myRank_( space.gridPart().comm().rank() ),
594 mySize_( space.gridPart().comm().size() ),
595 dirichlet_( dirichlet ),
596 space_( space ),
597 mapper_( mapper )
598 {
599 }
600 bool contains ( int dim, int codim ) const
601 {
602 return mapper_.contains( codim );
603 }
604
605 bool fixedSize ( int dim, int codim ) const
606 {
607 return false;
608 }
609
611 template< class MessageBuffer, class Entity >
612 inline void gather ( MessageBuffer &buffer,
613 const Entity &entity ) const
614 {
615 unsigned int localBlocks = mapper_.numEntityDofs( entity );
616 std::vector<std::size_t> globalBlockDofs(localBlocks);
617 mapper_.mapEntityDofs( entity, globalBlockDofs );
618 assert( size(entity) == globalBlockDofs.size()*blockSize );
619 for( unsigned int localBlock = 0 ; localBlock < globalBlockDofs.size(); ++ localBlock )
620 {
621 int global = globalBlockDofs[localBlock];
622 for (int r=0;r<blockSize;++r)
623 if (dirichlet_.dirichletBlocks_[ global ][r] )
624 buffer.write( 1 );
625 else
626 buffer.write( 0 );
627 }
628 }
629
634 template< class MessageBuffer, class EntityType >
635 inline void scatter ( MessageBuffer &buffer,
636 const EntityType &entity,
637 size_t n )
638 {
639 unsigned int localBlocks = mapper_.numEntityDofs( entity );
640 std::vector<std::size_t> globalBlockDofs(localBlocks);
641 mapper_.mapEntityDofs( entity, globalBlockDofs );
642 assert( n == globalBlockDofs.size()*blockSize );
643 assert( n == size(entity) );
644 for( unsigned int localBlock = 0 ; localBlock < globalBlockDofs.size(); ++ localBlock )
645 {
646 int global = globalBlockDofs[localBlock];
647 for (int r=0;r<blockSize;++r)
648 {
649 int val;
650 buffer.read(val);
651 if ( !dirichlet_.dirichletBlocks_[ global ][r] && val == 1)
652 dirichlet_.dirichletBlocks_[ global ][r] = true;
653 }
654 }
655 }
657 template< class Entity >
658 size_t size ( const Entity &entity ) const
659 {
660 return blockSize * mapper_.numEntityDofs( entity );
661 }
662};
663
664} // end namespace Dune
665#endif
discrete function space
Traits::ConstDofIteratorType ConstDofIteratorType
type of the const dof iterator
Definition: discretefunction.hh:628
ConstDofIteratorType dbegin() const
Obtain the constant iterator pointing to the first dof.
Definition: discretefunction.hh:761
ConstDofIteratorType dend() const
Obtain the constant iterator pointing to the last dof.
Definition: discretefunction.hh:773
DiscreteFunctionSpaceType::RangeType RangeType
type of range vector
Definition: discretefunction.hh:614
Traits::DofIteratorType DofIteratorType
type of the dof iterator
Definition: discretefunction.hh:626
Traits::FunctionSpaceType FunctionSpaceType
type of function space
Definition: discretefunctionspace.hh:194
A vector valued function space.
Definition: functionspace.hh:60
Definition: nonblockmapper.hh:284
forward declaration
Definition: discretefunction.hh:51
Default exception for dummy implementations.
Definition: exceptions.hh:357
#define DUNE_THROW(E,...)
Definition: exceptions.hh:314
@ ForwardCommunication
communicate as given in InterfaceType
Definition: gridenums.hh:171
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
constexpr std::bool_constant<((II==value)||...)> contains(std::integer_sequence< T, II... >, std::integral_constant< T, value >)
Checks whether or not a given sequence contains a value.
Definition: integersequence.hh:137
std::size_t fixedSize
The number of data items per index if it is fixed, 0 otherwise.
Definition: variablesizecommunicator.hh:264
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Sep 4, 22:38, 2025)