Xpetra Version of the Day
Loading...
Searching...
No Matches
Xpetra_IteratorOps.hpp
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
10#ifndef PACKAGES_XPETRA_SUP_UTILS_XPETRA_ITERATOROPS_HPP_
11#define PACKAGES_XPETRA_SUP_UTILS_XPETRA_ITERATOROPS_HPP_
12
13#include "Xpetra_ConfigDefs.hpp"
14
15#include "Xpetra_Matrix.hpp"
16#include "Xpetra_MatrixMatrix.hpp"
17#include "Xpetra_CrsMatrixWrap.hpp"
18
19#include "Xpetra_Map.hpp"
20#include "Xpetra_StridedMap.hpp"
21#include "Xpetra_StridedMapFactory.hpp"
22#include "Xpetra_MapExtractor.hpp"
23#include "Xpetra_MatrixFactory.hpp"
24#include "Xpetra_BlockedCrsMatrix.hpp"
25
26namespace Xpetra {
27
28// General implementation
29// Epetra variant throws
30template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
31void Jacobi(
32 Scalar omega,
37 bool call_FillComplete_on_result = true,
38 bool doOptimizeStorage = true,
39 const std::string& label = std::string(),
40 const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) {
41 typedef Scalar SC;
42 typedef LocalOrdinal LO;
43 typedef GlobalOrdinal GO;
44 typedef Node NO;
45
47 "XpetraExt::MatrixMatrix::Jacobi: row map of C is not same as row map of A")
49 "XpetraExt::MatrixMatrix::Jacobi: row map of C is not same as row map of B");
52
53 bool haveMultiplyDoFillComplete = call_FillComplete_on_result && doOptimizeStorage;
54
55 if (C.getRowMap()->lib() == Xpetra::UseEpetra) {
56#ifndef HAVE_XPETRA_EPETRAEXT
57 throw(Xpetra::Exceptions::RuntimeError("Xpetra::MatrixMatrix::Jacobi requires EpetraExt to be compiled."));
58#else
59 throw(Xpetra::Exceptions::RuntimeError("Xpetra::MatrixMatrix::Jacobi requires you to use an Epetra-compatible data type."));
60#endif
61 } else if (C.getRowMap()->lib() == Xpetra::UseTpetra) {
62#ifdef HAVE_XPETRA_TPETRA
63 const Tpetra::CrsMatrix<SC, LO, GO, NO>& tpA = Xpetra::Helpers<SC, LO, GO, NO>::Op2TpetraCrs(A);
64 const Tpetra::CrsMatrix<SC, LO, GO, NO>& tpB = Xpetra::Helpers<SC, LO, GO, NO>::Op2TpetraCrs(B);
65 Tpetra::CrsMatrix<SC, LO, GO, NO>& tpC = Xpetra::Helpers<SC, LO, GO, NO>::Op2NonConstTpetraCrs(C);
67 Tpetra::MatrixMatrix::Jacobi(omega, *tpD, tpA, tpB, tpC, haveMultiplyDoFillComplete, label, params);
68#else
69 throw(Xpetra::Exceptions::RuntimeError("Xpetra must be compiled with Tpetra."));
70#endif
71 }
72
73 if (call_FillComplete_on_result && !haveMultiplyDoFillComplete) {
75 fillParams->set("Optimize Storage", doOptimizeStorage);
76 C.fillComplete(B.getDomainMap(), B.getRangeMap(), fillParams);
77 }
78
79 // transfer striding information
80 RCP<Xpetra::Matrix<SC, LO, GO, NO> > rcpA = Teuchos::rcp_const_cast<Xpetra::Matrix<SC, LO, GO, NO> >(Teuchos::rcpFromRef(A));
81 RCP<Xpetra::Matrix<SC, LO, GO, NO> > rcpB = Teuchos::rcp_const_cast<Xpetra::Matrix<SC, LO, GO, NO> >(Teuchos::rcpFromRef(B));
82 C.CreateView("stridedMaps", rcpA, false, rcpB, false); // TODO use references instead of RCPs
83} // end Jacobi
84
85#if defined(HAVE_XPETRA_EPETRA) && !defined(XPETRA_EPETRA_NO_32BIT_GLOBAL_INDICES)
86template <>
92 bool call_FillComplete_on_result,
93 bool doOptimizeStorage,
94 const std::string& label,
96#endif
97
98#if defined(HAVE_XPETRA_EPETRA) && !defined(XPETRA_EPETRA_NO_64BIT_GLOBAL_INDICES)
99template <>
105 bool call_FillComplete_on_result,
106 bool doOptimizeStorage,
107 const std::string& label,
109#endif
110
118template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
120#undef XPETRA_ITERATOROPS_SHORT
122 public:
123 static RCP<Matrix>
124 Jacobi(SC omega, const Vector& Dinv, const Matrix& A, const Matrix& B, RCP<Matrix> C_in, Teuchos::FancyOStream& fos, const std::string& label, RCP<ParameterList>& params) {
127
128 RCP<Matrix> C = C_in;
129 if (C == Teuchos::null) {
131 } else {
132 C->resumeFill(); // why this is not done inside of Tpetra Jacobi?
133 fos << "Reuse C pattern" << std::endl;
134 }
135
136 Xpetra::Jacobi<Scalar, LocalOrdinal, GlobalOrdinal, Node>(omega, Dinv, A, B, *C, true, true, label, params);
137 C->CreateView("stridedMaps", rcpFromRef(A), false, rcpFromRef(B), false);
138
139 return C;
140 }
141};
142
143} // end namespace Xpetra
144
145#define XPETRA_ITERATOROPS_SHORT
146
147#endif /* PACKAGES_XPETRA_SUP_UTILS_XPETRA_ITERATOROPS_HPP_ */
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< const Tpetra::CrsMatrix< SC, LO, GO, NO > > Op2TpetraCrs(RCP< Matrix > Op)
Xpetra utility class containing iteration operators.
static RCP< Matrix > Jacobi(SC omega, const Vector &Dinv, const Matrix &A, const Matrix &B, RCP< Matrix > C_in, Teuchos::FancyOStream &fos, const std::string &label, RCP< ParameterList > &params)
static RCP< Matrix > Build(const RCP< const Map > &rowMap)
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)
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)
void Jacobi(Scalar omega, const Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Dinv, const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &C, bool call_FillComplete_on_result=true, bool doOptimizeStorage=true, const std::string &label=std::string(), const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)