This file is indexed.

/usr/include/deal.II/hp/dof_faces.h is in libdeal.ii-dev 8.1.0-4.

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
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
// ---------------------------------------------------------------------
// $Id: dof_faces.h 31932 2013-12-08 02:15:54Z heister $
//
// Copyright (C) 2006 - 2013 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 __deal2__hp_dof_faces_h
#define __deal2__hp_dof_faces_h


#include <deal.II/base/config.h>
#include <deal.II/hp/fe_collection.h>

#include <vector>

DEAL_II_NAMESPACE_OPEN

namespace hp
{
  template <int dim, int spacedim>
  class FECollection;
}


namespace internal
{
  namespace hp
  {
    /**
     * Store the indices of the degrees of freedom which are located on
     * objects of dimension @p structdim < dim, i.e., for faces or edges
     * of cells. This is opposed to the internal::hp::DoFLevels class
     * that stores the DoF indices on cells.
     *
     * The things we store here is very similar to what is stored in the
     * internal::DoFHandler::DoFObjects classes (see there for more
     * information, in particular on the layout of the class hierarchy,
     * and the use of file names).
     *
     * <h4>Offset computations</h4>
     *
     * For hp methods, not all cells may use the same finite element, and
     * it is consequently more complicated to determine where the DoF
     * indices for a given line, quad, or hex are stored. As described in
     * the documentation of the internal::DoFHandler::DoFLevel class, we
     * can compute the location of the first line DoF, for example, by
     * calculating the offset as <code>line_index *
     * dof_handler.get_fe().dofs_per_line</code>. This of course doesn't
     * work any more if different lines may have different numbers of
     * degrees of freedom associated with them. Consequently, rather than
     * using this simple multiplication, the dofs array has an associated
     * array dof_offsets. The data corresponding to a
     * line then starts at index <code>line_dof_offsets[line_index]</code>
     * within the <code>line_dofs</code> array.
     *
     *
     * <h4>Multiple data sets per object</h4>
     *
     * If two adjacent cells use different finite elements, then
     * the face that they share needs to store DoF indices for both
     * involved finite elements. While faces therefore have to have at
     * most two sets of DoF indices, it is easy to see that edges and
     * vertices can have as many sets of DoF indices associated with them
     * as there are adjacent cells.
     *
     * Consequently, for objects that have a lower dimensionality than
     * cells, we have to store a map from the finite element index to the
     * set of DoF indices associated. Since real sets are typically very
     * inefficient to store, and since most of the time we expect the
     * number of individual keys to be small (frequently, adjacent cells
     * will have the same finite element, and only a single entry will
     * exist in the map), what we do is instead to store a linked list. In
     * this format, the first entry starting at position
     * <code>lines.dofs[lines.dof_offsets[line_index]]</code> will denote
     * the finite element index of the set of DoF indices following; after
     * this set, we will store the finite element index of the second set
     * followed by the corresponding DoF indices; and so on. Finally, when
     * all finite element indices adjacent to this object have been
     * covered, we write a -1 to indicate the end of the list.
     *
     * Access to this kind of data, as well as the distinction between
     * cells and objects of lower dimensionality are encoded in the
     * accessor functions, DoFObjects::set_dof_index() and
     * DoFLevel::get_dof_index(). They are able to traverse this
     * list and pick out or set a DoF index given the finite element index
     * and its location within the set of DoFs corresponding to this
     * finite element.
     *
     *
     * @ingroup hp
     * @author Tobias Leicht, 2006
     */
    template <int structdim>
    class DoFIndicesOnFacesOrEdges
    {
    public:
      /**
       * Store the start index for
       * the degrees of freedom of each
       * object in the @p dofs array.
       *
       * The type we store is then obviously the type the @p dofs array
       * uses for indexing.
       */
      std::vector<unsigned int> dof_offsets;

      /**
       * Store the global indices of
       * the degrees of freedom. See
       * DoFLevel() for detailed
       * information.
       */
      std::vector<types::global_dof_index> dofs;

      /**
       * Set the global index of
       * the @p local_index-th
       * degree of freedom located
       * on the object with number @p
       * obj_index to the value
       * given by @p global_index. The @p
       * dof_handler argument is
       * used to access the finite
       * element that is to be used
       * to compute the location
       * where this data is stored.
       *
       * The third argument, @p
       * fe_index, denotes which of
       * the finite elements
       * associated with this
       * object we shall
       * access. Refer to the
       * general documentation of
       * the internal::hp::DoFLevel
       * class template for more
       * information.
       */
      template <int dim, int spacedim>
      void
      set_dof_index (const dealii::hp::DoFHandler<dim,spacedim> &dof_handler,
                     const unsigned int               obj_index,
                     const unsigned int               fe_index,
                     const unsigned int               local_index,
                     const types::global_dof_index    global_index,
                     const unsigned int               obj_level);

      /**
       * Return the global index of
       * the @p local_index-th
       * degree of freedom located
       * on the object with number @p
       * obj_index. The @p
       * dof_handler argument is
       * used to access the finite
       * element that is to be used
       * to compute the location
       * where this data is stored.
       *
       * The third argument, @p
       * fe_index, denotes which of
       * the finite elements
       * associated with this
       * object we shall
       * access. Refer to the
       * general documentation of
       * the internal::hp::DoFLevel
       * class template for more
       * information.
       */
      template <int dim, int spacedim>
      types::global_dof_index
      get_dof_index (const dealii::hp::DoFHandler<dim,spacedim> &dof_handler,
                     const unsigned int               obj_index,
                     const unsigned int               fe_index,
                     const unsigned int               local_index,
                     const unsigned int               obj_level) const;

      /**
       * Return the number of
       * finite elements that are
       * active on a given
       * object. If this is a cell,
       * the answer is of course
       * one. If it is a face, the
       * answer may be one or two,
       * depending on whether the
       * two adjacent cells use the
       * same finite element or
       * not. If it is an edge in
       * 3d, the possible return
       * value may be one or any
       * other value larger than
       * that.
       *
       * If the object is not part
       * of an active cell, then no
       * degrees of freedom have
       * been distributed and zero
       * is returned.
       */
      template <int dim, int spacedim>
      unsigned int
      n_active_fe_indices (const dealii::hp::DoFHandler<dim,spacedim> &dof_handler,
                           const unsigned int               obj_index) const;

      /**
       * Return the fe_index of the
       * n-th active finite element
       * on this object.
       */
      template <int dim, int spacedim>
      types::global_dof_index
      nth_active_fe_index (const dealii::hp::DoFHandler<dim,spacedim> &dof_handler,
                           const unsigned int               obj_level,
                           const unsigned int               obj_index,
                           const unsigned int               n) const;

      /**
       * Check whether a given
       * finite element index is
       * used on the present
       * object or not.
       */
      template <int dim, int spacedim>
      bool
      fe_index_is_active (const dealii::hp::DoFHandler<dim,spacedim> &dof_handler,
                          const unsigned int               obj_index,
                          const unsigned int               fe_index,
                          const unsigned int               obj_level) const;

      /**
       * Determine an estimate for the
       * memory consumption (in bytes)
       * of this object.
       */
      std::size_t memory_consumption () const;
    };



    /**
     * These classes are similar to the internal::hp::DoFLevel classes. We here store
     * information that is associated with faces, rather than cells, as this information is
     * independent of the hierarchical structure of cells, which are organized in levels. In 2D
     * we store information on degrees of freedom located on lines whereas in 3D we store
     * information on drefrees of freedom located on quads and lines. In 1D we do nothing, as
     * the faces of lines are vertices which are treated separately.
     *
     * Apart from the internal::hp::DoFObjects object containing the data to store
     * (degree of freedom indices) and all the access-functionality to this data, we do not
     * store any data or provide any functionality. However, we do implement a function to
     * determine an estimate of the memory consumption of the contained
     * internal::hp::DoFObjects object(s).
     *
     * The data contained isn't usually directly accessed. Rather, except for some access from
     * the DoFHandler class, access is usually through the DoFAccessor::set_dof_index() and
     * DoFAccessor::dof_index() functions or similar functions of derived classes that in turn
     * access the member variables using the DoFHandler::get_dof_index() and corresponding
     * setter functions. Knowledge of the actual data format is therefore encapsulated to the
     * present hierarchy of classes as well as the ::DoFHandler class.
     *
     * @ingroup dofs
     * @author Tobias Leicht, 2006
     */
    template<int dim>
    class DoFIndicesOnFaces;


    /**
     * Store the indices of degrees of freedom on faces in 1D. As these would be vertices, which
     * are treated separately, don't do anything.
     *
     * @ingroup hp
     * @author Tobias Leicht, 2006
     */
    template<>
    class DoFIndicesOnFaces<1>
    {
    public:
      /**
       * Determine an estimate for the
       * memory consumption (in bytes)
       * of this object.
       */
      std::size_t memory_consumption () const;
    };

    /**
     * Store the indices of degrees of freedom on faces in 2D, which are lines.
     *
     * @ingroup hp
     * @author Tobias Leicht, 2006
     */
    template<>
    class DoFIndicesOnFaces<2>
    {
    public:
      /**
       * Indices of DoFs on the lines that bound cells.
       */
      internal::hp::DoFIndicesOnFacesOrEdges<1> lines;

      /**
       * Determine an estimate for the
       * memory consumption (in bytes)
       * of this object.
       */
      std::size_t memory_consumption () const;
    };

    /**
     * Store the indices of degrees of freedom on faces in 3D, which are
     * quads, additionally also on lines.
     *
     * @ingroup hp
     * @author Tobias Leicht, 2006
     */
    template<>
    class DoFIndicesOnFaces<3>
    {
    public:
      /**
       * Indices of DoFs on the lines that form the edges of cells.
       */
      internal::hp::DoFIndicesOnFacesOrEdges<1> lines;

      /**
       * Indices of DoFs on the quads that bound cells.
       */
      internal::hp::DoFIndicesOnFacesOrEdges<2> quads;

      /**
       * Determine an estimate for the
       * memory consumption (in bytes)
       * of this object.
       */
      std::size_t memory_consumption () const;
    };


    // --------------------- inline and template functions ------------------

    template <int structdim>
    template <int dim, int spacedim>
    inline
    types::global_dof_index
    DoFIndicesOnFacesOrEdges<structdim>::
    get_dof_index (const dealii::hp::DoFHandler<dim,spacedim> &dof_handler,
                   const unsigned int                obj_index,
                   const unsigned int                fe_index,
                   const unsigned int                local_index,
                   const unsigned int                obj_level) const
    {
      Assert ((fe_index != dealii::hp::DoFHandler<dim,spacedim>::default_fe_index),
              ExcMessage ("You need to specify a FE index when working "
                          "with hp DoFHandlers"));
      Assert (&dof_handler != 0,
              ExcMessage ("No DoFHandler is specified for this iterator"));
      Assert (&dof_handler.get_fe() != 0,
              ExcMessage ("No finite element collection is associated with "
                          "this DoFHandler"));
      Assert (fe_index < dof_handler.get_fe().size(),
              ExcIndexRange (fe_index, 0, dof_handler.get_fe().size()));
      Assert (local_index <
              dof_handler.get_fe()[fe_index].template n_dofs_per_object<structdim>(),
              ExcIndexRange(local_index, 0,
                            dof_handler.get_fe()[fe_index]
                            .template n_dofs_per_object<structdim>()));
      Assert (obj_index < dof_offsets.size(),
              ExcIndexRange (obj_index, 0, dof_offsets.size()));

      // make sure we are on an
      // object for which DoFs have
      // been allocated at all
      Assert (dof_offsets[obj_index] != numbers::invalid_unsigned_int,
              ExcMessage ("You are trying to access degree of freedom "
                          "information for an object on which no such "
                          "information is available"));

      Assert (structdim<dim, ExcMessage ("This object can not be used for cells."));

      // there may be multiple finite elements associated with
      // this object. hop along the list of index sets until we
      // find the one with the correct fe_index, and then poke
      // into that part. trigger an exception if we can't find a
      // set for this particular fe_index
      const types::global_dof_index starting_offset = dof_offsets[obj_index];
      const types::global_dof_index *pointer        = &dofs[starting_offset];
      while (true)
        {
          Assert (*pointer != numbers::invalid_dof_index,
                  ExcInternalError());
          if (*pointer == fe_index)
            return *(pointer + 1 + local_index);
          else
            pointer += static_cast<types::global_dof_index>(
                         dof_handler.get_fe()[*pointer]
                         .template n_dofs_per_object<structdim>() + 1);
        }
    }



    template <int structdim>
    template <int dim, int spacedim>
    inline
    void
    DoFIndicesOnFacesOrEdges<structdim>::
    set_dof_index (const dealii::hp::DoFHandler<dim,spacedim> &dof_handler,
                   const unsigned int                obj_index,
                   const unsigned int                fe_index,
                   const unsigned int                local_index,
                   const types::global_dof_index     global_index,
                   const unsigned int                obj_level)
    {
      Assert ((fe_index != dealii::hp::DoFHandler<dim,spacedim>::default_fe_index),
              ExcMessage ("You need to specify a FE index when working "
                          "with hp DoFHandlers"));
      Assert (&dof_handler != 0,
              ExcMessage ("No DoFHandler is specified for this iterator"));
      Assert (&dof_handler.get_fe() != 0,
              ExcMessage ("No finite element collection is associated with "
                          "this DoFHandler"));
      Assert (fe_index < dof_handler.get_fe().size(),
              ExcIndexRange (fe_index, 0, dof_handler.get_fe().size()));
      Assert (local_index <
              dof_handler.get_fe()[fe_index].template n_dofs_per_object<structdim>(),
              ExcIndexRange(local_index, 0,
                            dof_handler.get_fe()[fe_index]
                            .template n_dofs_per_object<structdim>()));
      Assert (obj_index < dof_offsets.size(),
              ExcIndexRange (obj_index, 0, dof_offsets.size()));

      // make sure we are on an
      // object for which DoFs have
      // been allocated at all
      Assert (dof_offsets[obj_index] != numbers::invalid_unsigned_int,
              ExcMessage ("You are trying to access degree of freedom "
                          "information for an object on which no such "
                          "information is available"));

      Assert (structdim<dim, ExcMessage ("This object can not be used for cells."));

      // there may be multiple finite elements associated with
      // this object.  hop along the list of index sets until we
      // find the one with the correct fe_index, and then poke
      // into that part. trigger an exception if we can't find a
      // set for this particular fe_index
      const types::global_dof_index starting_offset = dof_offsets[obj_index];
      types::global_dof_index      *pointer         = &dofs[starting_offset];
      while (true)
        {
          Assert (*pointer != numbers::invalid_dof_index,
                  ExcInternalError());
          if (*pointer == fe_index)
            {
              *(pointer + 1 + local_index) = global_index;
              return;
            }
          else
            pointer += dof_handler.get_fe()[*pointer]
                       .template n_dofs_per_object<structdim>() + 1;
        }
    }



    template <int structdim>
    template <int dim, int spacedim>
    inline
    unsigned int
    DoFIndicesOnFacesOrEdges<structdim>::
    n_active_fe_indices (const dealii::hp::DoFHandler<dim,spacedim> &dof_handler,
                         const unsigned int                obj_index) const
    {
      Assert (&dof_handler != 0,
              ExcMessage ("No DoFHandler is specified for this iterator"));
      Assert (&dof_handler.get_fe() != 0,
              ExcMessage ("No finite element collection is associated with "
                          "this DoFHandler"));
      Assert (obj_index < dof_offsets.size(),
              ExcIndexRange (obj_index, 0, dof_offsets.size()));

      // make sure we are on an
      // object for which DoFs have
      // been allocated at all
      if (dof_offsets[obj_index] == numbers::invalid_unsigned_int)
        return 0;

      Assert (structdim<dim, ExcMessage ("This object can not be used for cells."));

      // there may be multiple finite elements associated with this
      // object. hop along the list of index sets until we find the
      // one with the correct fe_index, and then poke into that
      // part. trigger an exception if we can't find a set for this
      // particular fe_index
      const unsigned int starting_offset = dof_offsets[obj_index];
      const types::global_dof_index *pointer        = &dofs[starting_offset];
      unsigned int counter = 0;
      while (true)
        {
          if (*pointer == numbers::invalid_dof_index)
            // end of list reached
            return counter;
          else
            {
              ++counter;
              pointer += dof_handler.get_fe()[*pointer]
                         .template n_dofs_per_object<structdim>() + 1;
            }
        }
    }



    template <int structdim>
    template <int dim, int spacedim>
    inline
    types::global_dof_index
    DoFIndicesOnFacesOrEdges<structdim>::
    nth_active_fe_index (const dealii::hp::DoFHandler<dim,spacedim> &dof_handler,
                         const unsigned int                obj_level,
                         const unsigned int                obj_index,
                         const unsigned int                n) const
    {
      Assert (&dof_handler != 0,
              ExcMessage ("No DoFHandler is specified for this iterator"));
      Assert (&dof_handler.get_fe() != 0,
              ExcMessage ("No finite element collection is associated with "
                          "this DoFHandler"));
      Assert (obj_index < dof_offsets.size(),
              ExcIndexRange (obj_index, 0, dof_offsets.size()));

      // make sure we are on an
      // object for which DoFs have
      // been allocated at all
      Assert (dof_offsets[obj_index] != numbers::invalid_unsigned_int,
              ExcMessage ("You are trying to access degree of freedom "
                          "information for an object on which no such "
                          "information is available"));

      Assert (structdim<dim, ExcMessage ("This object can not be used for cells."));

      Assert (n < n_active_fe_indices(dof_handler, obj_index),
              ExcIndexRange (n, 0,
                             n_active_fe_indices(dof_handler, obj_index)));

      // there may be multiple finite elements associated with
      // this object. hop along the list of index sets until we
      // find the one with the correct fe_index, and then poke
      // into that part. trigger an exception if we can't find a
      // set for this particular fe_index
      const unsigned int starting_offset = dof_offsets[obj_index];
      const types::global_dof_index *pointer = &dofs[starting_offset];
      unsigned int counter = 0;
      while (true)
        {
          Assert (*pointer != numbers::invalid_dof_index,
                  ExcInternalError());

          const unsigned int fe_index = *pointer;

          Assert (fe_index < dof_handler.get_fe().size(),
                  ExcInternalError());

          if (counter == n)
            return fe_index;

          ++counter;
          pointer += dof_handler.get_fe()[fe_index]
                     .template n_dofs_per_object<structdim>() + 1;
        }
    }



    template <int structdim>
    template <int dim, int spacedim>
    inline
    bool
    DoFIndicesOnFacesOrEdges<structdim>::
    fe_index_is_active (const dealii::hp::DoFHandler<dim,spacedim> &dof_handler,
                        const unsigned int                obj_index,
                        const unsigned int                fe_index,
                        const unsigned int                obj_level) const
    {
      Assert (&dof_handler != 0,
              ExcMessage ("No DoFHandler is specified for this iterator"));
      Assert (&dof_handler.get_fe() != 0,
              ExcMessage ("No finite element collection is associated with "
                          "this DoFHandler"));
      Assert (obj_index < dof_offsets.size(),
              ExcIndexRange (obj_index, 0, static_cast<unsigned int>(dof_offsets.size())));
      Assert ((fe_index != dealii::hp::DoFHandler<dim,spacedim>::default_fe_index),
              ExcMessage ("You need to specify a FE index when working "
                          "with hp DoFHandlers"));
      Assert (fe_index < dof_handler.get_fe().size(),
              ExcIndexRange (fe_index, 0, dof_handler.get_fe().size()));

      // make sure we are on an
      // object for which DoFs have
      // been allocated at all
      Assert (dof_offsets[obj_index] != numbers::invalid_unsigned_int,
              ExcMessage ("You are trying to access degree of freedom "
                          "information for an object on which no such "
                          "information is available"));

      Assert (structdim<dim, ExcMessage ("This object can not be used for cells."));

      // there may be multiple finite elements associated with
      // this object. hop along the list of index sets until we
      // find the one with the correct fe_index, and then poke
      // into that part. trigger an exception if we can't find a
      // set for this particular fe_index
      const types::global_dof_index starting_offset = dof_offsets[obj_index];
      const types::global_dof_index *pointer = &dofs[starting_offset];
      while (true)
        {
          if (*pointer == numbers::invalid_dof_index)
            // end of list reached
            return false;
          else if (*pointer == fe_index)
            return true;
          else
            pointer += static_cast<types::global_dof_index>(
                         dof_handler.get_fe()[*pointer]
                         .template n_dofs_per_object<structdim>()+1);
        }
    }


  } // namespace hp

} // namespace internal

DEAL_II_NAMESPACE_CLOSE

#endif