3 #ifndef DUNE_DUAL_Q1_LOCALBASIS_HH
4 #define DUNE_DUAL_Q1_LOCALBASIS_HH
6 #include <dune/common/fvector.hh>
7 #include <dune/common/fmatrix.hh>
22 template<
class D,
class R,
int dim>
29 void setCoefficients(
const Dune::array<Dune::FieldVector<R, (1<<dim)> ,(1<<dim)>& coefficients)
31 coefficients_ = coefficients;
42 std::vector<typename Traits::RangeType>& out)
const
45 std::vector<typename Traits::RangeType> q1Values(
size());
47 for (
size_t i=0; i<
size(); i++) {
51 for (
int j=0; j<dim; j++)
53 q1Values[i] *= (i & (1<<j)) ? in[j] : 1-in[j];
59 for (
size_t i=0; i<
size(); i++)
62 for (
size_t i=0; i<
size(); i++)
63 for (
size_t j=0; j<
size(); j++)
64 out[i] += coefficients_[i][j]*q1Values[j];
72 std::vector<typename Traits::JacobianType>& out)
const
75 std::vector<typename Traits::JacobianType> q1Jacs(
size());
78 for (
size_t i=0; i<
size(); i++) {
81 for (
int j=0; j<dim; j++) {
85 q1Jacs[i][0][j] = (i & (1<<j)) ? 1 : -1;
87 for (
int k=0; k<dim; k++) {
91 q1Jacs[i][0][j] *= (i & (1<<k)) ? in[k] : 1-in[k];
101 for (
size_t i=0; i<
size(); i++)
104 for (
size_t i=0; i<
size(); i++)
105 for (
size_t j=0; j<
size(); j++)
106 out[i].axpy(coefficients_[i][j],q1Jacs[j]);
117 Dune::array<Dune::FieldVector<R, (1<<dim)> ,(1<<dim)> coefficients_;
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:39
LocalBasisTraits< D, dim, Dune::FieldVector< D, dim >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, dim > > Traits
Definition: dualq1localbasis.hh:27
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: dualq1localbasis.hh:71
unsigned int order() const
Polynomial order of the shape functions.
Definition: dualq1localbasis.hh:111
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: dualq1localbasis.hh:41
void setCoefficients(const Dune::array< Dune::FieldVector< R,(1<< dim)>,(1<< dim)> &coefficients)
Definition: dualq1localbasis.hh:29
Dual Lagrange shape functions of order 1 on the reference cube.
Definition: dualq1localbasis.hh:23
unsigned int size() const
number of shape functions
Definition: dualq1localbasis.hh:35
D DomainType
domain type
Definition: localbasis.hh:51