This file is indexed.

/usr/include/deal.II/lac/solver_bicgstab.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
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
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
// ---------------------------------------------------------------------
//
// Copyright (C) 1998 - 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__solver_bicgstab_h
#define dealii__solver_bicgstab_h


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

DEAL_II_NAMESPACE_OPEN

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

/**
 * Bicgstab algorithm by van der Vorst.
 *
 * 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.
 *
 * The Bicgstab-method has two additional parameters: the first is a boolean,
 * deciding whether to compute the actual residual in each step (@p true) or
 * to use the length of the computed orthogonal residual (@p false). Note that
 * computing the residual causes a third matrix-vector-multiplication, though
 * no additional preconditioning, in each step. The reason for doing this is,
 * that the size of the orthogonalized residual computed during the iteration
 * may be larger by orders of magnitude than the true residual. This is due to
 * numerical instabilities related to badly conditioned matrices. Since this
 * instability results in a bad stopping criterion, the default for this
 * parameter is @p true. Whenever the user knows that the estimated residual
 * works reasonably as well, the flag should be set to @p false in order to
 * increase the performance of the solver.
 *
 * The second parameter is the size of a breakdown criterion. It is difficult
 * to find a general good criterion, so if things do not work for you, try to
 * change this value.
 *
 *
 * <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.
 *
 */
template <typename VectorType = Vector<double> >
class SolverBicgstab : public Solver<VectorType>
{
public:
  /**
   * There are two possibilities to compute the residual: one is an estimate
   * using the computed value @p tau. The other is exact computation using
   * another matrix vector multiplication. This increases the costs of the
   * algorithm, so it is should be set to false whenever the problem allows
   * it.
   *
   * Bicgstab is susceptible to breakdowns, so we need a parameter telling us,
   * which numbers are considered zero.
   */
  struct AdditionalData
  {
    /**
     * Constructor.
     *
     * The default is to perform an exact residual computation and breakdown
     * parameter 1e-10.
     */
    explicit
    AdditionalData(const bool   exact_residual = true,
                   const double breakdown      = 1.e-10)
      : exact_residual(exact_residual),
        breakdown(breakdown)
    {}
    /**
     * Flag for exact computation of residual.
     */
    bool exact_residual;
    /**
     * Breakdown threshold.
     */
    double breakdown;
  };

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

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

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

  /**
   * Solve primal problem only.
   */
  template<typename MatrixType, typename PreconditionerType>
  void
  solve (const MatrixType         &A,
         VectorType               &x,
         const VectorType         &b,
         const PreconditionerType &precondition);

protected:
  /**
   * Computation of the stopping criterion.
   */
  template <typename MatrixType>
  double criterion (const MatrixType &A, const VectorType &x, const VectorType &b);

  /**
   * 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;

  /**
   * Auxiliary vector.
   */
  VectorType *Vx;
  /**
   * Auxiliary vector.
   */
  VectorType *Vr;
  /**
   * Auxiliary vector.
   */
  VectorType *Vrbar;
  /**
   * Auxiliary vector.
   */
  VectorType *Vp;
  /**
   * Auxiliary vector.
   */
  VectorType *Vy;
  /**
   * Auxiliary vector.
   */
  VectorType *Vz;
  /**
   * Auxiliary vector.
   */
  VectorType *Vt;
  /**
   * Auxiliary vector.
   */
  VectorType *Vv;
  /**
   * Right hand side vector.
   */
  const VectorType *Vb;

  /**
   * Auxiliary value.
   */
  double alpha;
  /**
   * Auxiliary value.
   */
  double beta;
  /**
   * Auxiliary value.
   */
  double omega;
  /**
   * Auxiliary value.
   */
  double rho;
  /**
   * Auxiliary value.
   */
  double rhobar;

  /**
   * Current iteration step.
   */
  unsigned int step;

  /**
   * Residual.
   */
  double res;

  /**
   * Additional parameters.
   */
  AdditionalData additional_data;

private:
  /**
   * Everything before the iteration loop.
   */
  template <typename MatrixType>
  SolverControl::State start(const MatrixType &A);

  /**
   * A structure returned by the iterate() function representing what it found
   * is happening during the iteration.
   */
  struct IterationResult
  {
    bool                 breakdown;
    SolverControl::State state;
    unsigned int         last_step;
    double               last_residual;

    IterationResult (const bool breakdown,
                     const SolverControl::State state,
                     const unsigned int         last_step,
                     const double               last_residual);
  };

  /**
   * The iteration loop itself. The function returns a structure indicating
   * what happened in this function.
   */
  template<typename MatrixType, typename PreconditionerType>
  IterationResult
  iterate(const MatrixType         &A,
          const PreconditionerType &precondition);
};

/*@}*/
/*-------------------------Inline functions -------------------------------*/

#ifndef DOXYGEN


template<typename VectorType>
SolverBicgstab<VectorType>::IterationResult::IterationResult
(const bool                 breakdown,
 const SolverControl::State state,
 const unsigned int         last_step,
 const double               last_residual)
  :
  breakdown (breakdown),
  state (state),
  last_step (last_step),
  last_residual (last_residual)
{}


template<typename VectorType>
SolverBicgstab<VectorType>::SolverBicgstab (SolverControl            &cn,
                                            VectorMemory<VectorType> &mem,
                                            const AdditionalData     &data)
  :
  Solver<VectorType>(cn,mem),
  additional_data(data)
{}



template<typename VectorType>
SolverBicgstab<VectorType>::SolverBicgstab (SolverControl        &cn,
                                            const AdditionalData &data)
  :
  Solver<VectorType>(cn),
  additional_data(data)
{}



template<typename VectorType>
SolverBicgstab<VectorType>::~SolverBicgstab ()
{}



template <typename VectorType>
template <typename MatrixType>
double
SolverBicgstab<VectorType>::criterion (const MatrixType &A,
                                       const VectorType &x,
                                       const VectorType &b)
{
  A.vmult(*Vt, x);
  Vt->add(-1.,b);
  res = Vt->l2_norm();

  return res;
}



template <typename VectorType >
template <typename MatrixType>
SolverControl::State
SolverBicgstab<VectorType>::start (const MatrixType &A)
{
  A.vmult(*Vr, *Vx);
  Vr->sadd(-1.,1.,*Vb);
  res = Vr->l2_norm();

  return this->iteration_status(step, res, *Vx);
}



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



template<typename VectorType>
template<typename MatrixType, typename PreconditionerType>
typename SolverBicgstab<VectorType>::IterationResult
SolverBicgstab<VectorType>::iterate(const MatrixType         &A,
                                    const PreconditionerType &precondition)
{
//TODO:[GK] Implement "use the length of the computed orthogonal residual" in the BiCGStab method.
  SolverControl::State state = SolverControl::iterate;
  alpha = omega = rho = 1.;

  VectorType &r = *Vr;
  VectorType &rbar = *Vrbar;
  VectorType &p = *Vp;
  VectorType &y = *Vy;
  VectorType &z = *Vz;
  VectorType &t = *Vt;
  VectorType &v = *Vv;

  rbar = r;
  bool startup = true;

  do
    {
      ++step;

      rhobar = r*rbar;
      beta   = rhobar * alpha / (rho * omega);
      rho    = rhobar;
      if (startup == true)
        {
          p = r;
          startup = false;
        }
      else
        {
          p.sadd(beta, 1., r);
          p.add(-beta*omega, v);
        }

      precondition.vmult(y,p);
      A.vmult(v,y);
      rhobar = rbar * v;

      alpha = rho/rhobar;

//TODO:[?] Find better breakdown criterion

      if (std::fabs(alpha) > 1.e10)
        return IterationResult(true, state, step, res);

      res = std::sqrt(r.add_and_dot(-alpha, v, r));

      // check for early success, see the lac/bicgstab_early testcase as to
      // why this is necessary
      //
      // note: the vector *Vx we pass to the iteration_status signal here is only
      // the current approximation, not the one we will return with,
      // which will be x=*Vx + alpha*y
      if (this->iteration_status(step, res, *Vx) == SolverControl::success)
        {
          Vx->add(alpha, y);
          print_vectors(step, *Vx, r, y);
          return IterationResult(false, SolverControl::success, step, res);
        }

      precondition.vmult(z,r);
      A.vmult(t,z);
      rhobar = t*r;
      omega = rhobar/(t*t);
      Vx->add(alpha, y, omega, z);

      if (additional_data.exact_residual)
        {
          r.add(-omega, t);
          res = criterion(A, *Vx, *Vb);
        }
      else
        res = std::sqrt(r.add_and_dot(-omega, t, r));

      state = this->iteration_status(step, res, *Vx);
      print_vectors(step, *Vx, r, y);
    }
  while (state == SolverControl::iterate);
  return IterationResult(false, state, step, res);
}


template<typename VectorType>
template<typename MatrixType, typename PreconditionerType>
void
SolverBicgstab<VectorType>::solve(const MatrixType         &A,
                                  VectorType               &x,
                                  const VectorType         &b,
                                  const PreconditionerType &precondition)
{
  deallog.push("Bicgstab");
  Vr    = this->memory.alloc();
  Vr->reinit(x, true);
  Vrbar = this->memory.alloc();
  Vrbar->reinit(x, true);
  Vp    = this->memory.alloc();
  Vp->reinit(x, true);
  Vy    = this->memory.alloc();
  Vy->reinit(x, true);
  Vz    = this->memory.alloc();
  Vz->reinit(x, true);
  Vt    = this->memory.alloc();
  Vt->reinit(x, true);
  Vv    = this->memory.alloc();
  Vv->reinit(x, true);

  Vx = &x;
  Vb = &b;

  step = 0;

  IterationResult state(false,SolverControl::failure,0,0);

  // iterate while the inner iteration returns a breakdown
  do
    {
      if (step != 0)
        deallog << "Restart step " << step << std::endl;
      if (start(A) == SolverControl::success)
        {
          state.state = SolverControl::success;
          break;
        }
      state = iterate(A, precondition);
    }
  while (state.breakdown == true);

  this->memory.free(Vr);
  this->memory.free(Vrbar);
  this->memory.free(Vp);
  this->memory.free(Vy);
  this->memory.free(Vz);
  this->memory.free(Vt);
  this->memory.free(Vv);

  deallog.pop();

  // in case of failure: throw exception
  AssertThrow(state.state == SolverControl::success,
              SolverControl::NoConvergence (state.last_step,
                                            state.last_residual));
  // otherwise exit as normal
}

#endif // DOXYGEN

DEAL_II_NAMESPACE_CLOSE

#endif