DUNE PDELab (git)

bdm1simplex2dfem.hh
1// -*- tab-width: 4; indent-tabs-mode: nil -*-
2#ifndef DUNE_PDELAB_FINITEELEMENTMAP_BDM1SIMPLEX2DFEM_HH
3#define DUNE_PDELAB_FINITEELEMENTMAP_BDM1SIMPLEX2DFEM_HH
4
5#include <vector>
6#include <dune/localfunctions/brezzidouglasmarini/brezzidouglasmarini1simplex2d.hh>
7#include "finiteelementmap.hh"
8
9namespace Dune {
10 namespace PDELab {
11
14 template<typename GV, typename D, typename R>
17 LocalFiniteElementMapTraits< Dune::BDM1Simplex2DLocalFiniteElement<D,R> >,
18 BDM1Simplex2DLocalFiniteElementMap<GV,D,R> >
19 {
21 typedef typename GV::IndexSet IndexSet;
22
23 public:
26
28 static constexpr int dimension = GV::dimension;
29
32 : gv(gv_), is(gv_.indexSet()), orient(gv_.size(0))
33 {
34 // create all variants
35 for (int i = 0; i < 8; i++)
36 {
37 variant[i] = FE(i);
38 }
39
40 // compute orientation for all elements
41 //--------------------------------------------
42 // loop once over the grid
43 for (const auto& cell : elements(gv)) {
44 unsigned int myId = is.template index<0>(cell);
45 orient[myId] = 0;
46
47 for (const auto& intersection : intersections(gv,cell)) {
48 if (intersection.neighbor()
49 && is.template index<0>(intersection.outside()) > myId)
50 {
51 orient[myId] |= 1 << intersection.indexInInside();
52 }
53 }
54 }
55 }
56
58 template<class EntityType>
59 const typename Traits::FiniteElementType& find(const EntityType& e) const
60 {
61 return variant[orient[is.index(e)]];
62 }
63
64 static constexpr bool fixedSize()
65 {
66 return true;
67 }
68
69 static constexpr bool hasDOFs(int codim)
70 {
71 return codim == 1;
72 }
73
74 static constexpr std::size_t size(GeometryType gt)
75 {
76 switch (gt.dim())
77 {
78 case 1:
79 return 2;
80 default:
81 return 0;
82 }
83 }
84
85 static constexpr std::size_t maxLocalSize()
86 {
87 return 6;
88 }
89
90 private:
91 GV gv;
92 FE variant[8];
93 const IndexSet& is;
94 std::vector<unsigned char> orient;
95 };
96 } // end namespace PDELab
97} // end namespace Dune
98
99#endif // DUNE_PDELAB_FINITEELEMENTMAP_BDM1SIMPLEX2DFEM_HH
First order Brezzi-Douglas-Marini shape functions on triangles.
Definition: brezzidouglasmarini1simplex2d.hh:28
Definition: bdm1simplex2dfem.hh:19
BDM1Simplex2DLocalFiniteElementMap(const GV &gv_)
Use when Imp has a standard constructor.
Definition: bdm1simplex2dfem.hh:31
const Traits::FiniteElementType & find(const EntityType &e) const
get local basis functions for entity
Definition: bdm1simplex2dfem.hh:59
static constexpr int dimension
The dimension of the finite elements returned by this map.
Definition: bdm1simplex2dfem.hh:28
LocalFiniteElementMapTraits< FE > Traits
export type of the signature
Definition: bdm1simplex2dfem.hh:25
interface for a finite element map
Definition: finiteelementmap.hh:43
bool gt(const T &first, const T &second, typename EpsilonType< T >::Type epsilon)
test if first greater than second
Definition: float_cmp.cc:158
Dune namespace.
Definition: alignedallocator.hh:13
T FiniteElementType
Type of finite element from local functions.
Definition: finiteelementmap.hh:30
collect types exported by a finite element map
Definition: finiteelementmap.hh:38
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)