dune-localfunctions  2.4.1
prismp1localbasis.hh
Go to the documentation of this file.
1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 #ifndef DUNE_PRISM_P1_LOCALBASIS_HH
4 #define DUNE_PRISM_P1_LOCALBASIS_HH
5 
6 #include <dune/common/fmatrix.hh>
7 
9 
10 namespace Dune
11 {
22  template<class D, class R>
24  {
25  public:
27  typedef LocalBasisTraits<D,3,Dune::FieldVector<D,3>,R,1,Dune::FieldVector<R,1>,
28  Dune::FieldMatrix<R,1,3> > Traits;
29 
31  unsigned int size () const
32  {
33  return 6;
34  }
35 
37  inline void evaluateFunction (const typename Traits::DomainType& in,
38  std::vector<typename Traits::RangeType>& out) const
39  {
40  out.resize(6);
41  out[0] = (1.0-in[0]-in[1])*(1.0-in[2]);
42  out[1] = in[0]*(1-in[2]);
43  out[2] = in[1]*(1-in[2]);
44  out[3] = in[2]*(1.0-in[0]-in[1]);
45  out[4] = in[0]*in[2];
46  out[5] = in[1]*in[2];
47  }
48 
50  inline void
51  evaluateJacobian (const typename Traits::DomainType& in, // position
52  std::vector<typename Traits::JacobianType>& out) const // return value
53  {
54  out.resize(6);
55  out[0][0][0] = in[2]-1; out[0][0][1] = in[2]-1; out[0][0][2] = in[0]+in[1]-1; // basis function 0
56  out[1][0][0] = 1-in[2]; out[1][0][1] = 0; out[1][0][2] = -in[0]; // basis function 1
57  out[2][0][0] = 0; out[2][0][1] = 1-in[2]; out[2][0][2] = -in[1]; // basis function 2
58  out[3][0][0] = -in[2]; out[3][0][1] = -in[2]; out[3][0][2] = 1-in[0]-in[1]; // basis function 3
59  out[4][0][0] = in[2]; out[4][0][1] = 0; out[4][0][2] = in[0]; // basis function 4
60  out[5][0][0] = 0; out[5][0][1] = in[2]; out[5][0][2] = in[1]; // basis function 5
61  }
62 
64  unsigned int order () const
65  {
66  return 1;
67  }
68  };
69 }
70 #endif
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: prismp1localbasis.hh:51
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: prismp1localbasis.hh:37
unsigned int size() const
number of shape functions
Definition: prismp1localbasis.hh:31
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:14
unsigned int order() const
Polynomial order of the shape functions.
Definition: prismp1localbasis.hh:64
Linear Lagrange shape functions on the prism.
Definition: prismp1localbasis.hh:23
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:37
D DomainType
domain type
Definition: localbasis.hh:49
LocalBasisTraits< D, 3, Dune::FieldVector< D, 3 >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, 3 > > Traits
export type traits for function signature
Definition: prismp1localbasis.hh:28