This file is indexed.

/usr/include/deal.II/lac/solver_minres.h is in libdeal.ii-dev 8.5.1-3.

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
// ---------------------------------------------------------------------
//
// Copyright (C) 2000 - 2017 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__solver_minres_h
#define dealii__solver_minres_h


#include <deal.II/base/config.h>
#include <deal.II/lac/solver.h>
#include <deal.II/lac/solver_control.h>
#include <deal.II/base/logstream.h>
#include <cmath>
#include <deal.II/base/subscriptor.h>
#include <deal.II/base/signaling_nan.h>

DEAL_II_NAMESPACE_OPEN

/*!@addtogroup Solvers */
/*@{*/

/**
 * Minimal residual method for symmetric matrices.
 *
 * For the requirements on matrices and vectors in order to work with this
 * class, see the documentation of the Solver base class.
 *
 * Like all other solver classes, this class has a local structure called @p
 * AdditionalData which is used to pass additional parameters to the solver,
 * like damping parameters or the number of temporary vectors. We use this
 * additional structure instead of passing these values directly to the
 * constructor because this makes the use of the @p SolverSelector and other
 * classes much easier and guarantees that these will continue to work even if
 * number or type of the additional parameters for a certain solver changes.
 *
 * However, since the MinRes method does not need additional data, the
 * respective structure is empty and does not offer any functionality. The
 * constructor has a default argument, so you may call it without the
 * additional parameter.
 *
 * The preconditioner has to be positive definite and symmetric
 *
 * The algorithm is taken from the Master thesis of Astrid Battermann with
 * some changes. The full text can be found at
 * http://scholar.lib.vt.edu/theses/public/etd-12164379662151/etd-title.html
 *
 *
 * <h3>Observing the progress of linear solver iterations</h3>
 *
 * The solve() function of this class uses the mechanism described in the
 * Solver base class to determine convergence. This mechanism can also be used
 * to observe the progress of the iteration.
 *
 *
 * @author Thomas Richter, 2000, Luca Heltai, 2006
 */
template <class VectorType = Vector<double> >
class SolverMinRes : public Solver<VectorType>
{
public:
  /**
   * Standardized data struct to pipe additional data to the solver. This
   * solver does not need additional data yet.
   */
  struct AdditionalData
  {
  };

  /**
   * Constructor.
   */
  SolverMinRes (SolverControl            &cn,
                VectorMemory<VectorType> &mem,
                const AdditionalData     &data=AdditionalData());

  /**
   * Constructor. Use an object of type GrowingVectorMemory as a default to
   * allocate memory.
   */
  SolverMinRes (SolverControl        &cn,
                const AdditionalData &data=AdditionalData());

  /**
   * Virtual destructor.
   */
  virtual ~SolverMinRes ();

  /**
   * Solve the linear system $Ax=b$ for x.
   */
  template<typename MatrixType, typename PreconditionerType>
  void
  solve (const MatrixType         &A,
         VectorType               &x,
         const VectorType         &b,
         const PreconditionerType &precondition);

  /**
   * @addtogroup Exceptions
   * @{
   */

  /**
   * Exception
   */
  DeclException0 (ExcPreconditionerNotDefinite);
  //@}

protected:
  /**
   * Implementation of the computation of the norm of the residual.
   */
  virtual double criterion();
  /**
   * Interface for derived class. This function gets the current iteration
   * vector, the residual and the update vector in each step. It can be used
   * for a graphical output of the convergence history.
   */
  virtual void print_vectors(const unsigned int step,
                             const VectorType   &x,
                             const VectorType   &r,
                             const VectorType   &d) const;

  /**
   * Temporary vectors, allocated through the @p VectorMemory object at the
   * start of the actual solution process and deallocated at the end.
   */
  VectorType *Vu0, *Vu1, *Vu2;
  VectorType *Vm0, *Vm1, *Vm2;
  VectorType *Vv;

  /**
   * Within the iteration loop, the square of the residual vector is stored in
   * this variable. The function @p criterion uses this variable to compute
   * the convergence value, which in this class is the norm of the residual
   * vector and thus the square root of the @p res2 value.
   */
  double res2;
};

/*@}*/
/*------------------------- Implementation ----------------------------*/

#ifndef DOXYGEN

template<class VectorType>
SolverMinRes<VectorType>::SolverMinRes (SolverControl            &cn,
                                        VectorMemory<VectorType> &mem,
                                        const AdditionalData &)
  :
  Solver<VectorType>(cn,mem)
{}



template<class VectorType>
SolverMinRes<VectorType>::SolverMinRes (SolverControl        &cn,
                                        const AdditionalData &)
  :
  Solver<VectorType>(cn),
  Vu0(NULL),
  Vu1(NULL),
  Vu2(NULL),
  Vm0(NULL),
  Vm1(NULL),
  Vm2(NULL),
  Vv(NULL),
  res2(numbers::signaling_nan<double>())
{}


template<class VectorType>
SolverMinRes<VectorType>::~SolverMinRes ()
{}



template<class VectorType>
double
SolverMinRes<VectorType>::criterion()
{
  return res2;
}


template<class VectorType>
void
SolverMinRes<VectorType>::print_vectors(const unsigned int,
                                        const VectorType &,
                                        const VectorType &,
                                        const VectorType &) const
{}



template<class VectorType>
template<typename MatrixType, typename PreconditionerType>
void
SolverMinRes<VectorType>::solve (const MatrixType         &A,
                                 VectorType               &x,
                                 const VectorType         &b,
                                 const PreconditionerType &precondition)
{
  deallog.push("minres");

  // Memory allocation
  Vu0  = this->memory.alloc();
  Vu1  = this->memory.alloc();
  Vu2  = this->memory.alloc();
  Vv   = this->memory.alloc();
  Vm0  = this->memory.alloc();
  Vm1  = this->memory.alloc();
  Vm2  = this->memory.alloc();
  // define some aliases for simpler access
  typedef VectorType *vecptr;
  vecptr u[3] = {Vu0, Vu1, Vu2};
  vecptr m[3] = {Vm0, Vm1, Vm2};
  VectorType &v   = *Vv;
  // resize the vectors, but do not set
  // the values since they'd be overwritten
  // soon anyway.
  u[0]->reinit(b,true);
  u[1]->reinit(b,true);
  u[2]->reinit(b,true);
  m[0]->reinit(b,true);
  m[1]->reinit(b,true);
  m[2]->reinit(b,true);
  v.reinit(b,true);

  // some values needed
  double delta[3] = { 0, 0, 0 };
  double f[2] = { 0, 0 };
  double e[2] = { 0, 0 };

  double r_l2 = 0;
  double r0   = 0;
  double tau  = 0;
  double c    = 0;
  double s    = 0;
  double d_   = 0;

  // The iteration step.
  unsigned int j = 1;


  // Start of the solution process
  A.vmult(*m[0],x);
  *u[1] = b;
  *u[1] -= *m[0];
  // Precondition is applied.
  // The preconditioner has to be
  // positive definite and symmetric

  // M v = u[1]
  precondition.vmult (v,*u[1]);

  delta[1] = v * (*u[1]);
  // Preconditioner positive
  Assert (delta[1]>=0, ExcPreconditionerNotDefinite());

  r0 = std::sqrt(delta[1]);
  r_l2 = r0;


  u[0]->reinit(b);
  delta[0] = 1.;
  m[0]->reinit(b);
  m[1]->reinit(b);
  m[2]->reinit(b);

  SolverControl::State conv = this->iteration_status(0, r_l2, x);
  while (conv==SolverControl::iterate)
    {
      if (delta[1]!=0)
        v *= 1./std::sqrt(delta[1]);
      else
        v.reinit(b);

      A.vmult(*u[2],v);
      u[2]->add (-std::sqrt(delta[1]/delta[0]), *u[0]);

      const double gamma = *u[2] * v;
      u[2]->add (-gamma / std::sqrt(delta[1]), *u[1]);
      *m[0] = v;

      // precondition: solve M v = u[2]
      // Preconditioner has to be positive
      // definite and symmetric.
      precondition.vmult(v,*u[2]);

      delta[2] = v * (*u[2]);

      Assert (delta[2]>=0, ExcPreconditionerNotDefinite());

      if (j==1)
        {
          d_ = gamma;
          e[1] = std::sqrt(delta[2]);
        }
      if (j>1)
        {
          d_ = s * e[0] - c * gamma;
          e[0] = c * e[0] + s * gamma;
          f[1] = s * std::sqrt(delta[2]);
          e[1] = -c * std::sqrt(delta[2]);
        }

      const double d = std::sqrt (d_*d_ + delta[2]);

      if (j>1)
        tau *= s / c;
      c = d_ / d;
      tau *= c;

      s = std::sqrt(delta[2]) / d;

      if (j==1)
        tau = r0 * c;

      m[0]->add (-e[0], *m[1]);
      if (j>1)
        m[0]->add (-f[0], *m[2]);
      *m[0] *= 1./d;
      x.add (tau, *m[0]);
      r_l2 *= std::fabs(s);

      conv = this->iteration_status(j, r_l2, x);

      // next iteration step
      ++j;
      // All vectors have to be shifted
      // one iteration step.
      // This should be changed one time.
      //
      // the previous code was like this:
      //   m[2] = m[1];
      //   m[1] = m[0];
      // but it can be made more efficient,
      // since the value of m[0] is no more
      // needed in the next iteration
      swap (*m[2], *m[1]);
      swap (*m[1], *m[0]);

      // likewise, but reverse direction:
      //   u[0] = u[1];
      //   u[1] = u[2];
      swap (*u[0], *u[1]);
      swap (*u[1], *u[2]);

      // these are scalars, so need
      // to bother
      f[0] = f[1];
      e[0] = e[1];
      delta[0] = delta[1];
      delta[1] = delta[2];
    }

  // Deallocation of Memory
  this->memory.free(Vu0);
  this->memory.free(Vu1);
  this->memory.free(Vu2);
  this->memory.free(Vv);
  this->memory.free(Vm0);
  this->memory.free(Vm1);
  this->memory.free(Vm2);
  // Output
  deallog.pop ();

  // in case of failure: throw exception
  AssertThrow(conv == SolverControl::success,
              SolverControl::NoConvergence (j, r_l2));

  // otherwise exit as normal
}

#endif // DOXYGEN

DEAL_II_NAMESPACE_CLOSE

#endif