/usr/include/dune/localfunctions/lagrange/pk3d/pk3dlocalbasis.hh is in libdune-localfunctions-dev 2.2.1-2.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 | #ifndef DUNE_PK3DLOCALBASIS_HH
#define DUNE_PK3DLOCALBASIS_HH
#include <dune/common/fmatrix.hh>
#include <dune/localfunctions/common/localbasis.hh>
namespace Dune
{
/**@ingroup LocalBasisImplementation
\brief Lagrange shape functions of arbitrary order on the reference tetrahedron.
Lagrange shape functions of arbitrary order have the property that
\f$\hat\phi^i(x_j) = \delta_{i,j}\f$ for certain points \f$x_j\f$.
\tparam D Type to represent the field in the domain.
\tparam R Type to represent the field in the range.
\tparam k Polynomial order.
\nosubgrouping
*/
template<class D, class R, unsigned int k>
class Pk3DLocalBasis
{
public:
enum {N = (k+1)*(k+2)*(k+3)/6};
enum {O = k};
typedef LocalBasisTraits<D,3,Dune::FieldVector<D,3>,R,1,Dune::FieldVector<R,1>,
Dune::FieldMatrix<R,1,3> > Traits;
//! \brief Standard constructor
Pk3DLocalBasis () {}
//! \brief number of shape functions
unsigned int size () const
{
return N;
}
//! \brief Evaluate all shape functions
inline void evaluateFunction (const typename Traits::DomainType& x,
std::vector<typename Traits::RangeType>& out) const
{
out.resize(N);
typename Traits::DomainType kx = x;
kx *= k;
unsigned int n = 0;
unsigned int i[4];
R factor[4];
for (i[2] = 0; i[2] <= k; ++i[2])
{
factor[2] = 1.0;
for (unsigned int j = 0; j < i[2]; ++j)
factor[2] *= (kx[2]-j) / (i[2]-j);
for (i[1] = 0; i[1] <= k - i[2]; ++i[1])
{
factor[1] = 1.0;
for (unsigned int j = 0; j < i[1]; ++j)
factor[1] *= (kx[1]-j) / (i[1]-j);
for (i[0] = 0; i[0] <= k - i[1] - i[2]; ++i[0])
{
factor[0] = 1.0;
for (unsigned int j = 0; j < i[0]; ++j)
factor[0] *= (kx[0]-j) / (i[0]-j);
i[3] = k - i[0] - i[1] - i[2];
D kx3 = k - kx[0] - kx[1] - kx[2];
factor[3] = 1.0;
for (unsigned int j = 0; j < i[3]; ++j)
factor[3] *= (kx3-j) / (i[3]-j);
out[n++] = factor[0] * factor[1] * factor[2] * factor[3];
}
}
}
}
//! \brief Evaluate Jacobian of all shape functions
inline void
evaluateJacobian (const typename Traits::DomainType& x, // position
std::vector<typename Traits::JacobianType>& out) const // return value
{
out.resize(N);
typename Traits::DomainType kx = x;
kx *= k;
unsigned int n = 0;
unsigned int i[4];
R factor[4];
for (i[2] = 0; i[2] <= k; ++i[2])
{
factor[2] = 1.0;
for (unsigned int j = 0; j < i[2]; ++j)
factor[2] *= (kx[2]-j) / (i[2]-j);
for (i[1] = 0; i[1] <= k - i[2]; ++i[1])
{
factor[1] = 1.0;
for (unsigned int j = 0; j < i[1]; ++j)
factor[1] *= (kx[1]-j) / (i[1]-j);
for (i[0] = 0; i[0] <= k - i[1] - i[2]; ++i[0])
{
factor[0] = 1.0;
for (unsigned int j = 0; j < i[0]; ++j)
factor[0] *= (kx[0]-j) / (i[0]-j);
i[3] = k - i[0] - i[1] - i[2];
D kx3 = k - kx[0] - kx[1] - kx[2];
R sum3 = 0.0;
factor[3] = 1.0;
for (unsigned int j = 0; j < i[3]; ++j)
factor[3] /= i[3] - j;
R prod_all = factor[0] * factor[1] * factor[2] * factor[3];
for (unsigned int j = 0; j < i[3]; ++j)
{
R prod = prod_all;
for (unsigned int l = 0; l < i[3]; ++l)
if (j == l)
prod *= -R(k);
else
prod *= kx3 - l;
sum3 += prod;
}
for (unsigned int j = 0; j < i[3]; ++j)
factor[3] *= kx3 - j;
for (unsigned int m = 0; m < 3; ++m)
{
out[n][0][m] = sum3;
for (unsigned int j = 0; j < i[m]; ++j)
{
R prod = factor[3];
for (unsigned int p = 0; p < 3; ++p)
{
if (m == p)
for (unsigned int l = 0; l < i[p]; ++l)
if (j == l)
prod *= R(k) / (i[p]-l);
else
prod *= (kx[p]-l) / (i[p]-l);
else
prod *= factor[p];
}
out[n][0][m] += prod;
}
}
n++;
}
}
}
}
//! \brief Polynomial order of the shape functions
unsigned int order () const
{
return k;
}
};
//Specialization for k=0
template<class D, class R>
class Pk3DLocalBasis<D,R,0>
{
public:
typedef LocalBasisTraits<D,3,Dune::FieldVector<D,3>,R,1,Dune::FieldVector<R,1>,
Dune::FieldMatrix<R,1,3> > Traits;
/** \brief Export the number of degrees of freedom */
enum {N = 1};
/** \brief Export the element order */
enum {O = 0};
unsigned int size () const
{
return 1;
}
inline void evaluateFunction (const typename Traits::DomainType& in,
std::vector<typename Traits::RangeType>& out) const
{
out.resize(1);
out[0] = 1;
}
// evaluate derivative of a single component
inline void
evaluateJacobian (const typename Traits::DomainType& in, // position
std::vector<typename Traits::JacobianType>& out) const // return value
{
out.resize(1);
out[0][0][0] = 0;
out[0][0][1] = 0;
out[0][0][2] = 0;
}
// local interpolation of a function
template<typename E, typename F, typename C>
void interpolate (const E& e, const F& f, std::vector<C>& out) const
{
typename Traits::DomainType x;
typename Traits::RangeType y;
x[0] = 1.0/4.0;
x[1] = 1.0/4.0;
x[2] = 1.0/4.0;
f.eval_local(e,x,y);
out[0] = y;
}
unsigned int order () const
{
return 0;
}
};
}
#endif
|