10 #ifndef Thyra_BlockedTriangularLinearOpWithSolveFactory_hpp
11 #define Thyra_BlockedTriangularLinearOpWithSolveFactory_hpp
13 #include "Thyra_LinearOpWithSolveBase.hpp"
14 #include "Thyra_DefaultBlockedLinearOp.hpp"
15 #include "Thyra_DefaultBlockedTriangularLinearOpWithSolve.hpp"
16 #include "Thyra_LinearOpSourceBase.hpp"
61 template <
class Scalar>
129 const std::string &precFactoryName);
138 std::string *precFactoryName);
142 virtual RCP<LinearOpWithSolveBase<Scalar> >
createOp()
const;
200 template <
class Scalar>
201 RCP<BlockedTriangularLinearOpWithSolveFactory<Scalar> >
213 template <
class Scalar>
214 RCP<BlockedTriangularLinearOpWithSolveFactory<Scalar> >
224 template <
class Scalar>
228 : lowsf_(lowsf.size())
230 for (
Ordinal i = 0; i < lowsf.size(); ++i) {
234 lowsf_[i].initialize(lowsf[i]);
238 template <
class Scalar>
242 : lowsf_(lowsf.size())
244 for (
Ordinal i = 0; i < lowsf.size(); ++i) {
248 lowsf_[i].initialize(lowsf[i]);
252 template <
class Scalar>
253 Array<RCP<LinearOpWithSolveFactoryBase<Scalar> > >
256 Array<RCP<LinearOpWithSolveFactoryBase<Scalar> > > lowsf(lowsf_.size());
257 for (
Ordinal i = 0; i < lowsf_.size(); ++i) {
258 lowsf[i] = lowsf_[i].getNonconstObj();
263 template <
class Scalar>
264 Array<RCP<const LinearOpWithSolveFactoryBase<Scalar> > >
267 Array<RCP<const LinearOpWithSolveFactoryBase<Scalar> > > lowsf(lowsf_.size());
268 for (
Ordinal i = 0; i < lowsf_.size(); ++i) {
269 lowsf[i] = lowsf_[i].getConstObj();
276 template <
class Scalar>
280 std::ostringstream oss;
282 for (
Ordinal i = 0; i < lowsf_.size(); ++i) {
284 if (!
is_null(lowsf_[i].getConstObj()))
285 oss << lowsf_[i].getConstObj()->description();
297 template <
class Scalar>
299 RCP<ParameterList>
const ¶mList)
301 for (
Ordinal i = 0; i < lowsf_.size(); ++i) {
302 lowsf_[i].getNonconstObj()->setParameterList(paramList);
306 template <
class Scalar>
310 return lowsf_[0].getNonconstObj()->getNonconstParameterList();
313 template <
class Scalar>
317 RCP<ParameterList> pl;
318 for (
Ordinal i = 0; i < lowsf_.size(); ++i) {
319 pl = lowsf_[i].getNonconstObj()->unsetParameterList();
324 template <
class Scalar>
325 RCP<const ParameterList>
328 return lowsf_[0].getConstObj()->getParameterList();
331 template <
class Scalar>
332 RCP<const ParameterList>
335 return lowsf_[0].getConstObj()->getValidParameters();
340 template <
class Scalar>
342 Scalar>::acceptsPreconditionerFactory()
const
347 template <
class Scalar>
355 true, std::logic_error,
356 "Error, we don't support a preconditioner factory!");
359 template <
class Scalar>
360 RCP<PreconditionerFactoryBase<Scalar> >
364 return Teuchos::null;
367 template <
class Scalar>
375 true, std::logic_error,
376 "Error, we don't support a preconditioner factory!");
379 template <
class Scalar>
388 template <
class Scalar>
389 RCP<LinearOpWithSolveBase<Scalar> >
392 return defaultBlockedTriangularLinearOpWithSolve<Scalar>();
395 template <
class Scalar>
403 using Teuchos::rcp_dynamic_cast;
410 for (
Ordinal i = 0; i < lowsf_.size(); ++i) {
411 lowsf_[i].getConstObj()->setOStream(this->getOStream());
412 lowsf_[i].getConstObj()->setVerbLevel(this->getVerbLevel());
417 const RCP<const PBLOB> blo =
418 rcp_dynamic_cast<
const PBLOB>(fwdOpSrc->getOp().assert_not_null());
424 DBTLOWS &btlows =
dyn_cast<DBTLOWS>(*Op);
429 const bool firstTime =
is_null(btlows.range());
434 btlows.beginBlockFill(blo->productRange(), blo->productDomain());
436 const int N = blo->productRange()->numBlocks();
437 for (
int k = 0; k < N; ++k) {
438 const RCP<const LinearOpBase<Scalar> > fwdOp_k =
439 blo->getBlock(k, k).assert_not_null();
443 btlows.setNonconstLOWSBlock(
444 k, k, linearOpWithSolve<Scalar>(*lowsf_[k].getConstObj(), fwdOp_k));
450 RCP<LinearOpWithSolveBase<Scalar> > invOp_k =
451 btlows.getNonconstLOWSBlock(k, k).assert_not_null();
452 Thyra::initializeOp<Scalar>(*lowsf_[k].getConstObj(), fwdOp_k,
459 if (firstTime) btlows.endBlockFill();
464 btlows.setBlocks(blo);
467 btlows.setOStream(this->getOStream());
468 btlows.setVerbLevel(this->getVerbLevel());
471 template <
class Scalar>
480 template <
class Scalar>
490 using Teuchos::rcp_dynamic_cast;
491 using Teuchos::rcp_implicit_cast;
494 DBTLOWS &btlowsOp =
dyn_cast<DBTLOWS>(*Op);
496 const RCP<const LinearOpBase<Scalar> > fwdOp = btlowsOp.getBlocks();
498 *fwdOpSrc = defaultLinearOpSource<Scalar>(fwdOp);
500 *fwdOpSrc = Teuchos::null;
502 if (prec) *prec = Teuchos::null;
503 if (approxFwdOpSrc) *approxFwdOpSrc = Teuchos::null;
506 template <
class Scalar>
520 template <
class Scalar>
530 true, std::logic_error,
531 "Error, we don't support an external preconditioner!");
534 template <
class Scalar>
544 true, std::logic_error,
545 "Error, we don't support an external preconditioner!");
550 template <
class Scalar>
552 Scalar>::informUpdatedVerbosityState()
const
554 for (
Ordinal i = 0; i < lowsf_.size(); ++i) {
555 lowsf_[i].getConstObj()->setVerbLevel(this->getVerbLevel());
556 lowsf_[i].getConstObj()->setOStream(this->getOStream());
562 #endif // Thyra_BlockedTriangularLinearOpWithSolveFactory_hpp
bool is_null(const boost::shared_ptr< T > &p)
virtual RCP< PreconditionerFactoryBase< Scalar > > getPreconditionerFactory() const
Returns null .
RCP< const ParameterList > getParameterList() const
void informUpdatedVerbosityState() const
Overridden from Teuchos::VerboseObjectBase.
Teuchos::ConstNonconstObjectContainer< LinearOpWithSolveFactoryBase< Scalar > > LOWSF_t
BlockedTriangularLinearOpWithSolveFactory()
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
RCP< ParameterList > getNonconstParameterList()
T_To & dyn_cast(T_From &from)
RCP< BlockedTriangularLinearOpWithSolveFactory< Scalar > > blockedTriangularLinearOpWithSolveFactory(const Array< RCP< const LinearOpWithSolveFactoryBase< Scalar > > > &lowsf)
Nonmember constructor.
RCP< ParameterList > unsetParameterList()
virtual bool isCompatible(const LinearOpSourceBase< Scalar > &fwdOpSrc) const
virtual void uninitializeOp(LinearOpWithSolveBase< Scalar > *Op, RCP< const LinearOpSourceBase< Scalar > > *fwdOpSrc, RCP< const PreconditionerBase< Scalar > > *prec, RCP< const LinearOpSourceBase< Scalar > > *approxFwdOpSrc, ESupportSolveUse *supportSolveUse) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
virtual void initializeOp(const RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, LinearOpWithSolveBase< Scalar > *Op, const ESupportSolveUse supportSolveUse) const
RCP< const ParameterList > getValidParameters() const
virtual bool acceptsPreconditionerFactory() const
returns false.
virtual RCP< LinearOpWithSolveBase< Scalar > > createOp() const
virtual void initializePreconditionedOp(const RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, const RCP< const PreconditionerBase< Scalar > > &prec, LinearOpWithSolveBase< Scalar > *Op, const ESupportSolveUse supportSolveUse) const
virtual std::string description() const
virtual void unsetPreconditionerFactory(RCP< PreconditionerFactoryBase< Scalar > > *precFactory, std::string *precFactoryName)
Throws exception.
virtual void setPreconditionerFactory(const RCP< PreconditionerFactoryBase< Scalar > > &precFactory, const std::string &precFactoryName)
Throws exception.
#define TEUCHOS_UNREACHABLE_RETURN(dummyReturnVal)
virtual void initializeAndReuseOp(const RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, LinearOpWithSolveBase< Scalar > *Op) const
Array< RCP< LinearOpWithSolveFactoryBase< Scalar > > > getUnderlyingLOWSF()
virtual bool supportsPreconditionerInputType(const EPreconditionerInputType precOpType) const
void setParameterList(RCP< ParameterList > const ¶mList)
virtual void initializeApproxPreconditionedOp(const RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, const RCP< const LinearOpSourceBase< Scalar > > &approxFwdOpSrc, LinearOpWithSolveBase< Scalar > *Op, const ESupportSolveUse supportSolveUse) const
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
RCP< BlockedTriangularLinearOpWithSolveFactory< Scalar > > blockedTriangularLinearOpWithSolveFactory(const Array< RCP< LinearOpWithSolveFactoryBase< Scalar > > > &lowsf)
Nonmember constructor.
Implicit subclass that takes a blocked triangular LOWB object and turns it into a LOWSB object...
std::string description() const