This file is indexed.

/usr/include/dune/localfunctions/lagrange/pk2d/pk2dlocalbasis.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
#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> > 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 *= 1.0/(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 *= -1.0/(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 *= 1.0/(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 *= -1.0/(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 Polynomial order of the shape functions
  unsigned int order () const
  {
	return k;
  }

private:
  R pos[k+1]; // positions on the interval
};

}
#endif