dune-pdelab  2.4-dev
bdm1cube2dfem.hh
Go to the documentation of this file.
1 // -*- tab-width: 4; indent-tabs-mode: nil -*-
2 #ifndef DUNE_PDELAB_FINITEELEMENTMAP_BDM1CUBE2DFEM_HH
3 #define DUNE_PDELAB_FINITEELEMENTMAP_BDM1CUBE2DFEM_HH
4 
5 #include <vector>
6 #include <dune/localfunctions/brezzidouglasmarini/brezzidouglasmarini1cube2d.hh>
7 #include "finiteelementmap.hh"
8 
9 namespace Dune {
10  namespace PDELab {
11 
14  template<typename GV, typename D, typename R>
17  LocalFiniteElementMapTraits<Dune::BDM1Cube2DLocalFiniteElement<D,R> >,
18  BDM1Cube2DLocalFiniteElementMap<GV,D,R> >
19  {
20  typedef Dune::BDM1Cube2DLocalFiniteElement<D,R> FE;
21  typedef typename GV::IndexSet IndexSet;
22 
23  public:
26 
29  : gv(gv_), is(gv_.indexSet()), orient(gv_.size(0))
30  {
31  // create all variants
32  for (int i = 0; i < 16; i++)
33  {
34  variant[i] = FE(i);
35  }
36 
37  // compute orientation for all elements
38  //--------------------------------------------
39  // loop once over the grid
40  for (const auto& cell : elements(gv)) {
41  unsigned int myId = is.template index<0>(cell);
42  orient[myId] = 0;
43 
44  for (const auto& intersection : intersections(gv,cell)) {
45  if (intersection.neighbor()
46  && is.template index<0>(intersection.outside()) > myId)
47  {
48  orient[myId] |= 1 << intersection.indexInInside();
49  }
50  }
51  }
52  }
53 
55  template<class EntityType>
56  const typename Traits::FiniteElementType& find(const EntityType& e) const
57  {
58  return variant[orient[is.index(e)]];
59  }
60 
61  bool fixedSize() const
62  {
63  return true;
64  }
65 
66  std::size_t size(GeometryType gt) const
67  {
68  switch (gt.dim())
69  {
70  case 1:
71  return 2;
72  default:
73  return 0;
74  }
75  }
76 
77  std::size_t maxLocalSize() const
78  {
79  return 8;
80  }
81 
82  private:
83  GV gv;
84  FE variant[16];
85  const IndexSet& is;
86  std::vector<unsigned char> orient;
87  };
88  } // end namespace PDELab
89 } // end namespace Dune
90 
91 #endif // DUNE_PDELAB_FINITEELEMENTMAP_BDM1CUBE2DFEM_HH
bool fixedSize() const
Definition: bdm1cube2dfem.hh:61
std::size_t size(GeometryType gt) const
Definition: bdm1cube2dfem.hh:66
LocalFiniteElementMapTraits< FE > Traits
export type of the signature
Definition: bdm1cube2dfem.hh:25
const E & e
Definition: interpolate.hh:172
collect types exported by a finite element map
Definition: finiteelementmap.hh:38
Definition: bdm1cube2dfem.hh:15
const Traits::FiniteElementType & find(const EntityType &e) const
get local basis functions for entity
Definition: bdm1cube2dfem.hh:56
T FiniteElementType
Type of finite element from local functions.
Definition: finiteelementmap.hh:30
Definition: adaptivity.hh:27
interface for a finite element map
Definition: finiteelementmap.hh:42
std::size_t maxLocalSize() const
Definition: bdm1cube2dfem.hh:77
BDM1Cube2DLocalFiniteElementMap(const GV &gv_)
Use when Imp has a standard constructor.
Definition: bdm1cube2dfem.hh:28