All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Xpetra_MatrixFactory.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 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 
47 // WARNING: This code is experimental. Backwards compatibility should not be expected.
48 
49 #ifndef XPETRA_MATRIXFACTORY_HPP
50 #define XPETRA_MATRIXFACTORY_HPP
51 
52 #include "Xpetra_ConfigDefs.hpp"
54 #include "Xpetra_Matrix.hpp"
55 #include "Xpetra_CrsMatrixWrap.hpp"
57 #include "Xpetra_Map.hpp"
58 #include "Xpetra_Vector.hpp"
59 #include "Xpetra_Exceptions.hpp"
60 
61 namespace Xpetra {
62 
63  /*template <class Scalar = Matrix<>::scalar_type,
64  class LocalOrdinal = typename Matrix<Scalar>::local_ordinal_type,
65  class GlobalOrdinal =
66  typename Matrix<Scalar, LocalOrdinal>::global_ordinal_type,
67  class Node =
68  typename Matrix<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>*/
69  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
70  class MatrixFactory2 {
71 #undef XPETRA_MATRIXFACTORY2_SHORT
72 #include "Xpetra_UseShortNames.hpp"
73 
74  public:
76  RCP<const CrsMatrixWrap> oldOp = Teuchos::rcp_dynamic_cast<const CrsMatrixWrap>(A);
77  if (oldOp == Teuchos::null)
78  throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
79 
80  RCP<const CrsMatrix> oldCrsOp = oldOp->getCrsMatrix();
81 
82  UnderlyingLib lib = A->getRowMap()->lib();
83 
85  "Not Epetra or Tpetra matrix");
86 
87 #ifdef HAVE_XPETRA_EPETRA
88  if (lib == UseEpetra) {
89  // NOTE: The proper Epetra conversion in Xpetra_MatrixFactory.cpp
90  throw Exceptions::RuntimeError("Xpetra::BuildCopy(): matrix templates are incompatible with Epetra");
91  }
92 #endif
93 
94 #ifdef HAVE_XPETRA_TPETRA
95  if (lib == UseTpetra) {
96  // Underlying matrix is Tpetra
97  RCP<const TpetraCrsMatrix> oldTCrsOp = Teuchos::rcp_dynamic_cast<const TpetraCrsMatrix>(oldCrsOp);
98 
99  if (oldTCrsOp != Teuchos::null) {
100  RCP<TpetraCrsMatrix> newTCrsOp(new TpetraCrsMatrix(*oldTCrsOp));
102  newOp->SetFixedBlockSize(A->GetFixedBlockSize());
103 
104  return newOp;
105  } else {
106  throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::TpetraCrsMatrix failed");
107  }
108  }
109 #endif
110 
111  return Teuchos::null;
112  }
113  };
114 #define XPETRA_MATRIXFACTORY2_SHORT
115 
116  //template<>
117  //class MatrixFactory2<double,int,int,typename Xpetra::Matrix<double, int, int>::node_type> {
118  template<class Node>
119  class MatrixFactory2<double,int,int,Node> {
120  typedef double Scalar;
121  typedef int LocalOrdinal;
122  typedef int GlobalOrdinal;
123  //typedef Matrix<double, int, GlobalOrdinal>::node_type Node;
124 #undef XPETRA_MATRIXFACTORY2_SHORT
125 #include "Xpetra_UseShortNames.hpp"
126  public:
128  RCP<const CrsMatrixWrap> oldOp = Teuchos::rcp_dynamic_cast<const CrsMatrixWrap>(A);
129  if (oldOp == Teuchos::null)
130  throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
131 
132  RCP<const CrsMatrix> oldCrsOp = oldOp->getCrsMatrix();
133 
134 #ifdef HAVE_XPETRA_EPETRA
135 #ifndef XPETRA_EPETRA_NO_32BIT_GLOBAL_INDICES
136  RCP<const EpetraCrsMatrixT<GlobalOrdinal,Node> > oldECrsOp = Teuchos::rcp_dynamic_cast<const EpetraCrsMatrixT<GlobalOrdinal,Node> >(oldCrsOp);
137  if (oldECrsOp != Teuchos::null) {
138  // Underlying matrix is Epetra
139  RCP<CrsMatrix> newECrsOp(new EpetraCrsMatrixT<GlobalOrdinal,Node>(*oldECrsOp));
140  RCP<CrsMatrixWrap> newOp (new CrsMatrixWrap (newECrsOp));
141  newOp->SetFixedBlockSize(A->GetFixedBlockSize());
142  return newOp;
143  }
144 #endif
145 #endif
146 
147 #ifdef HAVE_XPETRA_TPETRA
148  // Underlying matrix is Tpetra
149  RCP<const TpetraCrsMatrix> oldTCrsOp = Teuchos::rcp_dynamic_cast<const TpetraCrsMatrix>(oldCrsOp);
150  if (oldTCrsOp != Teuchos::null) {
151  RCP<CrsMatrix> newTCrsOp(new TpetraCrsMatrix(*oldTCrsOp));
152  RCP<CrsMatrixWrap> newOp (new CrsMatrixWrap(newTCrsOp));
153  newOp->SetFixedBlockSize(A->GetFixedBlockSize());
154  return newOp;
155  }
156  return Teuchos::null;
157 #else
158  throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::EpetraCrsMatrix or Xpetra::TpetraCrsMatrix failed");
159  TEUCHOS_UNREACHABLE_RETURN(Teuchos::null); // make compiler happy
160 #endif
161 
162  } //BuildCopy
163  };
164 
165 #define XPETRA_MATRIXFACTORY2_SHORT
166 
167 #ifdef HAVE_XPETRA_INT_LONG_LONG
168  //template<>
169  //class MatrixFactory2<double,int,long long,typename Xpetra::Matrix<double, int, long long>::node_type> {
170  template<class Node>
171  class MatrixFactory2<double, int, long long, Node> {
172  typedef double Scalar;
173  typedef int LocalOrdinal;
174  typedef long long GlobalOrdinal;
175  //typedef Matrix<double, int, GlobalOrdinal>::node_type Node;
176 #undef XPETRA_MATRIXFACTORY2_SHORT
177 #include "Xpetra_UseShortNames.hpp"
178  public:
179  static RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > BuildCopy(const RCP<const Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > A) {
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");
183 
184  RCP<const CrsMatrix> oldCrsOp = oldOp->getCrsMatrix();
185 
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) {
190  // Underlying matrix is Epetra
191  RCP<CrsMatrix> newECrsOp(new EpetraCrsMatrixT<GlobalOrdinal,Node>(*oldECrsOp));
192  RCP<CrsMatrixWrap> newOp (new CrsMatrixWrap (newECrsOp));
193  newOp->SetFixedBlockSize(A->GetFixedBlockSize());
194  return newOp;
195  }
196 #endif
197 #endif
198 
199 #ifdef HAVE_XPETRA_TPETRA
200  // Underlying matrix is 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());
206  return newOp;
207  }
208 #else
209  throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::EpetraCrsMatrix or Xpetra::TpetraCrsMatrix failed");
210 #endif
211 
212  return Teuchos::null; // make compiler happy
213  }
214  };
215 #endif // HAVE_XPETRA_INT_LONG_LONG
216 
217 #define XPETRA_MATRIXFACTORY2_SHORT
218 
219 
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,
224  class Node =
225  typename Matrix<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
226  class MatrixFactory {
227 #undef XPETRA_MATRIXFACTORY_SHORT
228 #include "Xpetra_UseShortNames.hpp"
229 
230  private:
233 
234  public:
235 
237  static RCP<Matrix> Build(const RCP<const Map>& rowMap, size_t maxNumEntriesPerRow, Xpetra::ProfileType pftype = Xpetra::DynamicProfile) {
238  return rcp(new CrsMatrixWrap(rowMap, maxNumEntriesPerRow, pftype));
239  }
240 
242  static RCP<Matrix> Build(const RCP<const Map>& rowMap, const RCP<const Map>& colMap, size_t maxNumEntriesPerRow, Xpetra::ProfileType pftype = Xpetra::DynamicProfile) {
243  return rcp(new CrsMatrixWrap(rowMap, colMap, maxNumEntriesPerRow, pftype));
244  }
245 
247  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) {
248  return rcp(new CrsMatrixWrap(rowMap, colMap, NumEntriesPerRowToAlloc, pftype));
249  }
250 
251 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
252  static RCP<Matrix> Build (
254  const Teuchos::RCP<const Map>& rowMap,
255  const Teuchos::RCP<const Map>& colMap,
257  const Teuchos::RCP<Teuchos::ParameterList>& params = null) {
258  XPETRA_MONITOR("MatrixFactory::Build");
259  return rcp(new CrsMatrixWrap(rowMap, colMap, lclMatrix, params));
260  }
262  static RCP<Matrix> Build (
264  const Teuchos::RCP<const Map>& rowMap,
265  const Teuchos::RCP<const Map>& colMap,
266  const Teuchos::RCP<const Map>& domainMap = Teuchos::null,
267  const Teuchos::RCP<const Map>& rangeMap = Teuchos::null,
268  const Teuchos::RCP<Teuchos::ParameterList>& params = null) {
269  XPETRA_MONITOR("MatrixFactory::Build");
270  return rcp(new CrsMatrixWrap(lclMatrix, rowMap, colMap, domainMap, rangeMap, params));
271  }
272 #endif
273 
275  static RCP<Matrix> Build(const RCP<const Map> &rowMap, const ArrayRCP<const size_t> &NumEntriesPerRowToAlloc, ProfileType pftype = Xpetra::DynamicProfile) {
276  return rcp( new CrsMatrixWrap(rowMap, NumEntriesPerRowToAlloc, pftype) );
277  }
278 
280  static RCP<Matrix> Build(const RCP<const CrsGraph>& graph, const RCP<ParameterList>& paramList = Teuchos::null) {
281  return rcp(new CrsMatrixWrap(graph, paramList));
282  }
283 
285  static RCP<Matrix> Build(const RCP<const Vector>& diagonal) {
286  Teuchos::ArrayRCP<const Scalar> vals = diagonal->getData(0);
287  LocalOrdinal NumMyElements = diagonal->getMap()->getNodeNumElements();
288  Teuchos::ArrayView<const GlobalOrdinal> MyGlobalElements = diagonal->getMap()->getNodeElementList();
289 
291 
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]) );
296  }
297  mtx->fillComplete();
298  return mtx;
299  }
300 
302  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>& params = Teuchos::null) {
303  RCP<const CrsMatrixWrap> crsOp = Teuchos::rcp_dynamic_cast<const CrsMatrixWrap>(sourceMatrix);
304  if (crsOp == Teuchos::null)
305  throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
306 
307  RCP<CrsMatrix> originalCrs = crsOp->getCrsMatrix();
308  RCP<CrsMatrix> newCrs = CrsMatrixFactory::Build(originalCrs, importer, domainMap, rangeMap, params);
309  if (newCrs->hasMatrix())
310  return rcp(new CrsMatrixWrap(newCrs));
311  else
312  return Teuchos::null;
313  }
314 
316  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>& params) {
317  RCP<const CrsMatrixWrap> crsOp = Teuchos::rcp_dynamic_cast<const CrsMatrixWrap>(sourceMatrix);
318  if (crsOp == Teuchos::null)
319  throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
320 
321  RCP<CrsMatrix> originalCrs = crsOp->getCrsMatrix();
322  return rcp(new CrsMatrixWrap(CrsMatrixFactory::Build(originalCrs, exporter, domainMap, rangeMap, params)));
323  }
324 
326  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>& params) {
327  RCP<const CrsMatrixWrap> crsOp = Teuchos::rcp_dynamic_cast<const CrsMatrixWrap>(sourceMatrix);
328  if (crsOp == Teuchos::null)
329  throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
330 
331  RCP<CrsMatrix> originalCrs = crsOp->getCrsMatrix();
332  RCP<CrsMatrix> newCrs = CrsMatrixFactory::Build(originalCrs, RowImporter, Teuchos::rcpFromRef(DomainImporter), domainMap, rangeMap, params);
333  if (newCrs->hasMatrix())
334  return rcp(new CrsMatrixWrap(newCrs));
335  else
336  return Teuchos::null;
337  }
338 
340  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>& params = Teuchos::null) {
341  RCP<const CrsMatrixWrap> crsOp = Teuchos::rcp_dynamic_cast<const CrsMatrixWrap>(sourceMatrix);
342  if (crsOp == Teuchos::null)
343  throw Exceptions::BadCast("Cast from Xpetra::Matrix to Xpetra::CrsMatrixWrap failed");
344 
345  RCP<CrsMatrix> originalCrs = crsOp->getCrsMatrix();
346  RCP<CrsMatrix> newCrs = CrsMatrixFactory::Build(originalCrs, RowExporter, Teuchos::rcpFromRef(DomainExporter), domainMap, rangeMap, params);
347  if (newCrs->hasMatrix())
348  return rcp(new CrsMatrixWrap(newCrs));
349  else
350  return Teuchos::null;
351  }
352 
353 
357  RCP<const BlockedCrsMatrix> input = Teuchos::rcp_dynamic_cast<const BlockedCrsMatrix>(A);
358  if(input == Teuchos::null)
360 
361  // deep copy of MapExtractors (and underlying maps)
362  RCP<const MapExtractor> rgMapExt = Teuchos::rcp(new MapExtractor(*(input->getRangeMapExtractor())));
363  RCP<const MapExtractor> doMapExt = Teuchos::rcp(new MapExtractor(*(input->getDomainMapExtractor())));
364 
365  // create new BlockedCrsMatrix object
366  RCP<BlockedCrsMatrix> bop = Teuchos::rcp(new BlockedCrsMatrix(rgMapExt, doMapExt, input->getNodeMaxNumRowEntries()));
367 
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) {
371  // make a deep copy of the matrix
372  // This is a recursive call to this function
373  RCP<Matrix> mat =
375  bop->setMatrix(r,c,mat);
376  }
377  }
378 
379  if(input->isFillComplete())
380  bop->fillComplete();
381  return bop;
382  }
383  };
384 #define XPETRA_MATRIXFACTORY_SHORT
385 
386 }
387 
388 #define XPETRA_MATRIXFACTORY_SHORT
389 #define XPETRA_MATRIXFACTORY2_SHORT
390 #endif
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 > &params=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 > &paramList=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 > &params)
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 > &params=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 > &params)
Constructor to create a Matrix using a fusedImport-style construction. The originalMatrix must be a X...