This file is indexed.

/usr/include/dune/istl/paamg/transfer.hh is in libdune-istl-dev 2.4.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
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#ifndef DUNE_AMGTRANSFER_HH
#define DUNE_AMGTRANSFER_HH

#include <dune/istl/bvector.hh>
#include <dune/istl/matrixredistribute.hh>
#include <dune/istl/paamg/pinfo.hh>
#include <dune/istl/owneroverlapcopy.hh>
#include <dune/istl/paamg/aggregates.hh>
#include <dune/common/exceptions.hh>
#include <dune/common/unused.hh>

namespace Dune
{
  namespace Amg
  {

    /**
     * @addtogroup ISTL_PAAMG
     *
     * @{
     */

    /** @file
     * @author Markus Blatt
     * @brief Prolongation and restriction for amg.
     */
    template<class V1, class V2, class T>
    class Transfer
    {

    public:
      typedef V1 Vertex;
      typedef V2 Vector;

      template<typename T1, typename R>
      static void prolongateVector(const AggregatesMap<Vertex>& aggregates, Vector& coarse, Vector& fine,
                                   Vector& fineRedist,T1 damp, R& redistributor=R());

      template<typename T1, typename R>
      static void prolongateVector(const AggregatesMap<Vertex>& aggregates, Vector& coarse, Vector& fine,
                                   T1 damp);

      static void restrictVector(const AggregatesMap<Vertex>& aggregates, Vector& coarse, const Vector& fine,
                                 T& comm);
    };

    template<class V,class V1>
    class Transfer<V,V1, SequentialInformation>
    {
    public:
      typedef V Vertex;
      typedef V1 Vector;
      typedef RedistributeInformation<SequentialInformation> Redist;
      template<typename T1>
      static void prolongateVector(const AggregatesMap<Vertex>& aggregates, Vector& coarse, Vector& fine,
                                   Vector& fineRedist, T1 damp,
                                   const SequentialInformation& comm=SequentialInformation(),
                                   const Redist& redist=Redist());
      template<typename T1>
      static void prolongateVector(const AggregatesMap<Vertex>& aggregates, Vector& coarse, Vector& fine,
                                   T1 damp,
                                   const SequentialInformation& comm=SequentialInformation());


      static void restrictVector(const AggregatesMap<Vertex>& aggregates, Vector& coarse, const Vector& fine,
                                 const SequentialInformation& comm);
    };

#if HAVE_MPI

    template<class V,class V1, class T1, class T2>
    class Transfer<V,V1,OwnerOverlapCopyCommunication<T1,T2> >
    {
    public:
      typedef V Vertex;
      typedef V1 Vector;
      typedef RedistributeInformation<OwnerOverlapCopyCommunication<T1,T2> > Redist;
      template<typename T3>
      static void prolongateVector(const AggregatesMap<Vertex>& aggregates, Vector& coarse, Vector& fine,
                                   Vector& fineRedist, T3 damp, OwnerOverlapCopyCommunication<T1,T2>& comm,
                                   const Redist& redist);
      template<typename T3>
      static void prolongateVector(const AggregatesMap<Vertex>& aggregates, Vector& coarse, Vector& fine,
                                   T3 damp, OwnerOverlapCopyCommunication<T1,T2>& comm);

      static void restrictVector(const AggregatesMap<Vertex>& aggregates, Vector& coarse, const Vector& fine,
                                 OwnerOverlapCopyCommunication<T1,T2>& comm);
    };

#endif

    template<class V, class V1>
    template<typename T>
    inline void
    Transfer<V,V1,SequentialInformation>::prolongateVector(const AggregatesMap<Vertex>& aggregates,
                                                           Vector& coarse, Vector& fine, Vector& fineRedist,
                                                           T damp,
                                                           const SequentialInformation& comm,
                                                           const Redist& redist)
    {
      DUNE_UNUSED_PARAMETER(fineRedist);
      DUNE_UNUSED_PARAMETER(comm);
      DUNE_UNUSED_PARAMETER(redist);
      prolongateVector(aggregates, coarse, fine, damp);
    }
    template<class V, class V1>
    template<typename T>
    inline void
    Transfer<V,V1,SequentialInformation>::prolongateVector(const AggregatesMap<Vertex>& aggregates,
                                                           Vector& coarse, Vector& fine,
                                                           T damp,
                                                           const SequentialInformation& comm)
    {
      DUNE_UNUSED_PARAMETER(comm);
      typedef typename Vector::iterator Iterator;

      Iterator end = coarse.end();
      Iterator begin= coarse.begin();
      for(; begin!=end; ++begin)
        *begin*=damp;
      end=fine.end();
      begin=fine.begin();

      for(Iterator block=begin; block != end; ++block) {
        std::ptrdiff_t index=block-begin;
        const Vertex& vertex = aggregates[index];
        if(vertex != AggregatesMap<Vertex>::ISOLATED)
          *block += coarse[aggregates[index]];
      }
    }

    template<class V, class V1>
    inline void
    Transfer<V,V1,SequentialInformation>::restrictVector(const AggregatesMap<Vertex>& aggregates,
                                                         Vector& coarse,
                                                         const Vector& fine,
                                                         const SequentialInformation& comm)
    {
      DUNE_UNUSED_PARAMETER(comm);
      // Set coarse vector to zero
      coarse=0;

      typedef typename Vector::const_iterator Iterator;
      Iterator end = fine.end();
      Iterator begin=fine.begin();

      for(Iterator block=begin; block != end; ++block) {
        const Vertex& vertex = aggregates[block-begin];
        if(vertex != AggregatesMap<Vertex>::ISOLATED)
          coarse[vertex] += *block;
      }
    }

#if HAVE_MPI
    template<class V, class V1, class T1, class T2>
    template<typename T3>
    inline void Transfer<V,V1,OwnerOverlapCopyCommunication<T1,T2> >::prolongateVector(const AggregatesMap<Vertex>& aggregates,
                                                                                       Vector& coarse, Vector& fine,
                                                                                       Vector& fineRedist, T3 damp,
                                                                                       OwnerOverlapCopyCommunication<T1,T2>& comm,
                                                                                       const Redist& redist)
    {
      if(fineRedist.size()>0)
        // we operated on the coarse level
        Transfer<V,V1,SequentialInformation>::prolongateVector(aggregates, coarse, fineRedist, damp);

      // TODO This could be accomplished with one communication, too!
      redist.redistributeBackward(fine, fineRedist);
      comm.copyOwnerToAll(fine,fine);
    }

    template<class V, class V1, class T1, class T2>
    template<typename T3>
    inline void Transfer<V,V1,OwnerOverlapCopyCommunication<T1,T2> >::prolongateVector(const AggregatesMap<Vertex>& aggregates,
                                                                                       Vector& coarse, Vector& fine,
                                                                                       T3 damp,
                                                                                       OwnerOverlapCopyCommunication<T1,T2>& comm)
    {
      DUNE_UNUSED_PARAMETER(comm);
      Transfer<V,V1,SequentialInformation>::prolongateVector(aggregates, coarse, fine, damp);
    }
    template<class V, class V1, class T1, class T2>
    inline void Transfer<V,V1,OwnerOverlapCopyCommunication<T1,T2> >::restrictVector(const AggregatesMap<Vertex>& aggregates,
                                                                                     Vector& coarse, const Vector& fine,
                                                                                     OwnerOverlapCopyCommunication<T1,T2>& comm)
    {
      Transfer<V,V1,SequentialInformation>::restrictVector(aggregates, coarse, fine, SequentialInformation());
      // We need this here to avoid it in the smoothers on the coarse level.
      // There (in the preconditioner d is const.
      comm.project(coarse);
    }
#endif
    /** @} */
  }    // namspace Amg
}     // namspace Dune
#endif