DUNE-FEM (unstable)

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