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