MueLu Version of the Day
MueLu_BrickAggregationFactory_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_BRICKAGGREGATIONFACTORY_DECL_HPP_
47#define MUELU_BRICKAGGREGATIONFACTORY_DECL_HPP_
48
49#include "MueLu_ConfigDefs.hpp"
50
51#include <Xpetra_Import_fwd.hpp>
52#include <Xpetra_ImportFactory_fwd.hpp>
53#include <Xpetra_Map_fwd.hpp>
54#include <Xpetra_MapFactory_fwd.hpp>
55#include <Xpetra_Matrix_fwd.hpp>
56#include <Xpetra_MultiVector_fwd.hpp>
57#include <Xpetra_MultiVectorFactory_fwd.hpp>
58
61
62#include "MueLu_Level_fwd.hpp"
64#include "MueLu_Exceptions.hpp"
66
76namespace MueLu {
77
78 template<class Scalar = DefaultScalar,
81 class Node = DefaultNode>
83#undef MUELU_BRICKAGGREGATIONFACTORY_SHORT
85 private:
86 typedef Teuchos::ScalarTraits<Scalar> STS;
87
88 // Comparator for doubles
89 // Generally, the coordinates for coarser levels would come out of averaging of fine level coordinates
90 // It is possible that the result of the averaging differs slightly between clusters, as we might have
91 // 3x2 and 2x2 cluster which would result in averaging 6 and 4 y-coordinates respectively, leading to
92 // slightly different results.
93 // Therefore, we hardcode a constant so that close points are considered the same.
94 class compare {
95 public:
96 bool operator()(const Scalar& x, const Scalar& y) const {
97 if (STS::magnitude(x - y) < 1e-14)
98 return false;
99 return STS::real(x) < STS::real(y);
100 }
101 };
102 typedef std::map<Scalar,GlobalOrdinal,compare> container;
103
104 public:
106
107
109 BrickAggregationFactory() : nDim_(-1), nx_(-1), ny_(-1), nz_(-1), bx_(-1), by_(-1), bz_(-1) { };
110
113
114 RCP<const ParameterList> GetValidParameterList() const;
115
117
118 // Options shared by all aggregation algorithms
119
121
122
123 void DeclareInput(Level &currentLevel) const;
124
126
128
129
131 void Build(Level &currentLevel) const;
132
134
135 private:
136 void Setup(const RCP<const Teuchos::Comm<int> >& comm, const RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO> >& coords, const RCP<const Map>& map) const;
137 RCP<container> Construct1DMap(const RCP<const Teuchos::Comm<int> >& comm, const ArrayRCP<const typename Teuchos::ScalarTraits<Scalar>::magnitudeType>& x) const;
138
139 bool isDirichlet(LocalOrdinal LID) const;
140 bool isRoot (LocalOrdinal LID) const;
143
144 void getIJK(LocalOrdinal LID, int &i, int &j, int &k) const;
145 void getAggIJK(LocalOrdinal LID, int &i, int &j, int &k) const;
146
147 mutable
148 int nDim_;
149 mutable
150 RCP<container> xMap_, yMap_, zMap_;
151 mutable
152 ArrayRCP<const typename Teuchos::ScalarTraits<Scalar>::magnitudeType> x_, y_, z_;
153 mutable
154 int nx_, ny_, nz_;
155 mutable
156 int bx_, by_, bz_;
157 mutable
159 mutable
161
162 mutable
163 std::map<GlobalOrdinal,GlobalOrdinal> revMap_;
164 }; // class BrickAggregationFactory
165
166}
167
168#define MUELU_BRICKAGGREGATIONFACTORY_SHORT
169#endif /* MUELU_BRICKAGGREGATIONFACTORY_DECL_HPP_ */
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
bool operator()(const Scalar &x, const Scalar &y) const
void Setup(const RCP< const Teuchos::Comm< int > > &comm, const RCP< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::magnitudeType, LO, GO, NO > > &coords, const RCP< const Map > &map) const
ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::magnitudeType > x_
GlobalOrdinal getRoot(LocalOrdinal LID) const
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
std::map< GlobalOrdinal, GlobalOrdinal > revMap_
ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::magnitudeType > y_
ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::magnitudeType > z_
GlobalOrdinal getAggGID(LocalOrdinal LID) const
void DeclareInput(Level &currentLevel) const
Input.
void Build(Level &currentLevel) const
Build aggregates.
void getAggIJK(LocalOrdinal LID, int &i, int &j, int &k) const
RCP< container > Construct1DMap(const RCP< const Teuchos::Comm< int > > &comm, const ArrayRCP< const typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &x) const
void getIJK(LocalOrdinal LID, int &i, int &j, int &k) const
std::map< Scalar, GlobalOrdinal, compare > container
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
Base class for factories that use one level (currentLevel).
Namespace for MueLu classes and methods.
KokkosClassic::DefaultNode::DefaultNodeType DefaultNode