This file is indexed.

/usr/include/deal.II/fe/fe_poly_face.templates.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
// ---------------------------------------------------------------------
//
// Copyright (C) 2009 - 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.
//
// ---------------------------------------------------------------------


#include <deal.II/base/qprojector.h>
#include <deal.II/base/polynomial_space.h>
#include <deal.II/fe/fe_values.h>
#include <deal.II/fe/fe_poly_face.h>


DEAL_II_NAMESPACE_OPEN

template <class PolynomialType, int dim, int spacedim>
FE_PolyFace<PolynomialType,dim,spacedim>::FE_PolyFace (
  const PolynomialType &poly_space,
  const FiniteElementData<dim> &fe_data,
  const std::vector<bool> &restriction_is_additive_flags):
  FiniteElement<dim,spacedim> (fe_data,
                               restriction_is_additive_flags,
                               std::vector<ComponentMask> (1, ComponentMask(1,true))),
  poly_space(poly_space)
{
  AssertDimension(dim, PolynomialType::dimension+1);
}


template <class PolynomialType, int dim, int spacedim>
unsigned int
FE_PolyFace<PolynomialType,dim,spacedim>::get_degree () const
{
  return this->degree;
}


//---------------------------------------------------------------------------
// Auxiliary functions
//---------------------------------------------------------------------------


template <class PolynomialType, int dim, int spacedim>
UpdateFlags
FE_PolyFace<PolynomialType,dim,spacedim>::requires_update_flags (const UpdateFlags flags) const
{
  UpdateFlags out = flags & update_values;
  if (flags & update_gradients)
    out |= update_gradients | update_covariant_transformation;
  if (flags & update_hessians)
    out |= update_hessians | update_covariant_transformation;
  if (flags & update_cell_normal_vectors)
    out |= update_cell_normal_vectors | update_JxW_values;

  return out;
}


//---------------------------------------------------------------------------
// Fill data of FEValues
//---------------------------------------------------------------------------
template <class PolynomialType, int dim, int spacedim>
void
FE_PolyFace<PolynomialType,dim,spacedim>::
fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &,
                const CellSimilarity::Similarity                                     ,
                const Quadrature<dim> &,
                const Mapping<dim,spacedim> &,
                const typename Mapping<dim,spacedim>::InternalDataBase &,
                const dealii::internal::FEValues::MappingRelatedData<dim, spacedim> &,
                const typename FiniteElement<dim,spacedim>::InternalDataBase &,
                dealii::internal::FEValues::FiniteElementRelatedData<dim, spacedim> &) const
{
  // Do nothing, since we do not have
  // values in the interior
}



template <class PolynomialType, int dim, int spacedim>
void
FE_PolyFace<PolynomialType,dim,spacedim>::
fill_fe_face_values (const typename Triangulation<dim,spacedim>::cell_iterator &,
                     const unsigned int                                                   face_no,
                     const Quadrature<dim-1>                                             &quadrature,
                     const Mapping<dim,spacedim> &,
                     const typename Mapping<dim,spacedim>::InternalDataBase &,
                     const dealii::internal::FEValues::MappingRelatedData<dim, spacedim> &,
                     const typename FiniteElement<dim,spacedim>::InternalDataBase        &fe_internal,
                     dealii::internal::FEValues::FiniteElementRelatedData<dim, spacedim> &output_data) const
{
  // convert data object to internal
  // data for this class. fails with
  // an exception if that is not
  // possible
  Assert (dynamic_cast<const InternalData *> (&fe_internal) != 0, ExcInternalError());
  const InternalData &fe_data = static_cast<const InternalData &> (fe_internal);

  if (fe_data.update_each & update_values)
    for (unsigned int i=0; i<quadrature.size(); ++i)
      {
        for (unsigned int k=0; k<this->dofs_per_cell; ++k)
          output_data.shape_values(k,i) = 0.;
        switch (dim)
          {
          case 3:
          {
            // Fill data for quad shape functions
            if (this->dofs_per_quad !=0)
              {
                const unsigned int foffset = this->first_quad_index + this->dofs_per_quad * face_no;
                for (unsigned int k=0; k<this->dofs_per_quad; ++k)
                  output_data.shape_values(foffset+k,i) = fe_data.shape_values[k+this->first_face_quad_index][i];
              }

            // fall through...
          }

          case 2:
          {
            // Fill data for line shape functions
            if (this->dofs_per_line != 0)
              {
                const unsigned int foffset = this->first_line_index;
                for (unsigned int line=0; line<GeometryInfo<dim>::lines_per_face; ++line)
                  {
                    for (unsigned int k=0; k<this->dofs_per_line; ++k)
                      output_data.shape_values(foffset+GeometryInfo<dim>::face_to_cell_lines(face_no, line)*this->dofs_per_line+k,i)
                        = fe_data.shape_values[k+(line*this->dofs_per_line)+this->first_face_line_index][i];
                  }
              }

            // fall through...
          }

          case 1:
          {
            // Fill data for vertex shape functions
            if (this->dofs_per_vertex != 0)
              for (unsigned int lvertex=0; lvertex<GeometryInfo<dim>::vertices_per_face; ++lvertex)
                output_data.shape_values(GeometryInfo<dim>::face_to_cell_vertices(face_no, lvertex),i)
                  = fe_data.shape_values[lvertex][i];
            break;
          }
          }
      }
}


template <class PolynomialType, int dim, int spacedim>
void
FE_PolyFace<PolynomialType,dim,spacedim>::
fill_fe_subface_values (const typename Triangulation<dim,spacedim>::cell_iterator &,
                        const unsigned int                                                   face_no,
                        const unsigned int                                                   sub_no,
                        const Quadrature<dim-1>                                             &quadrature,
                        const Mapping<dim,spacedim> &,
                        const typename Mapping<dim,spacedim>::InternalDataBase &,
                        const dealii::internal::FEValues::MappingRelatedData<dim, spacedim> &,
                        const typename FiniteElement<dim,spacedim>::InternalDataBase        &fe_internal,
                        dealii::internal::FEValues::FiniteElementRelatedData<dim, spacedim> &output_data) const
{
  // convert data object to internal
  // data for this class. fails with
  // an exception if that is not
  // possible
  Assert (dynamic_cast<const InternalData *> (&fe_internal) != 0, ExcInternalError());
  const InternalData &fe_data = static_cast<const InternalData &> (fe_internal);

  const unsigned int foffset = fe_data.shape_values.size() * face_no;
  const unsigned int offset = sub_no*quadrature.size();

  if (fe_data.update_each & update_values)
    {
      for (unsigned int k=0; k<this->dofs_per_cell; ++k)
        for (unsigned int i=0; i<quadrature.size(); ++i)
          output_data.shape_values(k,i) = 0.;
      for (unsigned int k=0; k<fe_data.shape_values.size(); ++k)
        for (unsigned int i=0; i<quadrature.size(); ++i)
          output_data.shape_values(foffset+k,i) = fe_data.shape_values[k][i+offset];
    }

  Assert (!(fe_data.update_each & update_gradients), ExcNotImplemented());
  Assert (!(fe_data.update_each & update_hessians), ExcNotImplemented());
}

DEAL_II_NAMESPACE_CLOSE