This file is indexed.

/usr/include/deal.II/lac/petsc_matrix_free.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
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
// ---------------------------------------------------------------------
//
// Copyright (C) 2012 - 2015 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__petsc_matrix_free_h
#define dealii__petsc_matrix_free_h


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

#ifdef DEAL_II_WITH_PETSC

#  include <deal.II/lac/exceptions.h>
#  include <deal.II/lac/petsc_matrix_base.h>
#  include <deal.II/lac/petsc_vector.h>

DEAL_II_NAMESPACE_OPEN



namespace PETScWrappers
{
  /**
   * Implementation of a parallel matrix class based on PETSc
   * <tt>MatShell</tt> matrix-type. This base class implements only the
   * interface to the PETSc matrix object, while all the functionality is
   * contained in the matrix-vector multiplication which must be reimplemented
   * in derived classes.
   *
   * This interface is an addition to the dealii::MatrixFree class to realize
   * user-defined matrix-classes together with PETSc solvers and
   * functionalities. See also the documentation of dealii::MatrixFree class
   * and step-37 and step-48.
   *
   * Similar to other matrix classes in namespaces PETScWrappers and
   * PETScWrappers::MPI, the MatrixFree class provides the usual matrix-vector
   * multiplication <tt>vmult(VectorBase &dst, const VectorBase &src)</tt>
   * which is pure virtual and must be reimplemented in derived classes.
   * Besides the usual interface, this class has a matrix-vector
   * multiplication <tt>vmult(Vec  &dst, const Vec  &src)</tt> taking PETSc
   * Vec objects, which will be called by <tt>matrix_free_mult(Mat A, Vec src,
   * Vec dst)</tt> registered as matrix-vector multiplication of this PETSc
   * matrix object. The default implementation of the vmult function in the
   * base class translates the given PETSc <tt>Vec*</tt> vectors into a
   * deal.II vector, calls the usual vmult function with the usual interface
   * and converts the result back to PETSc <tt>Vec*</tt>. This could be made
   * much more efficient in derived classes without allocating new memory.
   *
   * @ingroup PETScWrappers
   * @ingroup Matrix1
   * @author Wolfgang Bangerth, Martin Steigemann, 2012
   */
  class MatrixFree : public MatrixBase
  {
  public:

    /**
     * Default constructor. Create an empty matrix object.
     */
    MatrixFree ();

    /**
     * Create a matrix object of dimensions @p m times @p n with communication
     * happening over the provided @p communicator.
     *
     * For the meaning of the @p local_rows and @p local_columns parameters,
     * see the PETScWrappers::MPI::SparseMatrix class documentation.
     *
     * As other PETSc matrices, also the the matrix-free object needs to have
     * a size and to perform matrix vector multiplications efficiently in
     * parallel also @p local_rows and @p local_columns. But in contrast to
     * PETSc::SparseMatrix classes a PETSc matrix-free object does not need
     * any estimation of non_zero entries and has no option
     * <tt>is_symmetric</tt>.
     */
    MatrixFree (const MPI_Comm     &communicator,
                const unsigned int  m,
                const unsigned int  n,
                const unsigned int  local_rows,
                const unsigned int  local_columns);

    /**
     * Create a matrix object of dimensions @p m times @p n with communication
     * happening over the provided @p communicator.
     *
     * As other PETSc matrices, also the the matrix-free object needs to have
     * a size and to perform matrix vector multiplications efficiently in
     * parallel also @p local_rows and @p local_columns. But in contrast to
     * PETSc::SparseMatrix classes a PETSc matrix-free object does not need
     * any estimation of non_zero entries and has no option
     * <tt>is_symmetric</tt>.
     */
    MatrixFree (const MPI_Comm     &communicator,
                const unsigned int  m,
                const unsigned int  n,
                const std::vector<unsigned int> &local_rows_per_process,
                const std::vector<unsigned int> &local_columns_per_process,
                const unsigned int  this_process);

    /**
     * Constructor for the serial case: Same function as
     * <tt>MatrixFree()</tt>, see above, with <tt>communicator =
     * MPI_COMM_WORLD</tt>.
     */
    MatrixFree (const unsigned int  m,
                const unsigned int  n,
                const unsigned int  local_rows,
                const unsigned int  local_columns);

    /**
     * Constructor for the serial case: Same function as
     * <tt>MatrixFree()</tt>, see above, with <tt>communicator =
     * MPI_COMM_WORLD</tt>.
     */
    MatrixFree (const unsigned int  m,
                const unsigned int  n,
                const std::vector<unsigned int> &local_rows_per_process,
                const std::vector<unsigned int> &local_columns_per_process,
                const unsigned int  this_process);

    /**
     * Throw away the present matrix and generate one that has the same
     * properties as if it were created by the constructor of this class with
     * the same argument list as the present function.
     */
    void reinit (const MPI_Comm     &communicator,
                 const unsigned int  m,
                 const unsigned int  n,
                 const unsigned int  local_rows,
                 const unsigned int  local_columns);

    /**
     * Throw away the present matrix and generate one that has the same
     * properties as if it were created by the constructor of this class with
     * the same argument list as the present function.
     */
    void reinit (const MPI_Comm     &communicator,
                 const unsigned int  m,
                 const unsigned int  n,
                 const std::vector<unsigned int> &local_rows_per_process,
                 const std::vector<unsigned int> &local_columns_per_process,
                 const unsigned int  this_process);

    /**
     * Calls the @p reinit() function above with <tt>communicator =
     * MPI_COMM_WORLD</tt>.
     */
    void reinit (const unsigned int  m,
                 const unsigned int  n,
                 const unsigned int  local_rows,
                 const unsigned int  local_columns);

    /**
     * Calls the @p reinit() function above with <tt>communicator =
     * MPI_COMM_WORLD</tt>.
     */
    void reinit (const unsigned int  m,
                 const unsigned int  n,
                 const std::vector<unsigned int> &local_rows_per_process,
                 const std::vector<unsigned int> &local_columns_per_process,
                 const unsigned int  this_process);

    /**
     * Release all memory and return to a state just like after having called
     * the default constructor.
     */
    void clear ();

    /**
     * Return a reference to the MPI communicator object in use with this
     * matrix.
     */
    const MPI_Comm &get_mpi_communicator () const;

    /**
     * Matrix-vector multiplication: let <i>dst = M*src</i> with <i>M</i>
     * being this matrix.
     *
     * Source and destination must not be the same vector.
     *
     * Note that if the current object represents a parallel distributed
     * matrix (of type PETScWrappers::MPI::SparseMatrix), then both vectors
     * have to be distributed vectors as well. Conversely, if the matrix is
     * not distributed, then neither of the vectors may be.
     */
    virtual
    void vmult (VectorBase       &dst,
                const VectorBase &src) const = 0;

    /**
     * Matrix-vector multiplication: let <i>dst = M<sup>T</sup>*src</i> with
     * <i>M</i> being this matrix. This function does the same as @p vmult()
     * but takes the transposed matrix.
     *
     * Source and destination must not be the same vector.
     *
     * Note that if the current object represents a parallel distributed
     * matrix then both vectors have to be distributed vectors as well.
     * Conversely, if the matrix is not distributed, then neither of the
     * vectors may be.
     */
    virtual
    void Tvmult (VectorBase       &dst,
                 const VectorBase &src) const = 0;

    /**
     * Adding Matrix-vector multiplication. Add <i>M*src</i> on <i>dst</i>
     * with <i>M</i> being this matrix.
     *
     * Source and destination must not be the same vector.
     *
     * Note that if the current object represents a parallel distributed
     * matrix then both vectors have to be distributed vectors as well.
     * Conversely, if the matrix is not distributed, then neither of the
     * vectors may be.
     */
    virtual
    void vmult_add (VectorBase       &dst,
                    const VectorBase &src) const = 0;

    /**
     * Adding Matrix-vector multiplication. Add <i>M<sup>T</sup>*src</i> to
     * <i>dst</i> with <i>M</i> being this matrix. This function does the same
     * as @p vmult_add() but takes the transposed matrix.
     *
     * Source and destination must not be the same vector.
     *
     * Note that if the current object represents a parallel distributed
     * matrix then both vectors have to be distributed vectors as well.
     * Conversely, if the matrix is not distributed, then neither of the
     * vectors may be.
     */
    virtual
    void Tvmult_add (VectorBase       &dst,
                     const VectorBase &src) const = 0;

    /**
     * The matrix-vector multiplication called by @p matrix_free_mult(). This
     * function can be reimplemented in derived classes for efficiency. The
     * default implementation copies the given vectors into
     * PETScWrappers::*::Vector and calls <tt>vmult(VectorBase &dst, const
     * VectorBase &src)</tt> which is purely virtual and must be reimplemented
     * in derived classes.
     */
    virtual
    void vmult (Vec  &dst, const Vec  &src) const;

  private:

    /**
     * Copy of the communicator object to be used for this parallel matrix-
     * free object.
     */
    MPI_Comm  communicator;

    /**
     * Callback-function registered as the matrix-vector multiplication of
     * this matrix-free object called by PETSc routines. This function must be
     * static and takes a PETSc matrix @p A, and vectors @p src and @p dst,
     * where <i>dst = A*src</i>
     *
     * Source and destination must not be the same vector.
     *
     * This function calls <tt>vmult(Vec &dst, const Vec &src)</tt> which
     * should be reimplemented in derived classes.
     */
    static int matrix_free_mult (Mat  A, Vec  src, Vec  dst);

    /**
     * Do the actual work for the respective @p reinit() function and the
     * matching constructor, i.e. create a matrix object. Getting rid of the
     * previous matrix is left to the caller.
     */
    void do_reinit (const unsigned int  m,
                    const unsigned int  n,
                    const unsigned int  local_rows,
                    const unsigned int  local_columns);
  };



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

  inline
  const MPI_Comm &
  MatrixFree::get_mpi_communicator () const
  {
    return communicator;
  }
}

DEAL_II_NAMESPACE_CLOSE

#endif // DEAL_II_WITH_PETSC


/*----------------------------   petsc_matrix_free.h     ---------------------------*/

#endif
/*----------------------------   petsc_matrix_free.h     ---------------------------*/