MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_LowPrecisionFactory_def.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_LOWPRECISIONFACTORY_DEF_HPP
47 #define MUELU_LOWPRECISIONFACTORY_DEF_HPP
48 
49 #include <Xpetra_Matrix.hpp>
50 #include <Xpetra_Operator.hpp>
53 
55 
56 #include "MueLu_Level.hpp"
57 #include "MueLu_Monitor.hpp"
58 
59 namespace MueLu {
60 
61 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
63  RCP<ParameterList> validParamList = rcp(new ParameterList());
64 
65  validParamList->set<std::string>("matrix key", "A", "");
66  validParamList->set<RCP<const FactoryBase> >("R", Teuchos::null, "Generating factory of the matrix A to be converted to lower precision");
67  validParamList->set<RCP<const FactoryBase> >("A", Teuchos::null, "Generating factory of the matrix A to be converted to lower precision");
68  validParamList->set<RCP<const FactoryBase> >("P", Teuchos::null, "Generating factory of the matrix A to be converted to lower precision");
69 
70  return validParamList;
71 }
72 
73 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
75  const ParameterList& pL = GetParameterList();
76  std::string matrixKey = pL.get<std::string>("matrix key");
77  Input(currentLevel, matrixKey);
78 }
79 
80 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
83 
84  const ParameterList& pL = GetParameterList();
85  std::string matrixKey = pL.get<std::string>("matrix key");
86 
87  FactoryMonitor m(*this, "Converting " + matrixKey + " to half precision", currentLevel);
88 
89  RCP<Matrix> A = Get<RCP<Matrix> >(currentLevel, matrixKey);
90 
91  GetOStream(Warnings) << "Matrix not converted to half precision. This only works for Tpetra and when both Scalar and HalfScalar have been instantiated." << std::endl;
92  Set(currentLevel, matrixKey, A);
93 }
94 
95 #if defined(HAVE_TPETRA_INST_DOUBLE) && defined(HAVE_TPETRA_INST_FLOAT)
96 template <class LocalOrdinal, class GlobalOrdinal, class Node>
98  RCP<ParameterList> validParamList = rcp(new ParameterList());
99 
100  validParamList->set<std::string>("matrix key", "A", "");
101  validParamList->set<RCP<const FactoryBase> >("R", Teuchos::null, "Generating factory of the matrix A to be converted to lower precision");
102  validParamList->set<RCP<const FactoryBase> >("A", Teuchos::null, "Generating factory of the matrix A to be converted to lower precision");
103  validParamList->set<RCP<const FactoryBase> >("P", Teuchos::null, "Generating factory of the matrix A to be converted to lower precision");
104 
105  return validParamList;
106 }
107 
108 template <class LocalOrdinal, class GlobalOrdinal, class Node>
110  const ParameterList& pL = GetParameterList();
111  std::string matrixKey = pL.get<std::string>("matrix key");
112  Input(currentLevel, matrixKey);
113 }
114 
115 template <class LocalOrdinal, class GlobalOrdinal, class Node>
118  using HalfScalar = typename Teuchos::ScalarTraits<Scalar>::halfPrecision;
119 
120  const ParameterList& pL = GetParameterList();
121  std::string matrixKey = pL.get<std::string>("matrix key");
122 
123  FactoryMonitor m(*this, "Converting " + matrixKey + " to half precision", currentLevel);
124 
125  RCP<Matrix> A = Get<RCP<Matrix> >(currentLevel, matrixKey);
126 
127  if ((A->getRowMap()->lib() == Xpetra::UseTpetra) && std::is_same<Scalar, double>::value) {
128  auto tpA = rcp_dynamic_cast<TpetraCrsMatrix>(rcp_dynamic_cast<CrsMatrixWrap>(A)->getCrsMatrix(), true)->getTpetra_CrsMatrix();
129  auto tpLowA = tpA->template convert<HalfScalar>();
131  auto xpTpLowOpA = rcp(new TpetraOperator(tpLowOpA));
132  auto xpLowOpA = rcp_dynamic_cast<Operator>(xpTpLowOpA);
133  Set(currentLevel, matrixKey, xpLowOpA);
134  return;
135  }
136 
137  GetOStream(Warnings) << "Matrix not converted to half precision. This only works for Tpetra and when both Scalar and HalfScalar have been instantiated." << std::endl;
138  Set(currentLevel, matrixKey, A);
139 }
140 #endif
141 
142 #if defined(HAVE_TPETRA_INST_COMPLEX_DOUBLE) && defined(HAVE_TPETRA_INST_COMPLEX_FLOAT)
143 template <class LocalOrdinal, class GlobalOrdinal, class Node>
144 RCP<const ParameterList> LowPrecisionFactory<std::complex<double>, LocalOrdinal, GlobalOrdinal, Node>::GetValidParameterList() const {
145  RCP<ParameterList> validParamList = rcp(new ParameterList());
146 
147  validParamList->set<std::string>("matrix key", "A", "");
148  validParamList->set<RCP<const FactoryBase> >("R", Teuchos::null, "Generating factory of the matrix A to be converted to lower precision");
149  validParamList->set<RCP<const FactoryBase> >("A", Teuchos::null, "Generating factory of the matrix A to be converted to lower precision");
150  validParamList->set<RCP<const FactoryBase> >("P", Teuchos::null, "Generating factory of the matrix A to be converted to lower precision");
151 
152  return validParamList;
153 }
154 
155 template <class LocalOrdinal, class GlobalOrdinal, class Node>
156 void LowPrecisionFactory<std::complex<double>, LocalOrdinal, GlobalOrdinal, Node>::DeclareInput(Level& currentLevel) const {
157  const ParameterList& pL = GetParameterList();
158  std::string matrixKey = pL.get<std::string>("matrix key");
159  Input(currentLevel, matrixKey);
160 }
161 
162 template <class LocalOrdinal, class GlobalOrdinal, class Node>
163 void LowPrecisionFactory<std::complex<double>, LocalOrdinal, GlobalOrdinal, Node>::Build(Level& currentLevel) const {
165  using HalfScalar = typename Teuchos::ScalarTraits<Scalar>::halfPrecision;
166 
167  const ParameterList& pL = GetParameterList();
168  std::string matrixKey = pL.get<std::string>("matrix key");
169 
170  FactoryMonitor m(*this, "Converting " + matrixKey + " to half precision", currentLevel);
171 
172  RCP<Matrix> A = Get<RCP<Matrix> >(currentLevel, matrixKey);
173 
174  if ((A->getRowMap()->lib() == Xpetra::UseTpetra) && std::is_same<Scalar, std::complex<double> >::value) {
175  auto tpA = rcp_dynamic_cast<TpetraCrsMatrix>(rcp_dynamic_cast<CrsMatrixWrap>(A)->getCrsMatrix(), true)->getTpetra_CrsMatrix();
176  auto tpLowA = tpA->template convert<HalfScalar>();
178  auto xpTpLowOpA = rcp(new TpetraOperator(tpLowOpA));
179  auto xpLowOpA = rcp_dynamic_cast<Operator>(xpTpLowOpA);
180  Set(currentLevel, matrixKey, xpLowOpA);
181  return;
182  }
183 
184  GetOStream(Warnings) << "Matrix not converted to half precision. This only works for Tpetra and when both Scalar and HalfScalar have been instantiated." << std::endl;
185  Set(currentLevel, matrixKey, A);
186 }
187 #endif
188 
189 } // namespace MueLu
190 
191 #endif // MUELU_LOWPRECISIONFACTORY_DEF_HPP
MueLu::DefaultLocalOrdinal LocalOrdinal
T & get(const std::string &name, T def_value)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Timer to be used in factories. Similar to Monitor but with additional timers.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
MueLu::DefaultNode Node
void Build(Level &currentLevel) const
Build method.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
void DeclareInput(Level &currentLevel) const
Input.
Print all warning messages.