MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_Utilities_def.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// MueLu: A package for multigrid based preconditioning
4//
5// Copyright 2012 NTESS and the MueLu contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef MUELU_UTILITIES_DEF_HPP
11#define MUELU_UTILITIES_DEF_HPP
12
13#include <Teuchos_DefaultComm.hpp>
14#include <Teuchos_ParameterList.hpp>
15#include <iostream>
16
17#include "MueLu_ConfigDefs.hpp"
18#include "Xpetra_TpetraRowMatrix.hpp"
19
20#include <MatrixMarket_Tpetra.hpp>
21#include <Tpetra_RowMatrixTransposer.hpp>
22#include <TpetraExt_MatrixMatrix.hpp>
23#include <Xpetra_TpetraMultiVector.hpp>
24#include <Xpetra_TpetraOperator.hpp>
25#include <Xpetra_TpetraCrsMatrix.hpp>
26#include <Xpetra_TpetraBlockCrsMatrix.hpp>
27
28#include <Xpetra_BlockedCrsMatrix.hpp>
29//#include <Xpetra_DefaultPlatform.hpp>
30#include <Xpetra_IO.hpp>
31#include <Xpetra_Map.hpp>
32#include <Xpetra_MapFactory.hpp>
33#include <Xpetra_Matrix.hpp>
34#include <Xpetra_MatrixFactory.hpp>
35#include <Xpetra_MultiVector.hpp>
36#include <Xpetra_MultiVectorFactory.hpp>
37#include <Xpetra_Operator.hpp>
38#include <Xpetra_Vector.hpp>
39#include <Xpetra_VectorFactory.hpp>
40
41#include <Xpetra_MatrixMatrix.hpp>
42
44
45namespace MueLu {
46
47template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
48RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
50 Transpose(Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& Op, bool /* optimizeTranspose */, const std::string& label, const Teuchos::RCP<Teuchos::ParameterList>& params) {
51 std::string TorE = "tpetra";
52
53 if (TorE == "tpetra") {
54 using Helpers = Xpetra::Helpers<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
55 /***************************************************************/
56 if (Helpers::isTpetraCrs(Op)) {
57 const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& tpetraOp = toTpetra(Op);
58
59 RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> A;
60 Tpetra::RowMatrixTransposer<Scalar, LocalOrdinal, GlobalOrdinal, Node> transposer(rcpFromRef(tpetraOp), label); // more than meets the eye
61
62 {
63 using Teuchos::ParameterList;
64 using Teuchos::rcp;
65 RCP<ParameterList> transposeParams = params.is_null() ? rcp(new ParameterList) : rcp(new ParameterList(*params));
66 transposeParams->set("sort", false);
67 A = transposer.createTranspose(transposeParams);
68 }
69
70 RCP<Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> AA = rcp(new Xpetra::TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(A));
71 RCP<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> AAA = rcp_implicit_cast<Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(AA);
72 RCP<Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> AAAA = rcp(new Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>(AAA));
73 if (!AAAA->isFillComplete())
74 AAAA->fillComplete(Op.getRangeMap(), Op.getDomainMap());
75
76 if (Op.IsView("stridedMaps"))
77 AAAA->CreateView("stridedMaps", Teuchos::rcpFromRef(Op), true /*doTranspose*/);
78
79 return AAAA;
80 } else if (Helpers::isTpetraBlockCrs(Op)) {
81 using XMatrix = Xpetra::Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
82 using XCrsMatrix = Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
83 using XCrsMatrixWrap = Xpetra::CrsMatrixWrap<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
84 using BCRS = Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
85 // using CRS = Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>;
86 const BCRS& tpetraOp = toTpetraBlock(Op);
87
88 RCP<BCRS> At;
89 {
90 Tpetra::BlockCrsMatrixTransposer<Scalar, LocalOrdinal, GlobalOrdinal, Node> transposer(rcpFromRef(tpetraOp), label);
91
92 using Teuchos::ParameterList;
93 using Teuchos::rcp;
94 RCP<ParameterList> transposeParams = params.is_null() ? rcp(new ParameterList) : rcp(new ParameterList(*params));
95 transposeParams->set("sort", false);
96 At = transposer.createTranspose(transposeParams);
97 }
98
99 RCP<Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> AA = rcp(new Xpetra::TpetraBlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(At));
100 RCP<XCrsMatrix> AAA = rcp_implicit_cast<XCrsMatrix>(AA);
101 RCP<XMatrix> AAAA = rcp(new XCrsMatrixWrap(AAA));
102
103 if (Op.IsView("stridedMaps"))
104 AAAA->CreateView("stridedMaps", Teuchos::rcpFromRef(Op), true /*doTranspose*/);
105
106 return AAAA;
107 } else {
108 throw Exceptions::RuntimeError("Utilities::Transpose failed, perhaps because matrix is not a Crs matrix");
109 }
110 } // if
111
112 // Epetra case
113 std::cout << "Utilities::Transpose() not implemented for Epetra" << std::endl;
114 return Teuchos::null;
115
116} // Transpose
117
118template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
119RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
121 RealValuedToScalarMultiVector(RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::coordinateType, LocalOrdinal, GlobalOrdinal, Node>> X) {
122 RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Xscalar;
123#if (defined(HAVE_TPETRA_INST_COMPLEX_DOUBLE) || defined(HAVE_TPETRA_INST_COMPLEX_FLOAT))
124 using range_type = Kokkos::RangePolicy<LocalOrdinal, typename Node::execution_space>;
125
126 // Need to cast the real-valued multivector to Scalar=complex
127 if ((typeid(Scalar).name() == typeid(std::complex<double>).name()) ||
128 (typeid(Scalar).name() == typeid(std::complex<float>).name())) {
129 size_t numVecs = X->getNumVectors();
130 Xscalar = Xpetra::MultiVectorFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(X->getMap(), numVecs);
131 auto XVec = X->getLocalViewDevice(Tpetra::Access::ReadOnly);
132 auto XVecScalar = Xscalar->getLocalViewDevice(Tpetra::Access::ReadWrite);
133
134 Kokkos::parallel_for(
135 "MueLu:Utils::RealValuedToScalarMultiVector", range_type(0, X->getLocalLength()),
136 KOKKOS_LAMBDA(const size_t i) {
137 for (size_t j = 0; j < numVecs; j++)
138 XVecScalar(i, j) = XVec(i, j);
139 });
140 } else
141#endif
142 Xscalar = rcp_dynamic_cast<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>(X);
143 return Xscalar;
144}
145
146template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
147RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LocalOrdinal, GlobalOrdinal, Node>>
149 ExtractCoordinatesFromParameterList(ParameterList& paramList) {
150 RCP<Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LocalOrdinal, GlobalOrdinal, Node>> coordinates = Teuchos::null;
151
152 // check whether coordinates are contained in parameter list
153 if (paramList.isParameter("Coordinates") == false)
154 return coordinates;
155
156 // define Tpetra::MultiVector type with Scalar=float only if
157 // * ETI is turned off, since then the compiler will instantiate it automatically OR
158 // * Tpetra is instantiated on Scalar=float
159#if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT)
160 typedef Tpetra::MultiVector<float, LocalOrdinal, GlobalOrdinal, Node> tfMV;
161 RCP<tfMV> floatCoords = Teuchos::null;
162#endif
163
164 // define Tpetra::MultiVector type with Scalar=double only if
165 // * ETI is turned off, since then the compiler will instantiate it automatically OR
166 // * Tpetra is instantiated on Scalar=double
167#if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_DOUBLE)
168 typedef Tpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::magnitudeType, LocalOrdinal, GlobalOrdinal, Node> tdMV;
169 RCP<tdMV> doubleCoords = Teuchos::null;
170 if (paramList.isType<RCP<tdMV>>("Coordinates")) {
171 // Coordinates are stored as a double vector
172 doubleCoords = paramList.get<RCP<tdMV>>("Coordinates");
173 paramList.remove("Coordinates");
174 }
175#if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT)
176 else if (paramList.isType<RCP<tfMV>>("Coordinates")) {
177 // check if coordinates are stored as a float vector
178 floatCoords = paramList.get<RCP<tfMV>>("Coordinates");
179 paramList.remove("Coordinates");
180 doubleCoords = rcp(new tdMV(floatCoords->getMap(), floatCoords->getNumVectors()));
181 deep_copy(*doubleCoords, *floatCoords);
182 }
183#endif
184 // We have the coordinates in a Tpetra double vector
185 if (doubleCoords != Teuchos::null) {
186 // rcp(new Xpetra::TpetraMultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>(Vtpetra));
187 coordinates = Xpetra::toXpetra(doubleCoords);
188 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(coordinates));
189 TEUCHOS_TEST_FOR_EXCEPT(doubleCoords->getNumVectors() != coordinates->getNumVectors());
190 }
191#else
192 // coordinates usually are stored as double vector
193 // Tpetra is not instantiated on scalar=double
194 throw Exceptions::RuntimeError("ExtractCoordinatesFromParameterList: The coordinates vector in parameter list is expected to be a Tpetra multivector with SC=double or float.");
195#endif
196
197 // check for Xpetra coordinates vector
198 if (paramList.isType<decltype(coordinates)>("Coordinates")) {
199 coordinates = paramList.get<decltype(coordinates)>("Coordinates");
200 }
201
202 return coordinates;
203} // ExtractCoordinatesFromParameterList
204
205} // namespace MueLu
206
207#define MUELU_UTILITIES_SHORT
208#endif // MUELU_UTILITIES_DEF_HPP
209
210// LocalWords: LocalOrdinal
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultScalar Scalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
Exception throws to report errors in the internal logical of the program.
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Transpose(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, bool optimizeTranspose=false, const std::string &label=std::string(), const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
static RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > RealValuedToScalarMultiVector(RCP< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::coordinateType, LocalOrdinal, GlobalOrdinal, Node > > X)
static RCP< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::magnitudeType, LocalOrdinal, GlobalOrdinal, Node > > ExtractCoordinatesFromParameterList(ParameterList &paramList)
Namespace for MueLu classes and methods.