46#ifndef THYRA_MUELU_MAXWELL1_PRECONDITIONER_FACTORY_DEF_HPP
47#define THYRA_MUELU_MAXWELL1_PRECONDITIONER_FACTORY_DEF_HPP
51#include <Xpetra_CrsMatrixWrap.hpp>
52#include <Xpetra_CrsMatrix.hpp>
53#include <Xpetra_Matrix.hpp>
54#include <Xpetra_ThyraUtils.hpp>
55#include <MueLu_Maxwell1.hpp>
56#include <Xpetra_TpetraHalfPrecisionOperator.hpp>
59#if defined(HAVE_MUELU_STRATIMIKOS) && defined(HAVE_MUELU_THYRA)
62#if((defined(HAVE_TPETRA_INST_DOUBLE) && defined(HAVE_TPETRA_INST_FLOAT) && !defined(HAVE_TPETRA_INST_COMPLEX_DOUBLE) && !defined(HAVE_TPETRA_INST_COMPLEX_FLOAT)) || \
63 (!defined(HAVE_TPETRA_INST_DOUBLE) && !defined(HAVE_TPETRA_INST_FLOAT) && defined(HAVE_TPETRA_INST_COMPLEX_DOUBLE) && defined(HAVE_TPETRA_INST_COMPLEX_FLOAT)) || \
64 (defined(HAVE_TPETRA_INST_DOUBLE) && defined(HAVE_TPETRA_INST_FLOAT) && defined(HAVE_TPETRA_INST_COMPLEX_DOUBLE) && defined(HAVE_TPETRA_INST_COMPLEX_FLOAT)))
65# define MUELU_CAN_USE_MIXED_PRECISION
72 using Teuchos::ParameterList;
73 using Teuchos::rcp_dynamic_cast;
74 using Teuchos::rcp_const_cast;
78 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
79 MueLuMaxwell1PreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::MueLuMaxwell1PreconditionerFactory() :
80 paramList_(rcp(new ParameterList()))
85 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
86 bool MueLuMaxwell1PreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isCompatible(
const LinearOpSourceBase<Scalar>& fwdOpSrc)
const {
87 const RCP<const LinearOpBase<Scalar> > fwdOp = fwdOpSrc.getOp();
89 if (Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isTpetra(fwdOp))
return true;
91#ifdef HAVE_MUELU_EPETRA
92 if (Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isEpetra(fwdOp))
return true;
99 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
100 RCP<PreconditionerBase<Scalar> > MueLuMaxwell1PreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::createPrec()
const {
101 return Teuchos::rcp(
new DefaultPreconditioner<Scalar>);
104 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
105 void MueLuMaxwell1PreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::
106 initializePrec(
const RCP<
const LinearOpSourceBase<Scalar> >& fwdOpSrc, PreconditionerBase<Scalar>* prec,
const ESupportSolveUse supportSolveUse)
const {
109 typedef Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node> XpOp;
110 typedef Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpThyUtils;
111 typedef Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpMat;
112 typedef Thyra::LinearOpBase<Scalar> ThyLinOpBase;
114#if defined(MUELU_CAN_USE_MIXED_PRECISION)
115 typedef Xpetra::TpetraHalfPrecisionOperator<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpHalfPrecOp;
116 typedef Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpMV;
117 typedef typename XpHalfPrecOp::HalfScalar HalfScalar;
118 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType Magnitude;
119 typedef typename Teuchos::ScalarTraits<Magnitude>::halfPrecision HalfMagnitude;
120 typedef Xpetra::MultiVector<HalfScalar,LocalOrdinal,GlobalOrdinal,Node> XphMV;
121 typedef Xpetra::MultiVector<Magnitude,LocalOrdinal,GlobalOrdinal,Node > XpmMV;
122 typedef Xpetra::MultiVector<HalfMagnitude,LocalOrdinal,GlobalOrdinal,Node > XphmMV;
123 typedef Xpetra::Matrix<HalfScalar,LocalOrdinal,GlobalOrdinal,Node> XphMat;
125 Teuchos::TimeMonitor tM(*Teuchos::TimeMonitor::getNewTimer(std::string(
"ThyraMueLuMaxwell1::initializePrec")));
128 TEUCHOS_ASSERT(Teuchos::nonnull(fwdOpSrc));
129 TEUCHOS_ASSERT(this->isCompatible(*fwdOpSrc));
130 TEUCHOS_ASSERT(prec);
133 ParameterList paramList = *paramList_;
136 const RCP<const ThyLinOpBase> fwdOp = fwdOpSrc->getOp();
137 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(fwdOp));
140 bool bIsEpetra = XpThyUtils::isEpetra(fwdOp);
141 bool bIsTpetra = XpThyUtils::isTpetra(fwdOp);
142 TEUCHOS_TEST_FOR_EXCEPT((bIsEpetra ==
true && bIsTpetra ==
true));
146 RCP<XpMat> A = XpThyUtils::toXpetra(Teuchos::rcp_const_cast<ThyLinOpBase>(fwdOp));
147 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(A));
150 const Teuchos::Ptr<DefaultPreconditioner<Scalar> > defaultPrec = Teuchos::ptr(
dynamic_cast<DefaultPreconditioner<Scalar> *
>(prec));
151 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(defaultPrec));
154 RCP<ThyLinOpBase> thyra_precOp = Teuchos::null;
155 thyra_precOp = rcp_dynamic_cast<Thyra::LinearOpBase<Scalar> >(defaultPrec->getNonconstUnspecifiedPrecOp(),
true);
160 const bool startingOver = (thyra_precOp.is_null() || !paramList.isParameter(
"Maxwell1: enable reuse") || !paramList.get<
bool>(
"Maxwell1: enable reuse"));
161 const bool useHalfPrecision = paramList.get<
bool>(
"half precision",
false) && bIsTpetra;
164 if (startingOver ==
true) {
167 std::list<std::string> convertXpetra = {
"Coordinates",
"Nullspace",
"Kn",
"D0"};
168 for (
auto it = convertXpetra.begin(); it != convertXpetra.end(); ++it)
169 Converters<Scalar,LocalOrdinal,GlobalOrdinal,Node>::replaceWithXpetra(paramList,*it);
171 std::list<std::string> sublists = {
"maxwell1: 11list",
"maxwell1: 22list"};
172 for (
auto itSublist = sublists.begin(); itSublist != sublists.end(); ++itSublist)
173 if (paramList.isSublist(*itSublist)) {
174 ParameterList& sublist = paramList.sublist(*itSublist);
175 for (
int lvlNo=0; lvlNo < 10; ++lvlNo) {
176 if (sublist.isSublist(
"level " + std::to_string(lvlNo) +
" user data")) {
177 ParameterList& lvlList = sublist.sublist(
"level " + std::to_string(lvlNo) +
" user data");
178 std::list<std::string> convertKeys;
179 for (
auto it = lvlList.begin(); it != lvlList.end(); ++it)
180 convertKeys.push_back(lvlList.name(it));
181 for (
auto it = convertKeys.begin(); it != convertKeys.end(); ++it)
182 Converters<Scalar,LocalOrdinal,GlobalOrdinal,Node>::replaceWithXpetra(lvlList,*it);
187 ParameterList& sublist = paramList.sublist(
"maxwell1: 11list");
188 if (sublist.isParameter(
"D0")) {
189 Converters<Scalar,LocalOrdinal,GlobalOrdinal,Node>::replaceWithXpetra(sublist,
"D0");
192 paramList.set<
bool>(
"Maxwell1: use as preconditioner",
true);
193 if (useHalfPrecision) {
194#if defined(MUELU_CAN_USE_MIXED_PRECISION)
197 RCP<XphMat> halfA = Xpetra::convertToHalfPrecision(A);
198 if (paramList.isType<RCP<XpmMV> >(
"Coordinates")) {
199 RCP<XpmMV> coords = paramList.get<RCP<XpmMV> >(
"Coordinates");
200 paramList.remove(
"Coordinates");
201 RCP<XphmMV> halfCoords = Xpetra::convertToHalfPrecision(coords);
202 paramList.set(
"Coordinates",halfCoords);
204 if (paramList.isType<RCP<XpMV> >(
"Nullspace")) {
205 RCP<XpMV> nullspace = paramList.get<RCP<XpMV> >(
"Nullspace");
206 paramList.remove(
"Nullspace");
207 RCP<XphMV> halfNullspace = Xpetra::convertToHalfPrecision(nullspace);
208 paramList.set(
"Nullspace",halfNullspace);
210 std::list<std::string> convertMat = {
"Kn",
"D0"};
211 for (
auto it = convertMat.begin(); it != convertMat.end(); ++it) {
212 if (paramList.isType<RCP<XpMat> >(*it)) {
213 RCP<XpMat> M = paramList.get<RCP<XpMat> >(*it);
214 paramList.remove(*it);
215 RCP<XphMat> halfM = Xpetra::convertToHalfPrecision(M);
216 paramList.set(*it,halfM);
222 xpPrecOp = rcp(
new XpHalfPrecOp(halfPrec));
224 TEUCHOS_TEST_FOR_EXCEPT(
true);
230 xpPrecOp = rcp_dynamic_cast<XpOp>(preconditioner);
235 RCP<ThyXpOp> thyXpOp = rcp_dynamic_cast<ThyXpOp>(thyra_precOp,
true);
236 RCP<XpOp> xpOp = thyXpOp->getXpetraOperator();
237#if defined(MUELU_CAN_USE_MIXED_PRECISION)
238 RCP<XpHalfPrecOp> xpHalfPrecOp = rcp_dynamic_cast<XpHalfPrecOp>(xpOp);
239 if (!xpHalfPrecOp.is_null()) {
240 RCP<MueLu::Maxwell1<HalfScalar,LocalOrdinal,GlobalOrdinal,Node> > preconditioner = rcp_dynamic_cast<MueLu::Maxwell1<HalfScalar,LocalOrdinal,GlobalOrdinal,Node>>(xpHalfPrecOp->GetHalfPrecisionOperator(),
true);
241 RCP<XphMat> halfA = Xpetra::convertToHalfPrecision(A);
242 preconditioner->resetMatrix(halfA);
243 xpPrecOp = rcp_dynamic_cast<XpOp>(preconditioner);
247 RCP<MueLu::Maxwell1<Scalar,LocalOrdinal,GlobalOrdinal,Node> > preconditioner = rcp_dynamic_cast<MueLu::Maxwell1<Scalar,LocalOrdinal,GlobalOrdinal,Node>>(xpOp,
true);
248 preconditioner->resetMatrix(A);
249 xpPrecOp = rcp_dynamic_cast<XpOp>(preconditioner);
254 RCP<const VectorSpaceBase<Scalar> > thyraRangeSpace = Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node>::toThyra(xpPrecOp->getRangeMap());
255 RCP<const VectorSpaceBase<Scalar> > thyraDomainSpace = Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node>::toThyra(xpPrecOp->getDomainMap());
258 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(thyraPrecOp));
260 defaultPrec->initializeUnspecified(thyraPrecOp);
264 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
265 void MueLuMaxwell1PreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::
266 uninitializePrec(PreconditionerBase<Scalar>* prec, RCP<
const LinearOpSourceBase<Scalar> >* fwdOp, ESupportSolveUse* supportSolveUse)
const {
267 TEUCHOS_ASSERT(prec);
270 const Teuchos::Ptr<DefaultPreconditioner<Scalar> > defaultPrec = Teuchos::ptr(
dynamic_cast<DefaultPreconditioner<Scalar> *
>(prec));
271 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(defaultPrec));
275 *fwdOp = Teuchos::null;
278 if (supportSolveUse) {
280 *supportSolveUse = Thyra::SUPPORT_SOLVE_UNSPECIFIED;
283 defaultPrec->uninitialize();
288 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
289 void MueLuMaxwell1PreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::setParameterList(RCP<ParameterList>
const& paramList) {
290 TEUCHOS_TEST_FOR_EXCEPT(Teuchos::is_null(paramList));
291 paramList_ = paramList;
294 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
295 RCP<ParameterList> MueLuMaxwell1PreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getNonconstParameterList() {
299 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
300 RCP<ParameterList> MueLuMaxwell1PreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::unsetParameterList() {
301 RCP<ParameterList> savedParamList = paramList_;
302 paramList_ = Teuchos::null;
303 return savedParamList;
306 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
307 RCP<const ParameterList> MueLuMaxwell1PreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getParameterList()
const {
311 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
312 RCP<const ParameterList> MueLuMaxwell1PreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getValidParameters()
const {
313 static RCP<const ParameterList> validPL;
315 if (Teuchos::is_null(validPL))
316 validPL = rcp(
new ParameterList());
322 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
323 std::string MueLuMaxwell1PreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::description()
const {
324 return "Thyra::MueLuMaxwell1PreconditionerFactory";
Preconditioner (wrapped as a Xpetra::Operator) for Maxwell's equations in curl-curl form.
Concrete Thyra::LinearOpBase subclass for Xpetra::Operator.
RCP< XpetraLinearOp< Scalar, LocalOrdinal, GlobalOrdinal, Node > > xpetraLinearOp(const RCP< const VectorSpaceBase< Scalar > > &rangeSpace, const RCP< const VectorSpaceBase< Scalar > > &domainSpace, const RCP< Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &xpetraOperator)
Nonmmeber constructor for XpetraLinearOp.