This file is indexed.

/usr/include/dune/localfunctions/lagrange/pk1d/pk1dlocalbasis.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
#ifndef DUNE_Pk1DLOCALBASIS_HH
#define DUNE_Pk1DLOCALBASIS_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 1D 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 Pk1DLocalBasis
    {
        public:

        /** \brief Export the number of degrees of freedom */
        enum {N = k+1};

        /** \brief Export the element order */
        enum {O = k};

        typedef LocalBasisTraits<D,
                                 1,
                                 Dune::FieldVector<D,1>,
                                 R,
                                 1,
                                 Dune::FieldVector<R,1>,
                                 Dune::FieldMatrix<R,1,1>
                                 > Traits;

        //! \brief Standard constructor
        Pk1DLocalBasis ()
        {
            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);

            for (unsigned int i=0; i<N; i++)
            {
                out[i] = 1.0;
                for (unsigned int alpha=0; alpha<i; alpha++) 
                    out[i] *= (x[0]-pos[alpha])/(pos[i]-pos[alpha]);
                for (unsigned int gamma=i+1; gamma<=k; gamma++) 
                    out[i] *= (x[0]-pos[gamma])/(pos[i]-pos[gamma]);
            }
        }

        //! \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);

            for (unsigned int i=0; i<=k; i++) {

                // x_0 derivative
                out[i][0][0] = 0.0;
                R factor=1.0;
                for (unsigned int a=0; a<i; a++)
                {
                    R product=factor;
                    for (unsigned int alpha=0; alpha<i; alpha++)
                        product *= (alpha==a) ? 1.0/(pos[i]-pos[alpha])
                                              : (x[0]-pos[alpha])/(pos[i]-pos[alpha]);
                    for (unsigned int gamma=i+1; gamma<=k; gamma++) 
                        product *= (pos[gamma]-x[0])/(pos[gamma]-pos[i]);
                    out[i][0][0] += product;
                }
                for (unsigned int c=i+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+1; gamma<=k; gamma++) 
                        product *= (gamma==c) ? -1.0/(pos[gamma]-pos[i])
                                              : (pos[gamma]-x[0])/(pos[gamma]-pos[i]);
                    out[i][0][0] += product;
                }
            }
            
        }

        //! \brief Polynomial order of the shape functions
        unsigned int order () const
        {
            return k;
        }

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

}
#endif