This file is indexed.

/usr/include/dune/localfunctions/refined/common/refinedsimplexlocalbasis.hh is in libdune-localfunctions-dev 2.3.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
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#ifndef DUNE_REFINED_SIMPLEX_LOCALBASIS_HH
#define DUNE_REFINED_SIMPLEX_LOCALBASIS_HH

/** \file
    \brief Contains a base class for LocalBasis classes based on uniform refinement
 */

#include <dune/common/fvector.hh>
#include <dune/common/exceptions.hh>
#include <dune/localfunctions/common/localbasis.hh>

namespace Dune
{
  template<class D, int dim>
  class RefinedSimplexLocalBasis
  {
  protected:
    RefinedSimplexLocalBasis()
    {
      DUNE_THROW(Dune::NotImplemented,"RefinedSimplexLocalBasis not implemented for dim > 3.");
    }
  };

  /**@ingroup LocalBasisImplementation
     \brief Base class for LocalBasis classes based on uniform refinement in 1D; provides numbering and local coordinates of subelements

     \tparam D Type to represent the field in the domain.

     \nosubgrouping
   */
  template<class D>
  class RefinedSimplexLocalBasis<D,1>
  {
  protected:

    /** \brief Protected default constructor so this class can only be instantiated as a base class. */
    RefinedSimplexLocalBasis() {}

    /** \brief Get the number of the subelement containing a given point.
     *
     * The subelements are ordered according to
     *
     *     0       1
     * |-------:-------|
     *
     * \param[in] global Coordinates in the reference element
     * \returns Number of the subtriangle containing <tt>global</tt>
     */
    static int getSubElement(const FieldVector<D,1>& global)
    {
      if (global[0] <= 0.5)
        return 0;
      else if (global[0] <= 1.0)
        return 1;

      DUNE_THROW(InvalidStateException, "no subelement defined");
    }

    /** \brief Get local coordinates in the subelement

       \param[in] global Coordinates in the reference element
       \param[out] subElement Number of the subelement containing <tt>global</tt>
       \param[out] local The local coordinates in the subelement
     */
    static void getSubElement(const FieldVector<D,1>& global,
                              int& subElement,
                              FieldVector<D,1>& local)
    {
      if (global[0] <= 0.5) {
        subElement = 0;
        local[0] = 2.0 * global[0];
        return;
      }

      subElement = 1;
      local[0] = 2.0 * global[0] - 1.0;
    }

  };


  /**@ingroup LocalBasisImplementation
     \brief Base class for LocalBasis classes based on uniform refinement in 2D; provides numbering and local coordinates of subelements

     Shape functions like these are necessary for hierarchical error estimators
     for certain nonlinear problems.

     \tparam D Type to represent the field in the domain.

     \nosubgrouping
   */
  template<class D>
  class RefinedSimplexLocalBasis<D,2>
  {
  protected:

    /** \brief Protected default constructor so this class can only be instantiated as a base class. */
    RefinedSimplexLocalBasis() {}

    /** \brief Get the number of the subtriangle containing a given point.
     *
     * The triangles are ordered according to
     * \verbatim
     |\
     |\|2\
     |\|--\
     |\|\3|\
     |\|0\|1\
       ------ \endverbatim
     *
     * \param[in] global Coordinates in the reference triangle
     * \returns Number of the subtriangle containing <tt>global</tt>
     */
    static int getSubElement(const FieldVector<D,2>& global)
    {
      if (global[0] + global[1] <= 0.5)
        return 0;
      else if (global[0] >= 0.5)
        return 1;
      else if (global[1] >= 0.5)
        return 2;

      return 3;
    }

    /** \brief Get local coordinates in the subtriangle

       \param[in] global Coordinates in the reference triangle
       \param[out] subElement Number of the subtriangle containing <tt>global</tt>
       \param[out] local The local coordinates in the subtriangle
     */
    static void getSubElement(const FieldVector<D,2>& global,
                              int& subElement,
                              FieldVector<D,2>& local)
    {
      if (global[0] + global[1] <= 0.5) {
        subElement = 0;
        local[0] = 2*global[0];
        local[1] = 2*global[1];
        return;
      } else if (global[0] >= 0.5) {
        subElement = 1;
        local[0] = 2*global[0]-1;
        local[1] = 2*global[1];
        return;
      } else if (global[1] >= 0.5) {
        subElement = 2;
        local[0] = 2*global[0];
        local[1] = 2*global[1]-1;
        return;
      }

      subElement = 3;
      local[0] = -2 * global[0] + 1;
      local[1] = -2 * global[1] + 1;

    }


  };

  /**@ingroup LocalBasisImplementation
     \brief Base class for LocalBasis classes based on uniform refinement in 3D; provides numbering and local coordinates of subelements

     Shape functions like these are necessary for hierarchical error estimators
     for certain nonlinear problems.

     \tparam D Type to represent the field in the domain.

     \nosubgrouping
   */
  template<class D>
  class RefinedSimplexLocalBasis<D,3>
  {
  protected:

    /** \brief Protected default constructor so this class can only be instantiated as a base class. */
    RefinedSimplexLocalBasis() {}

    /** \brief Get the number of the subsimplex containing a given point in the reference element
     *
     * Defining the following points in the reference simplex
     *
     * 0: (0.0, 0.0, 0.0)
     * 1: (1.0, 0.0, 0.0)
     * 2: (0.0, 1.0, 0.0)
     * 3: (0.0, 0.0, 1.0)
     * 4: (0.5, 0.0, 0.0)
     * 5: (0.5, 0.5, 0.0)
     * 6: (0.0, 0.5, 0.0)
     * 7: (0.0, 0.0, 0.5)
     * 8: (0.5, 0.0, 0.5)
     * 9: (0.0, 0.5, 0.5)
     *
     * The subsimplices are numbered according to
     *
     * 0: 0467  -
     * 1: 4158   |_ "cut off" vertices
     * 2: 6529   |
     * 3: 7893  -
     *
     * 4: 6487  -
     * 5: 4568   |_  octahedron partition
     * 6: 6897   |
     * 7: 6895  -
     *
     * \param[in] global Coordinates in the reference simplex
     * \returns Number of the subsimplex containing <tt>global</tt>
     */
    static int getSubElement(const FieldVector<D,3>& global)
    {
      if (global[0] + global[1] + global[2] <= 0.5)
        return 0;
      else if (global[0] >= 0.5)
        return 1;
      else if (global[1] >= 0.5)
        return 2;
      else if (global[2] >= 0.5)
        return 3;
      else if ((global[0] + global[1] <= 0.5)and (global[1] + global[2] <= 0.5))
        return 4;
      else if ((global[0] + global[1] >= 0.5)and (global[1] + global[2] <= 0.5))
        return 5;
      else if ((global[0] + global[1] <= 0.5)and (global[1] + global[2] >= 0.5))
        return 6;
      else if ((global[0] + global[1] >= 0.5)and (global[1] + global[2] >= 0.5))
        return 7;

      DUNE_THROW(InvalidStateException, "no subelement defined");

    }
    /** \brief Get local coordinates in the subsimplex

       \param[in] global Coordinates in the reference simplex
       \param[out] subElement Number of the subsimplex containing <tt>global</tt>
       \param[out] local The local coordinates in the subsimplex
     */
    static void getSubElement(const FieldVector<D,3>& global,
                              int& subElement,
                              FieldVector<D,3>& local)
    {
      if (global[0] + global[1] + global[2] <= 0.5) {
        subElement = 0;
        local = global;
        local *= 2.0;
        return;
      } else if (global[0] >= 0.5) {
        subElement = 1;
        local = global;
        local[0] -= 0.5;
        local *= 2.0;
        return;
      } else if (global[1] >= 0.5) {
        subElement = 2;
        local = global;
        local[1] -= 0.5;
        local *= 2.0;
        return;
      } else if (global[2] >= 0.5) {
        subElement = 3;
        local = global;
        local[2] -= 0.5;
        local *= 2.0;
        return;
      } else if ((global[0] + global[1] <= 0.5)and (global[1] + global[2] <= 0.5)) {
        subElement = 4;
        local[0] = 2.0 * global[1];
        local[1] = 2.0 * (0.5 - global[0] - global[1]);
        local[2] = 2.0 * (-0.5 + global[0] + global[1] + global[2]);
        //              Dune::FieldMatrix<double,3,3> A(0.0);
        //              A[0][1] =  2.0;
        //              A[1][0] = -2.0;
        //              A[1][1] = -2.0;
        //              A[2][0] =  2.0;
        //              A[2][1] =  2.0;
        //              A[2][2] =  2.0;
        //              A.mv(global,local);
        //              local[1] += 1.0;
        //              local[2] -= 1.0;
        return;
      } else if ((global[0] + global[1] >= 0.5)and (global[1] + global[2] <= 0.5)) {
        subElement = 5;
        local[0] = 2.0 * (0.5 - global[0]);
        local[1] = 2.0 * (0.5 - global[1] - global[2]);
        local[2] = 2.0 * global[2];
        //              Dune::FieldMatrix<double,3,3> A(0.0);
        //              A[0][0] = -2.0;
        //              A[1][1] = -2.0;
        //              A[1][2] = -2.0;
        //              A[2][2] =  2.0;
        //              A.mv(global,local);
        //              local[0] += 1.0;
        //              local[1] += 1.0;
        return;
      } else if ((global[0] + global[1] <= 0.5)and (global[1] + global[2] >= 0.5)) {
        subElement = 6;
        local[0] = 2.0 * (0.5 - global[0] - global[1]);
        local[1] = 2.0 * global[0];
        local[2] = 2.0 * (-0.5 + global[1] + global[2]);
        //              Dune::FieldMatrix<double,3,3> A(0.0);
        //              A[0][0] = -2.0;
        //              A[0][1] = -2.0;
        //              A[1][0] =  2.0;
        //              A[2][1] =  2.0;
        //              A[2][2] =  2.0;
        //              A.mv(global,local);
        //              local[0] += 1.0;
        //              local[2] -= 1.0;
        return;
      } else if ((global[0] + global[1] >= 0.5)and (global[1] + global[2] >= 0.5)) {
        subElement = 7;
        local[0] = 2.0 * (-0.5 + global[1] + global[2]);
        local[1] = 2.0 * (0.5 - global[1]);
        local[2] = 2.0 * (-0.5 + global[0] + global[1]);
        //              Dune::FieldMatrix<double,3,3> A(0.0);
        //              A[0][1] =  2.0;
        //              A[0][2] =  2.0;
        //              A[1][1] = -2.0;
        //              A[2][0] =  2.0;
        //              A[2][1] =  2.0;
        //              A.mv(global,local);
        //              local[0] -= 1.0;
        //              local[1] += 1.0;
        //              local[2] -= 1.0;
        return;
      }

      DUNE_THROW(InvalidStateException, "no subelement defined");

    }

  };


}

#endif