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