/usr/include/dune/localfunctions/lagrange/pk2d/pk2dlocalbasis.hh is in libdune-localfunctions-dev 2.4.1-1.
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 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 | // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#ifndef DUNE_PK2DLOCALBASIS_HH
#define DUNE_PK2DLOCALBASIS_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 triangle.
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 Pk2DLocalBasis
{
public:
/** \brief Export the number of degrees of freedom */
enum {N = (k+1)*(k+2)/2};
/** \brief Export the element order
OS: Surprising that we need to export this both statically and dynamically!
*/
enum {O = k};
typedef LocalBasisTraits<D,2,Dune::FieldVector<D,2>,R,1,Dune::FieldVector<R,1>,
Dune::FieldMatrix<R,1,2>, 2 > Traits;
//! \brief Standard constructor
Pk2DLocalBasis ()
{
for (unsigned int i=0; i<=k; i++)
pos[i] = (1.0*i)/std::max(k,(unsigned int)1);
}
//! \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);
// specialization for k==0, not clear whether that is needed
if (k==0) {
out[0] = 1;
return;
}
int n=0;
for (unsigned int j=0; j<=k; j++)
for (unsigned int i=0; i<=k-j; i++)
{
out[n] = 1.0;
for (unsigned int alpha=0; alpha<i; alpha++)
out[n] *= (x[0]-pos[alpha])/(pos[i]-pos[alpha]);
for (unsigned int beta=0; beta<j; beta++)
out[n] *= (x[1]-pos[beta])/(pos[j]-pos[beta]);
for (unsigned int gamma=i+j+1; gamma<=k; gamma++)
out[n] *= (pos[gamma]-x[0]-x[1])/(pos[gamma]-pos[i]-pos[j]);
n++;
}
}
//! \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);
// specialization for k==0, not clear whether that is needed
if (k==0) {
out[0][0][0] = 0; out[0][0][1] = 0;
return;
}
int n=0;
for (unsigned int j=0; j<=k; j++)
for (unsigned int i=0; i<=k-j; i++)
{
// x_0 derivative
out[n][0][0] = 0.0;
R factor=1.0;
for (unsigned int beta=0; beta<j; beta++)
factor *= (x[1]-pos[beta])/(pos[j]-pos[beta]);
for (unsigned int a=0; a<i; a++)
{
R product=factor;
for (unsigned int alpha=0; alpha<i; alpha++)
if (alpha==a)
product *= D(1)/(pos[i]-pos[alpha]);
else
product *= (x[0]-pos[alpha])/(pos[i]-pos[alpha]);
for (unsigned int gamma=i+j+1; gamma<=k; gamma++)
product *= (pos[gamma]-x[0]-x[1])/(pos[gamma]-pos[i]-pos[j]);
out[n][0][0] += product;
}
for (unsigned int c=i+j+1; c<=k; c++)
{
R product=factor;
for (unsigned int alpha=0; alpha<i; alpha++)
product *= (x[0]-pos[alpha])/(pos[i]-pos[alpha]);
for (unsigned int gamma=i+j+1; gamma<=k; gamma++)
if (gamma==c)
product *= -D(1)/(pos[gamma]-pos[i]-pos[j]);
else
product *= (pos[gamma]-x[0]-x[1])/(pos[gamma]-pos[i]-pos[j]);
out[n][0][0] += product;
}
// x_1 derivative
out[n][0][1] = 0.0;
factor = 1.0;
for (unsigned int alpha=0; alpha<i; alpha++)
factor *= (x[0]-pos[alpha])/(pos[i]-pos[alpha]);
for (unsigned int b=0; b<j; b++)
{
R product=factor;
for (unsigned int beta=0; beta<j; beta++)
if (beta==b)
product *= D(1)/(pos[j]-pos[beta]);
else
product *= (x[1]-pos[beta])/(pos[j]-pos[beta]);
for (unsigned int gamma=i+j+1; gamma<=k; gamma++)
product *= (pos[gamma]-x[0]-x[1])/(pos[gamma]-pos[i]-pos[j]);
out[n][0][1] += product;
}
for (unsigned int c=i+j+1; c<=k; c++)
{
R product=factor;
for (unsigned int beta=0; beta<j; beta++)
product *= (x[1]-pos[beta])/(pos[j]-pos[beta]);
for (unsigned int gamma=i+j+1; gamma<=k; gamma++)
if (gamma==c)
product *= -D(1)/(pos[gamma]-pos[i]-pos[j]);
else
product *= (pos[gamma]-x[0]-x[1])/(pos[gamma]-pos[i]-pos[j]);
out[n][0][1] += product;
}
n++;
}
}
//! \brief Evaluate higher derivatives of all shape functions
template<unsigned int dOrder> //order of derivative
inline void evaluate(const std::array<int,dOrder>& directions, //direction of derivative
const typename Traits::DomainType& in, //position
std::vector<typename Traits::RangeType>& out) const //return value
{
out.resize(N);
if (dOrder > Traits::diffOrder)
DUNE_THROW(NotImplemented, "Desired derivative order is not implemented");
if (dOrder==0)
evaluateFunction(in, out);
else if (dOrder==1)
{
int n=0;
for (unsigned int j=0; j<=k; j++)
for (unsigned int i=0; i<=k-j; i++, n++)
{
out[n] = 0.0;
for (unsigned int no1=0; no1 < k; no1++)
{
R factor = lagrangianFactorDerivative(directions[0], no1, i, j, in);
for (unsigned int no2=0; no2 < k; no2++)
if (no1 != no2)
factor *= lagrangianFactor(no2, i, j, in);
out[n] += factor;
}
}
}
else if (dOrder==2)
{
// specialization for k<2, not clear whether that is needed
if (k<2) {
std::fill(out.begin(), out.end(), 0.0);
return;
}
//f = prod_{i} f_i -> dxa dxb f = sum_{i} {dxa f_i sum_{k \neq i} dxb f_k prod_{l \neq k,i} f_l
int n=0;
for (unsigned int j=0; j<=k; j++)
for (unsigned int i=0; i<=k-j; i++, n++)
{
R res = 0.0;
for (unsigned int no1=0; no1 < k; no1++)
{
R factor1 = lagrangianFactorDerivative(directions[0], no1, i, j, in);
for (unsigned int no2=0; no2 < k; no2++)
{
if (no1 == no2)
continue;
R factor2 = factor1*lagrangianFactorDerivative(directions[1], no2, i, j, in);
for (unsigned int no3=0; no3 < k; no3++)
{
if (no3 == no1 || no3 == no2)
continue;
factor2 *= lagrangianFactor(no3, i, j, in);
}
res += factor2;
}
}
out[n] = res;
}
}
}
//! \brief Polynomial order of the shape functions
unsigned int order () const
{
return k;
}
private:
/** \brief Returns a single Lagrangian factor of l_ij evaluated at x */
typename Traits::RangeType lagrangianFactor(const int no, const int i, const int j, const typename Traits::DomainType& x) const
{
if ( no < i)
return (x[0]-pos[no])/(pos[i]-pos[no]);
if (no < i+j)
return (x[1]-pos[no-i])/(pos[j]-pos[no-i]);
return (pos[no+1]-x[0]-x[1])/(pos[no+1]-pos[i]-pos[j]);
}
/** \brief Returns the derivative of a single Lagrangian factor of l_ij evaluated at x
* \param direction Derive in x-direction if this is 0, otherwise derive in y direction
*/
typename Traits::RangeType lagrangianFactorDerivative(const int direction, const int no, const int i, const int j, const typename Traits::DomainType& x) const
{
if ( no < i)
return (direction == 0) ? 1.0/(pos[i]-pos[no]) : 0;
if (no < i+j)
return (direction == 0) ? 0: 1.0/(pos[j]-pos[no-i]);
return -1.0/(pos[no+1]-pos[i]-pos[j]);
}
D pos[k+1]; // positions on the interval
};
}
#endif
|