Xpetra Version of the Day
Loading...
Searching...
No Matches
Xpetra_IteratorOps.hpp
Go to the documentation of this file.
1// @HEADER
2//
3// ***********************************************************************
4//
5// Xpetra: A linear algebra interface package
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// Tobias Wiesner (tawiesn@sandia.gov)
43//
44// ***********************************************************************
45//
46// @HEADER
47#ifndef PACKAGES_XPETRA_SUP_UTILS_XPETRA_ITERATOROPS_HPP_
48#define PACKAGES_XPETRA_SUP_UTILS_XPETRA_ITERATOROPS_HPP_
49
50#include "Xpetra_ConfigDefs.hpp"
51
52#include "Xpetra_Matrix.hpp"
54#include "Xpetra_CrsMatrixWrap.hpp"
55
56#include "Xpetra_Map.hpp"
57#include "Xpetra_StridedMap.hpp"
58#include "Xpetra_StridedMapFactory.hpp"
59#include "Xpetra_MapExtractor.hpp"
62
63namespace Xpetra {
64
65 // General implementation
66 // Epetra variant throws
67 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
68 void Jacobi(
69 Scalar omega,
74 bool call_FillComplete_on_result = true,
75 bool doOptimizeStorage = true,
76 const std::string & label = std::string(),
77 const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) {
78 typedef Scalar SC;
79 typedef LocalOrdinal LO;
80 typedef GlobalOrdinal GO;
81 typedef Node NO;
82
84 "XpetraExt::MatrixMatrix::Jacobi: row map of C is not same as row map of A")
86 "XpetraExt::MatrixMatrix::Jacobi: row map of C is not same as row map of B");
89
90 bool haveMultiplyDoFillComplete = call_FillComplete_on_result && doOptimizeStorage;
91
92 if (C.getRowMap()->lib() == Xpetra::UseEpetra) {
93#ifndef HAVE_XPETRA_EPETRAEXT
94 throw(Xpetra::Exceptions::RuntimeError("Xpetra::MatrixMatrix::Jacobi requires EpetraExt to be compiled."));
95#else
96 throw(Xpetra::Exceptions::RuntimeError("Xpetra::MatrixMatrix::Jacobi requires you to use an Epetra-compatible data type."));
97#endif
98 } else if (C.getRowMap()->lib() == Xpetra::UseTpetra) {
99#ifdef HAVE_XPETRA_TPETRA
100 const Tpetra::CrsMatrix<SC,LO,GO,NO> & tpA = Xpetra::Helpers<SC,LO,GO,NO>::Op2TpetraCrs(A);
101 const Tpetra::CrsMatrix<SC,LO,GO,NO> & tpB = Xpetra::Helpers<SC,LO,GO,NO>::Op2TpetraCrs(B);
102 Tpetra::CrsMatrix<SC,LO,GO,NO> & tpC = Xpetra::Helpers<SC,LO,GO,NO>::Op2NonConstTpetraCrs(C);
103 const RCP<Tpetra::Vector<SC,LO,GO,NO> > & tpD = toTpetra(Dinv);
104 Tpetra::MatrixMatrix::Jacobi(omega, *tpD, tpA, tpB, tpC, haveMultiplyDoFillComplete, label, params);
105#else
106 throw(Xpetra::Exceptions::RuntimeError("Xpetra must be compiled with Tpetra."));
107#endif
108 }
109
110 if (call_FillComplete_on_result && !haveMultiplyDoFillComplete) {
112 fillParams->set("Optimize Storage", doOptimizeStorage);
113 C.fillComplete(B.getDomainMap(), B.getRangeMap(), fillParams);
114 }
115
116 // transfer striding information
117 RCP<Xpetra::Matrix<SC,LO,GO,NO> > rcpA = Teuchos::rcp_const_cast<Xpetra::Matrix<SC,LO,GO,NO> >(Teuchos::rcpFromRef(A));
118 RCP<Xpetra::Matrix<SC,LO,GO,NO> > rcpB = Teuchos::rcp_const_cast<Xpetra::Matrix<SC,LO,GO,NO> >(Teuchos::rcpFromRef(B));
119 C.CreateView("stridedMaps", rcpA, false, rcpB, false); // TODO use references instead of RCPs
120 } // end Jacobi
121
122#if defined(HAVE_XPETRA_EPETRA) && !defined(XPETRA_EPETRA_NO_32BIT_GLOBAL_INDICES)
123 template<>
124 void Jacobi<double,int,int,EpetraNode>(double omega,
129 bool call_FillComplete_on_result,
130 bool doOptimizeStorage,
131 const std::string & label,
133#endif
134
135#if defined(HAVE_XPETRA_EPETRA) && !defined(XPETRA_EPETRA_NO_64BIT_GLOBAL_INDICES)
136 template<>
142 bool call_FillComplete_on_result,
143 bool doOptimizeStorage,
144 const std::string & label,
146#endif
147
155 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
157#undef XPETRA_ITERATOROPS_SHORT
159 public:
160
161 static RCP<Matrix>
162 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) {
165
166 RCP<Matrix> C = C_in;
167 if (C == Teuchos::null) {
169 } else {
170 C->resumeFill(); // why this is not done inside of Tpetra Jacobi?
171 fos << "Reuse C pattern" << std::endl;
172 }
173
174
175 Xpetra::Jacobi<Scalar,LocalOrdinal,GlobalOrdinal,Node>(omega, Dinv, A, B, *C, true, true,label,params);
176 C->CreateView("stridedMaps", rcpFromRef(A), false, rcpFromRef(B), false);
177
178 return C;
179 }
180
181 };
182
183} // end namespace Xpetra
184
185#define XPETRA_ITERATOROPS_SHORT
186
187#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 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 RCP< const Map > & getRowMap() const
Returns the Map that describes the row distribution in this matrix.
virtual Teuchos::RCP< const Map > getRangeMap() const =0
The Map associated with the range of this operator, which must be compatible with Y....
virtual Teuchos::RCP< const Map > getDomainMap() const =0
The Map associated with the domain of this operator, which must be compatible with X....
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Xpetra namespace
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)