Dune Core Modules (2.5.2)

transfer.hh
Go to the documentation of this file.
1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 #ifndef DUNE_AMGTRANSFER_HH
4 #define DUNE_AMGTRANSFER_HH
5 
6 #include <dune/istl/bvector.hh>
8 #include <dune/istl/paamg/pinfo.hh>
12 #include <dune/common/unused.hh>
13 
14 namespace Dune
15 {
16  namespace Amg
17  {
18 
29  template<class V1, class V2, class T>
30  class Transfer
31  {
32 
33  public:
34  typedef V1 Vertex;
35  typedef V2 Vector;
36 
37  template<typename T1, typename R>
38  static void prolongateVector(const AggregatesMap<Vertex>& aggregates, Vector& coarse, Vector& fine,
39  Vector& fineRedist,T1 damp, R& redistributor=R());
40 
41  template<typename T1, typename R>
42  static void prolongateVector(const AggregatesMap<Vertex>& aggregates, Vector& coarse, Vector& fine,
43  T1 damp);
44 
45  static void restrictVector(const AggregatesMap<Vertex>& aggregates, Vector& coarse, const Vector& fine,
46  T& comm);
47  };
48 
49  template<class V,class V1>
50  class Transfer<V,V1, SequentialInformation>
51  {
52  public:
53  typedef V Vertex;
54  typedef V1 Vector;
55  typedef RedistributeInformation<SequentialInformation> Redist;
56  template<typename T1>
57  static void prolongateVector(const AggregatesMap<Vertex>& aggregates, Vector& coarse, Vector& fine,
58  Vector& fineRedist, T1 damp,
59  const SequentialInformation& comm=SequentialInformation(),
60  const Redist& redist=Redist());
61  template<typename T1>
62  static void prolongateVector(const AggregatesMap<Vertex>& aggregates, Vector& coarse, Vector& fine,
63  T1 damp,
64  const SequentialInformation& comm=SequentialInformation());
65 
66 
67  static void restrictVector(const AggregatesMap<Vertex>& aggregates, Vector& coarse, const Vector& fine,
68  const SequentialInformation& comm);
69  };
70 
71 #if HAVE_MPI
72 
73  template<class V,class V1, class T1, class T2>
74  class Transfer<V,V1,OwnerOverlapCopyCommunication<T1,T2> >
75  {
76  public:
77  typedef V Vertex;
78  typedef V1 Vector;
79  typedef RedistributeInformation<OwnerOverlapCopyCommunication<T1,T2> > Redist;
80  template<typename T3>
81  static void prolongateVector(const AggregatesMap<Vertex>& aggregates, Vector& coarse, Vector& fine,
82  Vector& fineRedist, T3 damp, OwnerOverlapCopyCommunication<T1,T2>& comm,
83  const Redist& redist);
84  template<typename T3>
85  static void prolongateVector(const AggregatesMap<Vertex>& aggregates, Vector& coarse, Vector& fine,
86  T3 damp, OwnerOverlapCopyCommunication<T1,T2>& comm);
87 
88  static void restrictVector(const AggregatesMap<Vertex>& aggregates, Vector& coarse, const Vector& fine,
89  OwnerOverlapCopyCommunication<T1,T2>& comm);
90  };
91 
92 #endif
93 
94  template<class V, class V1>
95  template<typename T>
96  inline void
97  Transfer<V,V1,SequentialInformation>::prolongateVector(const AggregatesMap<Vertex>& aggregates,
98  Vector& coarse, Vector& fine, Vector& fineRedist,
99  T damp,
100  const SequentialInformation& comm,
101  const Redist& redist)
102  {
103  DUNE_UNUSED_PARAMETER(fineRedist);
104  DUNE_UNUSED_PARAMETER(comm);
105  DUNE_UNUSED_PARAMETER(redist);
106  prolongateVector(aggregates, coarse, fine, damp);
107  }
108  template<class V, class V1>
109  template<typename T>
110  inline void
111  Transfer<V,V1,SequentialInformation>::prolongateVector(const AggregatesMap<Vertex>& aggregates,
112  Vector& coarse, Vector& fine,
113  T damp,
114  const SequentialInformation& comm)
115  {
116  DUNE_UNUSED_PARAMETER(comm);
117  typedef typename Vector::iterator Iterator;
118 
119  Iterator end = coarse.end();
120  Iterator begin= coarse.begin();
121  for(; begin!=end; ++begin)
122  *begin*=damp;
123  end=fine.end();
124  begin=fine.begin();
125 
126  for(Iterator block=begin; block != end; ++block) {
127  std::ptrdiff_t index=block-begin;
128  const Vertex& vertex = aggregates[index];
129  if(vertex != AggregatesMap<Vertex>::ISOLATED)
130  *block += coarse[aggregates[index]];
131  }
132  }
133 
134  template<class V, class V1>
135  inline void
136  Transfer<V,V1,SequentialInformation>::restrictVector(const AggregatesMap<Vertex>& aggregates,
137  Vector& coarse,
138  const Vector& fine,
139  const SequentialInformation& comm)
140  {
141  DUNE_UNUSED_PARAMETER(comm);
142  // Set coarse vector to zero
143  coarse=0;
144 
145  typedef typename Vector::const_iterator Iterator;
146  Iterator end = fine.end();
147  Iterator begin=fine.begin();
148 
149  for(Iterator block=begin; block != end; ++block) {
150  const Vertex& vertex = aggregates[block-begin];
151  if(vertex != AggregatesMap<Vertex>::ISOLATED)
152  coarse[vertex] += *block;
153  }
154  }
155 
156 #if HAVE_MPI
157  template<class V, class V1, class T1, class T2>
158  template<typename T3>
159  inline void Transfer<V,V1,OwnerOverlapCopyCommunication<T1,T2> >::prolongateVector(const AggregatesMap<Vertex>& aggregates,
160  Vector& coarse, Vector& fine,
161  Vector& fineRedist, T3 damp,
162  OwnerOverlapCopyCommunication<T1,T2>& comm,
163  const Redist& redist)
164  {
165  if(fineRedist.size()>0)
166  // we operated on the coarse level
167  Transfer<V,V1,SequentialInformation>::prolongateVector(aggregates, coarse, fineRedist, damp);
168 
169  // TODO This could be accomplished with one communication, too!
170  redist.redistributeBackward(fine, fineRedist);
171  comm.copyOwnerToAll(fine,fine);
172  }
173 
174  template<class V, class V1, class T1, class T2>
175  template<typename T3>
176  inline void Transfer<V,V1,OwnerOverlapCopyCommunication<T1,T2> >::prolongateVector(const AggregatesMap<Vertex>& aggregates,
177  Vector& coarse, Vector& fine,
178  T3 damp,
179  OwnerOverlapCopyCommunication<T1,T2>& comm)
180  {
181  DUNE_UNUSED_PARAMETER(comm);
182  Transfer<V,V1,SequentialInformation>::prolongateVector(aggregates, coarse, fine, damp);
183  }
184  template<class V, class V1, class T1, class T2>
185  inline void Transfer<V,V1,OwnerOverlapCopyCommunication<T1,T2> >::restrictVector(const AggregatesMap<Vertex>& aggregates,
186  Vector& coarse, const Vector& fine,
187  OwnerOverlapCopyCommunication<T1,T2>& comm)
188  {
189  Transfer<V,V1,SequentialInformation>::restrictVector(aggregates, coarse, fine, SequentialInformation());
190  // We need this here to avoid it in the smoothers on the coarse level.
191  // There (in the preconditioner d is const.
192  comm.project(coarse);
193  }
194 #endif
196  } // namspace Amg
197 } // namspace Dune
198 #endif
Provides classes for the Coloring process of AMG.
This file implements a vector space as a tensor product of a given vector space. The number of compon...
A few common exception classes.
static const Vertex ISOLATED
Identifier of isolated vertices.
Definition: aggregates.hh:554
Functionality for redistributing a sparse matrix.
Dune namespace.
Definition: alignment.hh:11
Classes providing communication interfaces for overlapping Schwarz methods.
Definition of the DUNE_UNUSED macro for the case that config.h is not available.
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentionally unused function parameters with.
Definition: unused.hh:18
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 9, 22:29, 2024)