MueLu  Version of the Day
MueLu_CoupledAggregationCommHelper_decl.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // MueLu: A package for multigrid based preconditioning
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact
39 // Jonathan Hu (jhu@sandia.gov)
40 // Andrey Prokopenko (aprokop@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 #ifndef MUELU_COUPLEDAGGREGATIONCOMMHELPER_DECL_HPP
47 #define MUELU_COUPLEDAGGREGATIONCOMMHELPER_DECL_HPP
48 
49 #include <Xpetra_Import_fwd.hpp>
50 #include <Xpetra_ImportFactory_fwd.hpp>
51 #include <Xpetra_Vector_fwd.hpp>
52 
53 #include "MueLu_ConfigDefs.hpp"
54 #include "MueLu_BaseClass.hpp"
56 
57 #include "MueLu_Aggregates.hpp"
58 
59 namespace MueLu {
60 
68  template <class LocalOrdinal = int, class GlobalOrdinal = LocalOrdinal, class Node = KokkosClassic::DefaultNode::DefaultNodeType>
70 
71  typedef double Scalar; // Scalar type only used for weight: always a double.
72 #undef MUELU_COUPLEDAGGREGATIONCOMMHELPER_SHORT
73 #include "MueLu_UseShortNames.hpp"
74 
75  public:
76 
78 
79 
81  CoupledAggregationCommHelper(const RCP<const Map> & uniqueMap, const RCP<const Map> & nonUniqueMap);
82 
85 
87 
99  void ArbitrateAndCommunicate(Vector &weights, Aggregates &aggregates, const bool perturb) const {
100  ArbitrateAndCommunicate(weights, *aggregates.GetProcWinner(), &*aggregates.GetVertex2AggId(), perturb);
101  }
102 
161  /*
162  Output:
163  @param weight \f$ weight[k] \Leftarrow \max(weight[k_1],\dots,weight[k_n]) \f$
164  where \f$ weight[k_j] \f$ live on different processors
165  but have the same GlobalId as weight[k] on this processor.
166 
167  @param procWinner procWinner[k] <-- MyPid associated with the
168  kj yielding the max in
169  Max(weight[k1],...,weight[kn]) .
170  See weight Output comments.
171  NOTE: If all input weight[kj]'s are zero,
172  then procWinner[k] is left untouched.
173 
174  @param companion If not null,
175  companion[k] <-- companion[kj] where
176  companion[kj] lives on processor procWinner[k].
177  and corresponds to the same GlobalId as k.
178  NOTE: If for a particlar GlobalId, no processor
179  has a value of procWinner that matches
180  its MyPid, the corresponding companion
181  is not altered.
182  */
183  void ArbitrateAndCommunicate(Vector &weight, LOVector &procWinner, LOVector *companion, const bool perturb) const; //ArbitrateAndCommunicate(Vector&, LOVector &, LOVector *, const bool) const
184 
205  void NonUnique2NonUnique(const Vector &source, Vector &dest, const Xpetra::CombineMode what) const;
206 
207  private:
208  RCP<const Import> import_;
209  mutable RCP<const Import> winnerImport_; //FIXME get rid of "mutable"
210  mutable RCP<Import> pushWinners_; //FIXME get rid of mutable
211  RCP<Vector> tempVec_;
212  mutable RCP<Vector> perturbWt_;
213  mutable RCP<Vector> postComm_;
214  mutable RCP<Vector> candidateWinners_;
215  mutable ArrayRCP<GO> myWinners_;
216  mutable int numMyWinners_;
217  mutable RCP<Map> winnerMap_;
218  mutable int numCalls_;
219  int myPID_;
220 
221  // uniqueMap A subset of weight.getMap() where each GlobalId
222  // has only one unique copy on one processor.
223  // Normally, weight.getMap() would have both locals
224  // and ghost elements while uniqueMap would just
225  // have the locals. It should be possible to
226  // remove this or make it an optional argument
227  // and use some existing Epetra/Tpetra capability to
228  // make a uniqueMap.
229  //
230  // import_ This corresponds precisely to
231  // Import import_(
232  // weight.getMap(), uniqueMap);
233  // This could also be eliminated and created
234  // here, but for efficiency user's must pass in.
235  //
236  };
237 
238 
239 }
240 
241 //JG:
242 // - procWinner is an array of proc ID -> LocalOrdinal
243 // - companion == aggregates.GetVertex2AggId() == local aggregate ID -> LocalOrdinal
244 
245 #define MUELU_COUPLEDAGGREGATIONCOMMHELPER_SHORT
246 #endif // MUELU_COUPLEDAGGREGATIONCOMMHELPER_DECL_HPP
const RCP< LOVector > & GetProcWinner() const
Returns constant vector that maps local node IDs to owning processor IDs.
Container class for aggregation information.
void NonUnique2NonUnique(const Vector &source, Vector &dest, const Xpetra::CombineMode what) const
Redistribute data in source to dest where both source and dest might have multiple copies of the same...
Namespace for MueLu classes and methods.
const RCP< LOVector > & GetVertex2AggId() const
Returns constant vector that maps local node IDs to local aggregates IDs.
void ArbitrateAndCommunicate(Vector &weights, Aggregates &aggregates, const bool perturb) const
This method assigns unknowns to aggregates.
Helper class for providing arbitrated communication across processors.
Base class for MueLu classes.
CoupledAggregationCommHelper(const RCP< const Map > &uniqueMap, const RCP< const Map > &nonUniqueMap)
Constructor.