Xpetra Version of the Day
Loading...
Searching...
No Matches
Xpetra_TripleMatrixMultiply_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#ifndef PACKAGES_XPETRA_SUP_UTILS_XPETRA_TRIPLEMATRIXMULTIPLY_DEF_HPP_
11#define PACKAGES_XPETRA_SUP_UTILS_XPETRA_TRIPLEMATRIXMULTIPLY_DEF_HPP_
12
14
15namespace Xpetra {
16
17template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
19 const Matrix& A, bool transposeA,
20 const Matrix& P, bool transposeP,
21 Matrix& Ac,
22 bool call_FillComplete_on_result,
23 bool doOptimizeStorage,
24 const std::string& label,
25 const RCP<ParameterList>& params) {
26 TEUCHOS_TEST_FOR_EXCEPTION(transposeR == false && Ac.getRowMap()->isSameAs(*R.getRowMap()) == false,
27 Exceptions::RuntimeError, "XpetraExt::TripleMatrixMultiply::MultiplyRAP: row map of Ac is not same as row map of R");
28 TEUCHOS_TEST_FOR_EXCEPTION(transposeR == true && Ac.getRowMap()->isSameAs(*R.getDomainMap()) == false,
29 Exceptions::RuntimeError, "XpetraExt::TripleMatrixMultiply::MultiplyRAP: row map of Ac is not same as domain map of R");
30
34
35 bool haveMultiplyDoFillComplete = call_FillComplete_on_result && doOptimizeStorage;
36
37 if (Ac.getRowMap()->lib() == Xpetra::UseTpetra) {
38 using helpers = Xpetra::Helpers<SC, LO, GO, NO>;
39 if (helpers::isTpetraCrs(R) && helpers::isTpetraCrs(A) && helpers::isTpetraCrs(P)) {
40 // All matrices are Crs
45
46 // 18Feb2013 JJH I'm reenabling the code that allows the matrix matrix multiply to do the fillComplete.
47 // Previously, Tpetra's matrix matrix multiply did not support fillComplete.
48 Tpetra::TripleMatrixMultiply::MultiplyRAP(tpR, transposeR, tpA, transposeA, tpP, transposeP, tpAc, haveMultiplyDoFillComplete, label, params);
49 } else if (helpers::isTpetraBlockCrs(R) && helpers::isTpetraBlockCrs(A) && helpers::isTpetraBlockCrs(P)) {
50 // All matrices are BlockCrs (except maybe Ac)
51 // FIXME: For the moment we're just going to clobber the innards of Ac, so no reuse. Once we have a reuse kernel,
52 // we'll need to think about refactoring BlockCrs so we can do something smarter here.
53
54 if (!A.getRowMap()->getComm()->getRank())
55 std::cout << "WARNING: Using inefficient BlockCrs Multiply Placeholder" << std::endl;
56
60 // Tpetra::BlockCrsMatrix<SC,LO,GO,NO> & tpAc = Xpetra::Helpers<SC,LO,GO,NO>::Op2NonConstTpetraBlockCrs(Ac);
61
66 // RCP<CRS> Accrs = Tpetra::convertToCrsMatrix(tpAc);
67
68 // FIXME: The lines below only works because we're assuming Ac is Point
69 RCP<CRS> Accrs = Teuchos::rcp(new CRS(Rcrs->getRowMap(), 0));
70 const bool do_fill_complete = true;
71 Tpetra::TripleMatrixMultiply::MultiplyRAP(*Rcrs, transposeR, *Acrs, transposeA, *Pcrs, transposeP, *Accrs, do_fill_complete, label, params);
72
73 // Temporary output matrix
77
78 // We can now cheat and replace the innards of Ac
79 RCP<Xpetra::CrsMatrixWrap<SC, LO, GO, NO> > Ac_w = Teuchos::rcp_dynamic_cast<Xpetra::CrsMatrixWrap<SC, LO, GO, NO> >(Teuchos::rcpFromRef(Ac));
80 Ac_w->replaceCrsMatrix(Ac_p);
81 } else {
82 // Mix and match
83 TEUCHOS_TEST_FOR_EXCEPTION(1, Exceptions::RuntimeError, "Mix-and-match Crs/BlockCrs Multiply not currently supported");
84 }
85 }
86
87 if (call_FillComplete_on_result && !haveMultiplyDoFillComplete) {
89 fillParams->set("Optimize Storage", doOptimizeStorage);
90 Ac.fillComplete((transposeP) ? P.getRangeMap() : P.getDomainMap(),
91 (transposeR) ? R.getDomainMap() : R.getRangeMap(),
92 fillParams);
93 }
94
95 // transfer striding information
96 RCP<const Map> domainMap = Teuchos::null;
97 RCP<const Map> rangeMap = Teuchos::null;
98
99 const std::string stridedViewLabel("stridedMaps");
100 const size_t blkSize = 1;
101 std::vector<size_t> stridingInfo(1, blkSize);
102 LocalOrdinal stridedBlockId = -1;
103
104 if (R.IsView(stridedViewLabel)) {
105 rangeMap = transposeR ? R.getColMap(stridedViewLabel) : R.getRowMap(stridedViewLabel);
106 } else {
107 rangeMap = transposeR ? R.getDomainMap() : R.getRangeMap();
108 rangeMap = StridedMapFactory::Build(rangeMap, stridingInfo, stridedBlockId);
109 }
110
111 if (P.IsView(stridedViewLabel)) {
112 domainMap = transposeP ? P.getRowMap(stridedViewLabel) : P.getColMap(stridedViewLabel);
113 } else {
114 domainMap = transposeP ? P.getRangeMap() : P.getDomainMap();
115 domainMap = StridedMapFactory::Build(domainMap, stridingInfo, stridedBlockId);
116 }
117 Ac.CreateView(stridedViewLabel, rangeMap, domainMap);
118
119} // end Multiply
120
121} // end namespace Xpetra
122
123#define XPETRA_TRIPLEMATRIXMULTIPLY_SHORT
124
125#endif /* PACKAGES_XPETRA_SUP_UTILS_XPETRA_TRIPLEMATRIXMULTIPLY_DEF_HPP_ */
Exception throws to report errors in the internal logical of the program.
static RCP< Tpetra::CrsMatrix< SC, LO, GO, NO > > Op2NonConstTpetraCrs(RCP< Matrix > Op)
static RCP< const Tpetra::BlockCrsMatrix< SC, LO, GO, NO > > Op2TpetraBlockCrs(RCP< Matrix > Op)
static RCP< const Tpetra::CrsMatrix< SC, LO, GO, NO > > Op2TpetraCrs(RCP< Matrix > Op)
Xpetra-specific matrix class.
void CreateView(viewLabel_t viewLabel, const RCP< const Map > &rowMap, const RCP< const Map > &colMap)
virtual const RCP< const Map > & getRowMap() const
Returns the Map that describes the row distribution in this matrix.
virtual bool isFillComplete() const =0
Returns true if fillComplete() has been called and the matrix is in compute mode.
virtual const RCP< const Map > & getColMap() const
Returns the Map that describes the column distribution in this matrix. This might be null until fillC...
virtual void fillComplete(const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const RCP< ParameterList > &params=null)=0
Signal that data entry is complete, specifying domain and range maps.
bool IsView(const viewLabel_t viewLabel) const
virtual LocalOrdinal GetStorageBlockSize() const =0
Returns the block size of the storage mechanism, which is usually 1, except for Tpetra::BlockCrsMatri...
virtual const Teuchos::RCP< const map_type > getDomainMap() const =0
The Map associated with the domain of this operator, which must be compatible with X....
virtual const Teuchos::RCP< const map_type > getRangeMap() const =0
The Map associated with the range of this operator, which must be compatible with Y....
static RCP< Xpetra::StridedMap< LocalOrdinal, GlobalOrdinal, Node > > Build(UnderlyingLib lib, global_size_t numGlobalElements, GlobalOrdinal indexBase, std::vector< size_t > &stridingInfo, const Teuchos::RCP< const Teuchos::Comm< int > > &comm, LocalOrdinal stridedBlockId=-1, GlobalOrdinal offset=0, LocalGlobal lg=Xpetra::GloballyDistributed)
Map constructor with Xpetra-defined contiguous uniform distribution.
static void MultiplyRAP(const Matrix &R, bool transposeR, const Matrix &A, bool transposeA, const Matrix &P, bool transposeP, Matrix &Ac, bool call_FillComplete_on_result=true, bool doOptimizeStorage=true, const std::string &label=std::string(), const RCP< ParameterList > &params=null)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void MultiplyRAP(const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &R, bool transposeR, const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, bool transposeA, const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &P, bool transposeP, CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Ac, bool call_FillComplete_on_result=true, const std::string &label=std::string(), const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Teuchos::RCP< CrsMatrix< Scalar, LO, GO, Node > > convertToCrsMatrix(const Tpetra::BlockCrsMatrix< Scalar, LO, GO, Node > &blockMatrix)
Teuchos::RCP< BlockCrsMatrix< Scalar, LO, GO, Node > > convertToBlockCrsMatrix(const Tpetra::CrsMatrix< Scalar, LO, GO, Node > &pointMatrix, const LO &blockSize, bool use_local_ID=true)