3 #ifndef DUNE_DUAL_Q1_LOCALBASIS_HH
4 #define DUNE_DUAL_Q1_LOCALBASIS_HH
8 #include <dune/common/fvector.hh>
9 #include <dune/common/fmatrix.hh>
24 template<
class D,
class R,
int dim>
31 void setCoefficients(
const std::array<Dune::FieldVector<R, (1<<dim)> ,(1<<dim)>& coefficients)
33 coefficients_ = coefficients;
44 std::vector<typename Traits::RangeType>& out)
const
47 std::vector<typename Traits::RangeType> q1Values(
size());
49 for (
size_t i=0; i<
size(); i++) {
53 for (
int j=0; j<dim; j++)
55 q1Values[i] *= (i & (1<<j)) ? in[j] : 1-in[j];
61 for (
size_t i=0; i<
size(); i++)
64 for (
size_t i=0; i<
size(); i++)
65 for (
size_t j=0; j<
size(); j++)
66 out[i] += coefficients_[i][j]*q1Values[j];
74 std::vector<typename Traits::JacobianType>& out)
const
77 std::vector<typename Traits::JacobianType> q1Jacs(
size());
80 for (
size_t i=0; i<
size(); i++) {
83 for (
int j=0; j<dim; j++) {
87 q1Jacs[i][0][j] = (i & (1<<j)) ? 1 : -1;
89 for (
int k=0; k<dim; k++) {
93 q1Jacs[i][0][j] *= (i & (1<<k)) ? in[k] : 1-in[k];
103 for (
size_t i=0; i<
size(); i++)
106 for (
size_t i=0; i<
size(); i++)
107 for (
size_t j=0; j<
size(); j++)
108 out[i].axpy(coefficients_[i][j],q1Jacs[j]);
119 std::array<Dune::FieldVector<R, (1<<dim)> ,(1<<dim)> coefficients_;
unsigned int order() const
Polynomial order of the shape functions.
Definition: dualq1localbasis.hh:113
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: dualq1localbasis.hh:43
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:14
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: dualq1localbasis.hh:73
unsigned int size() const
number of shape functions
Definition: dualq1localbasis.hh:37
void setCoefficients(const std::array< Dune::FieldVector< R,(1<< dim)>,(1<< dim)> &coefficients)
Definition: dualq1localbasis.hh:31
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:37
D DomainType
domain type
Definition: localbasis.hh:49
Dual Lagrange shape functions of order 1 on the reference cube.
Definition: dualq1localbasis.hh:25
LocalBasisTraits< D, dim, Dune::FieldVector< D, dim >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, dim > > Traits
Definition: dualq1localbasis.hh:29