This file is indexed.

/usr/include/deal.II/base/polynomials_rannacher_turek.h is in libdeal.ii-dev 8.4.2-2+b1.

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
// ---------------------------------------------------------------------
//
// Copyright (C) 2015 - 2016 by the deal.II authors
//
// This file is part of the deal.II library.
//
// The deal.II library is free software; you can use it, redistribute
// it, and/or modify it under the terms of the GNU Lesser General
// Public License as published by the Free Software Foundation; either
// version 2.1 of the License, or (at your option) any later version.
// The full text of the license can be found in the file LICENSE at
// the top level of the deal.II distribution.
//
// ---------------------------------------------------------------------


#ifndef dealii__polynomials_rannacher_turek_h
#define dealii__polynomials_rannacher_turek_h

#include <deal.II/base/point.h>
#include <deal.II/base/tensor.h>
#include <vector>

DEAL_II_NAMESPACE_OPEN


/**
 * Basis for polynomial space on the unit square used for lowest order
 * Rannacher Turek element.
 *
 * The i-th basis function is the dual basis element corresponding to the dof
 * which evaluates the function's mean value across the i-th face. The
 * numbering can be found in GeometryInfo.
 *
 * @ingroup Polynomials
 * @author Patrick Esser
 * @date 2015
 */
template <int dim>
class PolynomialsRannacherTurek
{
public:
  /**
   * Dimension we are working in.
   */
  static const unsigned int dimension = dim;

  /**
   * Constructor, checking that the basis is implemented in this dimension.
   */
  PolynomialsRannacherTurek();

  /**
   * Value of basis function @p i at @p p.
   */
  double compute_value(const unsigned int i,
                       const Point<dim> &p) const;

  /**
   * <tt>order</tt>-th of basis function @p i at @p p.
   *
   * Consider using compute() instead.
   */
  template <int order>
  Tensor<order,dim> compute_derivative (const unsigned int i,
                                        const Point<dim> &p) const;

  /**
   * Gradient of basis function @p i at @p p.
   */
  Tensor<1, dim> compute_grad(const unsigned int i,
                              const Point<dim> &p) const;

  /**
   * Gradient of gradient of basis function @p i at @p p.
   */
  Tensor<2, dim> compute_grad_grad(const unsigned int i,
                                   const Point<dim> &p) const;

  /**
   * Compute values and derivatives of all basis functions at @p unit_point.
   *
   * Size of the vectors must be either equal to the number of polynomials or
   * zero. A size of zero means that we are not computing the vector entries.
   */
  void compute(const Point<dim> &unit_point,
               std::vector<double> &values,
               std::vector<Tensor<1, dim> > &grads,
               std::vector<Tensor<2,dim> > &grad_grads,
               std::vector<Tensor<3,dim> > &third_derivatives,
               std::vector<Tensor<4,dim> > &fourth_derivatives) const;
};


namespace internal
{
  namespace PolynomialsRannacherTurek
  {
    template <int order, int dim>
    inline
    Tensor<order,dim>
    compute_derivative (const unsigned int,
                        const Point<dim> &)
    {
      Assert (dim == 2, ExcNotImplemented());
      return Tensor<order,dim>();
    }


    template <int order>
    inline
    Tensor<order,2>
    compute_derivative (const unsigned int i,
                        const Point<2> &p)
    {
      const unsigned int dim = 2;

      Tensor<order,dim> derivative;
      switch (order)
        {
        case 1:
        {
          Tensor<1,dim> &grad = *reinterpret_cast<Tensor<1,dim>*>(&derivative);
          if (i == 0)
            {
              grad[0] = -2.5 + 3*p(0);
              grad[1] = 1.5 - 3*p(1);
            }
          else if (i == 1)
            {
              grad[0] = -0.5 + 3.0*p(0);
              grad[1] = 1.5 - 3.0*p(1);
            }
          else if (i == 2)
            {
              grad[0] = 1.5 - 3.0*p(0);
              grad[1] = -2.5 + 3.0*p(1);
            }
          else if (i == 3)
            {
              grad[0] = 1.5 - 3.0*p(0);
              grad[1] = -0.5 + 3.0*p(1);
            }
          else
            {
              Assert(false, ExcNotImplemented());
            }
          return derivative;
        }
        case 2:
        {
          Tensor<2,dim> &grad_grad = *reinterpret_cast<Tensor<2,dim>*>(&derivative);
          if (i == 0)
            {
              grad_grad[0][0] = 3;
              grad_grad[0][1] = 0;
              grad_grad[1][0] = 0;
              grad_grad[1][1] = -3;
            }
          else if (i == 1)
            {
              grad_grad[0][0] = 3;
              grad_grad[0][1] = 0;
              grad_grad[1][0] = 0;
              grad_grad[1][1] = -3;
            }
          else if (i == 2)
            {
              grad_grad[0][0] = -3;
              grad_grad[0][1] = 0;
              grad_grad[1][0] = 0;
              grad_grad[1][1] = 3;
            }
          else if (i == 3)
            {
              grad_grad[0][0] = -3;
              grad_grad[0][1] = 0;
              grad_grad[1][0] = 0;
              grad_grad[1][1] = 3;
            }
          return derivative;
        }
        default:
        {
          // higher derivatives are all zero
          return Tensor<order,dim>();
        }
        }
    }
  }
}



// template functions
template <int dim>
template <int order>
Tensor<order,dim>
PolynomialsRannacherTurek<dim>::compute_derivative (const unsigned int i,
                                                    const Point<dim> &p) const
{
  return internal::PolynomialsRannacherTurek::compute_derivative<order> (i, p);
}


DEAL_II_NAMESPACE_CLOSE

#endif