2#ifndef DUNE_PDELAB_CONSTRAINTS_INTERIORNODE_HH
3#define DUNE_PDELAB_CONSTRAINTS_INTERIORNODE_HH
7#include <dune/grid/common/gridenums.hh>
8#include <dune/geometry/referenceelements.hh>
9#include <dune/localfunctions/common/interfaceswitch.hh>
10#include <dune/localfunctions/common/localkey.hh>
12#include <dune/typetree/typetree.hh>
26 std::vector<bool> interior;
28 enum{doBoundary=
false};
29 enum{doProcessor=
false};
30 enum{doSkeleton=
false};
40 template<
typename P,
typename EG,
typename LFS,
typename T>
41 void volume (
const P& param,
const EG& eg,
const LFS& lfs, T& trafo)
const
43 typedef typename EG::Entity
Entity;
47 typename T::RowType empty;
48 typedef typename LFS::Traits::SizeType size_type;
50 typename LFS::Traits::FiniteElementType
52 for (size_type i=0; i<lfs.size(); i++){
53 const LocalKey& key = FESwitch::coefficients(lfs.finiteElement()).localKey(i);
54 assert(key.
codim() == dim &&
"InteriorNodeConstraints only work for vertex DOFs");
55 assert(key.
index() == 0 &&
"InteriorNodeConstraints only work for P1 shape functions");
61 unsigned int idx = lfs.gridFunctionSpace().gridView().indexSet().subIndex(eg.entity(), local_idx, dim);
65 trafo[lfs.dofIndex(i)] = empty;
70 const std::vector<bool> & interiorNodes()
const
76 void updateInteriorNodes(
const GV & gv)
79 const int dim = GV::dimension;
80 typedef typename GV::Grid::ctype ctype;
82 interior.resize(gv.indexSet().size(dim));
83 for(
int i=0; i< interior.size(); i++)
87 for(
const auto& entity : elements(gv))
90 for (
const auto& intersection : intersections(gv,entity))
92 if (intersection.boundary())
95 unsigned int f = intersection.indexInInside();
98 assert(entity.geometry().type().isSimplex() &&
"InteriorNodeConstraints only work for simplicial meshes");
99 unsigned int sz = refelem.size(f,1, dim);
101 for (
unsigned int v = 0; v < sz; ++v)
103 unsigned int local_idx = refelem.subEntity (f,1, v,dim);
104 unsigned int idx = gv.indexSet().subIndex(entity, local_idx, dim);
105 interior[idx] =
false;
Wrapper class for entities.
Definition: entity.hh:64
@ dimension
Know the grid dimension.
Definition: entity.hh:109
Describe position of one degree of freedom.
Definition: localkey.hh:21
unsigned int index() const
Return offset within subentity.
Definition: localkey.hh:66
unsigned int codim() const
Return codim of associated entity.
Definition: localkey.hh:60
unsigned int subEntity() const
Return number of associated subentity.
Definition: localkey.hh:54
constraints all DOFs associated with interior vertices This allows to implement surface FEM using sta...
Definition: interiornode.hh:25
void volume(const P ¶m, const EG &eg, const LFS &lfs, T &trafo) const
volume constraints
Definition: interiornode.hh:41
Base class for leaf nodes in a dune-typetree.
Definition: leafnode.hh:26
Dune namespace.
Definition: alignedallocator.hh:11
Switch for uniform treatment of finite element with either the local or the global interface.
Definition: interfaceswitch.hh:28
static const ReferenceElement & simplex()
get simplex reference elements
Definition: referenceelements.hh:202