Loading [MathJax]/extensions/tex2jax.js

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
52namespace Dune {
53
54namespace Fem {
55class HasLocalFunction;
56}
57
58template < class Model, class DiscreteFunctionSpace, bool useIdentity >
59class DirichletConstraints
60{
61public:
62 enum Operation { set = 0, sub = 1, add = 2 };
63 typedef Model ModelType;
64 typedef DiscreteFunctionSpace DiscreteFunctionSpaceType;
65 typedef typename DiscreteFunctionSpaceType::DomainType DomainType;
66 typedef typename DiscreteFunctionSpaceType::RangeType RangeType;
67 typedef typename DiscreteFunctionSpaceType::JacobianRangeType JacobianRangeType;
68 typedef typename DiscreteFunctionSpaceType::HessianRangeType HessianRangeType;
69
71 typedef typename DiscreteFunctionSpaceType :: GridPartType GridPartType;
73 typedef typename DiscreteFunctionSpaceType :: GridType GridType;
74
75 // types for boundary treatment
76 // ----------------------------
77 static const int localBlockSize = DiscreteFunctionSpaceType :: localBlockSize ;
78 static const int dimRange = DiscreteFunctionSpaceType::FunctionSpaceType::dimRange;
79 static_assert( localBlockSize <= dimRange,
80 "local block size of the space must be less than or equal to the dimension of the range of the function space.");
81 typedef std::array<int,localBlockSize> DirichletBlock;
82 typedef std::vector< DirichletBlock > DirichletBlockVector;
83
84 class BoundaryWrapper
85 {
86 public:
87 typedef typename DiscreteFunctionSpaceType::EntityType EntityType;
89 typedef typename DiscreteFunctionSpaceType::DomainType DomainType;
90 typedef typename DiscreteFunctionSpace::RangeType RangeType;
91 typedef typename DiscreteFunctionSpace::JacobianRangeType JacobianRangeType;
92 typedef typename DiscreteFunctionSpace::HessianRangeType HessianRangeType;
93
94 private:
95 const ModelType& impl_;
96 const EntityType& entity_;
97 const int order_;
98 int bndId_;
99
100 public:
101 static const int dimRange = RangeType::dimension;
102 BoundaryWrapper( const ModelType& impl, const EntityType& entity, const int order, int bndId )
103 : impl_( impl ), entity_(entity), order_(order), bndId_(bndId) {}
104 bool valid() const { return true; }
105 const EntityType& entity() const { return entity_; }
106 const int order () const { return order_; }
107 template <class Point>
108 void evaluate( const Point& x, RangeType& ret ) const
109 { impl_.dirichlet(bndId_,Dune::Fem::coordinate(x),ret); }
110 template <class Point>
111 void jacobian( const Point& x, JacobianRangeType& ret ) const
112 { DUNE_THROW(Dune::NotImplemented,"rhs jacobian not implemented"); }
113 };
114
115 DirichletConstraints( ModelType &model, const DiscreteFunctionSpaceType& space )
116 : model_(model),
117 space_( space ),
118 dirichletBlocks_(),
119 // mark DoFs on the Dirichlet boundary
120 hasDirichletDofs_( false ),
121 sequence_( -1 )
122 {}
123
131 template < class DiscreteFunctionType >
132 void operator ()( const DiscreteFunctionType& u, DiscreteFunctionType& w ) const
133 {
134 updateDirichletDofs();
135
136 // if Dirichlet Dofs have been found, treat them
137 if( hasDirichletDofs_ )
138 {
139 typedef typename DiscreteFunctionType :: DofIteratorType DofIteratorType ;
140 typedef typename DiscreteFunctionType :: ConstDofIteratorType ConstDofIteratorType ;
141
142 ConstDofIteratorType uIt = u.dbegin();
143 DofIteratorType wIt = w.dbegin();
144
145 // loop over all blocks
146 const unsigned int blocks = space_.blockMapper().size();
147 for( unsigned int blockDof = 0; blockDof < blocks ; ++ blockDof )
148 {
149 for( int l = 0; l < localBlockSize ; ++ l, ++ wIt, ++ uIt )
150 {
151 if( dirichletBlocks_[ blockDof ][l] )
152 {
153 // copy dofs of the block
154 assert( uIt != u.dend() );
155 assert( wIt != w.dend() );
156 (*wIt) = (*uIt);
157 }
158 }
159 }
160 }
161 }
169 template < class DiscreteFunctionType >
170 void operator ()( const typename DiscreteFunctionType::RangeType& value, DiscreteFunctionType& w ) const
171 {
172 updateDirichletDofs();
173
174 // if Dirichlet Dofs have been found, treat them
175 if( hasDirichletDofs_ )
176 {
177 typedef typename DiscreteFunctionType :: DofIteratorType DofIteratorType ;
178
179 DofIteratorType wIt = w.dbegin();
180
181 // loop over all blocks
182 const unsigned int blocks = space_.blockMapper().size();
183 for( unsigned int blockDof = 0; blockDof < blocks ; ++ blockDof )
184 {
185 for( int l = 0; l < localBlockSize ; ++ l, ++ wIt )
186 {
187 if( dirichletBlocks_[ blockDof ][l] )
188 {
189 // copy dofs of the block
190 assert( wIt != w.dend() );
191 (*wIt) = value[l];
192 }
193 }
194 }
195 }
196 }
197
204 template < class DiscreteFunctionType >
205 void operator ()( DiscreteFunctionType& w, Operation op=Operation::set ) const
206 {
207 if (op != Operation::set)
208 {
210 return (*this)(u,w,op);
211 }
212
213 updateDirichletDofs();
214
215 if( hasDirichletDofs_ )
216 {
217 typedef typename DiscreteFunctionSpaceType :: IteratorType IteratorType;
218 typedef typename IteratorType :: Entity EntityType;
219
220 Dune::Fem::SetSelectedLocalContribution< DiscreteFunctionType > wLocal( w );
221 Dune::Fem::LocalInterpolation< DiscreteFunctionSpaceType > interpolation( space_ );
222 for( const EntityType &entity : space_ )
223 {
224 model_.init(entity);
225
226 // bind local contribution to entity
227 auto wGuard = Dune::Fem::bindGuard( wLocal, entity );
228 // bind interpolation to entity
229 auto iGuard = bindGuard( interpolation, entity );
230
231 // interpolate dirichlet dofs
232 dirichletDofTreatment( interpolation, wLocal );
233 }
234 }
235 }
236 template < class GridFunctionType, class DiscreteFunctionType,
237 typename = std::enable_if_t< std::is_base_of<Dune::Fem::HasLocalFunction, GridFunctionType>::value > >
238 void operator ()( const GridFunctionType &u,
239 DiscreteFunctionType& w, Operation op=Operation::set ) const
240 {
241 updateDirichletDofs();
242
243 if( hasDirichletDofs_ )
244 {
245 Dune::Fem::ConstLocalFunction< GridFunctionType > uLocal( u );
246 Dune::Fem::SetSelectedLocalContribution< DiscreteFunctionType > wLocal( w );
247 Dune::Fem::LocalInterpolation< DiscreteFunctionSpaceType > interpolation( space_ );
248
249 for( const auto &entity : space_ )
250 {
251 auto guard = Dune::Fem::bindGuard( std::tie(uLocal,wLocal), entity );
252 // bind interpolation to entity
253 auto iGuard = bindGuard( interpolation, entity );
254
255 // interpolate dirichlet dofs
256 if (op == Operation::sub || op == Operation::add)
257 model_.init(entity);
258 dirichletDofTreatment( interpolation, uLocal, wLocal, op );
259 }
260 }
261 }
262
271 template <class LinearOperator>
272 void applyToOperator( LinearOperator& linearOperator ) const
273 {
274 updateDirichletDofs();
275
276 typedef typename DiscreteFunctionSpaceType :: IteratorType IteratorType;
277 typedef typename IteratorType :: Entity EntityType;
278
279 // if Dirichlet Dofs have been found, treat them
280 if( hasDirichletDofs_ )
281 {
282 typedef typename DiscreteFunctionSpaceType :: BlockMapperType BlockMapperType;
284
285 std::vector<std::size_t> globalBlockDofs;
286 std::vector<std::size_t> globalDofs;
287
288 // storage for unit rows and auxiliary rows, should be unique and sorted
289 std::set<std::size_t> unitRows;
290 std::set<std::size_t> auxRows;
291
292 const auto& space = space_; // linearOperator.rangeSpace();
293 // get auxiliary dof structure (for parallel runs) /*@LST0S@*/
294 //const auto &auxiliaryDofs = space.auxiliaryDofs();
295
296 const IteratorType end = space_.end();
297 for( IteratorType it = space_.begin(); it != end; ++it )
298 {
299 const EntityType &entity = *it;
300
301 // get number of basis functions
302 const int localBlocks = space.blockMapper().numDofs( entity );
303
304 // map local to global dofs
305 globalBlockDofs.resize(localBlocks);
306 // obtain all DofBlocks for this element
307 space.blockMapper().map( entity, globalBlockDofs );
308
309 // obtain all non-blocked dofs
310 globalDofs.resize(localBlocks * localBlockSize);
311 mapper.map( entity, globalDofs );
312
313 // counter for all local dofs (i.e. localBlockDof * localBlockSize + ... )
314 int localDof = 0;
315 // iterate over face dofs and set unit row
316 for( int localBlockDof = 0 ; localBlockDof < localBlocks; ++ localBlockDof )
317 {
318 int global = globalBlockDofs[localBlockDof];
319 // all diagonals were set to one anyway independent of auxDof or not
320 // so don't need this distinction anymore
321 // std::set<std::size_t>& rows = auxiliaryDofs.contains( global ) ? auxRows : unitRows;
322 std::set<std::size_t>& rows = useIdentity ? unitRows : auxRows;
323 for( int l = 0; l < localBlockSize; ++ l, ++ localDof )
324 {
325 if( dirichletBlocks_[global][l] )
326 {
327 // push non-blocked dof
328 rows.insert( globalDofs[ localDof ] );
329 }
330 }
331 }
332 } // end for elements
333
334 // set unit and auxiliary rows at once
335 linearOperator.setUnitRows( unitRows, auxRows );
336 }
337 }
338
339 const DirichletBlockVector &dirichletBlocks() const
340 {
341 updateDirichletDofs();
342 return dirichletBlocks_;
343 }
344
345protected:
347 template< class LocalInterpolationType, class LocalFunctionType >
348 void dirichletDofTreatment( const LocalInterpolationType& interpolation, LocalFunctionType &wLocal) const
349 {
350 // get entity
351 const typename LocalFunctionType::EntityType &entity = wLocal.entity();
352
353 // get number of Lagrange Points
354 const int localBlocks = space_.blockMapper().numDofs( entity );
355
356 // map local to global BlockDofs
357 std::vector<std::size_t> globalBlockDofs(localBlocks);
358 space_.blockMapper().map( entity, globalBlockDofs );
359 std::vector<typename LocalFunctionType::RangeFieldType> values( localBlocks*localBlockSize );
360
361 int localDof = 0;
362
363 // iterate over face dofs and set unit row
364 for( int localBlock = 0 ; localBlock < localBlocks; ++ localBlock )
365 {
366 // store result to dof vector
367 int global = globalBlockDofs[localBlock];
368 for( int l = 0; l < localBlockSize ; ++ l, ++localDof )
369 {
370 if( dirichletBlocks_[ global ][l] )
371 {
372 interpolation( BoundaryWrapper(model_, entity, wLocal.order(), dirichletBlocks_[global][l]), values );
373 // store result
374 assert( (unsigned int)localDof < wLocal.size() );
375 wLocal[ localDof ] = values[ localDof ];
376 }
377 }
378 }
379 }
380
381 template< class LocalInterpolationType, class GridLocalFunctionType, class LocalFunctionType >
382 void dirichletDofTreatment( const LocalInterpolationType& interpolation,
383 const GridLocalFunctionType &uLocal,
384 LocalFunctionType &wLocal, Operation op ) const
385 {
386 // get entity
387 const typename LocalFunctionType::EntityType &entity = wLocal.entity();
388
389 // get number of Lagrange Points
390 const int localBlocks = space_.blockMapper().numDofs( entity );
391
392 typedef typename DiscreteFunctionSpaceType::BlockMapperType::GlobalKeyType GlobalKeyType;
393
394 // map local to global BlockDofs
395 std::vector< GlobalKeyType > globalBlockDofs(localBlocks);
396 space_.blockMapper().map(entity,globalBlockDofs);
397 std::vector<double> values( localBlocks*localBlockSize,0 );
398 std::vector<double> valuesModel( localBlocks*localBlockSize );
399
400 if constexpr ( useIdentity )
401 interpolation( uLocal, values );
402
403 int localDof = 0;
404
405 // iterate over face dofs and set unit row
406 for( int localBlock = 0 ; localBlock < localBlocks; ++ localBlock )
407 {
408 // store result to dof vector
409 int global = globalBlockDofs[localBlock];
410 for( int l = 0; l < localBlockSize ; ++ l, ++localDof )
411 {
412 if( dirichletBlocks_[ global ][l] )
413 {
414 if (op == Operation::sub)
415 {
416 interpolation(BoundaryWrapper(model_, entity, wLocal.order(), dirichletBlocks_[global][l]), valuesModel);
417 values[ localDof ] -= valuesModel[ localDof ];
418 }
419 else if (op == Operation::add)
420 {
421 interpolation(BoundaryWrapper(model_, entity, wLocal.order(), dirichletBlocks_[global][l]), valuesModel);
422 values[ localDof ] += valuesModel[ localDof ];
423 }
424 // store result
425 assert( (unsigned int)localDof < wLocal.size() );
426 wLocal[ localDof ] = values[ localDof ];
427 }
428 }
429 }
430 }
431
432protected:
433 // detect all DoFs on the Dirichlet boundary
434 void updateDirichletDofs() const
435 {
436 if( sequence_ != space_.sequence() )
437 {
438 // only start search if Dirichlet boundary is present
439 if( ! model_.hasDirichletBoundary() )
440 {
441 hasDirichletDofs_ = false ;
442 return ;
443 }
444
445 // resize flag vector with number of blocks and reset flags
446 const int blocks = space_.blockMapper().size() ;
447 dirichletBlocks_.resize( blocks );
448 for( int i=0; i<blocks; ++i )
449 dirichletBlocks_[ i ].fill(0);
450
451 typedef typename DiscreteFunctionSpaceType :: IteratorType IteratorType;
452 typedef typename IteratorType :: Entity EntityType;
453
454 bool hasDirichletBoundary = false;
455 const IteratorType end = space_.end();
456 for( IteratorType it = space_.begin(); it != end; ++it )
457 {
458 const EntityType &entity = *it;
459 // if entity has boundary intersections
460 if( entity.hasBoundaryIntersections() )
461 {
462 hasDirichletBoundary |= searchEntityDirichletDofs( entity, model_ );
463 }
464 }
465
466 // update sequence number
467 sequence_ = space_.sequence();
468 if( space_.gridPart().comm().size() > 1 )
469 {
470 try
471 {
472 DirichletBuilder handle( *this, space_ , space_.blockMapper() );
473 space_.gridPart().communicate
474 ( handle, GridPartType::indexSetInterfaceType, ForwardCommunication );
475 }
476 // catch possible exceptions here to have a clue where it happend
477 catch( const Exception &e )
478 {
479 std::cerr << e << std::endl;
480 std::cerr << "Exception thrown in: " << __FILE__ << " line:" << __LINE__ << std::endl;
481 abort();
482 }
483 hasDirichletDofs_ = space_.gridPart().grid().comm().max( hasDirichletBoundary );
484 }
485 else
486 {
487 hasDirichletDofs_ = hasDirichletBoundary;
488 }
489 }
490 }
491
492 // detect all DoFs on the Dirichlet boundary of the given entity
493 template< class EntityType >
494 bool searchEntityDirichletDofs( const EntityType &entity, ModelType& model ) const
495 {
496 typedef typename DiscreteFunctionSpaceType :: GridPartType GridPartType;
497
498 typedef typename GridPartType :: IntersectionIteratorType
499 IntersectionIteratorType;
500
501 const GridPartType &gridPart = space_.gridPart();
502
503 // default is false
504 bool hasDirichletBoundary = false;
505
506 //map local to global BlockDofs
507 std::vector<size_t> globalBlockDofs(space_.blockMapper().numDofs(entity));
508 space_.blockMapper().map(entity,globalBlockDofs);
509
510 std::vector<char> subEntityFilter(space_.blockMapper().numDofs(entity));
511
512 IntersectionIteratorType it = gridPart.ibegin( entity );
513 const IntersectionIteratorType endit = gridPart.iend( entity );
514 for( ; it != endit; ++it )
515 {
516 typedef typename IntersectionIteratorType :: Intersection IntersectionType;
517 const IntersectionType& intersection = *it;
518
519 // if intersection is with boundary, adjust data
520 if( intersection.boundary() )
521 {
522 // get dirichlet information from model
523 DirichletBlock block;
524 block.fill(0);
525 std::array<int,dimRange> dblock;
526
527 model.init(intersection);
528 const bool isDirichletIntersection = model.isDirichletIntersection( intersection, dblock );
529 if (isDirichletIntersection)
530 {
531 // get face number of boundary intersection
532 const int face = intersection.indexInInside();
533 space_.blockMapper().onSubEntity(entity,face,1,subEntityFilter);
534 for( unsigned int i=0;i<globalBlockDofs.size();++i)
535 {
536 assert( i < subEntityFilter.size());
537 if ( subEntityFilter[i]==0) continue;
538 // mark global DoF number
539 for(int k = 0; k < localBlockSize; ++k)
540 {
541 int comp = subEntityFilter[i]-1+k;
542 assert( comp < dimRange );
543 dirichletBlocks_[ globalBlockDofs[ i ] ][k] = dblock[comp];
544 }
545 // we have Dirichlet values
546 hasDirichletBoundary = true ;
547 }
548 }
549 }
550 }
551
552 return hasDirichletBoundary;
553 }
554
555 ModelType &model_;
556 const DiscreteFunctionSpaceType& space_;
557 mutable DirichletBlockVector dirichletBlocks_;
558 mutable bool hasDirichletDofs_ ;
559 mutable int sequence_ ;
560
561 class DirichletBuilder;
562};
563
564template < class Model, class Space, bool useIdentity >
565class DirichletConstraints< Model,Space, useIdentity > :: DirichletBuilder
566 : public CommDataHandleIF< DirichletBuilder, int >
567{
568public:
569 typedef Space SpaceType;
570 typedef typename SpaceType::BlockMapperType MapperType;
571
572 enum { nCodim = SpaceType :: GridType :: dimension + 1 };
573
574public:
575 typedef int DataType;
576
577 const int myRank_;
578 const int mySize_;
579
580 typedef DirichletConstraints< Model,Space, useIdentity > DirichletType;
581 const DirichletType &dirichlet_;
582
583 const SpaceType &space_;
584 const MapperType &mapper_;
585
586 static const int blockSize = SpaceType::localBlockSize;
587
588public:
589 DirichletBuilder( const DirichletType &dirichlet,
590 const SpaceType &space,
591 const MapperType &mapper )
592 : myRank_( space.gridPart().comm().rank() ),
593 mySize_( space.gridPart().comm().size() ),
594 dirichlet_( dirichlet ),
595 space_( space ),
596 mapper_( mapper )
597 {
598 }
599 bool contains ( int dim, int codim ) const
600 {
601 return mapper_.contains( codim );
602 }
603
604 bool fixedSize ( int dim, int codim ) const
605 {
606 return false;
607 }
608
610 template< class MessageBuffer, class Entity >
611 inline void gather ( MessageBuffer &buffer,
612 const Entity &entity ) const
613 {
614 unsigned int localBlocks = mapper_.numEntityDofs( entity );
615 std::vector<std::size_t> globalBlockDofs(localBlocks);
616 mapper_.mapEntityDofs( entity, globalBlockDofs );
617 assert( size(entity) == globalBlockDofs.size()*blockSize );
618 for( unsigned int localBlock = 0 ; localBlock < globalBlockDofs.size(); ++ localBlock )
619 {
620 int global = globalBlockDofs[localBlock];
621 for (int r=0;r<blockSize;++r)
622 if (dirichlet_.dirichletBlocks_[ global ][r] )
623 buffer.write( 1 );
624 else
625 buffer.write( 0 );
626 }
627 }
628
633 template< class MessageBuffer, class EntityType >
634 inline void scatter ( MessageBuffer &buffer,
635 const EntityType &entity,
636 size_t n )
637 {
638 unsigned int localBlocks = mapper_.numEntityDofs( entity );
639 std::vector<std::size_t> globalBlockDofs(localBlocks);
640 mapper_.mapEntityDofs( entity, globalBlockDofs );
641 assert( n == globalBlockDofs.size()*blockSize );
642 assert( n == size(entity) );
643 for( unsigned int localBlock = 0 ; localBlock < globalBlockDofs.size(); ++ localBlock )
644 {
645 int global = globalBlockDofs[localBlock];
646 for (int r=0;r<blockSize;++r)
647 {
648 int val;
649 buffer.read(val);
650 if ( !dirichlet_.dirichletBlocks_[ global ][r] && val == 1)
651 dirichlet_.dirichletBlocks_[ global ][r] = true;
652 }
653 }
654 }
656 template< class Entity >
657 size_t size ( const Entity &entity ) const
658 {
659 return blockSize * mapper_.numEntityDofs( entity );
660 }
661};
662
663} // end namespace Dune
664#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 (Apr 24, 22:34, 2025)