40#ifndef DUNE_DIRICHLETCONSTRAINTS_HH
41#define DUNE_DIRICHLETCONSTRAINTS_HH
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>
56class HasLocalFunction;
59template <
class Model,
class DiscreteFunctionSpace,
bool useIdentity >
60class DirichletConstraints
63 enum Operation { set = 0, sub = 1, add = 2 };
64 typedef Model ModelType;
66 typedef typename DiscreteFunctionSpaceType::DomainType DomainType;
67 typedef typename DiscreteFunctionSpaceType::RangeType RangeType;
68 typedef typename DiscreteFunctionSpaceType::JacobianRangeType JacobianRangeType;
69 typedef typename DiscreteFunctionSpaceType::HessianRangeType HessianRangeType;
72 typedef typename DiscreteFunctionSpaceType :: GridPartType GridPartType;
74 typedef typename DiscreteFunctionSpaceType :: GridType GridType;
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;
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;
96 const ModelType& impl_;
97 const EntityType& entity_;
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 Po
int>
109 void evaluate(
const Point& x, RangeType& ret )
const
110 { impl_.dirichlet(bndId_,Dune::Fem::coordinate(x),ret); }
111 template <
class Po
int>
112 void jacobian(
const Point& x, JacobianRangeType& ret )
const
116 DirichletConstraints( ModelType &model,
const DiscreteFunctionSpaceType& space )
121 hasDirichletDofs_( false ),
132 template <
class DiscreteFunctionType >
135 updateDirichletDofs();
138 if( hasDirichletDofs_ )
143 ConstDofIteratorType uIt = u.
dbegin();
144 DofIteratorType wIt = w.
dbegin();
147 const unsigned int blocks = space_.blockMapper().size();
148 for(
unsigned int blockDof = 0; blockDof < blocks ; ++ blockDof )
150 for(
int l = 0; l < localBlockSize ; ++ l, ++ wIt, ++ uIt )
152 if( dirichletBlocks_[ blockDof ][l] )
155 assert( uIt != u.
dend() );
156 assert( wIt != w.
dend() );
170 template <
class DiscreteFunctionType >
173 updateDirichletDofs();
176 if( hasDirichletDofs_ )
180 DofIteratorType wIt = w.
dbegin();
183 const unsigned int blocks = space_.blockMapper().size();
184 for(
unsigned int blockDof = 0; blockDof < blocks ; ++ blockDof )
186 for(
int l = 0; l < localBlockSize ; ++ l, ++ wIt )
188 if( dirichletBlocks_[ blockDof ][l] )
191 assert( wIt != w.
dend() );
205 template <
class DiscreteFunctionType >
208 if (op != Operation::set)
211 return (*
this)(u,w,op);
214 updateDirichletDofs();
216 if( hasDirichletDofs_ )
218 typedef typename DiscreteFunctionSpaceType :: IteratorType IteratorType;
219 typedef typename IteratorType :: Entity EntityType;
221 Dune::Fem::SetSelectedLocalContribution< DiscreteFunctionType > wLocal( w );
222 Dune::Fem::LocalInterpolation< DiscreteFunctionSpaceType > interpolation( space_ );
223 for(
const EntityType &entity : space_ )
228 auto wGuard = Dune::Fem::bindGuard( wLocal, entity );
230 auto iGuard = bindGuard( interpolation, entity );
233 dirichletDofTreatment( interpolation, wLocal );
238 typename = std::enable_if_t< std::is_base_of<Dune::Fem::HasLocalFunction, GridFunctionType>::value > >
239 void operator ()(
const GridFunctionType &u,
242 updateDirichletDofs();
244 if( hasDirichletDofs_ )
246 Dune::Fem::ConstLocalFunction< GridFunctionType > uLocal( u );
247 Dune::Fem::SetSelectedLocalContribution< DiscreteFunctionType > wLocal( w );
248 Dune::Fem::LocalInterpolation< DiscreteFunctionSpaceType > interpolation( space_ );
250 for(
const auto &entity : space_ )
252 auto guard = Dune::Fem::bindGuard( std::tie(uLocal,wLocal), entity );
254 auto iGuard = bindGuard( interpolation, entity );
257 if (op == Operation::sub || op == Operation::add)
259 dirichletDofTreatment( interpolation, uLocal, wLocal, op );
272 template <
class LinearOperator>
273 void applyToOperator( LinearOperator& linearOperator )
const
275 updateDirichletDofs();
277 typedef typename DiscreteFunctionSpaceType :: IteratorType IteratorType;
278 typedef typename IteratorType :: Entity EntityType;
281 if( hasDirichletDofs_ )
283 typedef typename DiscreteFunctionSpaceType :: BlockMapperType BlockMapperType;
286 std::vector<std::size_t> globalBlockDofs;
287 std::vector<std::size_t> globalDofs;
290 std::set<std::size_t> unitRows;
291 std::set<std::size_t> auxRows;
293 const auto& space = space_;
297 const IteratorType end = space_.end();
298 for( IteratorType it = space_.begin(); it != end; ++it )
300 const EntityType &entity = *it;
303 const int localBlocks = space.blockMapper().numDofs( entity );
306 globalBlockDofs.resize(localBlocks);
308 space.blockMapper().map( entity, globalBlockDofs );
311 globalDofs.resize(localBlocks * localBlockSize);
312 mapper.map( entity, globalDofs );
317 for(
int localBlockDof = 0 ; localBlockDof < localBlocks; ++ localBlockDof )
319 int global = globalBlockDofs[localBlockDof];
323 std::set<std::size_t>& rows = useIdentity ? unitRows : auxRows;
324 for(
int l = 0; l < localBlockSize; ++ l, ++ localDof )
326 if( dirichletBlocks_[global][l] )
329 rows.insert( globalDofs[ localDof ] );
336 linearOperator.setUnitRows( unitRows, auxRows );
340 const DirichletBlockVector &dirichletBlocks()
const
342 updateDirichletDofs();
343 return dirichletBlocks_;
348 template<
class LocalInterpolationType,
class LocalFunctionType >
349 void dirichletDofTreatment(
const LocalInterpolationType& interpolation, LocalFunctionType &wLocal)
const
352 const typename LocalFunctionType::EntityType &entity = wLocal.entity();
355 const int localBlocks = space_.blockMapper().numDofs( entity );
358 std::vector<std::size_t> globalBlockDofs(localBlocks);
359 space_.blockMapper().map( entity, globalBlockDofs );
360 std::vector<typename LocalFunctionType::RangeFieldType> values( localBlocks*localBlockSize );
365 for(
int localBlock = 0 ; localBlock < localBlocks; ++ localBlock )
368 int global = globalBlockDofs[localBlock];
369 for(
int l = 0; l < localBlockSize ; ++ l, ++localDof )
371 if( dirichletBlocks_[ global ][l] )
373 interpolation( BoundaryWrapper(model_, entity, wLocal.order(), dirichletBlocks_[global][l]), values );
375 assert( (
unsigned int)localDof < wLocal.size() );
376 wLocal[ localDof ] = values[ localDof ];
382 template<
class LocalInterpolationType,
class Gr
idLocalFunctionType,
class LocalFunctionType >
383 void dirichletDofTreatment(
const LocalInterpolationType& interpolation,
384 const GridLocalFunctionType &uLocal,
385 LocalFunctionType &wLocal, Operation op )
const
388 const typename LocalFunctionType::EntityType &entity = wLocal.entity();
391 const int localBlocks = space_.blockMapper().numDofs( entity );
393 typedef typename DiscreteFunctionSpaceType::BlockMapperType::GlobalKeyType GlobalKeyType;
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 );
401 if constexpr ( useIdentity )
402 interpolation( uLocal, values );
407 for(
int localBlock = 0 ; localBlock < localBlocks; ++ localBlock )
410 int global = globalBlockDofs[localBlock];
411 for(
int l = 0; l < localBlockSize ; ++ l, ++localDof )
413 if( dirichletBlocks_[ global ][l] )
415 if (op == Operation::sub)
417 interpolation(BoundaryWrapper(model_, entity, wLocal.order(), dirichletBlocks_[global][l]), valuesModel);
418 values[ localDof ] -= valuesModel[ localDof ];
420 else if (op == Operation::add)
422 interpolation(BoundaryWrapper(model_, entity, wLocal.order(), dirichletBlocks_[global][l]), valuesModel);
423 values[ localDof ] += valuesModel[ localDof ];
426 assert( (
unsigned int)localDof < wLocal.size() );
427 wLocal[ localDof ] = values[ localDof ];
435 void updateDirichletDofs()
const
437 if( sequence_ != space_.sequence() )
440 if( ! model_.hasDirichletBoundary() )
442 hasDirichletDofs_ = false ;
447 const int blocks = space_.blockMapper().size() ;
448 dirichletBlocks_.resize( blocks );
449 for(
int i=0; i<blocks; ++i )
450 dirichletBlocks_[ i ].fill(0);
452 typedef typename DiscreteFunctionSpaceType :: IteratorType IteratorType;
453 typedef typename IteratorType :: Entity EntityType;
455 bool hasDirichletBoundary =
false;
456 const IteratorType end = space_.end();
457 for( IteratorType it = space_.begin(); it != end; ++it )
459 const EntityType &entity = *it;
461 if( Dune::Fem::HasBoundaryIntersection<GridPartType>::apply(entity) )
463 hasDirichletBoundary |= searchEntityDirichletDofs( entity, model_ );
468 sequence_ = space_.sequence();
469 if( space_.gridPart().comm().size() > 1 )
473 DirichletBuilder handle( *
this, space_ , space_.blockMapper() );
474 space_.gridPart().communicate
478 catch(
const Exception &e )
480 std::cerr << e << std::endl;
481 std::cerr <<
"Exception thrown in: " << __FILE__ <<
" line:" << __LINE__ << std::endl;
484 hasDirichletDofs_ = space_.gridPart().grid().comm().max( hasDirichletBoundary );
488 hasDirichletDofs_ = hasDirichletBoundary;
494 template<
class EntityType >
495 bool searchEntityDirichletDofs(
const EntityType &entity, ModelType& model )
const
497 typedef typename DiscreteFunctionSpaceType :: GridPartType GridPartType;
499 typedef typename GridPartType :: IntersectionIteratorType
500 IntersectionIteratorType;
502 const GridPartType &gridPart = space_.gridPart();
505 bool hasDirichletBoundary =
false;
508 std::vector<size_t> globalBlockDofs(space_.blockMapper().numDofs(entity));
509 space_.blockMapper().map(entity,globalBlockDofs);
511 std::vector<char> subEntityFilter(space_.blockMapper().numDofs(entity));
513 IntersectionIteratorType it = gridPart.ibegin( entity );
514 const IntersectionIteratorType endit = gridPart.iend( entity );
515 for( ; it != endit; ++it )
517 typedef typename IntersectionIteratorType :: Intersection IntersectionType;
518 const IntersectionType& intersection = *it;
521 if( intersection.boundary() )
524 DirichletBlock block;
526 std::array<int,dimRange> dblock;
528 model.init(intersection);
529 const bool isDirichletIntersection = model.isDirichletIntersection( intersection, dblock );
530 if (isDirichletIntersection)
533 const int face = intersection.indexInInside();
534 space_.blockMapper().onSubEntity(entity,face,1,subEntityFilter);
535 for(
unsigned int i=0;i<globalBlockDofs.size();++i)
537 assert( i < subEntityFilter.size());
538 if ( subEntityFilter[i]==0)
continue;
540 for(
int k = 0; k < localBlockSize; ++k)
542 int comp = subEntityFilter[i]-1+k;
543 assert( comp < dimRange );
544 dirichletBlocks_[ globalBlockDofs[ i ] ][k] = dblock[comp];
547 hasDirichletBoundary = true ;
553 return hasDirichletBoundary;
557 const DiscreteFunctionSpaceType& space_;
558 mutable DirichletBlockVector dirichletBlocks_;
559 mutable bool hasDirichletDofs_ ;
560 mutable int sequence_ ;
562 class DirichletBuilder;
565template <
class Model,
class Space,
bool useIdentity >
566class DirichletConstraints< Model,Space, useIdentity > :: DirichletBuilder
567 :
public CommDataHandleIF< DirichletBuilder, int >
571 typedef typename SpaceType::BlockMapperType MapperType;
573 enum { nCodim = SpaceType :: GridType :: dimension + 1 };
576 typedef int DataType;
581 typedef DirichletConstraints< Model,Space, useIdentity > DirichletType;
582 const DirichletType &dirichlet_;
585 const MapperType &mapper_;
587 static const int blockSize = SpaceType::localBlockSize;
590 DirichletBuilder(
const DirichletType &dirichlet,
592 const MapperType &mapper )
593 : myRank_( space.gridPart().comm().rank() ),
594 mySize_( space.gridPart().comm().
size() ),
595 dirichlet_( dirichlet ),
600 bool contains (
int dim,
int codim )
const
602 return mapper_.contains( codim );
605 bool fixedSize (
int dim,
int codim )
const
611 template<
class MessageBuffer,
class Entity >
612 inline void gather ( MessageBuffer &buffer,
613 const Entity &entity )
const
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 )
621 int global = globalBlockDofs[localBlock];
622 for (
int r=0;r<blockSize;++r)
623 if (dirichlet_.dirichletBlocks_[ global ][r] )
634 template<
class MessageBuffer,
class EntityType >
635 inline void scatter ( MessageBuffer &buffer,
636 const EntityType &entity,
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 )
646 int global = globalBlockDofs[localBlock];
647 for (
int r=0;r<blockSize;++r)
651 if ( !dirichlet_.dirichletBlocks_[ global ][r] && val == 1)
652 dirichlet_.dirichletBlocks_[ global ][r] =
true;
657 template<
class Entity >
658 size_t size (
const Entity &entity )
const
660 return blockSize * mapper_.numEntityDofs( entity );
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
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