5#ifndef DUNE_DUAL_Q1_LOCALBASIS_HH
6#define DUNE_DUAL_Q1_LOCALBASIS_HH
11#include <dune/common/fvector.hh>
12#include <dune/common/fmatrix.hh>
28 template<
class D,
class R,
int dim>
35 void setCoefficients(
const std::array<Dune::FieldVector<R, (1<<dim)> ,(1<<dim)>& coefficients)
37 coefficients_ = coefficients;
48 std::vector<typename Traits::RangeType>& out)
const
51 std::vector<typename Traits::RangeType> q1Values(
size());
53 for (
size_t i=0; i<
size(); i++) {
57 for (
int j=0; j<dim; j++)
59 q1Values[i] *= (i & (1<<j)) ? in[j] : 1-in[j];
65 for (
size_t i=0; i<
size(); i++)
68 for (
size_t i=0; i<
size(); i++)
69 for (
size_t j=0; j<
size(); j++)
70 out[i] += coefficients_[i][j]*q1Values[j];
78 std::vector<typename Traits::JacobianType>& out)
const
81 std::vector<typename Traits::JacobianType> q1Jacs(
size());
84 for (
size_t i=0; i<
size(); i++) {
87 for (
int j=0; j<dim; j++) {
91 q1Jacs[i][0][j] = (i & (1<<j)) ? 1 : -1;
93 for (
int k=0; k<dim; k++) {
97 q1Jacs[i][0][j] *= (i & (1<<k)) ? in[k] : 1-in[k];
107 for (
size_t i=0; i<
size(); i++)
110 for (
size_t i=0; i<
size(); i++)
111 for (
size_t j=0; j<
size(); j++)
112 out[i].axpy(coefficients_[i][j],q1Jacs[j]);
119 std::vector<typename Traits::RangeType>& out)
const
121 auto totalOrder = std::accumulate(
order.begin(),
order.end(), 0);
122 if (totalOrder == 0) {
125 DUNE_THROW(NotImplemented,
"Desired derivative order is not implemented");
136 std::array<Dune::FieldVector<R, (1<<dim)> ,(1<<dim)> coefficients_;
Definition bdfmcube.hh:18
Type traits for LocalBasisVirtualInterface.
Definition common/localbasis.hh:35
D DomainType
domain type
Definition common/localbasis.hh:43
Dual Lagrange shape functions of order 1 on the reference cube.
Definition dualq1localbasis.hh:30
unsigned int size() const
number of shape functions
Definition dualq1localbasis.hh:41
unsigned int order() const
Polynomial order of the shape functions.
Definition dualq1localbasis.hh:130
void setCoefficients(const std::array< Dune::FieldVector< R,(1<< dim)>,(1<< dim)> &coefficients)
Definition dualq1localbasis.hh:35
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition dualq1localbasis.hh:47
LocalBasisTraits< D, dim, Dune::FieldVector< D, dim >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, dim > > Traits
Definition dualq1localbasis.hh:33
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition dualq1localbasis.hh:77
void partial(const std::array< unsigned int, dim > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate partial derivatives of all shape functions.
Definition dualq1localbasis.hh:117