Xpetra Version of the Day
Loading...
Searching...
No Matches
Xpetra_TpetraHalfPrecisionOperator.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 XPETRA_TPETRAHALFPRECISIONOPEARTOR_HPP
47#define XPETRA_TPETRAHALFPRECISIONOPEARTOR_HPP
48
49#include "Xpetra_ConfigDefs.hpp"
50
52
53#include <Tpetra_CrsMatrix.hpp>
55#include <Xpetra_MultiVector.hpp>
56#include <Xpetra_CrsMatrixWrap.hpp>
57#include <Xpetra_MultiVectorFactory.hpp>
58
59namespace Xpetra {
60
61 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
62 RCP<Xpetra::Matrix<typename Teuchos::ScalarTraits<Scalar>::halfPrecision, LocalOrdinal, GlobalOrdinal, Node> >
64#if defined(HAVE_TPETRA_INST_DOUBLE) && defined(HAVE_TPETRA_INST_FLOAT)
65 using HalfScalar = typename Teuchos::ScalarTraits<Scalar>::halfPrecision;
67
68 RCP<XpCrs> xpCrs = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(A, true)->getCrsMatrix();
69 auto tpCrs = Teuchos::rcp_dynamic_cast<Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(xpCrs, true)->getTpetra_CrsMatrix();
70 auto newTpCrs = tpCrs->template convert<HalfScalar>();
73
74 return newA;
75#else
77 "Xpetra::convertToHalfPrecision only available for Tpetra with SC=double and SC=float enabled");
79#endif
80 }
81
82 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
83 RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::halfPrecision, LocalOrdinal, GlobalOrdinal, Node> >
85#if defined(HAVE_TPETRA_INST_DOUBLE) && defined(HAVE_TPETRA_INST_FLOAT)
86 using HalfScalar = typename Teuchos::ScalarTraits<Scalar>::halfPrecision;
87
88 auto tpX = toTpetra(*X);
89 auto newTpX = tpX.template convert<HalfScalar>();
91
92 return newX;
93#else
95 "Xpetra::convertToHalfPrecision only available for Tpetra with SC=double and SC=float enabled");
97#endif
98 }
99
102 template <class Scalar,
103 class LocalOrdinal,
104 class GlobalOrdinal,
106 class TpetraHalfPrecisionOperator : public Xpetra::Operator<Scalar,LocalOrdinal,GlobalOrdinal,Node> {
107 public:
109
111
112
115 Allocate(1);
116 }
117
118 void Allocate(int numVecs) {
121 }
122
125
127
130 return Op_->getDomainMap();
131 }
132
135 return Op_->getRangeMap();
136 }
137
139
146 Scalar alpha = Teuchos::ScalarTraits<Scalar>::one(),
147 Scalar beta = Teuchos::ScalarTraits<Scalar>::one()) const{
149 Tpetra::deep_copy(*Teuchos::rcp_dynamic_cast<tMVHalf>(X_)->getTpetra_MultiVector(),
150 toTpetra(X));
151 Op_->apply(*X_, *Y_, mode, Teuchos::as<HalfScalar>(alpha), Teuchos::as<HalfScalar>(beta));
152 Tpetra::deep_copy(toTpetra(Y),
153 *Teuchos::rcp_dynamic_cast<tMVHalf>(Y_)->getTpetra_MultiVector());
154 }
155
157 bool hasTransposeApply() const { return false; }
158
164 R.update(STS::one(),B,STS::zero());
165 this->apply (X, R, Teuchos::NO_TRANS, -STS::one(), STS::one());
166 }
167
169
170
173
175
177
178 private:
181 };
182
183} // namespace
184
185#endif // XPETRA_TPETRAHALFPRECISIONOPEARTOR_HPP
Concrete implementation of Xpetra::Matrix.
Exception throws to report errors in the internal logical of the program.
Xpetra-specific matrix class.
static Teuchos::RCP< MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Build(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Constructor specifying the number of non-zeros for all rows.
virtual void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta)=0
Update multi-vector values with scaled values of A, this = beta*this + alpha*A.
Wraps an existing halfer precision Xpetra::Operator as a Xpetra::Operator.
void SetHalfPrecisionOperator(const RCP< Xpetra::Operator< HalfScalar, LocalOrdinal, GlobalOrdinal, Node > > &op)
Teuchos::ScalarTraits< Scalar >::halfPrecision HalfScalar
RCP< Xpetra::MultiVector< HalfScalar, LocalOrdinal, GlobalOrdinal, Node > > X_
RCP< Xpetra::MultiVector< HalfScalar, LocalOrdinal, GlobalOrdinal, Node > > Y_
void residual(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &R) const
Compute a residual R = B - (*this) * X.
Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Tpetra::Map object associated with the range of this TpetraOperator.
void apply(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::one()) const
Returns in Y the result of a Xpetra::TpetraOperator applied to a Xpetra::MultiVector X.
RCP< Xpetra::Operator< HalfScalar, LocalOrdinal, GlobalOrdinal, Node > > Op_
TpetraHalfPrecisionOperator(const RCP< Xpetra::Operator< HalfScalar, LocalOrdinal, GlobalOrdinal, Node > > &op)
Constructor.
RCP< Xpetra::Operator< HalfScalar, LocalOrdinal, GlobalOrdinal, Node > > GetHalfPrecisionOperator() const
Direct access to the underlying TpetraOperator.
Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Tpetra::Map object associated with the domain of this TpetraOperator.
bool hasTransposeApply() const
Indicates whether this TpetraOperator supports applying the adjoint TpetraOperator.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
#define TEUCHOS_UNREACHABLE_RETURN(dummyReturnVal)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Xpetra namespace
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
RCP< Xpetra::Matrix< typename Teuchos::ScalarTraits< Scalar >::halfPrecision, LocalOrdinal, GlobalOrdinal, Node > > convertToHalfPrecision(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A)