4 #ifndef DUNE_LOCALFUNCTIONS_WHITNEY_EDGES0_5_BASIS_HH
5 #define DUNE_LOCALFUNCTIONS_WHITNEY_EDGES0_5_BASIS_HH
10 #include <dune/common/fmatrix.hh>
11 #include <dune/common/fvector.hh>
33 template<
class Geometry,
class RF>
35 private EdgeS0_5Common<Geometry::mydimension, typename Geometry::ctype>
48 typedef FieldVector<RangeField, dimRange>
Range;
50 typedef FieldMatrix<RangeField, dimRange, dimDomainGlobal>
Jacobian;
70 std::vector<typename P1Basis::Traits::Jacobian> p1j;
72 std::vector<typename Traits::DomainField> edgel;
82 template<
typename VertexOrder>
93 for(std::size_t i = 0; i <
s; ++i) {
94 edgel[i] = (geo.corner(
refelem.subEntity(i,dim-1,0,dim))-
95 geo.corner(
refelem.subEntity(i,dim-1,1,dim))
97 const typename VertexOrder::iterator& edgeVertexOrder =
98 vertexOrder.begin(dim-1, i);
99 if(edgeVertexOrder[0] > edgeVertexOrder[1])
105 std::size_t
size ()
const {
return s; }
109 std::vector<typename Traits::Range>& out)
const
115 std::vector<typename P1LocalBasis::Traits::RangeType> p1v;
118 for(std::size_t i = 0; i <
s; i++) {
119 const std::size_t i0 =
refelem.subEntity(i,dim-1,0,dim);
120 const std::size_t i1 =
refelem.subEntity(i,dim-1,1,dim);
121 out[i].axpy( p1v[i0], p1j[i1][0]);
122 out[i].axpy(-p1v[i1], p1j[i0][0]);
129 std::vector<typename Traits::Jacobian>& out)
const
133 for(std::size_t i = 0; i <
s; i++) {
134 const std::size_t i0 =
refelem.subEntity(i,dim-1,0,dim);
135 const std::size_t i1 =
refelem.subEntity(i,dim-1,1,dim);
136 for(std::size_t j = 0; j < dim; j++)
137 for(std::size_t k = 0; k < dim; k++)
138 out[i][j][k] = edgel[i] *
139 (p1j[i0][0][k]*p1j[i1][0][j]-p1j[i1][0][k]*p1j[i0][0][j]);
144 std::size_t
order ()
const {
return 1; }
147 template<
class Geometry,
class RF>
148 const typename EdgeS0_5Basis<Geometry, RF>::P1LocalBasis&
149 EdgeS0_5Basis<Geometry, RF>::p1LocalBasis = P1LocalBasis();
153 #endif // DUNE_LOCALFUNCTIONS_WHITNEY_EDGES0_5_BASIS_HH