43 #include "Thyra_DefaultDiagonalLinearOpWithSolve.hpp" 46 #include "Teuchos_dyn_cast.hpp" 48 #include "Epetra_RowMatrix.h" 49 #include "Epetra_Vector.h" 50 #include "Epetra_Map.h" 57 const LinearOpSourceBase<double> &fwdOpSrc
60 using Teuchos::outArg;
61 RCP<const LinearOpBase<double> >
62 fwdOp = fwdOpSrc.getOp();
64 if( ! (eFwdOp = dynamic_cast<const EpetraLinearOpBase*>(&*fwdOp)) )
66 RCP<const Epetra_Operator> epetraFwdOp;
67 EOpTransp epetraFwdOpTransp;
71 outArg(epetraFwdOpApplyAs), outArg(epetraFwdOpAdjointSupport) );
72 if( !dynamic_cast<const Epetra_RowMatrix*>(&*epetraFwdOp) )
78 RCP<LinearOpWithSolveBase<double> >
81 return Teuchos::rcp(
new DefaultDiagonalLinearOpWithSolve<double>());
86 const RCP<
const LinearOpSourceBase<double> > &fwdOpSrc
87 ,LinearOpWithSolveBase<double> *Op
88 ,
const ESupportSolveUse supportSolveUse
91 using Teuchos::outArg;
92 TEUCHOS_TEST_FOR_EXCEPT(Op==NULL);
93 TEUCHOS_TEST_FOR_EXCEPT(fwdOpSrc.get()==NULL);
94 TEUCHOS_TEST_FOR_EXCEPT(fwdOpSrc->getOp().get()==NULL);
95 RCP<const LinearOpBase<double> > fwdOp = fwdOpSrc->getOp();
97 RCP<const Epetra_Operator> epetraFwdOp;
98 EOpTransp epetraFwdOpTransp;
102 outArg(epetraFwdOpApplyAs), outArg(epetraFwdOpAdjointSupport) );
103 const Epetra_RowMatrix &eRMOp =
104 Teuchos::dyn_cast<
const Epetra_RowMatrix>(*epetraFwdOp);
105 const Epetra_Map &map = eRMOp.OperatorDomainMap();
107 e_diag = Teuchos::rcp(
new Epetra_Vector(map));
108 eRMOp.ExtractDiagonalCopy(*e_diag);
109 RCP< const VectorSpaceBase<double> >
111 RCP< const VectorBase<double> >
113 Teuchos::set_extra_data<RCP<const LinearOpSourceBase<double> > >(
114 fwdOpSrc,
"Thyra::DiagonalEpetraLinearOpWithSolveFactory::fwdOpSrc",
115 Teuchos::inOutArg(diag)
117 Teuchos::dyn_cast< DefaultDiagonalLinearOpWithSolve<double> >(*Op).initialize(
118 Teuchos::rcp_implicit_cast<
const VectorBase<double> >(diag)
125 LinearOpWithSolveBase<double> *Op
126 ,RCP<
const LinearOpSourceBase<double> > *fwdOpSrc
127 ,RCP<
const PreconditionerBase<double> > *prec
128 ,RCP<
const LinearOpSourceBase<double> > *approxFwdOpSrc
129 ,ESupportSolveUse *supportSolveUse
132 using Teuchos::get_extra_data;
133 TEUCHOS_TEST_FOR_EXCEPT(Op==NULL);
134 DefaultDiagonalLinearOpWithSolve<double>
135 &diagOp = Teuchos::dyn_cast<DefaultDiagonalLinearOpWithSolve<double> >(*Op);
136 RCP< const VectorBase<double> >
137 diag = diagOp.getDiag();
141 get_extra_data<RCP<const LinearOpSourceBase<double> > >(
142 diag,
"Thyra::DiagonalEpetraLinearOpWithSolveFactory::fwdOpSrc" 147 *fwdOpSrc = Teuchos::null;
149 if(prec) *prec = Teuchos::null;
150 if(approxFwdOpSrc) *approxFwdOpSrc = Teuchos::null;
158 RCP<Teuchos::ParameterList>
const& paramList
163 RCP<Teuchos::ParameterList>
166 return Teuchos::null;
170 RCP<Teuchos::ParameterList>
173 return Teuchos::null;
177 RCP<const Teuchos::ParameterList>
180 return Teuchos::null;
184 RCP<const Teuchos::ParameterList>
187 return Teuchos::null;
void setParameterList(Teuchos::RCP< Teuchos::ParameterList > const ¶mList)
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
RCP< const VectorSpaceBase< double > > create_VectorSpace(const RCP< const Epetra_Map > &epetra_map)
Create an VectorSpaceBase object given an Epetra_Map object.
Abstract base class for all LinearOpBase objects that can return an Epetra_Operator view of themselve...
RCP< VectorBase< double > > create_Vector(const RCP< Epetra_Vector > &epetra_v, const RCP< const VectorSpaceBase< double > > &space)
Create a non-const VectorBase object from a non-const Epetra_Vector object.
Teuchos::RCP< const Teuchos::ParameterList > getParameterList() const
EAdjointEpetraOp
Determine if adjoints are supported on Epetra_Opeator or not.
virtual void getEpetraOpView(const Ptr< RCP< const Epetra_Operator > > &epetraOp, const Ptr< EOpTransp > &epetraOpTransp, const Ptr< EApplyEpetraOpAs > &epetraOpApplyAs, const Ptr< EAdjointEpetraOp > &epetraOpAdjointSupport) const =0
Return a smart pointer to a const Epetra_Operator view of this object and how the object is applied t...
Teuchos::RCP< Teuchos::ParameterList > getNonconstParameterList()
void initializeOp(const Teuchos::RCP< const LinearOpSourceBase< double > > &fwdOpSrc, LinearOpWithSolveBase< double > *Op, const ESupportSolveUse supportSolveUse) const
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList()
EApplyEpetraOpAs
Determine how the apply an Epetra_Operator as a linear operator.
bool isCompatible(const LinearOpSourceBase< double > &fwdOpSrc) const
void uninitializeOp(LinearOpWithSolveBase< double > *Op, Teuchos::RCP< const LinearOpSourceBase< double > > *fwdOpSrc, Teuchos::RCP< const PreconditionerBase< double > > *prec, Teuchos::RCP< const LinearOpSourceBase< double > > *approxFwdOpSrc, ESupportSolveUse *supportSolveUse) const
Teuchos::RCP< LinearOpWithSolveBase< double > > createOp() const