MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_LeftoverAggregationAlgorithm_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_LEFTOVERAGGREGATIONALGORITHM_DECL_HPP
47 #define MUELU_LEFTOVERAGGREGATIONALGORITHM_DECL_HPP
48 
49 #include <Xpetra_Vector_fwd.hpp>
50 #include <Xpetra_VectorFactory_fwd.hpp>
51 
52 #include "MueLu_ConfigDefs.hpp"
53 #include "MueLu_BaseClass.hpp"
55 
56 #include "MueLu_Aggregates_fwd.hpp"
57 #include "MueLu_GraphBase.hpp"
58 
60 
61 namespace MueLu {
62 
63  template<class LocalOrdinal = DefaultLocalOrdinal,
65  class Node = DefaultNode>
67 #undef MUELU_LEFTOVERAGGREGATIONALGORITHM_SHORT
69 
70  typedef GO global_size_t; //TODO
71  typedef LO my_size_t; //TODO
72 
73  public:
74 
76 
77 
80 
83 
85 
87 
88 
89  void SetMinNodesPerAggregate(int minNodesPerAggregate) { minNodesPerAggregate_ = minNodesPerAggregate; }
90  void SetPhase3AggCreation(double phase3AggCreation) { phase3AggCreation_ = phase3AggCreation; }
91 
92  double GetPhase3AggCreation() const { return phase3AggCreation_; }
94 
95  // TODO: Set/GetGraphName
97 
99 
100 
101 #define MUELU_NOSCORE -100 /* indicates that a quality score has not */
102  /* yet been assigned when determining to */
103  /* which existing aggregate a vertex */
104  /* should be assigned. */
105 
106 #define MUELU_DISTONE_VERTEX_WEIGHT 100 /* Weights associated with all */
107  /* vertices that have a direct connection */
108  /* to the aggregate root. */
109 
110 #define INCR_SCALING 3 /* Determines how much of a penalty should */
111  /* be deduced from a score during Phase 5 */
112  /* for each Phase 5 vertex already added */
113  /* to this aggregate. Specifically the */
114  /* penalty associated with aggregate y is */
115  /* max (INCR_SCALING*NNewVtx, */
116  /* UnpenalizedScore*(1- */
117  /* MUELU_PENALTYFACTOR)) */
118  /* where NNewVtx is the number of phase 5 */
119  /* vertices already assigned to y. */
120 
121 #define MUELU_PENALTYFACTOR .30 /* Determines maximum allowable */
122  /* percentage of a score that can be */
123  /* deducted from this score for having */
124  /* already enlargened an aggregate to */
125  /* which we are contemplated adding another*/
126  /* vertex. Should be between 0 and 1. */
127 
128 
292  void AggregateLeftovers(GraphBase const &graph, Aggregates &aggregates) const; //AggregateLeftovers
293 
295 
297 
298 
303  void RootCandidates(my_size_t nVertices, ArrayView<const LO> & vertex2AggId, GraphBase const &graph,
304  ArrayRCP<LO> &candidates, my_size_t &nCandidates, global_size_t &nCandidatesGlobal) const; //RootCandidates
305 
307  int RemoveSmallAggs(Aggregates& aggregates, int min_size,
308  RCP<Xpetra::Vector<double,LO,GO,NO> > & distWeights, const MueLu::CoupledAggregationCommHelper<LO,GO,NO> & myWidget) const; //RemoveSmallAggs
309 
311 
312  private:
315 
316  // JG TODO: rename variables:
317  // Adjacent-> adjacent
318  // homogenization of variables names :
319  // - colj and j
320  // - i and iNode
321  // - k->kNode
322  // - ...
323 
324  }; //class LeftoverAggregationAlgorithm
325 
326 } //namespace MueLu
327 
328 // graphName_("UC_CleanUp")
329 
330 #define MUELU_LEFTOVERAGGREGATIONALGORITHM_SHORT
331 #endif // MUELU_LEFTOVERAGGREGATIONALGORITHM_DECL_HPP
void RootCandidates(my_size_t nVertices, ArrayView< const LO > &vertex2AggId, GraphBase const &graph, ArrayRCP< LO > &candidates, my_size_t &nCandidates, global_size_t &nCandidatesGlobal) const
Build a list of candidate root nodes.
MueLu::DefaultLocalOrdinal LocalOrdinal
KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
void AggregateLeftovers(GraphBase const &graph, Aggregates &aggregates) const
Take a partially aggregated graph and complete the aggregation.
MueLu::DefaultNode Node
Helper class for providing arbitrated communication across processors.
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Base class for MueLu classes.
int RemoveSmallAggs(Aggregates &aggregates, int min_size, RCP< Xpetra::Vector< double, LO, GO, NO > > &distWeights, const MueLu::CoupledAggregationCommHelper< LO, GO, NO > &myWidget) const
Attempt to clean up aggregates that are too small.