This file is indexed.

/usr/include/dune/istl/paamg/twolevelmethod.hh is in libdune-istl-dev 2.5.1-1.

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
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#ifndef DUNE_ISTL_TWOLEVELMETHOD_HH
#define DUNE_ISTL_TWOLEVELMETHOD_HH

#include <tuple>

#include<dune/istl/operators.hh>
#include"amg.hh"
#include"galerkin.hh"
#include<dune/istl/solver.hh>

#include<dune/common/unused.hh>

/**
 * @addtogroup ISTL_PAAMG
 * @{
 * @file
 * @author Markus Blatt
 * @brief Algebraic twolevel methods.
 */
namespace Dune
{
namespace Amg
{

/**
 * @brief Abstract base class for transfer between levels and creation
 * of the coarse level system.
 *
 * @tparam FO The type of the linear operator of the finel level system. Has to be
 * derived from AssembledLinearOperator.
 * @tparam CO The type of the linear operator of the coarse level system. Has to be
 * derived from AssembledLinearOperator.
 */
template<class FO, class CO>
class LevelTransferPolicy
{
public:
  /**
   * @brief The linear operator of the finel level system. Has to be
   * derived from AssembledLinearOperator.
   */
  typedef FO FineOperatorType;
  /**
   * @brief The type of the range of the fine level operator.
   */
  typedef typename FineOperatorType::range_type FineRangeType;
  /**
   * @brief The type of the domain of the fine level operator.
   */
  typedef typename FineOperatorType::domain_type FineDomainType;
  /**
   * @brief The linear operator of the finel level system. Has to be
   * derived from AssembledLinearOperator.
   */
  typedef CO CoarseOperatorType;
  /**
   * @brief The type of the range of the coarse level operator.
   */
  typedef typename CoarseOperatorType::range_type CoarseRangeType;
  /**
   * @brief The type of the domain of the coarse level operator.
   */
  typedef typename CoarseOperatorType::domain_type CoarseDomainType;
  /**
   * @brief Get the coarse level operator.
   * @return A shared pointer to the coarse level system.
   */
  std::shared_ptr<CoarseOperatorType>& getCoarseLevelOperator()
  {
    return operator_;
  }
  /**
   * @brief Get the coarse level right hand side.
   * @return The coarse level right hand side.
   */
  CoarseRangeType& getCoarseLevelRhs()
  {
    return rhs_;
  }

  /**
   * @brief Get the coarse level left hand side.
   * @return The coarse level leftt hand side.
   */
  CoarseDomainType& getCoarseLevelLhs()
  {
    return lhs_;
  }
  /**
   * @brief Transfers the data to the coarse level.
   *
   * Restricts the residual to the right hand side of the
   * coarse level system and initialies the left hand side
   * of the coarse level system. These can afterwards be accessed
   * usinf getCoarseLevelRhs() and getCoarseLevelLhs().
   * @param fineDefect The current residual of the fine level system.
   */
  virtual void moveToCoarseLevel(const FineRangeType& fineRhs)=0;
  /**
   * @brief Updates the fine level linear system after the correction
   * of the coarse levels system.
   *
   * After returning from this function the coarse level correction
   * will have been added to fine level system.
   * @param[inout] fineLhs The left hand side of the fine level to update
   * with the coarse level correction.
   */
  virtual void moveToFineLevel(FineDomainType& fineLhs)=0;
  /**
   * @brief Algebraically creates the coarse level system.
   *
   * After returning from this function the coarse level operator
   * can be accessed using getCoarseLevelOperator().
   * @param fineOperator The operator of the fine level system.
   */
  virtual void createCoarseLevelSystem(const FineOperatorType& fineOperator)=0;

  /** @brief Clone the current object. */
  virtual LevelTransferPolicy* clone() const =0;

  /** @brief Destructor. */
  virtual ~LevelTransferPolicy(){}

  protected:
  /** @brief The coarse level rhs. */
  CoarseRangeType rhs_;
  /** @brief The coarse level lhs. */
  CoarseDomainType lhs_;
  /** @brief the coarse level linear operator. */
  std::shared_ptr<CoarseOperatorType> operator_;
};

/**
 * @brief A LeveTransferPolicy that used aggregation to construct the coarse level system.
 * @tparam O The type of the fine and coarse level operator.
 * @tparam C The criterion that describes the aggregation procedure.
 */
template<class O, class C>
class AggregationLevelTransferPolicy
  : public LevelTransferPolicy<O,O>
{
  typedef Dune::Amg::AggregatesMap<typename O::matrix_type::size_type> AggregatesMap;
public:
  typedef LevelTransferPolicy<O,O> FatherType;
  typedef C Criterion;
  typedef SequentialInformation ParallelInformation;

  AggregationLevelTransferPolicy(const Criterion& crit)
  : criterion_(crit)
  {}

  void createCoarseLevelSystem(const O& fineOperator)
  {
    prolongDamp_ = criterion_.getProlongationDampingFactor();
    GalerkinProduct<ParallelInformation> productBuilder;
    typedef typename Dune::Amg::MatrixGraph<const typename O::matrix_type> MatrixGraph;
    typedef typename Dune::Amg::PropertiesGraph<MatrixGraph,Dune::Amg::VertexProperties,
      Dune::Amg::EdgeProperties,Dune::IdentityMap,Dune::IdentityMap> PropertiesGraph;
    MatrixGraph mg(fineOperator.getmat());
    PropertiesGraph pg(mg,Dune::IdentityMap(),Dune::IdentityMap());
    typedef NegateSet<typename ParallelInformation::OwnerSet> OverlapFlags;

    aggregatesMap_.reset(new AggregatesMap(pg.maxVertex()+1));

    int noAggregates, isoAggregates, oneAggregates, skippedAggregates;

    std::tie(noAggregates, isoAggregates, oneAggregates, skippedAggregates) =
       aggregatesMap_->buildAggregates(fineOperator.getmat(), pg, criterion_, true);
     std::cout<<"no aggregates="<<noAggregates<<" iso="<<isoAggregates<<" one="<<oneAggregates<<" skipped="<<skippedAggregates<<std::endl;
    // misuse coarsener to renumber aggregates
    Dune::Amg::IndicesCoarsener<Dune::Amg::SequentialInformation,int> renumberer;
    typedef std::vector<bool>::iterator Iterator;
    typedef Dune::IteratorPropertyMap<Iterator, Dune::IdentityMap> VisitedMap;
    std::vector<bool> excluded(fineOperator.getmat().N(), false);
    VisitedMap vm(excluded.begin(), Dune::IdentityMap());
    ParallelInformation pinfo;
    std::size_t aggregates = renumberer.coarsen(pinfo, pg, vm,
                                                *aggregatesMap_, pinfo,
                                                noAggregates);
    std::vector<bool>& visited=excluded;

    typedef std::vector<bool>::iterator Iterator;

    for(Iterator iter= visited.begin(), end=visited.end();
        iter != end; ++iter)
          *iter=false;
    matrix_.reset(productBuilder.build(mg, vm,
                                       SequentialInformation(),
                                       *aggregatesMap_,
                                       aggregates,
                                       OverlapFlags()));
    productBuilder.calculate(fineOperator.getmat(), *aggregatesMap_, *matrix_, pinfo, OverlapFlags());
    this->lhs_.resize(this->matrix_->M());
    this->rhs_.resize(this->matrix_->N());
    this->operator_.reset(new O(*matrix_));
  }

  void moveToCoarseLevel(const typename FatherType::FineRangeType& fineRhs)
  {
    Transfer<std::size_t,typename FatherType::FineRangeType,ParallelInformation>
      ::restrictVector(*aggregatesMap_, this->rhs_, fineRhs, ParallelInformation());
    this->lhs_=0;
  }

  void moveToFineLevel(typename FatherType::FineDomainType& fineLhs)
  {
    Transfer<std::size_t,typename FatherType::FineRangeType,ParallelInformation>
      ::prolongateVector(*aggregatesMap_, this->lhs_, fineLhs,
                         prolongDamp_, ParallelInformation());
  }

  AggregationLevelTransferPolicy* clone() const
  {
    return new AggregationLevelTransferPolicy(*this);
  }

private:
  typename O::matrix_type::field_type prolongDamp_;
  std::shared_ptr<AggregatesMap> aggregatesMap_;
  Criterion criterion_;
  std::shared_ptr<typename O::matrix_type> matrix_;
};

/**
 * @brief A policy class for solving the coarse level system using one step of AMG.
 * @tparam O The type of the linear operator used.
 * @tparam S The type of the smoother used in AMG.
 * @tparam C The type of the crition used for the aggregation within AMG.
 */
template<class O, class S, class C>
class OneStepAMGCoarseSolverPolicy
{
public:
  /** @brief The type of the linear operator used. */
  typedef O Operator;
  /** @brief The type of the range and domain of the operator. */
  typedef typename O::range_type X;
  /** @brief The type of the crition used for the aggregation within AMG.*/
  typedef C Criterion;
  /** @brief The type of the smoother used in AMG. */
  typedef S Smoother;
  /** @brief The type of the arguments used for constructing the smoother. */
  typedef typename Dune::Amg::SmootherTraits<S>::Arguments SmootherArgs;
  /** @brief The type of the AMG construct on the coarse level.*/
  typedef AMG<Operator,X,Smoother> AMGType;
  /**
   * @brief Constructs the coarse solver policy.
   * @param args The arguments used for constructing the smoother.
   * @param c The crition used for the aggregation within AMG.
   */
  OneStepAMGCoarseSolverPolicy(const SmootherArgs& args, const Criterion& c)
    : smootherArgs_(args), criterion_(c)
  {}
  /** @brief Copy constructor. */
  OneStepAMGCoarseSolverPolicy(const OneStepAMGCoarseSolverPolicy& other)
  : coarseOperator_(other.coarseOperator_), smootherArgs_(other.smootherArgs_),
    criterion_(other.criterion_)
  {}
private:
  /**
   * @brief A wrapper that makes an inverse operator out of AMG.
   *
   * The operator will use one step of AMG to approximately solve
   * the coarse level system.
   */
  struct AMGInverseOperator : public InverseOperator<X,X>
  {
    AMGInverseOperator(const typename AMGType::Operator& op,
                       const Criterion& crit,
                       const typename AMGType::SmootherArgs& args)
      : amg_(op, crit,args), first_(true)
    {}

    void apply(X& x, X& b, double reduction, InverseOperatorResult& res)
    {
      DUNE_UNUSED_PARAMETER(reduction);
      DUNE_UNUSED_PARAMETER(res);
      if(first_)
      {
        amg_.pre(x,b);
        first_=false;
        x_=x;
      }
      amg_.apply(x,b);
    }

    void apply(X& x, X& b, InverseOperatorResult& res)
    {
      return apply(x,b,1e-8,res);
    }

    ~AMGInverseOperator()
    {
      if(!first_)
        amg_.post(x_);
    }
    AMGInverseOperator(const AMGInverseOperator& other)
    : x_(other.x_), amg_(other.amg_), first_(other.first_)
    {
    }
  private:
    X x_;
    AMGType amg_;
    bool first_;
  };

public:
  /** @brief The type of solver constructed for the coarse level. */
  typedef AMGInverseOperator CoarseLevelSolver;

  /**
    * @brief Constructs a coarse level solver.
    *
    * @param transferPolicy The policy describing the transfer between levels.
    * @return A pointer to the constructed coarse level solver.
    * @tparam P The type of the level transfer policy.
    */
  template<class P>
  CoarseLevelSolver* createCoarseLevelSolver(P& transferPolicy)
  {
    coarseOperator_=transferPolicy.getCoarseLevelOperator();
    AMGInverseOperator* inv = new AMGInverseOperator(*coarseOperator_,
                                                     criterion_,
                                                     smootherArgs_);

    return inv; //std::shared_ptr<InverseOperator<X,X> >(inv);

  }

private:
  /** @brief The coarse level operator. */
  std::shared_ptr<Operator> coarseOperator_;
  /** @brief The arguments used to construct the smoother. */
  SmootherArgs smootherArgs_;
  /** @brief The coarsening criterion. */
  Criterion criterion_;
};

/**
 * @tparam FO The type of the fine level linear operator.
 * @tparam CSP The type of the coarse level solver policy.
 * @tparam S The type of the fine level smoother used.
 */
template<class FO, class CSP, class S>
class TwoLevelMethod :
    public Preconditioner<typename FO::domain_type, typename FO::range_type>
{
public:
  /** @brief The type of the policy for constructing the coarse level solver. */
  typedef CSP CoarseLevelSolverPolicy;
  /** @brief The type of the coarse level solver. */
  typedef typename CoarseLevelSolverPolicy::CoarseLevelSolver CoarseLevelSolver;
  /**
   * @brief The linear operator of the finel level system. Has to be
   * derived from AssembledLinearOperator.
   */
  typedef FO FineOperatorType;
  /**
   * @brief The type of the range of the fine level operator.
   */
  typedef typename FineOperatorType::range_type FineRangeType;
  /**
   * @brief The type of the domain of the fine level operator.
   */
  typedef typename FineOperatorType::domain_type FineDomainType;
  /**
   * @brief The linear operator of the finel level system. Has to be
   * derived from AssembledLinearOperator.
   */
  typedef typename CSP::Operator CoarseOperatorType;
  /**
   * @brief The type of the range of the coarse level operator.
   */
  typedef typename CoarseOperatorType::range_type CoarseRangeType;
  /**
   * @brief The type of the domain of the coarse level operator.
   */
  typedef typename CoarseOperatorType::domain_type CoarseDomainType;
  /**
   * @brief The type of the fine level smoother.
   */
  typedef S SmootherType;
  // define the category
  enum {
    //! \brief The category the preconditioner is part of.
    category=SolverCategory::sequential
  };

  /**
   * @brief Constructs a two level method.
   *
   * @tparam CoarseSolverPolicy The policy for constructing the coarse
   * solver, e.g. OneStepAMGCoarseSolverPolicy
   * @param op The fine level operator.
   * @param smoother The fine level smoother.
   * @param policy The level transfer policy.
   * @param coarsePolicy The policy for constructing the coarse level solver.
   * @param preSteps The number of smoothing steps to apply before the coarse
   * level correction.
   * @param preSteps The number of smoothing steps to apply after the coarse
   * level correction.
   */
  TwoLevelMethod(const FineOperatorType& op,
                 std::shared_ptr<SmootherType> smoother,
                 const LevelTransferPolicy<FineOperatorType,
                                           CoarseOperatorType>& policy,
                 CoarseLevelSolverPolicy& coarsePolicy,
                 std::size_t preSteps=1, std::size_t postSteps=1)
    : operator_(&op), smoother_(smoother),
      preSteps_(preSteps), postSteps_(postSteps)
  {
    policy_ = policy.clone();
    policy_->createCoarseLevelSystem(*operator_);
    coarseSolver_=coarsePolicy.createCoarseLevelSolver(*policy_);
  }

  TwoLevelMethod(const TwoLevelMethod& other)
  : operator_(other.operator_), coarseSolver_(new CoarseLevelSolver(*other.coarseSolver_)),
    smoother_(other.smoother_), policy_(other.policy_->clone()),
    preSteps_(other.preSteps_), postSteps_(other.postSteps_)
  {}

  ~TwoLevelMethod()
  {
    // Each instance has its own policy.
    delete policy_;
    delete coarseSolver_;
  }

  void pre(FineDomainType& x, FineRangeType& b)
  {
    smoother_->pre(x,b);
  }

  void post(FineDomainType& x)
  {
    DUNE_UNUSED_PARAMETER(x);
  }

  void apply(FineDomainType& v, const FineRangeType& d)
  {
    FineDomainType u(v);
    FineRangeType rhs(d);
    LevelContext context;
    SequentialInformation info;
    context.pinfo=&info;
    context.lhs=&u;
    context.update=&v;
    context.smoother=smoother_;
    context.rhs=&rhs;
    context.matrix=operator_;
    // Presmoothing
    presmooth(context, preSteps_);
    //Coarse grid correction
    policy_->moveToCoarseLevel(*context.rhs);
    InverseOperatorResult res;
    coarseSolver_->apply(policy_->getCoarseLevelLhs(), policy_->getCoarseLevelRhs(), res);
    *context.lhs=0;
    policy_->moveToFineLevel(*context.lhs);
    *context.update += *context.lhs;
    // Postsmoothing
    postsmooth(context, postSteps_);

  }

private:
  /**
   * @brief Struct containing the level information.
   */
  struct LevelContext
  {
    /** @brief The type of the smoother used. */
    typedef S SmootherType;
    /** @brief A pointer to the smoother. */
    std::shared_ptr<SmootherType> smoother;
    /** @brief The left hand side passed to the and returned by the smoother. */
    FineDomainType* lhs;
    /*
     * @brief The right hand side holding the current residual.
     *
     * This is passed to the smoother as the right hand side.
     */
    FineRangeType* rhs;
    /**
     * @brief The total update calculated by the preconditioner.
     *
     * I.e. all update from smoothing and coarse grid correction summed up.
     */
    FineDomainType* update;
    /** @parallel information */
    SequentialInformation* pinfo;
    /**
     * @brief The matrix that we are solving.
     *
     * Needed to update the residual.
     */
    const FineOperatorType* matrix;
  };
  const FineOperatorType* operator_;
  /** @brief The coarse level solver. */
  CoarseLevelSolver* coarseSolver_;
  /** @brief The fine level smoother. */
  std::shared_ptr<S> smoother_;
  /** @brief Policy for prolongation, restriction, and coarse level system creation. */
  LevelTransferPolicy<FO,typename CSP::Operator>* policy_;
  /** @brief The number of presmoothing steps to apply. */
  std::size_t preSteps_;
  /** @brief The number of postsmoothing steps to apply. */
  std::size_t postSteps_;
};
}// end namespace Amg
}// end namespace Dune

/** @} */
#endif