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 
48 #ifndef PACKAGES_XPETRA_SUP_UTILS_XPETRA_ITERATOROPS_HPP_
49 #define PACKAGES_XPETRA_SUP_UTILS_XPETRA_ITERATOROPS_HPP_
50 
51 #include "Xpetra_ConfigDefs.hpp"
52 
53 #include "Xpetra_Matrix.hpp"
54 #include "Xpetra_MatrixMatrix.hpp"
55 #include "Xpetra_CrsMatrixWrap.hpp"
56 
57 #include "Xpetra_Map.hpp"
58 #include "Xpetra_StridedMap.hpp"
60 #include "Xpetra_MapExtractor.hpp"
61 #include "Xpetra_MatrixFactory.hpp"
63 
64 
65 
66 namespace Xpetra {
67 
68 
69 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
70 void Jacobi(
71  Scalar omega,
76  bool call_FillComplete_on_result = true,
77  bool doOptimizeStorage = true,
78  const std::string & label = std::string()) {
79 
80  if(C.getRowMap()->isSameAs(*A.getRowMap()) == false) {
81  std::string msg = "XpetraExt::MatrixMatrix::Jacobi: row map of C is not same as row map of A";
83  }
84  else if(C.getRowMap()->isSameAs(*B.getRowMap()) == false) {
85  std::string msg = "XpetraExt::MatrixMatrix::Jacobi: row map of C is not same as row map of B";
87  }
88 
89  if (!A.isFillComplete())
90  throw(Xpetra::Exceptions::RuntimeError("A is not fill-completed"));
91  if (!B.isFillComplete())
92  throw(Xpetra::Exceptions::RuntimeError("B is not fill-completed"));
93 
94  bool haveMultiplyDoFillComplete = call_FillComplete_on_result && doOptimizeStorage;
95 
96  if (C.getRowMap()->lib() == Xpetra::UseEpetra) {
97 #ifndef HAVE_XPETRA_EPETRAEXT
98  throw(Xpetra::Exceptions::RuntimeError("Xpetra::MatrixMatrix::Jacobi requires EpetraExt to be compiled."));
99 #else
100  throw(Xpetra::Exceptions::RuntimeError("Xpetra::MatrixMatrix::Jacobi requires you to use an Epetra-compatible data type."));
101 #endif
102  } else if (C.getRowMap()->lib() == Xpetra::UseTpetra) {
103 #ifdef HAVE_XPETRA_TPETRA
104  const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> & tpA = Xpetra::Helpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2TpetraCrs(A);
105  const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> & tpB = Xpetra::Helpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2TpetraCrs(B);
106  Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> & tpC = Xpetra::Helpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Op2NonConstTpetraCrs(C);
108  Tpetra::MatrixMatrix::Jacobi(omega,*tpD,tpA,tpB,tpC,haveMultiplyDoFillComplete,label);
109 #else
110  throw(Xpetra::Exceptions::RuntimeError("Xpetra must be compiled with Tpetra."));
111 #endif
112  }
113 
114  if(call_FillComplete_on_result && !haveMultiplyDoFillComplete) {
116  params->set("Optimize Storage",doOptimizeStorage);
117  C.fillComplete(B.getDomainMap(),B.getRangeMap(),params);
118  }
119 
120  // transfer striding information
123  C.CreateView("stridedMaps", rcpA, false, rcpB, false); // TODO use references instead of RCPs
124 } // end Jacobi
125 
126 
127 template <class GlobalOrdinal>
128 inline void JacobiT(
129  double omega,
134  bool call_FillComplete_on_result,
135  bool doOptimizeStorage,
136  const std::string & label) {
137 
138  typedef double SC;
139  typedef int LO;
140  typedef GlobalOrdinal GO;
142 
143  if(C.getRowMap()->isSameAs(*A.getRowMap()) == false) {
144  std::string msg = "XpetraExt::MatrixMatrix::Jacobi: row map of C is not same as row map of A";
146  }
147  else if(C.getRowMap()->isSameAs(*B.getRowMap()) == false) {
148  std::string msg = "XpetraExt::MatrixMatrix::Jacobi: row map of C is not same as row map of B";
150  }
151 
152  if (!A.isFillComplete())
153  throw(Xpetra::Exceptions::RuntimeError("A is not fill-completed"));
154  if (!B.isFillComplete())
155  throw(Xpetra::Exceptions::RuntimeError("B is not fill-completed"));
156 
157  bool haveMultiplyDoFillComplete = call_FillComplete_on_result && doOptimizeStorage;
158 
159  if (C.getRowMap()->lib() == Xpetra::UseEpetra) {
160 # ifndef HAVE_XPETRA_EPETRAEXT
161  throw(Xpetra::Exceptions::RuntimeError("Xpetra::IteratorOps::Jacobi requires EpetraExt to be compiled."));
162 #else
163  Epetra_CrsMatrix & epA = Xpetra::Helpers<SC,LO,GO,NO>::Op2NonConstEpetraCrs(A);
164  Epetra_CrsMatrix & epB = Xpetra::Helpers<SC,LO,GO,NO>::Op2NonConstEpetraCrs(B);
165  Epetra_CrsMatrix & epC = Xpetra::Helpers<SC,LO,GO,NO>::Op2NonConstEpetraCrs(C);
166  // const Epetra_Vector & epD = toEpetra(Dinv);
167  XPETRA_DYNAMIC_CAST(const EpetraVectorT<GO>, Dinv, epD, "Xpetra::IteratorOps::Jacobi() only accepts Xpetra::EpetraVector as input argument.");
168 
169  int i = EpetraExt::MatrixMatrix::Jacobi(omega,*epD.getEpetra_Vector(),epA,epB,epC,haveMultiplyDoFillComplete);
170  if (i != 0) {
171  std::ostringstream buf;
172  buf << i;
173  std::string msg = "EpetraExt::MatrixMatrix::Jacobi return value of " + buf.str();
174  throw(Exceptions::RuntimeError(msg));
175  }
176 #endif
177  } else if (C.getRowMap()->lib() == Xpetra::UseTpetra) {
178 #ifdef HAVE_XPETRA_TPETRA
179  const Tpetra::CrsMatrix<SC, LO, GO, NO> & tpA = Xpetra::Helpers<SC,LO,GO,NO>::Op2TpetraCrs(A);
180  const Tpetra::CrsMatrix<SC, LO, GO, NO> & tpB = Xpetra::Helpers<SC,LO,GO,NO>::Op2TpetraCrs(B);
181  Tpetra::CrsMatrix<SC, LO, GO, NO> & tpC = Xpetra::Helpers<SC,LO,GO,NO>::Op2NonConstTpetraCrs(C);
182  const RCP<Tpetra::Vector<SC, LO, GO, NO> > & tpD = toTpetra(Dinv);
183  Tpetra::MatrixMatrix::Jacobi(omega,*tpD,tpA,tpB,tpC,haveMultiplyDoFillComplete,label);
184 #else
185  throw(Xpetra::Exceptions::RuntimeError("Xpetra must be compiled with Tpetra."));
186 #endif
187  }
188 
189  if(call_FillComplete_on_result && !haveMultiplyDoFillComplete) {
191  params->set("Optimize Storage",doOptimizeStorage);
192  C.fillComplete(B.getDomainMap(),B.getRangeMap(),params);
193  }
194 
195  // transfer striding information
196  RCP<Xpetra::Matrix<SC, LO, GO, NO> > rcpA = Teuchos::rcp_const_cast<Xpetra::Matrix<SC, LO, GO, NO> >(Teuchos::rcpFromRef(A));
197  RCP<Xpetra::Matrix<SC, LO, GO, NO> > rcpB = Teuchos::rcp_const_cast<Xpetra::Matrix<SC, LO, GO, NO> >(Teuchos::rcpFromRef(B));
198  C.CreateView("stridedMaps", rcpA, false, rcpB, false); // TODO use references instead of RCPs
199 } // end Jacobi
200 
201 #ifndef XPETRA_EPETRA_NO_32BIT_GLOBAL_INDICES
202 inline void Jacobi(
203  double omega,
204  const Xpetra::Vector<double,int,int> & Dinv,
208  bool call_FillComplete_on_result,
209  bool doOptimizeStorage,
210  const std::string & label) {
211  JacobiT<int>(omega, Dinv, A, B, C, call_FillComplete_on_result, doOptimizeStorage, label);
212 }
213 #endif
214 
215 #ifdef HAVE_XPETRA_INT_LONG_LONG
216 inline void Jacobi(
217  double omega,
222  bool call_FillComplete_on_result,
223  bool doOptimizeStorage,
224  const std::string & label) {
225  JacobiT<long long>(omega, Dinv, A, B, C, call_FillComplete_on_result, doOptimizeStorage, label);
226 }
227 #endif // HAVE_XPETRA_INT_LONG_LONG
228 
236 template <class Scalar,
237 class LocalOrdinal = int,
238 class GlobalOrdinal = LocalOrdinal,
240 class IteratorOps {
241 
242 private:
243 #undef XPETRA_ITERATOROPS_SHORT
244 #include "Xpetra_UseShortNames.hpp"
245 
246 public:
247 
249  Jacobi(Scalar omega,
250  const Vector& Dinv,
251  const Matrix& A,
252  const Matrix& B,
253  RCP<Matrix> C_in,
255  const std::string & label) {
256  // Sanity checks
257  if (!A.isFillComplete())
258  throw Exceptions::RuntimeError("A is not fill-completed");
259  if (!B.isFillComplete())
260  throw Exceptions::RuntimeError("B is not fill-completed");
261 
262  // Default case: Xpetra Jacobi
263  RCP<Matrix> C = C_in;
264  if (C == Teuchos::null)
266 
267  Xpetra::Jacobi(omega, Dinv, A, B, *C, true,true,label);
268  C->CreateView("stridedMaps", rcpFromRef(A),false, rcpFromRef(B), false);
269  return C;
270  } //Jacobi
271 
272 
273 };
274 
275 
276 } // end namespace Xpetra
277 
278 #define XPETRA_ITERATOROPS_SHORT
279 
280 #endif /* PACKAGES_XPETRA_SUP_UTILS_XPETRA_ITERATOROPS_HPP_ */
void Jacobi(double omega, const Xpetra::Vector< double, int, int > &Dinv, const Xpetra::Matrix< double, int, int > &A, const Xpetra::Matrix< double, int, int > &B, Xpetra::Matrix< double, int, int > &C, bool call_FillComplete_on_result, bool doOptimizeStorage, const std::string &label)
GlobalOrdinal GO
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Xpetra namespace
Exception throws to report errors in the internal logical of the program.
LocalOrdinal LO
Xpetra utility class containing iteration operators.
static RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraCrs(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
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.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
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())
virtual bool isFillComplete() const =0
Returns true if fillComplete() has been called and the matrix is in compute mode. ...
void CreateView(viewLabel_t viewLabel, const RCP< const Map > &rowMap, const RCP< const Map > &colMap)
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Jacobi(Scalar omega, const Vector &Dinv, const Matrix &A, const Matrix &B, RCP< Matrix > C_in, Teuchos::FancyOStream &fos, const std::string &label)
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
void JacobiT(double omega, const Xpetra::Vector< double, int, GlobalOrdinal > &Dinv, const Xpetra::Matrix< double, int, GlobalOrdinal > &A, const Xpetra::Matrix< double, int, GlobalOrdinal > &B, Xpetra::Matrix< double, int, GlobalOrdinal > &C, bool call_FillComplete_on_result, bool doOptimizeStorage, const std::string &label)
static RCP< Epetra_CrsMatrix > Op2NonConstEpetraCrs(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
virtual Teuchos::RCP< const Map > getRangeMap() const =0
The Map associated with the range of this operator, which must be compatible with Y...
static RCP< Matrix > Build(const RCP< const Map > &rowMap, size_t maxNumEntriesPerRow, Xpetra::ProfileType pftype=Xpetra::DynamicProfile)
Constructor specifying the number of non-zeros for all rows.
static RCP< const Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2TpetraCrs(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
virtual const RCP< const Map > & getRowMap() const
Returns the Map that describes the row distribution in this matrix.
virtual Teuchos::RCP< const Map > getDomainMap() const =0
The Map associated with the domain of this operator, which must be compatible with X...
Xpetra-specific matrix class.