46 #ifndef THYRA_MUELU_REFMAXWELL_PRECONDITIONER_FACTORY_DEF_HPP
47 #define THYRA_MUELU_REFMAXWELL_PRECONDITIONER_FACTORY_DEF_HPP
51 #if defined(HAVE_MUELU_STRATIMIKOS) && defined(HAVE_MUELU_THYRA)
62 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
63 MueLuRefMaxwellPreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::MueLuRefMaxwellPreconditionerFactory() :
64 paramList_(
rcp(new ParameterList()))
69 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
70 bool MueLuRefMaxwellPreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isCompatible(
const LinearOpSourceBase<Scalar>& fwdOpSrc)
const {
71 const RCP<const LinearOpBase<Scalar> > fwdOp = fwdOpSrc.getOp();
73 #ifdef HAVE_MUELU_TPETRA
74 if (Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isTpetra(fwdOp))
return true;
77 if (Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isBlockedOperator(fwdOp))
return true;
83 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
84 RCP<PreconditionerBase<Scalar> > MueLuRefMaxwellPreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::createPrec()
const {
88 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
89 void MueLuRefMaxwellPreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::
90 initializePrec(
const RCP<
const LinearOpSourceBase<Scalar> >& fwdOpSrc, PreconditionerBase<Scalar>* prec,
const ESupportSolveUse supportSolveUse)
const {
91 using Teuchos::rcp_dynamic_cast;
94 typedef Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> XpMap;
95 typedef Xpetra::Operator<Scalar, LocalOrdinal, GlobalOrdinal, Node> XpOp;
96 typedef Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpThyUtils;
97 typedef Xpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpCrsMat;
98 typedef Xpetra::BlockedCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpBlockedCrsMat;
99 typedef Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpMat;
100 typedef Xpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> XpMultVec;
102 typedef Thyra::LinearOpBase<Scalar> ThyLinOpBase;
103 typedef Thyra::DiagonalLinearOpBase<Scalar> ThyDiagLinOpBase;
112 ParameterList paramList = *paramList_;
115 const RCP<const ThyLinOpBase> fwdOp = fwdOpSrc->getOp();
119 bool bIsEpetra = XpThyUtils::isEpetra(fwdOp);
120 bool bIsTpetra = XpThyUtils::isTpetra(fwdOp);
121 bool bIsBlocked = XpThyUtils::isBlockedOperator(fwdOp);
126 RCP<XpMat> A = Teuchos::null;
129 Teuchos::rcp_dynamic_cast<
const Thyra::BlockedLinearOpBase<Scalar> >(fwdOp);
137 RCP<const XpCrsMat > xpetraFwdCrsMat00 = XpThyUtils::toXpetra(b00);
141 RCP<XpCrsMat> xpetraFwdCrsMatNonConst00 = Teuchos::rcp_const_cast<XpCrsMat>(xpetraFwdCrsMat00);
145 RCP<XpMat> A00 =
rcp(
new Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>(xpetraFwdCrsMatNonConst00));
148 RCP<const XpMap> rowmap00 = A00->getRowMap();
149 RCP< const Teuchos::Comm< int > > comm = rowmap00->getComm();
152 RCP<XpBlockedCrsMat> bMat =
Teuchos::rcp(
new XpBlockedCrsMat(ThyBlockedOp, comm));
158 RCP<const XpCrsMat > xpetraFwdCrsMat = XpThyUtils::toXpetra(fwdOp);
162 RCP<XpCrsMat> xpetraFwdCrsMatNonConst = Teuchos::rcp_const_cast<XpCrsMat>(xpetraFwdCrsMat);
166 A =
rcp(
new Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>(xpetraFwdCrsMatNonConst));
175 RCP<ThyLinOpBase> thyra_precOp = Teuchos::null;
176 thyra_precOp = rcp_dynamic_cast<Thyra::LinearOpBase<Scalar> >(defaultPrec->getNonconstUnspecifiedPrecOp(),
true);
179 RCP<MueLu::RefMaxwell<Scalar,LocalOrdinal,GlobalOrdinal,Node> > preconditioner = Teuchos::null;
184 const bool startingOver = (thyra_precOp.is_null() || !paramList.isParameter(
"reuse: type") || paramList.get<std::string>(
"reuse: type") ==
"none");
186 if (startingOver ==
true) {
192 paramList.set<RCP<XpMultVecDouble> >(
"Coordinates", coordinates);
196 #ifdef HAVE_MUELU_TPETRA
198 typedef Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> tV;
199 typedef Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> tMV;
200 typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> TpCrsMat;
204 RCP<tMV> tpetra_nullspace = paramList.get<RCP<tMV> >(
"Nullspace");
205 paramList.remove(
"Nullspace");
206 RCP<XpMultVec> nullspace = MueLu::TpetraMultiVector_To_XpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>(tpetra_nullspace);
207 paramList.set<RCP<XpMultVec> >(
"Nullspace", nullspace);
211 if (paramList.isParameter(
"M1")) {
214 RCP<TpCrsMat> tM1 = paramList.get<RCP<TpCrsMat> >(
"M1");
215 paramList.remove(
"M1");
216 RCP<XpCrsMat> xM1 = rcp_dynamic_cast<XpCrsMat>(tM1,
true);
217 paramList.set<RCP<XpCrsMat> >(
"M1", xM1);
220 RCP<const ThyLinOpBase> thyM1 = paramList.get<RCP<const ThyLinOpBase> >(
"M1");
221 paramList.remove(
"M1");
222 RCP<const XpCrsMat> crsM1 = XpThyUtils::toXpetra(thyM1);
225 RCP<XpCrsMat> crsM1NonConst = Teuchos::rcp_const_cast<XpCrsMat>(crsM1);
228 RCP<XpMat> M1 =
rcp(
new Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>(crsM1NonConst));
229 paramList.set<RCP<XpMat> >(
"M1", M1);
237 if (paramList.isParameter(
"Ms")) {
240 RCP<TpCrsMat> tMs = paramList.get<RCP<TpCrsMat> >(
"Ms");
241 paramList.remove(
"Ms");
242 RCP<XpCrsMat> xMs = rcp_dynamic_cast<XpCrsMat>(tMs,
true);
243 paramList.set<RCP<XpCrsMat> >(
"Ms", xMs);
246 RCP<const ThyLinOpBase> thyMs = paramList.get<RCP<const ThyLinOpBase> >(
"Ms");
247 paramList.remove(
"Ms");
248 RCP<const XpCrsMat> crsMs = XpThyUtils::toXpetra(thyMs);
251 RCP<XpCrsMat> crsMsNonConst = Teuchos::rcp_const_cast<XpCrsMat>(crsMs);
254 RCP<XpMat> Ms =
rcp(
new Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>(crsMsNonConst));
255 paramList.set<RCP<XpMat> >(
"Ms", Ms);
262 if (paramList.isParameter(
"D0")) {
265 RCP<TpCrsMat> tD0 = paramList.get<RCP<TpCrsMat> >(
"D0");
266 paramList.remove(
"D0");
267 RCP<XpCrsMat> xD0 = rcp_dynamic_cast<XpCrsMat>(tD0,
true);
268 paramList.set<RCP<XpCrsMat> >(
"D0", xD0);
271 RCP<const ThyLinOpBase> thyD0 = paramList.get<RCP<const ThyLinOpBase> >(
"D0");
272 paramList.remove(
"D0");
273 RCP<const XpCrsMat> crsD0 = XpThyUtils::toXpetra(thyD0);
276 RCP<XpCrsMat> crsD0NonConst = Teuchos::rcp_const_cast<XpCrsMat>(crsD0);
279 RCP<XpMat> D0 =
rcp(
new Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>(crsD0NonConst));
280 paramList.set<RCP<XpMat> >(
"D0", D0);
288 if (paramList.isParameter(
"M0inv")) {
291 RCP<TpCrsMat> tM0inv = paramList.get<RCP<TpCrsMat> >(
"M0inv");
292 paramList.remove(
"M0inv");
293 RCP<XpCrsMat> xM0inv = rcp_dynamic_cast<XpCrsMat>(tM0inv,
true);
294 paramList.set<RCP<XpCrsMat> >(
"M0inv", xM0inv);
297 RCP<const ThyDiagLinOpBase> thyM0inv = paramList.get<RCP<const ThyDiagLinOpBase> >(
"M0inv");
298 paramList.remove(
"M0inv");
299 RCP<const Thyra::VectorBase<Scalar> > diag = thyM0inv->getDiag();
300 RCP<const tV> tDiag = Thyra::TpetraOperatorVectorExtraction<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getConstTpetraVector(diag);
301 RCP<XpMat> M0inv = Xpetra::MatrixFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Build(Xpetra::toXpetra(tDiag));
302 paramList.set<RCP<XpMat> >(
"M0inv", M0inv);
305 RCP<const ThyLinOpBase> thyM0inv = paramList.get<RCP<const ThyLinOpBase> >(
"M0inv");
306 paramList.remove(
"M0inv");
307 RCP<const XpCrsMat> crsM0inv = XpThyUtils::toXpetra(thyM0inv);
310 RCP<XpCrsMat> crsM0invNonConst = Teuchos::rcp_const_cast<XpCrsMat>(crsM0inv);
313 RCP<XpMat> M0inv =
rcp(
new Xpetra::CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>(crsM0invNonConst));
314 paramList.set<RCP<XpMat> >(
"M0inv", M0inv);
328 paramList.set<
bool>(
"refmaxwell: use as preconditioner",
true);
334 preconditioner->resetMatrix(A);
338 RCP<ThyLinOpBase > thyraPrecOp = Teuchos::null;
339 RCP<const VectorSpaceBase<Scalar> > thyraRangeSpace = Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node>::toThyra(preconditioner->getRangeMap());
340 RCP<const VectorSpaceBase<Scalar> > thyraDomainSpace = Xpetra::ThyraUtils<Scalar,LocalOrdinal,GlobalOrdinal,Node>::toThyra(preconditioner->getDomainMap());
342 RCP<XpOp> xpOp = Teuchos::rcp_dynamic_cast<XpOp>(preconditioner);
343 thyraPrecOp = Thyra::xpetraLinearOp<Scalar, LocalOrdinal, GlobalOrdinal, Node>(thyraRangeSpace, thyraDomainSpace,xpOp);
347 defaultPrec->initializeUnspecified(thyraPrecOp);
351 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
352 void MueLuRefMaxwellPreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::
353 uninitializePrec(PreconditionerBase<Scalar>* prec, RCP<
const LinearOpSourceBase<Scalar> >* fwdOp, ESupportSolveUse* supportSolveUse)
const {
362 *fwdOp = Teuchos::null;
365 if (supportSolveUse) {
367 *supportSolveUse = Thyra::SUPPORT_SOLVE_UNSPECIFIED;
370 defaultPrec->uninitialize();
375 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
376 void MueLuRefMaxwellPreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::setParameterList(RCP<ParameterList>
const& paramList) {
378 paramList_ = paramList;
381 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
382 RCP<ParameterList> MueLuRefMaxwellPreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getNonconstParameterList() {
386 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
387 RCP<ParameterList> MueLuRefMaxwellPreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::unsetParameterList() {
388 RCP<ParameterList> savedParamList = paramList_;
389 paramList_ = Teuchos::null;
390 return savedParamList;
393 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
394 RCP<const ParameterList> MueLuRefMaxwellPreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getParameterList()
const {
398 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
399 RCP<const ParameterList> MueLuRefMaxwellPreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getValidParameters()
const {
400 static RCP<const ParameterList> validPL;
403 validPL =
rcp(
new ParameterList());
409 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
410 std::string MueLuRefMaxwellPreconditionerFactory<Scalar,LocalOrdinal,GlobalOrdinal,Node>::description()
const {
411 return "Thyra::MueLuRefMaxwellPreconditionerFactory";
415 #endif // HAVE_MUELU_STRATIMIKOS
417 #endif // ifdef THYRA_MUELU_REFMAXWELL_PRECONDITIONER_FACTORY_DEF_HPP
MueLu::DefaultLocalOrdinal LocalOrdinal
static RCP< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::magnitudeType, LocalOrdinal, GlobalOrdinal, Node > > ExtractCoordinatesFromParameterList(ParameterList ¶mList)
bool nonnull(const std::shared_ptr< T > &p)
bool is_null(const std::shared_ptr< T > &p)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Preconditioner (wrapped as a Xpetra::Operator) for Maxwell's equations in curl-curl form...
static RCP< Time > getNewTimer(const std::string &name)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Exception throws to report errors in the internal logical of the program.
#define TEUCHOS_ASSERT(assertion_test)
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)