This file is indexed.

/usr/include/dune/localfunctions/common/virtualinterface.hh is in libdune-localfunctions-dev 2.5.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
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#ifndef DUNE_VIRTUALINTERFACE_HH
#define DUNE_VIRTUALINTERFACE_HH

#include <array>

#include <dune/common/function.hh>

#include <dune/geometry/type.hh>

#include <dune/localfunctions/common/localbasis.hh>
#include <dune/localfunctions/common/localkey.hh>
#include <dune/localfunctions/common/localfiniteelementtraits.hh>

namespace Dune
{

  // forward declaration needed by the helper traits
  template<class DomainType, class RangeType>
  class LocalInterpolationVirtualInterface;

  template<class T>
  class LocalBasisVirtualInterface;

  // -----------------------------------------------------------------
  // Helper traits classes
  // -----------------------------------------------------------------

  /**
   * @brief Construct LocalBasisTraits with one diff order lower
   *
   * @tparam T A LocalBasisTraits class
   */
  template<class T>
  struct LowerOrderLocalBasisTraits
  {
    //! The LocalBasisTraits with one order lower
    typedef LocalBasisTraits<
        typename T::DomainFieldType,
        T::dimDomain,
        typename T::DomainType,
        typename T::RangeFieldType,
        T::dimRange,
        typename T::RangeType,
        typename T::JacobianType,
        T::diffOrder-1> Traits;
  };

  /**
   * @brief Construct LocalBasisTraits with fixed diff order
   *
   * @tparam T A LocalBasisTraits class
   * @tparam order The differentiation order
   */
  template<class T, int order>
  struct FixedOrderLocalBasisTraits
  {
    //! The LocalBasisTraits specified order
    typedef LocalBasisTraits<
        typename T::DomainFieldType,
        T::dimDomain,
        typename T::DomainType,
        typename T::RangeFieldType,
        T::dimRange,
        typename T::RangeType,
        typename T::JacobianType,
        order> Traits;
  };

  /**
   * @brief Return a proper base class for functions to use with LocalInterpolation.
   *
   * @tparam FE A FiniteElement type
   */
  template<class FE>
  class LocalFiniteElementFunctionBase
  {
    typedef typename FE::Traits::LocalBasisType::Traits::DomainType DomainType;
    typedef typename FE::Traits::LocalBasisType::Traits::RangeType RangeType;

    typedef LocalInterpolationVirtualInterface<DomainType, RangeType> Interface;
    typedef typename FE::Traits::LocalInterpolationType Implementation;

  public:

    typedef VirtualFunction<DomainType, RangeType> VirtualFunctionBase;
    typedef Function<const DomainType&, RangeType&> FunctionBase;

    /** \brief Base class type for functions to use with LocalInterpolation
     *
     * This is the VirtualFunction interface class if FE implements the virtual
     * interface and Function base class otherwise.
     */
    typedef typename std::conditional<std::is_base_of<Interface, Implementation>::value, VirtualFunctionBase, FunctionBase>::type type;
  };



  // -----------------------------------------------------------------
  // Basis
  // -----------------------------------------------------------------

  // current versions of doxygen (<= 1.6.2) enter an infinite loop when parsing
  // the following class
#ifndef DOXYGEN
  /**
   * @brief virtual base class for a local basis
   *
   * Provides the local basis interface with pure virtual methods.
   * The class derives from the interface with one differentiation order lower.
   *
   * This class defines the interface using pure virtual methods.
   * In applications you should use the derived class
   * LocalBasisVirtualInterface that also
   * contains an evaluate with order as template parameter.
   *
   * This template method cannot be defined in the same
   * class as the virtual method. Otherwise name resolution fails.
   */
  template<class T>
  class LocalBasisVirtualInterfaceBase :
    public virtual LocalBasisVirtualInterface<typename LowerOrderLocalBasisTraits<T>::Traits>
  {
    typedef LocalBasisVirtualInterface<typename LowerOrderLocalBasisTraits<T>::Traits> BaseInterface;
  public:
    typedef T Traits;

    //! \todo Please doc me!
    virtual void evaluate (
      const typename std::template array<int,Traits::diffOrder>& directions,
      const typename Traits::DomainType& in,
      std::vector<typename Traits::RangeType>& out) const = 0;

    using BaseInterface::evaluate;
  };
#endif // DOXYGEN

  /**
   * @brief virtual base class for a local basis
   *
   * Provides the local basis interface with pure virtual methods.
   * This is the base interface with differentiation order 0.
   */
  template<class DF, int n, class D, class RF, int m, class R, class J>
  class LocalBasisVirtualInterfaceBase<LocalBasisTraits<DF,n,D,RF,m,R,J,0> >
  {
  public:
    typedef LocalBasisTraits<DF,n,D,RF,m,R,J,0> Traits;

    virtual ~LocalBasisVirtualInterfaceBase() {}

    //! \brief Number of shape functions
    virtual unsigned int size () const = 0;

    //! \brief Polynomial order of the shape functions
    virtual unsigned int order () const = 0;

    /** \brief Evaluate all basis function at given position
     *
     * Evaluates all shape functions at the given position and returns
     * these values in a vector.
     */
    virtual void evaluateFunction (const typename Traits::DomainType& in,
                                   std::vector<typename Traits::RangeType>& out) const = 0;

    /** \brief Evaluate jacobian of all shape functions at given position.
     *
     * out[i][j][k] is \f$\partial_k \hat\phi_j^i \f$, when \f$\hat\phi^i \f$ is the
     *		i'th shape function.
     *
     * \param [in]  in  The position where evaluated
     * \param [out] out The result
     */
    virtual void evaluateJacobian(const typename Traits::DomainType& in,         // position
                                  std::vector<typename Traits::JacobianType>& out) const = 0;

    /** \brief Evaluate partial derivatives of any order of all shape functions
     * \param order Order of the partial derivatives, in the classic multi-index notation
     * \param in Position where to evaluate the derivatives
     * \param[out] out Return value: the desired partial derivatives
     */
    virtual void partial(const std::array<unsigned int,n>& order,
                         const typename Traits::DomainType& in,
                         std::vector<typename Traits::RangeType>& out) const = 0;

    //! \todo Please doc me!
    virtual void evaluate (
      const typename std::template array<int,Traits::diffOrder>& directions,
      const typename Traits::DomainType& in,
      std::vector<typename Traits::RangeType>& out) const = 0;

  };

  /**
   * @brief virtual base class for a local basis
   *
   * Provides the local basis interface with pure virtual methods.
   * The class derives from the interface with one differentiation order lower.
   *
   * This class defines the interface using pure virtual methods.
   * It also contains the evaluate method with order as template parameter.
   */
  template<class T>
  class LocalBasisVirtualInterface :
    public virtual LocalBasisVirtualInterfaceBase<T>
  {
    typedef LocalBasisVirtualInterfaceBase<T> BaseInterface;
  public:
    typedef T Traits;

    //! \todo Please doc me!
    template <int k>
    void evaluate (
      const typename std::template array<int,k>& directions,
      const typename Traits::DomainType& in,
      std::vector<typename Traits::RangeType>& out) const
    {
      typedef LocalBasisVirtualInterfaceBase<typename FixedOrderLocalBasisTraits<T,k>::Traits > OrderKBaseInterface;
      const OrderKBaseInterface& asBase = *this;
      asBase.evaluate(directions, in, out);
    }

    using BaseInterface::size;
    using BaseInterface::order;
    using BaseInterface::partial;
    using BaseInterface::evaluateFunction;
    using BaseInterface::evaluateJacobian;
    /* Unfortunately, the intel compiler cannot use the different evaluate
     * methods with varying argument lists. :-( */
#ifndef __INTEL_COMPILER
    using BaseInterface::evaluate;
#endif
  };




  // -----------------------------------------------------------------
  // Interpolation
  // -----------------------------------------------------------------

  /**
   * @brief virtual base class for a local interpolation
   *
   * This class defines the interface using pure virtual methods.
   * In applications you should use the derived class
   * LocalInterpolationVirtualInterface that also
   * contains a interpolate method where the function type
   * is a template parameter.
   *
   * This template method cannot be defined in the same
   * class as the virtual method. Otherwise name resolution fails.
   */
  template<class DomainType, class RangeType>
  class LocalInterpolationVirtualInterfaceBase
  {
  public:

    //! type of virtual function to interpolate
    typedef Dune::VirtualFunction<DomainType, RangeType> FunctionType;

    //! type of the coefficient vector in the interpolate method
    typedef typename RangeType::field_type CoefficientType;

    virtual ~LocalInterpolationVirtualInterfaceBase() {}

    /** \brief determine coefficients interpolating a given function
     *
     * This is the pure virtual method taking a VirtualFunction.
     *
     * \param[in]  f   Function instance used to interpolate.
     * \param[out] out Resulting coefficients vector.
     */
    virtual void interpolate (const FunctionType& f, std::vector<CoefficientType>& out) const = 0;
  };

  /**
   * @brief virtual base class for a local interpolation
   *
   * This class defines the interface using pure virtual methods.
   * It also contains the interpolate method with function
   * type as template parameter.
   */
  template<class DomainType, class RangeType>
  class LocalInterpolationVirtualInterface
    : public LocalInterpolationVirtualInterfaceBase<DomainType, RangeType>
  {
  public:

    //! type of virtual function to interpolate
    typedef Dune::VirtualFunction<DomainType, RangeType> FunctionType;

    //! type of the coefficient vector in the interpolate method
    typedef typename RangeType::field_type CoefficientType;


    virtual ~LocalInterpolationVirtualInterface() {}

    // This method is only notet again for to make the documentation complete.

    /** \brief determine coefficients interpolating a given function
     *
     * This is the pure virtual method taking a VirtualFunction.
     *
     * \param[in]  f   Function instance used to interpolate.
     * \param[out] out Resulting coefficients vector.
     */
    virtual void interpolate (const FunctionType& f, std::vector<CoefficientType>& out) const = 0;

    //! \copydoc LocalInterpolationVirtualInterfaceBase::interpolate
    //! This uses the pure virtual method by wrapping the template argument into a VirtualFunction
    template<class F>
    void interpolate (const F& f, std::vector<CoefficientType>& out) const
    {
      const LocalInterpolationVirtualInterfaceBase<DomainType, RangeType>& asBase = *this;
      asBase.interpolate(VirtualFunctionWrapper<F>(f),out);
    }

    template<class F, class C>
    void interpolate (const F& f, std::vector<C>& out) const
    {
      std::vector<CoefficientType> outDummy;
      const LocalInterpolationVirtualInterfaceBase<DomainType, RangeType>& asBase = *this;
      asBase.interpolate(VirtualFunctionWrapper<F>(f),outDummy);
      out.resize(outDummy.size());
      for(typename std::vector<CoefficientType>::size_type i=0; i<outDummy.size(); ++i)
        out[i] = outDummy[i];
    }

  private:

    template <typename F>
    struct VirtualFunctionWrapper
      : public FunctionType
    {
    public:
      VirtualFunctionWrapper(const F &f)
        : f_(f)
      {}

      virtual ~VirtualFunctionWrapper() {}

      virtual void evaluate(const DomainType& x, RangeType& y) const
      {
        f_.evaluate(x,y);
      }

      const F &f_;
    };
  };



  // -----------------------------------------------------------------
  // Coefficients
  // -----------------------------------------------------------------

  /**
   * @brief virtual base class for local coefficients
   *
   * This class defines the interface using pure virtual methods.
   */
  class LocalCoefficientsVirtualInterface
  {
  public:

    virtual ~LocalCoefficientsVirtualInterface() {}

    //! number of coefficients
    virtual std::size_t size () const = 0;

    //! get i'th index
    const virtual LocalKey& localKey (std::size_t i) const = 0;

  };



  // -----------------------------------------------------------------
  // Finite Element
  // -----------------------------------------------------------------

  /**
   * @brief virtual base class for local finite elements with functions
   *
   * This class defines the same interface using pure virtual methods.
   * The class derives from the interface with one differentiation order lower.
   */
  template<class T>
  class LocalFiniteElementVirtualInterface
    : public virtual LocalFiniteElementVirtualInterface<typename LowerOrderLocalBasisTraits<T>::Traits >
  {
    typedef LocalFiniteElementVirtualInterface<typename LowerOrderLocalBasisTraits<T>::Traits > BaseInterface;

  public:
    typedef LocalFiniteElementTraits<
        LocalBasisVirtualInterface<T>,
        LocalCoefficientsVirtualInterface,
        LocalInterpolationVirtualInterface<
            typename T::DomainType,
            typename T::RangeType> > Traits;

    //! \copydoc LocalFiniteElementVirtualInterface::localBasis
    virtual const typename Traits::LocalBasisType& localBasis () const = 0;

    using BaseInterface::localBasis;
    using BaseInterface::localCoefficients;
    using BaseInterface::localInterpolation;
    using BaseInterface::type;

    virtual LocalFiniteElementVirtualInterface<T>* clone() const = 0;
  };


  /**
   * @brief virtual base class for local finite elements with functions
   *
   * This class defines the same interface using pure virtual methods.
   * This is the base interface with differentiation order 0.
   */
  template<class DF, int n, class D, class RF, int m, class R, class J>
  class LocalFiniteElementVirtualInterface<LocalBasisTraits<DF,n,D,RF,m,R,J,0> >
  {
    typedef LocalBasisTraits<DF,n,D,RF,m,R,J,0> T;

  public:
    typedef LocalFiniteElementTraits<
        LocalBasisVirtualInterface<T>,
        LocalCoefficientsVirtualInterface,
        LocalInterpolationVirtualInterface<
            typename T::DomainType,
            typename T::RangeType> > Traits;

    virtual ~LocalFiniteElementVirtualInterface() {}

    //! \copydoc LocalFiniteElementVirtualInterface::localBasis
    virtual const typename Traits::LocalBasisType& localBasis () const = 0;

    //! \copydoc LocalFiniteElementVirtualInterface::localCoefficients
    virtual const typename Traits::LocalCoefficientsType& localCoefficients () const = 0;

    //! \copydoc LocalFiniteElementVirtualInterface::localInterpolation
    virtual const typename Traits::LocalInterpolationType& localInterpolation () const = 0;

    //! \copydoc LocalFiniteElementVirtualInterface::size
    virtual unsigned int size () const = 0;

    //! \copydoc LocalFiniteElementVirtualInterface::type
    virtual const GeometryType type () const = 0;

    virtual LocalFiniteElementVirtualInterface<T>* clone() const = 0;
  };

}
#endif