MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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>
51 #include <Xpetra_Vector_fwd.hpp>
53 
54 #include "MueLu_ConfigDefs.hpp"
55 #include "MueLu_BaseClass.hpp"
57 
58 #include "MueLu_Aggregates.hpp"
59 
60 
61 namespace MueLu {
62 
70  template<class LocalOrdinal = DefaultLocalOrdinal,
72  class Node = DefaultNode>
75 
76 #undef MUELU_COUPLEDAGGREGATIONCOMMHELPER_SHORT
77 #include "MueLu_UseShortNames.hpp"
78 
79  public:
80 
82 
83 
85  CoupledAggregationCommHelper(const RCP<const Map> & uniqueMap, const RCP<const Map> & nonUniqueMap);
86 
89 
91 
103  void ArbitrateAndCommunicate(Vector &weights, Aggregates &aggregates, const bool perturb) const {
104  ArbitrateAndCommunicate(weights, *aggregates.GetProcWinner(), &*aggregates.GetVertex2AggId(), perturb);
105  }
106 
165  /*
166  Output:
167  @param weight \f$ weight[k] \Leftarrow \max(weight[k_1],\dots,weight[k_n]) \f$
168  where \f$ weight[k_j] \f$ live on different processors
169  but have the same GlobalId as weight[k] on this processor.
170 
171  @param procWinner procWinner[k] <-- MyPid associated with the
172  kj yielding the max in
173  Max(weight[k1],...,weight[kn]) .
174  See weight Output comments.
175  NOTE: If all input weight[kj]'s are zero,
176  then procWinner[k] is left untouched.
177 
178  @param companion If not null,
179  companion[k] <-- companion[kj] where
180  companion[kj] lives on processor procWinner[k].
181  and corresponds to the same GlobalId as k.
182  NOTE: If for a particlar GlobalId, no processor
183  has a value of procWinner that matches
184  its MyPid, the corresponding companion
185  is not altered.
186  */
187  void ArbitrateAndCommunicate(Vector &weight, LOVector &procWinner, LOMultiVector *companion, const bool perturb) const; //ArbitrateAndCommunicate(Vector&, LOVector &, LOMultiVector *, const bool) const
188 
209  void NonUnique2NonUnique(const Vector &source, Vector &dest, const Xpetra::CombineMode what) const;
210 
211  private:
213  mutable RCP<const Import> winnerImport_; //FIXME get rid of "mutable"
214  mutable RCP<Import> pushWinners_; //FIXME get rid of mutable
220  mutable int numMyWinners_;
222  mutable int numCalls_;
223  int myPID_;
224 
225  // uniqueMap A subset of weight.getMap() where each GlobalId
226  // has only one unique copy on one processor.
227  // Normally, weight.getMap() would have both locals
228  // and ghost elements while uniqueMap would just
229  // have the locals. It should be possible to
230  // remove this or make it an optional argument
231  // and use some existing Epetra/Tpetra capability to
232  // make a uniqueMap.
233  //
234  // import_ This corresponds precisely to
235  // Import import_(
236  // weight.getMap(), uniqueMap);
237  // This could also be eliminated and created
238  // here, but for efficiency user's must pass in.
239  //
240  };
241 
242 
243 }
244 
245 //JG:
246 // - procWinner is an array of proc ID -> LocalOrdinal
247 // - companion == aggregates.GetVertex2AggId() == local aggregate ID -> LocalOrdinal
248 
249 #define MUELU_COUPLEDAGGREGATIONCOMMHELPER_SHORT
250 #endif // MUELU_COUPLEDAGGREGATIONCOMMHELPER_DECL_HPP
MueLu::DefaultLocalOrdinal LocalOrdinal
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
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...
MueLu::DefaultNode Node
void ArbitrateAndCommunicate(Vector &weights, Aggregates &aggregates, const bool perturb) const
This method assigns unknowns to aggregates.
Helper class for providing arbitrated communication across processors.
Tpetra::Details::DefaultTypes::scalar_type DefaultScalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
const RCP< LOMultiVector > & GetVertex2AggId() const
Returns constant vector that maps local node IDs to local aggregates IDs.
Base class for MueLu classes.
CoupledAggregationCommHelper(const RCP< const Map > &uniqueMap, const RCP< const Map > &nonUniqueMap)
Constructor.