MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_Utilities_kokkos_decl.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 MUELU_UTILITIES_KOKKOS_DECL_HPP
47 #define MUELU_UTILITIES_KOKKOS_DECL_HPP
48 
49 #include "MueLu_ConfigDefs.hpp"
50 #if defined(HAVE_MUELU_KOKKOS_REFACTOR)
51 
52 #include <string>
53 
54 #include "Teuchos_ScalarTraits.hpp"
56 
57 #include "Xpetra_BlockedCrsMatrix_fwd.hpp"
58 #include "Xpetra_CrsMatrix_fwd.hpp"
59 #include "Xpetra_CrsMatrixWrap_fwd.hpp"
60 #include "Xpetra_ExportFactory.hpp"
61 #include "Xpetra_ImportFactory_fwd.hpp"
62 #include "Xpetra_MapFactory_fwd.hpp"
63 #include "Xpetra_Map_fwd.hpp"
64 #include "Xpetra_MatrixFactory_fwd.hpp"
65 #include "Xpetra_Matrix_fwd.hpp"
66 #include "Xpetra_MultiVectorFactory_fwd.hpp"
67 #include "Xpetra_MultiVector_fwd.hpp"
68 #include "Xpetra_Operator_fwd.hpp"
69 #include "Xpetra_VectorFactory_fwd.hpp"
70 #include "Xpetra_Vector_fwd.hpp"
71 
72 #include "Xpetra_IO.hpp"
73 
74 #include "Kokkos_ArithTraits.hpp"
75 
76 #ifdef HAVE_MUELU_EPETRA
77 #include "Epetra_MultiVector.h"
78 #include "Epetra_CrsMatrix.h"
79 #include "Xpetra_EpetraCrsMatrix_fwd.hpp"
80 #include "Xpetra_EpetraMultiVector_fwd.hpp"
81 #endif
82 
83 #include "MueLu_Exceptions.hpp"
84 #include "MueLu_Utilities.hpp"
85 #include "MueLu_UtilitiesBase.hpp"
86 
87 #ifdef HAVE_MUELU_TPETRA
88 #include "Tpetra_CrsMatrix.hpp"
89 #include "Tpetra_Map.hpp"
90 #include "Tpetra_MultiVector.hpp"
91 #include "Xpetra_TpetraCrsMatrix_fwd.hpp"
92 #include "Xpetra_TpetraMultiVector_fwd.hpp"
93 #endif
94 
95 
96 namespace MueLu {
97 
105  template <class Scalar,
108  class Node = DefaultNode>
109  class Utilities_kokkos : public MueLu::UtilitiesBase<Scalar,LocalOrdinal,GlobalOrdinal,Node> {
110 #undef MUELU_UTILITIES_KOKKOS_SHORT
111 #include "MueLu_UseShortNames.hpp"
112 
113  public:
114  using TST = Teuchos::ScalarTraits<SC>;
115  using Magnitude = typename TST::magnitudeType;
116  using RealValuedMultiVector = Xpetra::MultiVector<Magnitude,LO,GO,NO>;
117 
118 #ifdef HAVE_MUELU_EPETRA
119  // @{
121  static RCP<const Epetra_MultiVector> MV2EpetraMV(RCP<MultiVector> const vec) { return Utilities::MV2EpetraMV(vec); }
122  static RCP< Epetra_MultiVector> MV2NonConstEpetraMV(RCP<MultiVector> vec) { return Utilities::MV2NonConstEpetraMV(vec); }
123 
124  static const Epetra_MultiVector& MV2EpetraMV(const MultiVector& vec) { return Utilities::MV2EpetraMV(vec); }
125  static Epetra_MultiVector& MV2NonConstEpetraMV(MultiVector& vec) { return Utilities::MV2NonConstEpetraMV(vec); }
126 
127  static RCP<const Epetra_CrsMatrix> Op2EpetraCrs(RCP<const Matrix> Op) { return Utilities::Op2EpetraCrs(Op); }
128  static RCP< Epetra_CrsMatrix> Op2NonConstEpetraCrs(RCP<Matrix> Op) { return Utilities::Op2NonConstEpetraCrs(Op); }
129 
130  static const Epetra_CrsMatrix& Op2EpetraCrs(const Matrix& Op) { return Utilities::Op2EpetraCrs(Op); }
131  static Epetra_CrsMatrix& Op2NonConstEpetraCrs(Matrix& Op) { return Utilities::Op2NonConstEpetraCrs(Op); }
132 
133  static const Epetra_Map& Map2EpetraMap(const Map& map) { return Utilities::Map2EpetraMap(map); }
134  // @}
135 #endif
136 
137 #ifdef HAVE_MUELU_TPETRA
138  // @{
140  static RCP<const Tpetra::MultiVector<SC,LO,GO,NO> > MV2TpetraMV(RCP<MultiVector> const vec) { return Utilities::MV2TpetraMV(vec); }
141  static RCP< Tpetra::MultiVector<SC,LO,GO,NO> > MV2NonConstTpetraMV(RCP<MultiVector> vec) { return Utilities::MV2NonConstTpetraMV(vec); }
142  static RCP< Tpetra::MultiVector<SC,LO,GO,NO> > MV2NonConstTpetraMV2(MultiVector& vec) { return Utilities::MV2NonConstTpetraMV2(vec); }
143 
144  static const Tpetra::MultiVector<SC,LO,GO,NO>& MV2TpetraMV(const MultiVector& vec) { return Utilities::MV2TpetraMV(vec); }
145  static Tpetra::MultiVector<SC,LO,GO,NO>& MV2NonConstTpetraMV(MultiVector& vec) { return Utilities::MV2NonConstTpetraMV(vec); }
146 
147  static RCP<const Tpetra::CrsMatrix<SC,LO,GO,NO> > Op2TpetraCrs(RCP<const Matrix> Op) { return Utilities::Op2TpetraCrs(Op); }
148  static RCP< Tpetra::CrsMatrix<SC,LO,GO,NO> > Op2NonConstTpetraCrs(RCP<Matrix> Op) { return Utilities::Op2NonConstTpetraCrs(Op); }
149 
150  static const Tpetra::CrsMatrix<SC,LO,GO,NO>& Op2TpetraCrs(const Matrix& Op) { return Utilities::Op2TpetraCrs(Op); }
151  static Tpetra::CrsMatrix<SC,LO,GO,NO>& Op2NonConstTpetraCrs(Matrix& Op) { return Utilities::Op2NonConstTpetraCrs(Op); }
152 
153  static RCP<const Tpetra::RowMatrix<SC,LO,GO,NO> > Op2TpetraRow(RCP<const Matrix> Op) { return Utilities::Op2TpetraRow(Op); }
154  static RCP< Tpetra::RowMatrix<SC,LO,GO,NO> > Op2NonConstTpetraRow(RCP<Matrix> Op) { return Utilities::Op2NonConstTpetraRow(Op); }
155 
156  static const RCP<const Tpetra::Map<LO, GO, NO> > Map2TpetraMap(const Map& map) { return Utilities::Map2TpetraMap(map); }
157 #endif
158 
159  static RCP<Xpetra::Matrix<SC,LO,GO,NO> > Crs2Op(RCP<CrsMatrix> Op) { return Utilities::Crs2Op(Op); }
160 
167  static Teuchos::ArrayRCP<SC> GetMatrixDiagonal(const Matrix& A); // FIXME
168 
175  static RCP<Vector> GetMatrixDiagonalInverse(const Matrix& A, Magnitude tol = TST::eps()*100); // FIXME
176 
177 
178 
186  static RCP<Vector> GetMatrixOverlappedDiagonal(const Matrix& A); // FIXME
187 
196 
197  // TODO: should NOT return an Array. Definition must be changed to:
198  // - ArrayRCP<> ResidualNorm(Matrix const &Op, MultiVector const &X, MultiVector const &RHS)
199  // or
200  // - void ResidualNorm(Matrix const &Op, MultiVector const &X, MultiVector const &RHS, Array &)
201  static Teuchos::Array<Magnitude> ResidualNorm(const Operator& Op, const MultiVector& X, const MultiVector& RHS) {
202  return Utilities::ResidualNorm(Op, X, RHS);
203  }
204 
205  static RCP<MultiVector> Residual(const Operator& Op, const MultiVector& X, const MultiVector& RHS) {
206  return Utilities::Residual(Op, X, RHS);
207  }
208 
209  static void PauseForDebugger();
210 
226  static SC PowerMethod(const Matrix& A, bool scaleByDiag = true,
227  LO niters = 10, Magnitude tolerance = 1e-2, bool verbose = false, unsigned int seed = 123) {
228  return Utilities::PowerMethod(A, scaleByDiag, niters, tolerance, verbose, seed);
229  }
230 
231  static void MyOldScaleMatrix(Matrix& Op, const Teuchos::ArrayRCP<const SC>& scalingVector, bool doInverse = true,
232  bool doFillComplete = true, bool doOptimizeStorage = true); // FIXME
233 
234  static void MyOldScaleMatrix_Tpetra(Matrix& Op, const Teuchos::ArrayRCP<SC>& scalingVector,
235  bool doFillComplete, bool doOptimizeStorage); // FIXME
236 
237  static void MyOldScaleMatrix_Epetra(Matrix& Op, const Teuchos::ArrayRCP<SC>& scalingVector,
238  bool doFillComplete, bool doOptimizeStorage); // FIXME
239 
240  static RCP<Teuchos::FancyOStream> MakeFancy(std::ostream& os) { return Utilities::MakeFancy(os); }
241 
249  static Kokkos::View<bool*, typename NO::device_type> DetectDirichletRows(const Matrix& A, const Magnitude& tol = Teuchos::ScalarTraits<SC>::zero(), const bool count_twos_as_dirichlet=false);
250 
262 
263 
264  static void ZeroDirichletRows(RCP<Matrix>& A, const Kokkos::View<const bool*, typename NO::device_type>& dirichletRows, SC replaceWith=Teuchos::ScalarTraits<SC>::zero());
265 
266  static void ZeroDirichletRows(RCP<MultiVector>& X, const Kokkos::View<const bool*, typename NO::device_type>& dirichletRows, SC replaceWith=Teuchos::ScalarTraits<SC>::zero());
267 
268  static void ZeroDirichletCols(RCP<Matrix>& A, const Kokkos::View<const bool*, typename NO::device_type>& dirichletCols, SC replaceWith=Teuchos::ScalarTraits<SC>::zero());
269 
270  static RCP<MultiVector> RealValuedToScalarMultiVector(RCP<RealValuedMultiVector> X);
271 
281  static void SetRandomSeed(const Teuchos::Comm<int> &comm) { return Utilities::SetRandomSeed(comm); }
282 
288  static RCP<Matrix> Transpose(Matrix& Op, bool optimizeTranspose = false, const std::string & label = std::string()) {
289  return Utilities::Transpose(Op, optimizeTranspose, label);
290  }
291 
292  static RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LocalOrdinal,GlobalOrdinal,Node> > ExtractCoordinatesFromParameterList(ParameterList& paramList) {
294  }
295 
296 
300  static RCP<Xpetra::Vector<LocalOrdinal,LocalOrdinal,GlobalOrdinal,Node> > ReverseCuthillMcKee(const Matrix &Op);
301 
305  static RCP<Xpetra::Vector<LocalOrdinal,LocalOrdinal,GlobalOrdinal,Node> > CuthillMcKee(const Matrix &Op);
306 
307  static void ApplyOAZToMatrixRows(RCP<Matrix>& A, const Kokkos::View<const bool*, typename Node::device_type>& dirichletRows);
308 
309  }; // class Utils
310 
311 
321  template <class Node>
322  class Utilities_kokkos<double,int,int,Node> : public UtilitiesBase<double,int,int,Node> {
323  public:
324  typedef double Scalar;
325  typedef int LocalOrdinal;
326  typedef int GlobalOrdinal;
327  typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType Magnitude;
328  typedef Xpetra::MultiVector<Magnitude,LocalOrdinal,GlobalOrdinal,Node> RealValuedMultiVector;
329 
330  private:
331 #undef MUELU_UTILITIES_KOKKOS_SHORT
332 #include "MueLu_UseShortNames.hpp"
333 
334  public:
335 
336 #ifdef HAVE_MUELU_EPETRA
337  // @{
339  static RCP<const Epetra_MultiVector> MV2EpetraMV(RCP<MultiVector> const vec) { return Utilities::MV2EpetraMV(vec); }
340  static RCP< Epetra_MultiVector> MV2NonConstEpetraMV(RCP<MultiVector> vec) { return Utilities::MV2NonConstEpetraMV(vec); }
341 
342  static const Epetra_MultiVector& MV2EpetraMV(const MultiVector& vec) { return Utilities::MV2EpetraMV(vec); }
343  static Epetra_MultiVector& MV2NonConstEpetraMV(MultiVector& vec) { return Utilities::MV2NonConstEpetraMV(vec); }
344 
345  static RCP<const Epetra_CrsMatrix> Op2EpetraCrs(RCP<const Matrix> Op) { return Utilities::Op2EpetraCrs(Op); }
346  static RCP< Epetra_CrsMatrix> Op2NonConstEpetraCrs(RCP<Matrix> Op) { return Utilities::Op2NonConstEpetraCrs(Op); }
347 
348  static const Epetra_CrsMatrix& Op2EpetraCrs(const Matrix& Op) { return Utilities::Op2EpetraCrs(Op); }
349  static Epetra_CrsMatrix& Op2NonConstEpetraCrs(Matrix& Op) { return Utilities::Op2NonConstEpetraCrs(Op); }
350 
351  static const Epetra_Map& Map2EpetraMap(const Map& map) { return Utilities::Map2EpetraMap(map); }
352  // @}
353 #endif
354 
355 #ifdef HAVE_MUELU_TPETRA
356  // @{
358  static RCP<const Tpetra::MultiVector<SC,LO,GO,NO> > MV2TpetraMV(RCP<MultiVector> const vec) { return Utilities::MV2TpetraMV(vec); }
359  static RCP< Tpetra::MultiVector<SC,LO,GO,NO> > MV2NonConstTpetraMV(RCP<MultiVector> vec) { return Utilities::MV2NonConstTpetraMV(vec); }
360  static RCP< Tpetra::MultiVector<SC,LO,GO,NO> > MV2NonConstTpetraMV2(MultiVector& vec) { return Utilities::MV2NonConstTpetraMV2(vec); }
361 
362  static const Tpetra::MultiVector<SC,LO,GO,NO>& MV2TpetraMV(const MultiVector& vec) { return Utilities::MV2TpetraMV(vec); }
363  static Tpetra::MultiVector<SC,LO,GO,NO>& MV2NonConstTpetraMV(MultiVector& vec) { return Utilities::MV2NonConstTpetraMV(vec); }
364 
365  static RCP<const Tpetra::CrsMatrix<SC,LO,GO,NO> > Op2TpetraCrs(RCP<const Matrix> Op) { return Utilities::Op2TpetraCrs(Op); }
366  static RCP< Tpetra::CrsMatrix<SC,LO,GO,NO> > Op2NonConstTpetraCrs(RCP<Matrix> Op) { return Utilities::Op2NonConstTpetraCrs(Op); }
367 
368  static const Tpetra::CrsMatrix<SC,LO,GO,NO>& Op2TpetraCrs(const Matrix& Op) { return Utilities::Op2TpetraCrs(Op); }
369  static Tpetra::CrsMatrix<SC,LO,GO,NO>& Op2NonConstTpetraCrs(Matrix& Op) { return Utilities::Op2NonConstTpetraCrs(Op); }
370 
371  static RCP<const Tpetra::RowMatrix<SC,LO,GO,NO> > Op2TpetraRow(RCP<const Matrix> Op) { return Utilities::Op2TpetraRow(Op); }
372  static RCP< Tpetra::RowMatrix<SC,LO,GO,NO> > Op2NonConstTpetraRow(RCP<Matrix> Op) { return Utilities::Op2NonConstTpetraRow(Op); }
373 
374  static const RCP<const Tpetra::Map<LO, GO, NO> > Map2TpetraMap(const Map& map) { return Utilities::Map2TpetraMap(map); }
375 #endif
376  static RCP<Xpetra::Matrix<SC,LO,GO,NO> > Crs2Op(RCP<CrsMatrix> Op) { return Utilities::Crs2Op(Op); }
377 
378  static ArrayRCP<SC> GetMatrixDiagonal(const Matrix& A) {
380  }
381  static RCP<Vector> GetMatrixDiagonalInverse(const Matrix& A, Magnitude tol = Teuchos::ScalarTraits<SC>::eps()*100) {
383  }
384  static ArrayRCP<SC> GetLumpedMatrixDiagonal(const Matrix& A) {
386  }
387  static RCP<Vector> GetLumpedMatrixDiagonal(RCP<const Matrix > A) {
389  }
390  static RCP<Vector> GetMatrixOverlappedDiagonal(const Matrix& A) {
392  }
393  static RCP<Vector> GetInverse(RCP<const Vector> v, Magnitude tol = Teuchos::ScalarTraits<SC>::eps()*100, SC tolReplacement = Teuchos::ScalarTraits<SC>::zero()) {
394  return UtilitiesBase::GetInverse(v,tol,tolReplacement);
395  }
396  static Teuchos::Array<Magnitude> ResidualNorm(const Operator& Op, const MultiVector& X, const MultiVector& RHS) {
397  return UtilitiesBase::ResidualNorm(Op,X,RHS);
398  }
399  static RCP<MultiVector> Residual(const Operator& Op, const MultiVector& X, const MultiVector& RHS) {
400  return UtilitiesBase::Residual(Op,X,RHS);
401  }
402  static void PauseForDebugger() {
404  }
405  static RCP<Teuchos::FancyOStream> MakeFancy(std::ostream& os) {
406  return UtilitiesBase::MakeFancy(os);
407  }
408  static void SetRandomSeed(const Teuchos::Comm<int> &comm) {
410  }
411 
412  // todo: move this to UtilitiesBase::kokkos
413  static Kokkos::View<bool*, typename Node::device_type> DetectDirichletRows(const Matrix& A, const Magnitude& tol = Teuchos::ScalarTraits<SC>::zero(), const bool count_twos_as_dirichlet=false);
414 
416 
417  static void ZeroDirichletRows(RCP<Matrix>& A, const Kokkos::View<const bool*, typename Node::device_type>& dirichletRows, SC replaceWith=Teuchos::ScalarTraits<SC>::zero());
418 
419  static void ZeroDirichletRows(RCP<MultiVector>& X, const Kokkos::View<const bool*, typename Node::device_type>& dirichletRows, SC replaceWith=Teuchos::ScalarTraits<SC>::zero());
420 
421  static void ZeroDirichletCols(RCP<Matrix>& A, const Kokkos::View<const bool*, typename Node::device_type>& dirichletCols, SC replaceWith=Teuchos::ScalarTraits<SC>::zero());
422 
423  static RCP<MultiVector> RealValuedToScalarMultiVector(RCP<RealValuedMultiVector> X);
424 
425  static Scalar PowerMethod(const Matrix& A, bool scaleByDiag = true, LO niters = 10, Magnitude tolerance = 1e-2, bool verbose = false, unsigned int seed = 123) {
426  return UtilitiesBase::PowerMethod(A,scaleByDiag,niters,tolerance,verbose,seed);
427  }
428 
429  static void MyOldScaleMatrix(Matrix& Op, const Teuchos::ArrayRCP<const SC>& scalingVector, bool doInverse = true,
430  bool doFillComplete = true, bool doOptimizeStorage = true) {
432  Teuchos::ArrayRCP<SC> sv(scalingVector.size());
433  if (doInverse) {
434  for (int i = 0; i < scalingVector.size(); ++i)
435  sv[i] = one / scalingVector[i];
436  } else {
437  for (int i = 0; i < scalingVector.size(); ++i)
438  sv[i] = scalingVector[i];
439  }
440 
441  switch (Op.getRowMap()->lib()) {
442  case Xpetra::UseTpetra:
443  MyOldScaleMatrix_Tpetra(Op, sv, doFillComplete, doOptimizeStorage);
444  break;
445 
446  case Xpetra::UseEpetra:
447  MyOldScaleMatrix_Epetra(Op, sv, doFillComplete, doOptimizeStorage);
448  break;
449 
450  default:
451  throw Exceptions::RuntimeError("Only Epetra and Tpetra matrices can be scaled.");
452 #ifndef __NVCC__ //prevent nvcc warning
453  break;
454 #endif
455  }
456  }
457 
458  // TODO This is the <double,int,int> specialization
459  static void MyOldScaleMatrix_Tpetra(Matrix& Op, const Teuchos::ArrayRCP<Scalar>& scalingVector,
460  bool doFillComplete, bool doOptimizeStorage) {
461  #ifdef HAVE_MUELU_TPETRA
462  #ifdef HAVE_MUELU_TPETRA_INST_INT_INT
463  try {
464  Tpetra::CrsMatrix<SC,LO,GO,NO>& tpOp = Op2NonConstTpetraCrs(Op);
465 
466  const RCP<const Tpetra::Map<LO,GO,NO> > rowMap = tpOp.getRowMap();
467  const RCP<const Tpetra::Map<LO,GO,NO> > domainMap = tpOp.getDomainMap();
468  const RCP<const Tpetra::Map<LO,GO,NO> > rangeMap = tpOp.getRangeMap();
469 
470  size_t maxRowSize = tpOp.getNodeMaxNumRowEntries();
471  if (maxRowSize == Teuchos::as<size_t>(-1)) // hasn't been determined yet
472  maxRowSize = 20;
473 
474  std::vector<SC> scaledVals(maxRowSize);
475  if (tpOp.isFillComplete())
476  tpOp.resumeFill();
477 
478  if (Op.isLocallyIndexed() == true) {
481 
482  for (size_t i = 0; i < rowMap->getNodeNumElements(); ++i) {
483  tpOp.getLocalRowView(i, cols, vals);
484  size_t nnz = tpOp.getNumEntriesInLocalRow(i);
485  if (nnz > maxRowSize) {
486  maxRowSize = nnz;
487  scaledVals.resize(maxRowSize);
488  }
489  for (size_t j = 0; j < nnz; ++j)
490  scaledVals[j] = vals[j]*scalingVector[i];
491 
492  if (nnz > 0) {
493  Teuchos::ArrayView<const SC> valview(&scaledVals[0], nnz);
494  tpOp.replaceLocalValues(i, cols, valview);
495  }
496  } //for (size_t i=0; ...
497 
498  } else {
501 
502  for (size_t i = 0; i < rowMap->getNodeNumElements(); ++i) {
503  GO gid = rowMap->getGlobalElement(i);
504  tpOp.getGlobalRowView(gid, cols, vals);
505  size_t nnz = tpOp.getNumEntriesInGlobalRow(gid);
506  if (nnz > maxRowSize) {
507  maxRowSize = nnz;
508  scaledVals.resize(maxRowSize);
509  }
510  // FIXME FIXME FIXME FIXME FIXME FIXME
511  for (size_t j = 0; j < nnz; ++j)
512  scaledVals[j] = vals[j]*scalingVector[i]; //FIXME i or gid?
513 
514  if (nnz > 0) {
515  Teuchos::ArrayView<const SC> valview(&scaledVals[0], nnz);
516  tpOp.replaceGlobalValues(gid, cols, valview);
517  }
518  } //for (size_t i=0; ...
519  }
520 
521  if (doFillComplete) {
522  if (domainMap == Teuchos::null || rangeMap == Teuchos::null)
523  throw Exceptions::RuntimeError("In Utilities::Scaling: cannot fillComplete because the domain and/or range map hasn't been defined");
524 
525  RCP<Teuchos::ParameterList> params = rcp(new Teuchos::ParameterList());
526  params->set("Optimize Storage", doOptimizeStorage);
527  params->set("No Nonlocal Changes", true);
528  Op.fillComplete(Op.getDomainMap(), Op.getRangeMap(), params);
529  }
530  } catch(...) {
531  throw Exceptions::RuntimeError("Only Tpetra::CrsMatrix types can be scaled (Err.1)");
532  }
533 #else
534  throw Exceptions::RuntimeError("Matrix scaling is not possible because Tpetra has not been compiled with support for LO=GO=int.");
535 #endif
536 #else
537  throw Exceptions::RuntimeError("Matrix scaling is not possible because Tpetra has not been enabled.");
538 #endif
539  }
540 
541  static void MyOldScaleMatrix_Epetra(Matrix& Op, const Teuchos::ArrayRCP<Scalar>& scalingVector, bool doFillComplete, bool doOptimizeStorage) {
542 #ifdef HAVE_MUELU_EPETRA
543  try {
544  //const Epetra_CrsMatrix& epOp = Utilities<double,int,int>::Op2NonConstEpetraCrs(Op);
545  const Epetra_CrsMatrix& epOp = Op2NonConstEpetraCrs(Op);
546 
547  Epetra_Map const &rowMap = epOp.RowMap();
548  int nnz;
549  double *vals;
550  int *cols;
551 
552  for (int i = 0; i < rowMap.NumMyElements(); ++i) {
553  epOp.ExtractMyRowView(i, nnz, vals, cols);
554  for (int j = 0; j < nnz; ++j)
555  vals[j] *= scalingVector[i];
556  }
557 
558  } catch (...){
559  throw Exceptions::RuntimeError("Only Epetra_CrsMatrix types can be scaled");
560  }
561 #else
562  throw Exceptions::RuntimeError("Matrix scaling is not possible because Epetra has not been enabled.");
563 #endif // HAVE_MUELU_EPETRA
564  }
565 
571  static RCP<Matrix> Transpose(Matrix& Op, bool /* optimizeTranspose */ = false, const std::string & label = std::string(),const Teuchos::RCP<Teuchos::ParameterList> &params=Teuchos::null) {
572  switch (Op.getRowMap()->lib()) {
573  case Xpetra::UseTpetra:
574  {
575 #ifdef HAVE_MUELU_TPETRA
576 #ifdef HAVE_MUELU_TPETRA_INST_INT_INT
577  try {
578  const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& tpetraOp = Utilities::Op2TpetraCrs(Op);
579 
580  // Compute the transpose A of the Tpetra matrix tpetraOp.
581  RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > A;
582  Tpetra::RowMatrixTransposer<SC,LO,GO,NO> transposer(rcpFromRef(tpetraOp),label);
583  {
585  using Teuchos::rcp;
586  RCP<ParameterList> transposeParams = params.is_null () ?
587  rcp (new ParameterList) :
588  rcp (new ParameterList (*params));
589  transposeParams->set ("sort", false);
590  A = transposer.createTranspose(transposeParams);
591  }
592 
593  RCP<Xpetra::TpetraCrsMatrix<SC,LO,GO,NO> > AA = rcp(new Xpetra::TpetraCrsMatrix<SC,LO,GO,NO>(A));
594  RCP<CrsMatrix> AAA = rcp_implicit_cast<CrsMatrix>(AA);
595  RCP<CrsMatrixWrap> AAAA = rcp(new CrsMatrixWrap(AAA));
596 
597  return AAAA;
598  }
599  catch (std::exception& e) {
600  std::cout << "threw exception '" << e.what() << "'" << std::endl;
601  throw Exceptions::RuntimeError("Utilities::Transpose failed, perhaps because matrix is not a Crs matrix");
602  }
603 #else
604  throw Exceptions::RuntimeError("Utilities::Transpose: Tpetra is not compiled with LO=GO=int. Add TPETRA_INST_INT_INT:BOOL=ON to your configuration!");
605 #endif
606 #else
607  throw Exceptions::RuntimeError("Utilities::Transpose: Tpetra is not compiled!");
608 #endif
609 #ifndef __NVCC__ //prevent nvcc warning
610  break;
611 #endif
612  }
613  case Xpetra::UseEpetra:
614  {
615 #if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_EPETRAEXT)
616  Teuchos::TimeMonitor tm(*Teuchos::TimeMonitor::getNewTimer("ZZ Entire Transpose"));
617  // Epetra case
620  Epetra_CrsMatrix * A = dynamic_cast<Epetra_CrsMatrix*>(&transposer(epetraOp));
621  transposer.ReleaseTranspose(); // So we can keep A in Muelu...
622 
623  RCP<Epetra_CrsMatrix> rcpA(A);
624  RCP<Xpetra::EpetraCrsMatrixT<GO,NO> > AA = rcp(new Xpetra::EpetraCrsMatrixT<GO,NO> (rcpA));
625  RCP<CrsMatrix> AAA = rcp_implicit_cast<CrsMatrix>(AA);
626  RCP<CrsMatrixWrap> AAAA = rcp( new CrsMatrixWrap(AAA));
627  AAAA->fillComplete(Op.getRangeMap(), Op.getDomainMap());
628 
629  return AAAA;
630 #else
631  throw Exceptions::RuntimeError("Epetra (Err. 2)");
632 #endif
633 #ifndef __NVCC__ //prevent nvcc warning
634  break;
635 #endif
636  }
637  default:
638  throw Exceptions::RuntimeError("Only Epetra and Tpetra matrices can be transposed.");
639 #ifndef __NVCC__ //prevent nvcc warning
640  break;
641 #endif
642  }
643 
644 #ifndef __NVCC__ //prevent nvcc warning
645  return Teuchos::null;
646 #endif
647  }
648 
651  static RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO> > ExtractCoordinatesFromParameterList(ParameterList& paramList) {
652  RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO> > coordinates = Teuchos::null;
653 
654  // check whether coordinates are contained in parameter list
655  if(paramList.isParameter ("Coordinates") == false)
656  return coordinates;
657 
658 #if defined(HAVE_MUELU_TPETRA)
659 #if ( defined(EPETRA_HAVE_OMP) && defined(HAVE_TPETRA_INST_OPENMP) && defined(HAVE_TPETRA_INST_INT_INT)) || \
660  (!defined(EPETRA_HAVE_OMP) && defined(HAVE_TPETRA_INST_SERIAL) && defined(HAVE_TPETRA_INST_INT_INT))
661 
662  // define Tpetra::MultiVector type with Scalar=float only if
663  // * ETI is turned off, since then the compiler will instantiate it automatically OR
664  // * Tpetra is instantiated on Scalar=float
665 #if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT)
666  typedef Tpetra::MultiVector<float, LO,GO,NO> tfMV;
667  RCP<tfMV> floatCoords = Teuchos::null;
668 #endif
669 
670  // define Tpetra::MultiVector type with Scalar=double only if
671  // * ETI is turned off, since then the compiler will instantiate it automatically OR
672  // * Tpetra is instantiated on Scalar=double
673  typedef Tpetra::MultiVector<SC,LO,GO,NO> tdMV;
674  RCP<tdMV> doubleCoords = Teuchos::null;
675  if (paramList.isType<RCP<tdMV> >("Coordinates")) {
676  // Coordinates are stored as a double vector
677  doubleCoords = paramList.get<RCP<tdMV> >("Coordinates");
678  paramList.remove("Coordinates");
679  }
680 #if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT)
681  else if (paramList.isType<RCP<tfMV> >("Coordinates")) {
682  // check if coordinates are stored as a float vector
683  floatCoords = paramList.get<RCP<tfMV> >("Coordinates");
684  paramList.remove("Coordinates");
685  doubleCoords = rcp(new tdMV(floatCoords->getMap(), floatCoords->getNumVectors()));
686  deep_copy(*doubleCoords, *floatCoords);
687  }
688 #endif
689  // We have the coordinates in a Tpetra double vector
690  if(doubleCoords != Teuchos::null) {
691  coordinates = Teuchos::rcp(new Xpetra::TpetraMultiVector<SC,LO,GO,NO>(doubleCoords));
692  TEUCHOS_TEST_FOR_EXCEPT(doubleCoords->getNumVectors() != coordinates->getNumVectors());
693  }
694 #endif // Tpetra instantiated on GO=int and EpetraNode
695 #endif // endif HAVE_TPETRA
696 
697 #if defined(HAVE_MUELU_EPETRA)
698  RCP<Epetra_MultiVector> doubleEpCoords;
699  if (paramList.isType<RCP<Epetra_MultiVector> >("Coordinates")) {
700  doubleEpCoords = paramList.get<RCP<Epetra_MultiVector> >("Coordinates");
701  paramList.remove("Coordinates");
702  RCP<Xpetra::EpetraMultiVectorT<GO,NO> > epCoordinates = Teuchos::rcp(new Xpetra::EpetraMultiVectorT<GO,NO>(doubleEpCoords));
703  coordinates = rcp_dynamic_cast<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType,LO,GO,NO> >(epCoordinates);
704  TEUCHOS_TEST_FOR_EXCEPT(doubleEpCoords->NumVectors() != Teuchos::as<int>(coordinates->getNumVectors()));
705  }
706 #endif
708  return coordinates;
709  }
710 
714  static RCP<Xpetra::Vector<LocalOrdinal,LocalOrdinal,GlobalOrdinal,Node> > ReverseCuthillMcKee(const Matrix &Op);
715 
719  static RCP<Xpetra::Vector<LocalOrdinal,LocalOrdinal,GlobalOrdinal,Node> > CuthillMcKee(const Matrix &Op);
720 
721  static void ApplyOAZToMatrixRows(RCP<Matrix>& A, const Kokkos::View<const bool*, typename Node::device_type>& dirichletRows);
722 
723  }; // class Utilities (specialization SC=double LO=GO=int)
724 
725 
726 
727  // Useful Kokkos conversions
728  template < class View, unsigned AppendValue >
729  struct AppendTrait {
730  // static_assert(false, "Error: NOT a Kokkos::View");
731  };
732 
733  template < class MT, unsigned T >
734  struct CombineMemoryTraits {
735  // empty
736  };
737 
738  template < unsigned U, unsigned T>
739  struct CombineMemoryTraits<Kokkos::MemoryTraits<U>, T> {
740  typedef Kokkos::MemoryTraits<U|T> type;
741  };
742 
743  template < class DataType, unsigned T, class... Pack >
744  struct AppendTrait< Kokkos::View< DataType, Pack... >, T> {
745  typedef Kokkos::View< DataType, Pack... > view_type;
747  };
748 
749 } //namespace MueLu
750 
751 #define MUELU_UTILITIES_KOKKOS_SHORT
752 
753 #endif // #if defined(HAVE_MUELU_KOKKOS_REFACTOR)
754 
755 #endif // MUELU_UTILITIES_KOKKOS_DECL_HPP
static void SetRandomSeed(const Teuchos::Comm< int > &comm)
static const RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > Map2TpetraMap(const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &map)
MueLu::DefaultLocalOrdinal LocalOrdinal
KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
static RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2TpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > const vec)
Helper utility to pull out the underlying Tpetra objects from an Xpetra object.
static RCP< const Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2TpetraCrs(RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static const Epetra_Map & Map2EpetraMap(const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &map)
Teuchos::RCP< Xpetra::Vector< LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node > > ReverseCuthillMcKee(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op)
static RCP< Epetra_CrsMatrix > Op2NonConstEpetraCrs(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static Scalar PowerMethod(const Matrix &A, bool scaleByDiag=true, LocalOrdinal niters=10, Magnitude tolerance=1e-2, bool verbose=false, unsigned int seed=123)
Simple transpose for Tpetra::CrsMatrix types.
static Teuchos::RCP< Vector > GetInverse(Teuchos::RCP< const Vector > v, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps()*100, Scalar tolReplacement=Teuchos::ScalarTraits< Scalar >::zero())
Return vector containing inverse of input vector.
static RCP< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::magnitudeType, LocalOrdinal, GlobalOrdinal, Node > > ExtractCoordinatesFromParameterList(ParameterList &paramList)
bool is_null(const std::shared_ptr< T > &p)
static RCP< const Epetra_MultiVector > MV2EpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > const vec)
Helper utility to pull out the underlying Epetra objects from an Xpetra object.
void ZeroDirichletCols(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, const Kokkos::View< const bool *, typename Node::device_type > &dirichletCols, Scalar replaceWith)
void ZeroDirichletRows(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, const Kokkos::View< const bool *, typename Node::device_type > &dirichletRows, Scalar replaceWith)
Teuchos::RCP< Xpetra::Vector< LocalOrdinal, LocalOrdinal, GlobalOrdinal, Node > > CuthillMcKee(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op)
void ApplyOAZToMatrixRows(Teuchos::RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, const Kokkos::View< const bool *, typename Node::device_type > &dirichletRows)
size_type size() const
MueLu::DefaultNode Node
static RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2NonConstTpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec)
static RCP< Time > getNewTimer(const std::string &name)
int ExtractMyRowView(int MyRow, int &NumEntries, double *&Values, int *&Indices) const
const Epetra_Map & RowMap() const
int NumMyElements() const
static RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Residual(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &RHS)
Kokkos::View< bool *, typename NO::device_type > DetectDirichletRows(const Xpetra::Matrix< SC, LO, GO, NO > &A, const typename Teuchos::ScalarTraits< SC >::magnitudeType &tol, const bool count_twos_as_dirichlet)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
static RCP< Epetra_MultiVector > MV2NonConstEpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > vec)
static RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraCrs(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2NonConstTpetraMV2(Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec)
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Transpose(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, bool optimizeTranspose=false, const std::string &label=std::string(), const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
static RCP< const Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2TpetraRow(RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static Teuchos::ArrayRCP< Scalar > GetLumpedMatrixDiagonal(const Matrix &A)
Extract Matrix Diagonal of lumped matrix.
static Teuchos::Array< Magnitude > ResidualNorm(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Crs2Op(RCP< Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static RCP< const Epetra_CrsMatrix > Op2EpetraCrs(RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
static Teuchos::Array< Magnitude > ResidualNorm(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &RHS)
static RCP< Teuchos::FancyOStream > MakeFancy(std::ostream &os)
static RCP< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op2NonConstTpetraRow(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Op)
static RCP< Vector > GetMatrixOverlappedDiagonal(const Matrix &A)
Extract Overlapped Matrix Diagonal.
Kokkos::View< bool *, typename NO::device_type > DetectDirichletCols(const Xpetra::Matrix< SC, LO, GO, NO > &A, const Kokkos::View< const bool *, typename NO::device_type > &dirichletRows)
static RCP< Teuchos::FancyOStream > MakeFancy(std::ostream &os)
static Scalar PowerMethod(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, bool scaleByDiag=true, LocalOrdinal niters=10, Magnitude tolerance=1e-2, bool verbose=false, unsigned int seed=123)
static RCP< Vector > GetMatrixDiagonalInverse(const Matrix &A, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps()*100)
Extract Matrix Diagonal.
static void SetRandomSeed(const Teuchos::Comm< int > &comm)
Set seed for random number generator.
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
void deep_copy(const DynRankView< DT, DP...> &dst, typename ViewTraits< DT, DP...>::const_value_type &value, typename std::enable_if< std::is_same< typename ViewTraits< DT, DP...>::specialize, void >::value >::type *=nullptr)
static Teuchos::ArrayRCP< Scalar > GetMatrixDiagonal(const Matrix &A)
Extract Matrix Diagonal.