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