KASKADE 7 development version
bddcSpaceTransfer.hh
Go to the documentation of this file.
1
2/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
3/* */
4/* This file is part of the library KASKADE 7 */
5/* https://www.zib.de/research/projects/kaskade7-finite-element-toolbox */
6/* */
7/* Copyright (C) 2024-2024 Zuse Institute Berlin */
8/* */
9/* KASKADE 7 is distributed under the terms of the ZIB Academic License. */
10/* see $KASKADE/academic.txt */
11/* */
12/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
13
14#ifndef BDDCSPACETRANSFER_HH
15#define BDDCSPACETRANSFER_HH
16
17#include "dune/common/fvector.hh"
18#include <dune/istl/bvector.hh>
19
20#include <vector>
21#include <map>
22namespace Kaskade::BDDC
23{
24
34 using Interfaces = std::vector<std::pair<int,std::vector<int>>>;
35
36 // ----------------------------------------------------------------------------------------------
37
42 template <int m, class Scalar>
44 {
45 public:
49
53 void setInterfaces(int id, Interfaces const& ifs);
54
61 void setInterfaceAveragingWeights(int s, XVector const& rho);
62
69 void getInterfaceAveragingWeights(int s, XVector& rho) const;
70
71 protected:
72
73 // Returns the interface-local indices of dofs shared with subdomain s.
74 std::vector<int> const& localIdx(int s) const;
75
76 int id; // own subdomain number
77
78 // We think in interface dofs and indices. We therefore maintain a mapping of
79 // interface dof indices to subdomain dof indices (globalIndex), and a
80 // list of interface indices (instead of subdomain indices) for each interface group
81 // (localIfs).
82 std::vector<int> globalIndex;
84
85 // This contains for each shared dof the own averaging weight w_i as well as the
86 // sum of all weights (i.e. the averaging denominator).
89 };
90
91 // ----------------------------------------------------------------------------------------------
92
118 template <int m, class Scalar=double, class TransmissionScalar=Scalar>
119 class SpaceTransfer: public SpaceTransferBase<m,Scalar>
120 {
122
123 public:
127
129
133 using TransmissionType = std::vector<TEntry>;
134
159 void preProlongate(Vector const& u);
160
165 void getInterfaceValues(int s, TransmissionType& values) const;
166
171 void setInterfaceValues(int s, TransmissionType const& values);
172
179 void postProlongate(Vector const& u, Vector& du) const;
180
203 void preRestrict(Vector const& r);
204
206
208
214
220 void setInActiveIdx(int gIdx);
221
222 void get_globalIndex(std::set<int> & globalIndex_notActive)
223 {
224 std::map<int, int>::iterator it;
225 for (it = globalIndex_active.begin(); it != globalIndex_active.end(); it++)
226 {
227 if(not it->second)
228 {
229 globalIndex_notActive.insert(it->first);
230 }
231 }
232 }
233
234 private:
235 using Base::localIdx;
236 using Base::globalIndex;
237 using Base::avgWeightSum;
238 using Base::myAvgWeight;
239
240 // We store the interface values in transmission message format to be sent to
241 // neighboring subdomains on demand. On restriction, this contains residual updates.
242 // On prolongation, this contains weighted correction values.
243 TransmissionType localIfValues;
244
245 // On prolongation, we need to accumulate the weighted interface values received from
246 // neighbors.
247 TransmissionType ifValues;
248
249 // We keep track of the local residual that we already have sent to neighbors (which may
250 // be inexact due to compression or mixed precision), such that we can send only the
251 // updates.
252 XVector localPublishedResidual;
253
254 // Correspondingly, we maintain the summed interface residual, which is updated with the
255 // (possibly inexact) residual updates we receive from the neighbors.
256 XVector residualSum;
257
258 std::map<int,int> globalIndex_active;
259 };
260
261
262}
263
264#endif
Base class for BDDC SpaceTransfer implementations.
std::vector< int > const & localIdx(int s) const
void setInterfaceAveragingWeights(int s, XVector const &rho)
Receives the averaging weights from the given subdomain.
void getInterfaceAveragingWeights(int s, XVector &rho) const
Provides the averaging weights to the given subdomain.
void setInterfaces(int id, Interfaces const &ifs)
Provides the interface dof list.
A trivial BDDC space transfer.
void postProlongate(Vector const &u, Vector &du) const
Computes difference of averaged interface values to own interface values and imbues them into the sub...
void getInterfaceResidualUpdate(int s, TransmissionType &dres) const
Dune::BlockVector< XEntry > XVector
void setInterfaceResidualUpdate(int s, TransmissionType const &dres)
void get_globalIndex(std::set< int > &globalIndex_notActive)
void preProlongate(Vector const &u)
Prepares exchange of interface values.
void postRestrict(Vector &r)
void preRestrict(Vector const &r)
void setInterfaceValues(int s, TransmissionType const &values)
Receives the weighted (and possibly compressed )interface dof values of an adjacent subdomain.
std::vector< TEntry > TransmissionType
the data type that is used for actual data exchange between subdomains
void getInterfaceValues(int s, TransmissionType &values) const
Provides the weighted own interface dof values , possibly compressed.
std::vector< std::pair< int, std::vector< int > > > Interfaces
A data structure describing the interfaces of a subdomain.