Xpetra Version of the Day
Loading...
Searching...
No Matches
Xpetra_IteratorOps.cpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Xpetra: A linear algebra interface package
4//
5// Copyright 2012 NTESS and the Xpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
11
12namespace Xpetra {
13
14#if defined(HAVE_XPETRA_EPETRA) && !defined(XPETRA_EPETRA_NO_32BIT_GLOBAL_INDICES)
15template <>
21 bool call_FillComplete_on_result,
22 bool doOptimizeStorage,
23 const std::string& label,
25 typedef double SC;
26 typedef int LO;
27 typedef int GO;
28 typedef EpetraNode NO;
29
31 "XpetraExt::MatrixMatrix::Jacobi: row map of C is not same as row map of A")
33 "XpetraExt::MatrixMatrix::Jacobi: row map of C is not same as row map of B");
36
37 bool haveMultiplyDoFillComplete = call_FillComplete_on_result && doOptimizeStorage;
38
39 if (C.getRowMap()->lib() == Xpetra::UseEpetra) {
40#ifndef HAVE_XPETRA_EPETRAEXT
41 throw(Xpetra::Exceptions::RuntimeError("Xpetra::IteratorOps::Jacobi requires EpetraExt to be compiled."));
42#else
46 // FIXME
47 XPETRA_DYNAMIC_CAST(const EpetraVectorT<GO XPETRA_COMMA NO>, Dinv, epD, "Xpetra::IteratorOps::Jacobi() only accepts Xpetra::EpetraVector as input argument.");
48
49 int i = EpetraExt::MatrixMatrix::Jacobi(omega, *epD.getEpetra_Vector(), epA, epB, epC, haveMultiplyDoFillComplete);
50 if (haveMultiplyDoFillComplete) {
51 // Due to Epetra wrapper intricacies, we need to explicitly call
52 // fillComplete on Xpetra matrix here. Specifically, EpetraCrsMatrix
53 // only keeps an internal variable to check whether we are in resumed
54 // state or not, but never touches the underlying Epetra object. As
55 // such, we need to explicitly update the state of Xpetra matrix to
56 // that of Epetra one afterwords
57 C.fillComplete();
58 }
59
60 if (i != 0) {
61 std::ostringstream buf;
62 buf << i;
63 std::string msg = "EpetraExt::MatrixMatrix::Jacobi return value of " + buf.str();
64 throw(Exceptions::RuntimeError(msg));
65 }
66#endif
67 } else if (C.getRowMap()->lib() == Xpetra::UseTpetra) {
68#ifdef HAVE_XPETRA_TPETRA
69#if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
70 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
71 throw(Xpetra::Exceptions::RuntimeError("Xpetra must be compiled with Tpetra GO=<double,int,int> enabled."));
72#else
73 const Tpetra::CrsMatrix<SC, LO, GO, NO>& tpA = Xpetra::Helpers<SC, LO, GO, NO>::Op2TpetraCrs(A);
74 const Tpetra::CrsMatrix<SC, LO, GO, NO>& tpB = Xpetra::Helpers<SC, LO, GO, NO>::Op2TpetraCrs(B);
75 Tpetra::CrsMatrix<SC, LO, GO, NO>& tpC = Xpetra::Helpers<SC, LO, GO, NO>::Op2NonConstTpetraCrs(C);
77 Tpetra::MatrixMatrix::Jacobi(omega, *tpD, tpA, tpB, tpC, haveMultiplyDoFillComplete, label, params);
78#endif
79#else
80 throw(Xpetra::Exceptions::RuntimeError("Xpetra must be compiled with Tpetra."));
81#endif
82 }
83
84 if (call_FillComplete_on_result && !haveMultiplyDoFillComplete) {
86 ppp->set("Optimize Storage", doOptimizeStorage);
87 C.fillComplete(B.getDomainMap(), B.getRangeMap(), ppp);
88 }
89
90 // transfer striding information
91 Teuchos::RCP<Xpetra::Matrix<SC, LO, GO, NO> > rcpA = Teuchos::rcp_const_cast<Xpetra::Matrix<SC, LO, GO, NO> >(Teuchos::rcpFromRef(A));
92 Teuchos::RCP<Xpetra::Matrix<SC, LO, GO, NO> > rcpB = Teuchos::rcp_const_cast<Xpetra::Matrix<SC, LO, GO, NO> >(Teuchos::rcpFromRef(B));
93 C.CreateView("stridedMaps", rcpA, false, rcpB, false); // TODO use references instead of RCPs
94}
95#endif
96
97#if defined(HAVE_XPETRA_EPETRA) && !defined(XPETRA_EPETRA_NO_64BIT_GLOBAL_INDICES)
98template <>
104 bool call_FillComplete_on_result,
105 bool doOptimizeStorage,
106 const std::string& label,
108 typedef double SC;
109 typedef int LO;
110 typedef long long GO;
111 typedef EpetraNode NO;
112
114 "XpetraExt::MatrixMatrix::Jacobi: row map of C is not same as row map of A")
116 "XpetraExt::MatrixMatrix::Jacobi: row map of C is not same as row map of B");
119
120 bool haveMultiplyDoFillComplete = call_FillComplete_on_result && doOptimizeStorage;
121
122 if (C.getRowMap()->lib() == Xpetra::UseEpetra) {
123#ifndef HAVE_XPETRA_EPETRAEXT
124 throw(Xpetra::Exceptions::RuntimeError("Xpetra::IteratorOps::Jacobi requires EpetraExt to be compiled."));
125#else
129 // FIXME
130 XPETRA_DYNAMIC_CAST(const EpetraVectorT<GO XPETRA_COMMA NO>, Dinv, epD, "Xpetra::IteratorOps::Jacobi() only accepts Xpetra::EpetraVector as input argument.");
131
132 int i = EpetraExt::MatrixMatrix::Jacobi(omega, *epD.getEpetra_Vector(), epA, epB, epC, haveMultiplyDoFillComplete);
133 if (haveMultiplyDoFillComplete) {
134 // Due to Epetra wrapper intricacies, we need to explicitly call
135 // fillComplete on Xpetra matrix here. Specifically, EpetraCrsMatrix
136 // only keeps an internal variable to check whether we are in resumed
137 // state or not, but never touches the underlying Epetra object. As
138 // such, we need to explicitly update the state of Xpetra matrix to
139 // that of Epetra one afterwords
140 C.fillComplete();
141 }
142
143 if (i != 0) {
144 std::ostringstream buf;
145 buf << i;
146 std::string msg = "EpetraExt::MatrixMatrix::Jacobi return value of " + buf.str();
147 throw(Exceptions::RuntimeError(msg));
148 }
149#endif
150 } else if (C.getRowMap()->lib() == Xpetra::UseTpetra) {
151#ifdef HAVE_XPETRA_TPETRA
152#if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))) || \
153 (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))))
154 throw(Xpetra::Exceptions::RuntimeError("Xpetra must be compiled with Tpetra GO=<double,int,long long> enabled."));
155#else
156 const Tpetra::CrsMatrix<SC, LO, GO, NO>& tpA = Xpetra::Helpers<SC, LO, GO, NO>::Op2TpetraCrs(A);
157 const Tpetra::CrsMatrix<SC, LO, GO, NO>& tpB = Xpetra::Helpers<SC, LO, GO, NO>::Op2TpetraCrs(B);
158 Tpetra::CrsMatrix<SC, LO, GO, NO>& tpC = Xpetra::Helpers<SC, LO, GO, NO>::Op2NonConstTpetraCrs(C);
159 const RCP<Tpetra::Vector<SC, LO, GO, NO> >& tpD = toTpetra(Dinv);
160 Tpetra::MatrixMatrix::Jacobi(omega, *tpD, tpA, tpB, tpC, haveMultiplyDoFillComplete, label, params);
161#endif
162#else
163 throw(Xpetra::Exceptions::RuntimeError("Xpetra must be compiled with Tpetra."));
164#endif
165 }
166
167 if (call_FillComplete_on_result && !haveMultiplyDoFillComplete) {
169 ppp->set("Optimize Storage", doOptimizeStorage);
170 C.fillComplete(B.getDomainMap(), B.getRangeMap(), ppp);
171 }
172
173 // transfer striding information
174 Teuchos::RCP<Xpetra::Matrix<SC, LO, GO, NO> > rcpA = Teuchos::rcp_const_cast<Xpetra::Matrix<SC, LO, GO, NO> >(Teuchos::rcpFromRef(A));
175 Teuchos::RCP<Xpetra::Matrix<SC, LO, GO, NO> > rcpB = Teuchos::rcp_const_cast<Xpetra::Matrix<SC, LO, GO, NO> >(Teuchos::rcpFromRef(B));
176 C.CreateView("stridedMaps", rcpA, false, rcpB, false); // TODO use references instead of RCPs
177}
178#endif
179
180} // namespace Xpetra
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
LocalOrdinal LO
GlobalOrdinal GO
Exception throws to report errors in the internal logical of the program.
static RCP< Tpetra::CrsMatrix< SC, LO, GO, NO > > Op2NonConstTpetraCrs(RCP< Matrix > Op)
static RCP< Epetra_CrsMatrix > Op2NonConstEpetraCrs(RCP< Matrix > Op)
static RCP< const Tpetra::CrsMatrix< SC, LO, GO, NO > > Op2TpetraCrs(RCP< Matrix > Op)
Xpetra-specific matrix class.
void CreateView(viewLabel_t viewLabel, const RCP< const Map > &rowMap, const RCP< const Map > &colMap)
virtual const RCP< const Map > & getRowMap() const
Returns the Map that describes the row distribution in this matrix.
virtual bool isFillComplete() const =0
Returns true if fillComplete() has been called and the matrix is in compute mode.
virtual void fillComplete(const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const RCP< ParameterList > &params=null)=0
Signal that data entry is complete, specifying domain and range maps.
virtual const Teuchos::RCP< const map_type > getDomainMap() const =0
The Map associated with the domain of this operator, which must be compatible with X....
virtual const Teuchos::RCP< const map_type > getRangeMap() const =0
The Map associated with the range of this operator, which must be compatible with Y....
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Tpetra::KokkosCompat::KokkosSerialWrapperNode EpetraNode
void Jacobi< double, int, long long, EpetraNode >(double omega, const Xpetra::Vector< double, int, long long, EpetraNode > &Dinv, const Xpetra::Matrix< double, int, long long, EpetraNode > &A, const Xpetra::Matrix< double, int, long long, EpetraNode > &B, Xpetra::Matrix< double, int, long long, EpetraNode > &C, bool call_FillComplete_on_result, bool doOptimizeStorage, const std::string &label, const Teuchos::RCP< Teuchos::ParameterList > &params)
void Jacobi< double, int, int, EpetraNode >(double omega, const Xpetra::Vector< double, int, int, EpetraNode > &Dinv, const Xpetra::Matrix< double, int, int, EpetraNode > &A, const Xpetra::Matrix< double, int, int, EpetraNode > &B, Xpetra::Matrix< double, int, int, EpetraNode > &C, bool call_FillComplete_on_result, bool doOptimizeStorage, const std::string &label, const Teuchos::RCP< Teuchos::ParameterList > &params)
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)