3 #ifndef DUNE_Q1_LOCALBASIS_HH 4 #define DUNE_Q1_LOCALBASIS_HH 8 #include <dune/common/fmatrix.hh> 25 template<
class D,
class R,
int dim>
40 std::vector<typename Traits::RangeType>& out)
const 44 for (
size_t i=0; i<
size(); i++) {
48 for (
int j=0; j<dim; j++)
50 out[i] *= (i & (1<<j)) ? in[j] : 1-in[j];
59 std::vector<typename Traits::JacobianType>& out)
const 64 for (
size_t i=0; i<
size(); i++) {
67 for (
int j=0; j<dim; j++) {
71 out[i][0][j] = (i & (1<<j)) ? 1 : -1;
73 for (
int k=0; k<dim; k++) {
77 out[i][0][j] *= (i & (1<<k)) ? in[k] : 1-in[k];
94 std::vector<typename Traits::RangeType>& out)
const 96 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
98 if (totalOrder == 0) {
101 else if (totalOrder == 1) {
104 auto direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
105 if (direction >= dim) {
106 DUNE_THROW(RangeError,
"Direction of partial derivative not found!");
110 for (std::size_t i = 0; i <
size(); ++i) {
114 out[i] = (i & (1<<direction)) ? 1 : -1;
116 for (
int k = 0; k < dim; ++k) {
119 out[i] *= (i & (1<<k)) ? in[k] : 1-in[k];
126 DUNE_THROW(NotImplemented,
"To be implemented!");
LocalBasisTraits< D, dim, Dune::FieldVector< D, dim >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, dim > > Traits
Definition: q1localbasis.hh:30
unsigned int size() const
number of shape functions
Definition: q1localbasis.hh:33
unsigned int order() const
Polynomial order of the shape functions.
Definition: q1localbasis.hh:131
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:15
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 any order of all shape functions.
Definition: q1localbasis.hh:92
D DomainType
domain type
Definition: localbasis.hh:49
Lagrange shape functions of order 1 on the reference cube.
Definition: q1localbasis.hh:26
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: q1localbasis.hh:58
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: q1localbasis.hh:39
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:37