Xpetra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Xpetra_MatrixFactory2_decl.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 // WARNING: This code is experimental. Backwards compatibility should not be expected.
11 
12 #ifndef XPETRA_MATRIXFACTORY2_DECL_HPP
13 #define XPETRA_MATRIXFACTORY2_DECL_HPP
14 
15 #include "Xpetra_ConfigDefs.hpp"
17 #include "Xpetra_Matrix.hpp"
18 #include "Xpetra_CrsMatrixWrap.hpp"
20 #include "Xpetra_Map.hpp"
21 #include "Xpetra_BlockedMap.hpp"
22 #include "Xpetra_Vector.hpp"
23 #include "Xpetra_BlockedVector.hpp"
24 #include "Xpetra_Exceptions.hpp"
25 
26 namespace Xpetra {
27 
28 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
30 #undef XPETRA_MATRIXFACTORY2_SHORT
31 #include "Xpetra_UseShortNames.hpp"
32 
33  public:
34  static RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > BuildCopy(const RCP<const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > A, bool setFixedBlockSize = true);
35 };
36 #define XPETRA_MATRIXFACTORY2_SHORT
37 
38 // template<>
39 // class MatrixFactory2<double,int,int,typename Xpetra::Matrix<double, int, int>::node_type> {
40 template <class Node>
41 class MatrixFactory2<double, int, int, Node> {
42  typedef double Scalar;
43  typedef int LocalOrdinal;
44  typedef int GlobalOrdinal;
45  // typedef Matrix<double, int, GlobalOrdinal>::node_type Node;
46 #undef XPETRA_MATRIXFACTORY2_SHORT
47 #include "Xpetra_UseShortNames.hpp"
48  public:
49  static RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > BuildCopy(const RCP<const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > A, bool setFixedBlockSize = true) {
50  RCP<const CrsMatrixWrap> oldOp = Teuchos::rcp_dynamic_cast<const CrsMatrixWrap>(A);
51  if (oldOp == Teuchos::null)
52  throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
53 
54  RCP<const CrsMatrix> oldCrsOp = oldOp->getCrsMatrix();
55 
56 #ifdef HAVE_XPETRA_EPETRA
57 #ifndef XPETRA_EPETRA_NO_32BIT_GLOBAL_INDICES
58  RCP<const EpetraCrsMatrixT<GlobalOrdinal, Node> > oldECrsOp = Teuchos::rcp_dynamic_cast<const EpetraCrsMatrixT<GlobalOrdinal, Node> >(oldCrsOp);
59  if (oldECrsOp != Teuchos::null) {
60  // Underlying matrix is Epetra
61  RCP<CrsMatrix> newECrsOp(new EpetraCrsMatrixT<GlobalOrdinal, Node>(*oldECrsOp));
62  RCP<CrsMatrixWrap> newOp(new CrsMatrixWrap(newECrsOp));
63  if (setFixedBlockSize)
64  newOp->SetFixedBlockSize(A->GetFixedBlockSize());
65  return newOp;
66  }
67 #endif
68 #endif
69 
70 #ifdef HAVE_XPETRA_TPETRA
71  // Underlying matrix is Tpetra
72  RCP<const TpetraCrsMatrix> oldTCrsOp = Teuchos::rcp_dynamic_cast<const TpetraCrsMatrix>(oldCrsOp);
73  if (oldTCrsOp != Teuchos::null) {
74  RCP<CrsMatrix> newTCrsOp(new TpetraCrsMatrix(*oldTCrsOp));
75  RCP<CrsMatrixWrap> newOp(new CrsMatrixWrap(newTCrsOp));
76  if (setFixedBlockSize)
77  newOp->SetFixedBlockSize(A->GetFixedBlockSize());
78  return newOp;
79  }
80  return Teuchos::null;
81 #else
82  throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::EpetraCrsMatrix or Xpetra::TpetraCrsMatrix failed");
83  TEUCHOS_UNREACHABLE_RETURN(Teuchos::null); // make compiler happy
84 #endif
85 
86  } // BuildCopy
87 };
88 
89 #define XPETRA_MATRIXFACTORY2_SHORT
90 
91 #ifdef HAVE_XPETRA_INT_LONG_LONG
92 template <class Node>
93 class MatrixFactory2<double, int, long long, Node> {
94  typedef double Scalar;
95  typedef int LocalOrdinal;
96  typedef long long GlobalOrdinal;
97  // typedef Matrix<double, int, GlobalOrdinal>::node_type Node;
98 #undef XPETRA_MATRIXFACTORY2_SHORT
99 #include "Xpetra_UseShortNames.hpp"
100  public:
101  static RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > BuildCopy(const RCP<const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > A, bool setFixedBlockSize = true) {
102  RCP<const CrsMatrixWrap> oldOp = Teuchos::rcp_dynamic_cast<const CrsMatrixWrap>(A);
103  if (oldOp == Teuchos::null)
104  throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
105 
106  RCP<const CrsMatrix> oldCrsOp = oldOp->getCrsMatrix();
107 
108 #ifdef HAVE_XPETRA_EPETRA
109 #ifndef XPETRA_EPETRA_NO_64BIT_GLOBAL_INDICES
110  RCP<const EpetraCrsMatrixT<GlobalOrdinal, Node> > oldECrsOp = Teuchos::rcp_dynamic_cast<const EpetraCrsMatrixT<GlobalOrdinal, Node> >(oldCrsOp);
111  if (oldECrsOp != Teuchos::null) {
112  // Underlying matrix is Epetra
113  RCP<CrsMatrix> newECrsOp(new EpetraCrsMatrixT<GlobalOrdinal, Node>(*oldECrsOp));
114  RCP<CrsMatrixWrap> newOp(new CrsMatrixWrap(newECrsOp));
115  if (setFixedBlockSize)
116  newOp->SetFixedBlockSize(A->GetFixedBlockSize());
117  return newOp;
118  }
119 #endif
120 #endif
121 
122 #ifdef HAVE_XPETRA_TPETRA
123  // Underlying matrix is Tpetra
124  RCP<const TpetraCrsMatrix> oldTCrsOp = Teuchos::rcp_dynamic_cast<const TpetraCrsMatrix>(oldCrsOp);
125  if (oldTCrsOp != Teuchos::null) {
126  RCP<CrsMatrix> newTCrsOp(new TpetraCrsMatrix(*oldTCrsOp));
127  RCP<CrsMatrixWrap> newOp(new CrsMatrixWrap(newTCrsOp));
128  if (setFixedBlockSize)
129  newOp->SetFixedBlockSize(A->GetFixedBlockSize());
130  return newOp;
131  }
132 #else
133  throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::EpetraCrsMatrix or Xpetra::TpetraCrsMatrix failed");
134 #endif
135 
136  return Teuchos::null; // make compiler happy
137  }
138 };
139 #endif // HAVE_XPETRA_INT_LONG_LONG
140 
141 #define XPETRA_MATRIXFACTORY2_SHORT
142 
143 } // namespace Xpetra
144 
145 #define XPETRA_MATRIXFACTORY2_SHORT
146 
147 #endif // ifndef XPETRA_MATRIXFACTORY2_DECL_HPP
Exception indicating invalid cast attempted.
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > BuildCopy(const RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > A, bool setFixedBlockSize=true)
Concrete implementation of Xpetra::Matrix.
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > BuildCopy(const RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > A, bool setFixedBlockSize=true)
Xpetra-specific matrix class.