49 #ifndef XPETRA_MATRIXFACTORY_HPP
50 #define XPETRA_MATRIXFACTORY_HPP
69 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
70 class MatrixFactory2 {
71 #undef XPETRA_MATRIXFACTORY2_SHORT
77 if (oldOp == Teuchos::null)
85 "Not Epetra or Tpetra matrix");
87 #ifdef HAVE_XPETRA_EPETRA
94 #ifdef HAVE_XPETRA_TPETRA
99 if (oldTCrsOp != Teuchos::null) {
102 newOp->SetFixedBlockSize(A->GetFixedBlockSize());
111 return Teuchos::null;
114 #define XPETRA_MATRIXFACTORY2_SHORT
124 #undef XPETRA_MATRIXFACTORY2_SHORT
129 if (oldOp == Teuchos::null)
134 #ifdef HAVE_XPETRA_EPETRA
135 #ifndef XPETRA_EPETRA_NO_32BIT_GLOBAL_INDICES
137 if (oldECrsOp != Teuchos::null) {
141 newOp->SetFixedBlockSize(A->GetFixedBlockSize());
147 #ifdef HAVE_XPETRA_TPETRA
150 if (oldTCrsOp != Teuchos::null) {
153 newOp->SetFixedBlockSize(A->GetFixedBlockSize());
156 return Teuchos::null;
158 throw Exceptions::BadCast(
"Cast from Xpetra::Matrix to Xpetra::EpetraCrsMatrix or Xpetra::TpetraCrsMatrix failed");
165 #define XPETRA_MATRIXFACTORY2_SHORT
167 #ifdef HAVE_XPETRA_INT_LONG_LONG
171 class MatrixFactory2<double, int, long long, Node> {
172 typedef double Scalar;
173 typedef int LocalOrdinal;
174 typedef long long GlobalOrdinal;
176 #undef XPETRA_MATRIXFACTORY2_SHORT
180 RCP<const CrsMatrixWrap> oldOp = Teuchos::rcp_dynamic_cast<
const CrsMatrixWrap>(A);
181 if (oldOp == Teuchos::null)
182 throw Exceptions::BadCast(
"Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
184 RCP<const CrsMatrix> oldCrsOp = oldOp->getCrsMatrix();
186 #ifdef HAVE_XPETRA_EPETRA
187 #ifndef XPETRA_EPETRA_NO_64BIT_GLOBAL_INDICES
188 RCP<const EpetraCrsMatrixT<GlobalOrdinal,Node> > oldECrsOp = Teuchos::rcp_dynamic_cast<
const EpetraCrsMatrixT<GlobalOrdinal,Node> >(oldCrsOp);
189 if (oldECrsOp != Teuchos::null) {
191 RCP<CrsMatrix> newECrsOp(
new EpetraCrsMatrixT<GlobalOrdinal,Node>(*oldECrsOp));
192 RCP<CrsMatrixWrap> newOp (
new CrsMatrixWrap (newECrsOp));
193 newOp->SetFixedBlockSize(A->GetFixedBlockSize());
199 #ifdef HAVE_XPETRA_TPETRA
201 RCP<const TpetraCrsMatrix> oldTCrsOp = Teuchos::rcp_dynamic_cast<
const TpetraCrsMatrix>(oldCrsOp);
202 if (oldTCrsOp != Teuchos::null) {
203 RCP<CrsMatrix> newTCrsOp(
new TpetraCrsMatrix(*oldTCrsOp));
204 RCP<CrsMatrixWrap> newOp (
new CrsMatrixWrap(newTCrsOp));
205 newOp->SetFixedBlockSize(A->GetFixedBlockSize());
209 throw Exceptions::BadCast(
"Cast from Xpetra::Matrix to Xpetra::EpetraCrsMatrix or Xpetra::TpetraCrsMatrix failed");
212 return Teuchos::null;
215 #endif // HAVE_XPETRA_INT_LONG_LONG
217 #define XPETRA_MATRIXFACTORY2_SHORT
220 template <class Scalar = Matrix<>::scalar_type,
221 class LocalOrdinal =
typename Matrix<Scalar>::local_ordinal_type,
222 class GlobalOrdinal =
223 typename Matrix<Scalar, LocalOrdinal>::global_ordinal_type,
225 typename Matrix<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
226 class MatrixFactory {
227 #undef XPETRA_MATRIXFACTORY_SHORT
248 return rcp(
new CrsMatrixWrap(rowMap, colMap, NumEntriesPerRowToAlloc, pftype));
251 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
262 static RCP<Matrix>
Build (
270 return rcp(
new CrsMatrixWrap(lclMatrix, rowMap, colMap, domainMap, rangeMap, params));
287 LocalOrdinal NumMyElements = diagonal->getMap()->getNodeNumElements();
292 for (LocalOrdinal i = 0; i < NumMyElements; ++i) {
293 mtx->insertGlobalValues(MyGlobalElements[i],
294 Teuchos::tuple<GlobalOrdinal>(MyGlobalElements[i]),
295 Teuchos::tuple<Scalar>(vals[i]) );
304 if (crsOp == Teuchos::null)
309 if (newCrs->hasMatrix())
312 return Teuchos::null;
318 if (crsOp == Teuchos::null)
328 if (crsOp == Teuchos::null)
333 if (newCrs->hasMatrix())
336 return Teuchos::null;
342 if (crsOp == Teuchos::null)
347 if (newCrs->hasMatrix())
350 return Teuchos::null;
358 if(input == Teuchos::null)
368 for (
size_t r = 0; r < input->Rows(); ++r) {
369 for (
size_t c = 0; c < input->Cols(); ++c)
370 if(input->getMatrix(r,c) != Teuchos::null) {
375 bop->setMatrix(r,c,mat);
379 if(input->isFillComplete())
384 #define XPETRA_MATRIXFACTORY_SHORT
388 #define XPETRA_MATRIXFACTORY_SHORT
389 #define XPETRA_MATRIXFACTORY2_SHORT
static RCP< Matrix > Build(const RCP< const Vector > &diagonal)
Constructor for creating a diagonal Xpetra::Matrix using the entries of a given vector for the diagon...
static RCP< Matrix > Build(const RCP< const Matrix > &sourceMatrix, const Import &importer, const RCP< const Map > &domainMap=Teuchos::null, const RCP< const Map > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Constructor to create a Matrix using a fusedImport-style construction. The originalMatrix must be a X...
static RCP< Matrix > Build(const RCP< const CrsGraph > &graph, const RCP< ParameterList > ¶mList=Teuchos::null)
Constructor specifying graph.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
static RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Build(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, Xpetra::ProfileType pftype=Xpetra::DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null)
Constructor specifying fixed number of entries for each row.
Exception throws to report errors in the internal logical of the program.
static RCP< Matrix > Build(const RCP< const Map > &rowMap, const RCP< const Map > &colMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, Xpetra::ProfileType pftype=Xpetra::DynamicProfile)
Constructor specifying the (possibly different) number of entries per row and providing column map...
static RCP< Matrix > Build(const RCP< const Map > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=Xpetra::DynamicProfile)
Constructor specifying (possibly different) number of entries in each row.
Exception indicating invalid cast attempted.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > BuildCopy(const RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > A)
static RCP< Matrix > Build(const RCP< const Map > &rowMap, const RCP< const Map > &colMap, size_t maxNumEntriesPerRow, Xpetra::ProfileType pftype=Xpetra::DynamicProfile)
Constructor specifying the max number of non-zeros per row and providing column map.
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > BuildCopy(const RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > A)
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > BuildCopy(const RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > A)
static RCP< Matrix > Build(const RCP< const Matrix > &sourceMatrix, const Export &exporter, const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Constructor to create a Matrix using a fusedExport-style construction. The originalMatrix must be a X...
#define TEUCHOS_UNREACHABLE_RETURN(dummyReturnVal)
TypeTo as(const TypeFrom &t)
Concrete implementation of Xpetra::Matrix.
#define XPETRA_MONITOR(funcName)
static RCP< Matrix > Build(const RCP< const Matrix > &sourceMatrix, const Export &RowExporter, const Export &DomainExporter, const RCP< const Map > &domainMap=Teuchos::null, const RCP< const Map > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Constructor to create a Matrix using a fusedExport-style construction. The originalMatrix must be a X...
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.
Xpetra-specific matrix class.
MatrixFactory()
Private constructor. This is a static class.
static RCP< Matrix > Build(const RCP< const Matrix > &sourceMatrix, const Import &RowImporter, const Import &DomainImporter, const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Constructor to create a Matrix using a fusedImport-style construction. The originalMatrix must be a X...