Xpetra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Xpetra_MatrixFactory_def.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_MATRIXFACTORY_DEF_HPP
13 #define XPETRA_MATRIXFACTORY_DEF_HPP
14 
17 #include "Xpetra_BlockedCrsMatrix.hpp"
18 
19 namespace Xpetra {
20 
21 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
22 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> MatrixFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(const RCP<const Map>& rowMap) {
23  return rcp(new CrsMatrixWrap(rowMap));
24 }
25 
26 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
27 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> MatrixFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(const RCP<const Map>& rowMap, size_t maxNumEntriesPerRow) {
28  return rcp(new CrsMatrixWrap(rowMap, maxNumEntriesPerRow));
29 }
30 
31 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
32 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> MatrixFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(const RCP<const Map>& rowMap, const RCP<const Map>& colMap, size_t maxNumEntriesPerRow) {
33  return rcp(new CrsMatrixWrap(rowMap, colMap, maxNumEntriesPerRow));
34 }
35 
36 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
37 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> MatrixFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(const RCP<const Map>& rowMap, const RCP<const Map>& colMap, const ArrayRCP<const size_t>& NumEntriesPerRowToAlloc) {
38  return rcp(new CrsMatrixWrap(rowMap, colMap, NumEntriesPerRowToAlloc));
39 }
40 
41 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
42 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> MatrixFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(
43  const Teuchos::RCP<const Map>& rowMap,
44  const Teuchos::RCP<const Map>& colMap,
46  const Teuchos::RCP<Teuchos::ParameterList>& params) {
47  XPETRA_MONITOR("MatrixFactory::Build");
48  return rcp(new CrsMatrixWrap(rowMap, colMap, lclMatrix, params));
49 }
50 
51 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
52 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> MatrixFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(
54  const Teuchos::RCP<const Map>& rowMap,
55  const Teuchos::RCP<const Map>& colMap,
56  const Teuchos::RCP<const Map>& domainMap,
57  const Teuchos::RCP<const Map>& rangeMap,
58  const Teuchos::RCP<Teuchos::ParameterList>& params) {
59  XPETRA_MONITOR("MatrixFactory::Build");
60  return rcp(new CrsMatrixWrap(lclMatrix, rowMap, colMap, domainMap, rangeMap, params));
61 }
62 
63 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
64 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> MatrixFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(const RCP<const Map>& rowMap, const ArrayRCP<const size_t>& NumEntriesPerRowToAlloc) {
65  return rcp(new CrsMatrixWrap(rowMap, NumEntriesPerRowToAlloc));
66 }
67 
68 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
69 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> MatrixFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(const RCP<const CrsGraph>& graph, const RCP<ParameterList>& paramList) {
70  return rcp(new CrsMatrixWrap(graph, paramList));
71 }
72 
73 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
74 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> MatrixFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(const RCP<const CrsGraph>& graph,
76  const RCP<ParameterList>& paramList) {
77  return rcp(new CrsMatrixWrap(graph, values, paramList));
78 }
79 
80 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
81 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> MatrixFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(const RCP<const Vector>& diagonal) {
82  RCP<const BlockedVector> bdiagonal = Teuchos::rcp_dynamic_cast<const BlockedVector>(diagonal);
83  Teuchos::RCP<Matrix> mtx = Teuchos::null;
84 
85  if (bdiagonal == Teuchos::null) {
86  Teuchos::ArrayRCP<const Scalar> vals = diagonal->getData(0);
87  LocalOrdinal numMyElements = diagonal->getMap()->getLocalNumElements();
88  Teuchos::ArrayView<const GlobalOrdinal> myGlobalElements = diagonal->getMap()->getLocalElementList();
89 
90  mtx = Teuchos::rcp(new CrsMatrixWrap(diagonal->getMap(), 1));
91 
92  for (LocalOrdinal i = 0; i < numMyElements; ++i) {
93  mtx->insertGlobalValues(myGlobalElements[i],
94  Teuchos::tuple<GlobalOrdinal>(myGlobalElements[i]),
95  Teuchos::tuple<Scalar>(vals[i]));
96  }
97  mtx->fillComplete();
98  } else {
99  RCP<BlockedCrsMatrix> bop = Teuchos::rcp(new BlockedCrsMatrix(bdiagonal->getBlockedMap(), bdiagonal->getBlockedMap(), 1));
100 
101  for (size_t r = 0; r < bdiagonal->getBlockedMap()->getNumMaps(); ++r) {
102  if (!bdiagonal->getMultiVector(r).is_null()) {
103  const RCP<MultiVector> subvec = bdiagonal->getMultiVector(r);
104  bop->setMatrix(r, r, Build(subvec->getVector(0)));
105  }
106  }
107  bop->fillComplete();
108  mtx = BuildCopy(bop);
109  }
110 
111  return mtx;
112 }
113 
114 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
115 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> MatrixFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(const RCP<const Matrix>& sourceMatrix, const Import& importer, const RCP<const Map>& domainMap, const RCP<const Map>& rangeMap, const Teuchos::RCP<Teuchos::ParameterList>& params) {
116  RCP<const CrsMatrixWrap> crsOp = Teuchos::rcp_dynamic_cast<const CrsMatrixWrap>(sourceMatrix);
117  if (crsOp == Teuchos::null)
118  throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
119 
120  RCP<CrsMatrix> originalCrs = crsOp->getCrsMatrix();
121  RCP<CrsMatrix> newCrs = CrsMatrixFactory::Build(originalCrs, importer, domainMap, rangeMap, params);
122  if (newCrs->hasMatrix())
123  return rcp(new CrsMatrixWrap(newCrs));
124  else
125  return Teuchos::null;
126 }
127 
128 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
129 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> MatrixFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(const RCP<const Matrix>& sourceMatrix, const Export& exporter, const RCP<const Map>& domainMap, const RCP<const Map>& rangeMap, const Teuchos::RCP<Teuchos::ParameterList>& params) {
130  RCP<const CrsMatrixWrap> crsOp = Teuchos::rcp_dynamic_cast<const CrsMatrixWrap>(sourceMatrix);
131  if (crsOp == Teuchos::null)
132  throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
133 
134  RCP<CrsMatrix> originalCrs = crsOp->getCrsMatrix();
135  return rcp(new CrsMatrixWrap(CrsMatrixFactory::Build(originalCrs, exporter, domainMap, rangeMap, params)));
136 }
137 
138 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
139 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> MatrixFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::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>& params) {
140  RCP<const CrsMatrixWrap> crsOp = Teuchos::rcp_dynamic_cast<const CrsMatrixWrap>(sourceMatrix);
141  if (crsOp == Teuchos::null)
142  throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
143 
144  RCP<CrsMatrix> originalCrs = crsOp->getCrsMatrix();
145  RCP<CrsMatrix> newCrs = CrsMatrixFactory::Build(originalCrs, RowImporter, Teuchos::rcpFromRef(DomainImporter), domainMap, rangeMap, params);
146  if (newCrs->hasMatrix())
147  return rcp(new CrsMatrixWrap(newCrs));
148  else
149  return Teuchos::null;
150 }
151 
152 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
153 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> MatrixFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(const RCP<const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& sourceMatrix, const Xpetra::Export<LocalOrdinal, GlobalOrdinal, Node>& RowExporter, const Xpetra::Export<LocalOrdinal, GlobalOrdinal, Node>& DomainExporter, const RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>& domainMap, const RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>& rangeMap, const Teuchos::RCP<Teuchos::ParameterList>& params) {
154  RCP<const CrsMatrixWrap> crsOp = Teuchos::rcp_dynamic_cast<const CrsMatrixWrap>(sourceMatrix);
155  if (crsOp == Teuchos::null)
156  throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
157 
158  RCP<CrsMatrix> originalCrs = crsOp->getCrsMatrix();
159  RCP<CrsMatrix> newCrs = CrsMatrixFactory::Build(originalCrs, RowExporter, Teuchos::rcpFromRef(DomainExporter), domainMap, rangeMap, params);
160  if (newCrs->hasMatrix())
161  return rcp(new CrsMatrixWrap(newCrs));
162  else
163  return Teuchos::null;
164 }
165 
166 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
167 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> MatrixFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::BuildCopy(const RCP<const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> A, bool setFixedBlockSize) {
168  RCP<const BlockedCrsMatrix> input = Teuchos::rcp_dynamic_cast<const BlockedCrsMatrix>(A);
169  if (input == Teuchos::null)
171 
172  // deep copy of MapExtractors (and underlying maps)
173  RCP<const MapExtractor> rgMapExt = Teuchos::rcp(new MapExtractor(*(input->getRangeMapExtractor())));
174  RCP<const MapExtractor> doMapExt = Teuchos::rcp(new MapExtractor(*(input->getDomainMapExtractor())));
175 
176  // create new BlockedCrsMatrix object
177  RCP<BlockedCrsMatrix> bop = Teuchos::rcp(new BlockedCrsMatrix(rgMapExt, doMapExt, input->getLocalMaxNumRowEntries()));
178 
179  for (size_t r = 0; r < input->Rows(); ++r) {
180  for (size_t c = 0; c < input->Cols(); ++c)
181  if (input->getMatrix(r, c) != Teuchos::null) {
182  // make a deep copy of the matrix
183  // This is a recursive call to this function
184  RCP<Matrix> mat =
185  Xpetra::MatrixFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::BuildCopy(input->getMatrix(r, c), setFixedBlockSize);
186  bop->setMatrix(r, c, mat);
187  }
188  }
189 
190  if (input->isFillComplete())
191  bop->fillComplete();
192  return bop;
193 }
194 
195 } // namespace Xpetra
196 
197 #endif
static RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Build(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rowMap)
Constructor for empty matrix (intended use is an import/export target - can&#39;t insert entries directly...
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > BuildCopy(const RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > A, bool setFixedBlockSize=true)
Exception indicating invalid cast attempted.
static RCP< Matrix > Build(const RCP< const Map > &rowMap)
virtual Teuchos::ArrayRCP< const Scalar > getData(size_t j) const
Const view of the local values in a particular vector of this vector.
Concrete implementation of Xpetra::Matrix.
#define XPETRA_MONITOR(funcName)
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > BuildCopy(const RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > A, bool setFixedBlockSize=true)
KokkosSparse::CrsMatrix< impl_scalar_type, LocalOrdinal, execution_space, void, typename local_graph_type::size_type > local_matrix_type
The specialization of Kokkos::CrsMatrix that represents the part of the sparse matrix on each MPI pro...
Xpetra-specific matrix class.