This file is indexed.

/usr/include/deal.II/multigrid/mg_tools.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
// ---------------------------------------------------------------------
// $Id: mg_tools.h 30036 2013-07-18 16:55:32Z maier $
//
// Copyright (C) 2005 - 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__mg_tools_h
#define __deal2__mg_tools_h

#include <deal.II/base/config.h>
#include <deal.II/base/index_set.h>
#include <deal.II/dofs/dof_tools.h>
#include <deal.II/dofs/dof_handler.h>

#include <vector>
#include <set>


DEAL_II_NAMESPACE_OPEN

template <class Object> class MGLevelObject;
template <int dim, int spacedim> class DoFHandler;
template <typename number> class Vector;
template <typename number> class SparseMatrix;
template <typename number> class BlockVector;
template <typename number> class BlockSparseMatrix;
template <typename number> class FullMatrix;
template <typename number> class BlockSparseMatrix;

/* !@addtogroup mg */
/* @{ */

/**
 * This is a collection of functions operating on, and manipulating
 * the numbers of degrees of freedom in a multilevel triangulation. It
 * is similar in purpose and function to the @p DoFTools class, but
 * operates on levels of DoFHandler
 * objects. See there and the documentation of the member functions
 * for more information.
 *
 * @author Wolfgang Bangerth, Guido Kanschat, 1999 - 2005, 2012
 */
namespace MGTools
{
  /**
   * Compute row length vector for
   * multilevel methods.
   */
  template <int dim, int spacedim>
  void
  compute_row_length_vector(const DoFHandler<dim,spacedim> &dofs,
                            const unsigned int level,
                            std::vector<unsigned int> &row_lengths,
                            const DoFTools::Coupling flux_couplings = DoFTools::none);

  /**
   * Compute row length vector for
   * multilevel methods with
   * optimization for block
   * couplings.
   */
  template <int dim, int spacedim>
  void
  compute_row_length_vector(const DoFHandler<dim,spacedim> &dofs,
                            const unsigned int level,
                            std::vector<unsigned int> &row_lengths,
                            const Table<2,DoFTools::Coupling> &couplings,
                            const Table<2,DoFTools::Coupling> &flux_couplings);

  /**
   * Write the sparsity structure
   * of the matrix belonging to the
   * specified @p level. The sparsity pattern
   * is not compressed, so before
   * creating the actual matrix
   * you have to compress the
   * matrix yourself, using
   * <tt>SparseMatrixStruct::compress()</tt>.
   *
   * There is no need to consider
   * hanging nodes here, since only
   * one level is considered.
   */
  template <class DH, class SparsityPattern>
  void
  make_sparsity_pattern (const DH &dof_handler,
                         SparsityPattern         &sparsity,
                         const unsigned int       level);

  /**
   * Make a sparsity pattern including fluxes
   * of discontinuous Galerkin methods.
   * @ref make_sparsity_pattern
   * @ref DoFTools
   */
  template <int dim, class SparsityPattern, int spacedim>
  void
  make_flux_sparsity_pattern (const DoFHandler<dim,spacedim> &dof_handler,
                              SparsityPattern         &sparsity,
                              const unsigned int       level);

  /**
   * Create sparsity pattern for
   * the fluxes at refinement
   * edges. The matrix maps a
   * function of the fine level
   * space @p level to the coarser
   * space.
   *
   * make_flux_sparsity_pattern()
   */
  template <int dim, class SparsityPattern, int spacedim>
  void
  make_flux_sparsity_pattern_edge (const DoFHandler<dim,spacedim> &dof_handler,
                                   SparsityPattern         &sparsity,
                                   const unsigned int       level);
  /**
   * This function does the same as
   * the other with the same name,
   * but it gets two additional
   * coefficient matrices. A matrix
   * entry will only be generated
   * for two basis functions, if
   * there is a non-zero entry
   * linking their associated
   * components in the coefficient
   * matrix.
   *
   * There is one matrix for
   * couplings in a cell and one
   * for the couplings occurring in
   * fluxes.
   */
  template <int dim, class SparsityPattern, int spacedim>
  void
  make_flux_sparsity_pattern (const DoFHandler<dim,spacedim> &dof,
                              SparsityPattern       &sparsity,
                              const unsigned int level,
                              const Table<2,DoFTools::Coupling> &int_mask,
                              const Table<2,DoFTools::Coupling> &flux_mask);

  /**
   * Create sparsity pattern for
   * the fluxes at refinement
   * edges. The matrix maps a
   * function of the fine level
   * space @p level to the coarser
   * space. This is the version
   * restricting the pattern to the
   * elements actually needed.
   *
   * make_flux_sparsity_pattern()
   */
  template <int dim, class SparsityPattern, int spacedim>
  void
  make_flux_sparsity_pattern_edge (const DoFHandler<dim,spacedim> &dof_handler,
                                   SparsityPattern         &sparsity,
                                   const unsigned int       level,
                                   const Table<2,DoFTools::Coupling> &flux_mask);

  /**
   * Count the dofs block-wise
   * on each level.
   *
   * Result is a vector containing
   * for each level a vector
   * containing the number of dofs
   * for each block (access is
   * <tt>result[level][block]</tt>).
   */
  template <class DH>
  void
  count_dofs_per_block (const DH     &dof_handler,
                        std::vector<std::vector<types::global_dof_index> > &dofs_per_block,
                        std::vector<unsigned int>  target_block = std::vector<unsigned int>());

  /**
   * Count the dofs component-wise
   * on each level.
   *
   * Result is a vector containing
   * for each level a vector
   * containing the number of dofs
   * for each component (access is
   * <tt>result[level][component]</tt>).
   */
  template <int dim, int spacedim>
  void
  count_dofs_per_component (const DoFHandler<dim,spacedim> &mg_dof,
                            std::vector<std::vector<types::global_dof_index> > &result,
                            const bool only_once = false,
                            std::vector<unsigned int> target_component = std::vector<unsigned int>());

  /**
   * @deprecated Wrapper for the
   * other function with same name,
   * introduced for compatibility.
   */
  template <int dim, int spacedim>
  void
  count_dofs_per_component (const DoFHandler<dim,spacedim> &mg_dof,
                            std::vector<std::vector<types::global_dof_index> > &result,
                            std::vector<unsigned int> target_component) DEAL_II_DEPRECATED;

  /**
   * Generate a list of those
   * degrees of freedom at the
   * boundary which should be
   * eliminated from the matrix.
   *
   * This is the multilevel
   * equivalent of
   * VectorTools::interpolate_boundary_values,
   * but since the multilevel
   * method does not have its own
   * right hand side, the function
   * values are ignored.
   *
   * @arg <tt>boundary_indices</tt>
   * is a vector which on return
   * contains all indices of
   * boundary constraint degrees of
   * freedom for each level. Its
   * length has to match the number
   * of levels.
   */
  template <int dim, int spacedim>
  void
  make_boundary_list (const DoFHandler<dim,spacedim>      &mg_dof,
                      const typename FunctionMap<dim>::type &function_map,
                      std::vector<std::set<types::global_dof_index> > &boundary_indices,
                      const ComponentMask                   &component_mask = ComponentMask());

  /**
   * The same function as above, but return
   * an IndexSet rather than a
   * std::set<unsigned int> on each level.
   */
  template <int dim, int spacedim>
  void
  make_boundary_list (const DoFHandler<dim,spacedim>      &mg_dof,
                      const typename FunctionMap<dim>::type &function_map,
                      std::vector<IndexSet>                 &boundary_indices,
                      const ComponentMask               &component_mask = ComponentMask());

  /**
   * @deprecated
   */
  template <typename number>
  void
  apply_boundary_values (const std::set<types::global_dof_index> &boundary_dofs,
                         SparseMatrix<number> &matrix,
                         const bool preserve_symmetry,
                         const bool ignore_zeros = false) DEAL_II_DEPRECATED;

  /**
   * @deprecated
   */
  template <typename number>
  void
  apply_boundary_values (const std::set<types::global_dof_index> &boundary_dofs,
                         BlockSparseMatrix<number> &matrix,
                         const bool preserve_symmetry) DEAL_II_DEPRECATED;

  /**
   * For each level in a multigrid
   * hierarchy, produce a boolean
   * mask that indicates which of
   * the degrees of freedom are
   * along interfaces of this level
   * to cells that only exist on
   * coarser levels. The function
   * returns the subset of these
   * indices in the last argument
   * that are not only on interior
   * interfaces (i.e. between cells
   * of a given level and adjacent
   * coarser levels) but also on
   * the external boundary of the
   * domain.
   */
  template <int dim, int spacedim>
  void
  extract_inner_interface_dofs (const DoFHandler<dim,spacedim> &mg_dof_handler,
                                std::vector<std::vector<bool> >  &interface_dofs,
                                std::vector<std::vector<bool> >  &boundary_interface_dofs);

  /**
   * Does the same as the function above,
   * but fills only the interface_dofs.
   */
  template <int dim, int spacedim>
  void
  extract_inner_interface_dofs (const DoFHandler<dim,spacedim> &mg_dof_handler,
                                std::vector<std::vector<bool> >  &interface_dofs);

  template <int dim, int spacedim>
  void
  extract_non_interface_dofs (const DoFHandler<dim,spacedim> &mg_dof_handler,
                              std::vector<std::set<types::global_dof_index> > &non_interface_dofs);
}

/* @} */

DEAL_II_NAMESPACE_CLOSE

#endif