This file is indexed.

/usr/include/deal.II/fe/fe_poly_face.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
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
// ---------------------------------------------------------------------
//
// 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.
//
// ---------------------------------------------------------------------

#ifndef dealii__fe_poly_face_h
#define dealii__fe_poly_face_h


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


DEAL_II_NAMESPACE_OPEN

/*!@addtogroup febase */
/*@{*/

/**
 * @warning This class is not sufficiently tested yet!
 *
 * This class gives a unified framework for the implementation of
 * FiniteElement classes only located on faces of the mesh. They are based on
 * polynomial spaces like the TensorProductPolynomials or a PolynomialSpace
 * classes.
 *
 * Every class that implements the following functions can be used as template
 * parameter PolynomialType.
 *
 * @code
 * double compute_value (const unsigned int i,
 *                       const Point<dim> &p) const;
 * @endcode
 * Example classes are TensorProductPolynomials, PolynomialSpace or
 * PolynomialsP.
 *
 * This class is not a fully implemented FiniteElement class. Instead there
 * are several pure virtual functions declared in the FiniteElement class
 * which cannot be implemented by this class but are left for implementation
 * in derived classes.
 *
 * @author Guido Kanschat, 2009
 */
template <class PolynomialType, int dim=PolynomialType::dimension+1, int spacedim=dim>
class FE_PolyFace : public FiniteElement<dim,spacedim>
{
public:
  /**
   * Constructor.
   */
  FE_PolyFace (const PolynomialType &poly_space,
               const FiniteElementData<dim> &fe_data,
               const std::vector<bool> &restriction_is_additive_flags);

  /**
   * Return the polynomial degree of this finite element, i.e. the value
   * passed to the constructor.
   */
  unsigned int get_degree () const;

  // for documentation, see the FiniteElement base class
  virtual
  UpdateFlags
  requires_update_flags (const UpdateFlags update_flags) const;

protected:
  /*
   * NOTE: The following functions have their definitions inlined into the class declaration
   * because we otherwise run into a compiler error with MS Visual Studio.
   */


  virtual
  typename FiniteElement<dim,spacedim>::InternalDataBase *
  get_data (const UpdateFlags                                                    /*update_flags*/,
            const Mapping<dim,spacedim>                                         &/*mapping*/,
            const Quadrature<dim>                                               &/*quadrature*/,
            dealii::internal::FEValues::FiniteElementRelatedData<dim, spacedim> &/*output_data*/) const
  {
    InternalData *data = new InternalData;
    return data;
  }

  typename FiniteElement<dim,spacedim>::InternalDataBase *
  get_face_data(const UpdateFlags                                                    update_flags,
                const Mapping<dim,spacedim>                                         &/*mapping*/,
                const Quadrature<dim-1>                                             &quadrature,
                dealii::internal::FEValues::FiniteElementRelatedData<dim, spacedim> &/*output_data*/) const
  {
    // generate a new data object and
    // initialize some fields
    InternalData *data = new InternalData;
    data->update_each = requires_update_flags(update_flags);

    const unsigned int n_q_points = quadrature.size();

    // some scratch arrays
    std::vector<double> values(0);
    std::vector<Tensor<1,dim-1> > grads(0);
    std::vector<Tensor<2,dim-1> > grad_grads(0);
    std::vector<Tensor<3,dim-1> > empty_vector_of_3rd_order_tensors;
    std::vector<Tensor<4,dim-1> > empty_vector_of_4th_order_tensors;

    // initialize fields only if really
    // necessary. otherwise, don't
    // allocate memory
    if (data->update_each & update_values)
      {
        values.resize (poly_space.n());
        data->shape_values.resize (poly_space.n(),
                                   std::vector<double> (n_q_points));
        for (unsigned int i=0; i<n_q_points; ++i)
          {
            poly_space.compute(quadrature.point(i),
                               values, grads, grad_grads,
                               empty_vector_of_3rd_order_tensors,
                               empty_vector_of_4th_order_tensors);

            for (unsigned int k=0; k<poly_space.n(); ++k)
              data->shape_values[k][i] = values[k];
          }
      }
    // No derivatives of this element
    // are implemented.
    if (data->update_each & update_gradients || data->update_each & update_hessians)
      {
        Assert(false, ExcNotImplemented());
      }

    return data;
  }

  typename FiniteElement<dim,spacedim>::InternalDataBase *
  get_subface_data(const UpdateFlags                                                    update_flags,
                   const Mapping<dim,spacedim>                                         &mapping,
                   const Quadrature<dim-1>                                             &quadrature,
                   dealii::internal::FEValues::FiniteElementRelatedData<dim, spacedim> &output_data) const
  {
    return get_face_data(update_flags, mapping,
                         QProjector<dim - 1>::project_to_all_children(quadrature),
                         output_data);
  }

  virtual
  void
  fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator           &cell,
                  const CellSimilarity::Similarity                                     cell_similarity,
                  const Quadrature<dim>                                               &quadrature,
                  const Mapping<dim,spacedim>                                         &mapping,
                  const typename Mapping<dim,spacedim>::InternalDataBase              &mapping_internal,
                  const dealii::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
                  const typename FiniteElement<dim,spacedim>::InternalDataBase        &fe_internal,
                  dealii::internal::FEValues::FiniteElementRelatedData<dim, spacedim> &output_data) const;

  virtual
  void
  fill_fe_face_values (const typename Triangulation<dim,spacedim>::cell_iterator           &cell,
                       const unsigned int                                                   face_no,
                       const Quadrature<dim-1>                                             &quadrature,
                       const Mapping<dim,spacedim>                                         &mapping,
                       const typename Mapping<dim,spacedim>::InternalDataBase              &mapping_internal,
                       const dealii::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
                       const typename FiniteElement<dim,spacedim>::InternalDataBase        &fe_internal,
                       dealii::internal::FEValues::FiniteElementRelatedData<dim, spacedim> &output_data) const;

  virtual
  void
  fill_fe_subface_values (const typename Triangulation<dim,spacedim>::cell_iterator           &cell,
                          const unsigned int                                                   face_no,
                          const unsigned int                                                   sub_no,
                          const Quadrature<dim-1>                                             &quadrature,
                          const Mapping<dim,spacedim>                                         &mapping,
                          const typename Mapping<dim,spacedim>::InternalDataBase              &mapping_internal,
                          const dealii::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
                          const typename FiniteElement<dim,spacedim>::InternalDataBase        &fe_internal,
                          dealii::internal::FEValues::FiniteElementRelatedData<dim, spacedim> &output_data) const;

  /**
   * Fields of cell-independent data.
   *
   * For information about the general purpose of this class, see the
   * documentation of the base class.
   */
  class InternalData : public FiniteElement<dim,spacedim>::InternalDataBase
  {
  public:
    /**
     * Array with shape function values in quadrature points on one face.
     * There is one row for each shape function, containing values for each
     * quadrature point.
     *
     * In this array, we store the values of the shape function in the
     * quadrature points on one face of the unit cell. Since these values do
     * not change under transformation to the real cell, we only need to copy
     * them over when visiting a concrete cell.
     *
     * In particular, we can simply copy the same set of values to each of the
     * faces.
     */
    std::vector<std::vector<double> > shape_values;
  };

  /**
   * The polynomial space. Its type is given by the template parameter
   * PolynomialType.
   */
  PolynomialType poly_space;
};

/*@}*/

DEAL_II_NAMESPACE_CLOSE

#endif