6#ifndef DUNE_LOCALFUNCTIONS_WHITNEY_EDGES0_5_BASIS_HH 
    7#define DUNE_LOCALFUNCTIONS_WHITNEY_EDGES0_5_BASIS_HH 
   15#include <dune/localfunctions/common/localtoglobaladaptors.hh> 
   16#include <dune/localfunctions/lagrange/lagrangesimplex.hh> 
   17#include <dune/localfunctions/whitney/edges0.5/common.hh> 
   35  template<
class Geometry, 
class RF>
 
   37    private EdgeS0_5Common<Geometry::mydimension, typename Geometry::ctype>
 
   48      typedef RF RangeField;
 
   49      static const std::size_t dimRange = dimDomainLocal;
 
   56    typedef Dune::Impl::LagrangeSimplexLocalBasis<
typename Traits::DomainField,
 
   57        typename Traits::RangeField,
 
   58        Traits::dimDomainLocal,
 
   63    static const P1LocalBasis& p1LocalBasis;
 
   64    static const std::size_t dim = Traits::dimDomainLocal;
 
   71    std::vector<typename P1Basis::Traits::Jacobian> p1j;
 
   73    std::vector<typename Traits::DomainField> edgel;
 
   83    template<
typename VertexOrder>
 
   91      P1Basis(p1LocalBasis, geo).evaluateJacobian(xl, p1j);
 
   94      for(std::size_t i = 0; i < s; ++i) {
 
   95        edgel[i] = (geo.
corner(refelem.subEntity(i,dim-1,0,dim))-
 
   96                    geo.
corner(refelem.subEntity(i,dim-1,1,dim))
 
   98        const typename VertexOrder::iterator& edgeVertexOrder =
 
   99          vertexOrder.begin(dim-1, i);
 
  100        if(edgeVertexOrder[0] > edgeVertexOrder[1])
 
  106    std::size_t 
size ()
 const { 
return s; }
 
  110                          std::vector<typename Traits::Range>& out)
 const 
  116      std::vector<typename P1LocalBasis::Traits::RangeType> p1v;
 
  117      p1LocalBasis.evaluateFunction(xl, p1v);
 
  119      for(std::size_t i = 0; i < s; i++) {
 
  120        const std::size_t i0 = refelem.subEntity(i,dim-1,0,dim);
 
  121        const std::size_t i1 = refelem.subEntity(i,dim-1,1,dim);
 
  122        out[i].axpy( p1v[i0], p1j[i1][0]);
 
  123        out[i].axpy(-p1v[i1], p1j[i0][0]);
 
  130                          std::vector<typename Traits::Jacobian>& out)
 const 
  134      for(std::size_t i = 0; i < s; i++) {
 
  135        const std::size_t i0 = refelem.subEntity(i,dim-1,0,dim);
 
  136        const std::size_t i1 = refelem.subEntity(i,dim-1,1,dim);
 
  137        for(std::size_t j = 0; j < dim; j++)
 
  138          for(std::size_t k = 0; k < dim; k++)
 
  139            out[i][j][k] = edgel[i] *
 
  140                           (p1j[i0][0][k]*p1j[i1][0][j]-p1j[i1][0][k]*p1j[i0][0][j]);
 
  147                  std::vector<typename Traits::Range>& out) 
const       
  150      if (totalOrder == 0) {
 
  152      } 
else if (totalOrder==1) {
 
  153        auto const k = std::distance(
order.begin(), std::find(
order.begin(), 
order.end(), 1));
 
  156        for (std::size_t i = 0; i < s; i++)
 
  158          const std::size_t i0 = refelem.subEntity(i,dim-1,0,dim);
 
  159          const std::size_t i1 = refelem.subEntity(i,dim-1,1,dim);
 
  160          for(std::size_t j = 0; j < dim; j++)
 
  161            out[i][j] = edgel[i] *
 
  162              (p1j[i0][0][k]*p1j[i1][0][j] - p1j[i1][0][k]*p1j[i0][0][j]);
 
  170    std::size_t 
order ()
 const { 
return 1; }
 
  173  template<
class Geometry, 
class RF>
 
  174  const typename EdgeS0_5Basis<Geometry, RF>::P1LocalBasis&
 
  175  EdgeS0_5Basis<Geometry, RF>::p1LocalBasis = P1LocalBasis();
 
Basis for order 0.5 (lowest order) edge elements on simplices.
Definition: basis.hh:38
 
EdgeS0_5Basis(const Geometry &geo, const VertexOrder &vertexOrder)
Construct an EdgeS0_5Basis.
Definition: basis.hh:84
 
void evaluateJacobian(const typename Traits::DomainLocal &, std::vector< typename Traits::Jacobian > &out) const
Evaluate all Jacobians.
Definition: basis.hh:129
 
void evaluateFunction(const typename Traits::DomainLocal &xl, std::vector< typename Traits::Range > &out) const
Evaluate all shape functions.
Definition: basis.hh:109
 
std::size_t size() const
number of shape functions
Definition: basis.hh:106
 
void partial(const std::array< unsigned int, dim > &order, const typename Traits::DomainLocal &in, std::vector< typename Traits::Range > &out) const
Evaluate partial derivatives of all shape functions.
Definition: basis.hh:145
 
std::size_t order() const
Polynomial order of the shape functions.
Definition: basis.hh:170
 
A dense n x m matrix.
Definition: fmatrix.hh:117
 
vector space out of a tensor product of fields.
Definition: fvector.hh:97
 
Wrapper class for geometries.
Definition: geometry.hh:71
 
GlobalCoordinate corner(int i) const
Obtain a corner of the geometry.
Definition: geometry.hh:219
 
static constexpr int coorddimension
dimension of embedding coordinate system
Definition: geometry.hh:97
 
GridImp::ctype ctype
define type used for coordinates in grid module
Definition: geometry.hh:100
 
static constexpr int mydimension
geometry dimension
Definition: geometry.hh:94
 
Default exception for dummy implementations.
Definition: exceptions.hh:357
 
Convert a simple scalar local basis into a global basis.
Definition: localtoglobaladaptors.hh:65
 
Implements a matrix constructed from a given type representing a field and compile-time given number ...
 
Implements a vector constructed from a given type representing a field and a compile-time given size.
 
#define DUNE_THROW(E,...)
Definition: exceptions.hh:314
 
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:280
 
Dune namespace.
Definition: alignedallocator.hh:13
 
export type traits for function signature
Definition: basis.hh:41
 
Common base class for edge elements.
Definition: common.hh:23
 
RefElem refelem
The reference element for this edge element.
Definition: common.hh:30
 
std::size_t s
The number of base functions.
Definition: common.hh:38