Thyra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Thyra_DiagonalEpetraLinearOpWithSolveFactory.cpp
1 // @HEADER
2 // *****************************************************************************
3 // Thyra: Interfaces and Support for Abstract Numerical Algorithms
4 //
5 // Copyright 2004 NTESS and the Thyra contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #include "Thyra_DiagonalEpetraLinearOpWithSolveFactory.hpp"
11 #include "Thyra_DefaultDiagonalLinearOpWithSolve.hpp"
12 #include "Thyra_EpetraLinearOp.hpp"
14 #include "Teuchos_dyn_cast.hpp"
15 
16 #include "Epetra_RowMatrix.h"
17 #include "Epetra_Vector.h"
18 #include "Epetra_Map.h"
19 
20 
21 namespace Thyra {
22 
23 
25  const LinearOpSourceBase<double> &fwdOpSrc
26  ) const
27 {
28  using Teuchos::outArg;
30  fwdOp = fwdOpSrc.getOp();
31  const EpetraLinearOpBase *eFwdOp = NULL;
32  if( ! (eFwdOp = dynamic_cast<const EpetraLinearOpBase*>(&*fwdOp)) )
33  return false;
34  RCP<const Epetra_Operator> epetraFwdOp;
35  EOpTransp epetraFwdOpTransp;
36  EApplyEpetraOpAs epetraFwdOpApplyAs;
37  EAdjointEpetraOp epetraFwdOpAdjointSupport;
38  eFwdOp->getEpetraOpView(outArg(epetraFwdOp), outArg(epetraFwdOpTransp),
39  outArg(epetraFwdOpApplyAs), outArg(epetraFwdOpAdjointSupport) );
40  if( !dynamic_cast<const Epetra_RowMatrix*>(&*epetraFwdOp) )
41  return false;
42  return true;
43 }
44 
45 
48 {
50 }
51 
52 
54  const RCP<const LinearOpSourceBase<double> > &fwdOpSrc
56  ,const ESupportSolveUse /* supportSolveUse */
57  ) const
58 {
59  using Teuchos::outArg;
60  TEUCHOS_TEST_FOR_EXCEPT(Op==NULL);
61  TEUCHOS_TEST_FOR_EXCEPT(fwdOpSrc.get()==NULL);
62  TEUCHOS_TEST_FOR_EXCEPT(fwdOpSrc->getOp().get()==NULL);
63  RCP<const LinearOpBase<double> > fwdOp = fwdOpSrc->getOp();
64  const EpetraLinearOpBase &eFwdOp = Teuchos::dyn_cast<const EpetraLinearOpBase>(*fwdOp);
65  RCP<const Epetra_Operator> epetraFwdOp;
66  EOpTransp epetraFwdOpTransp;
67  EApplyEpetraOpAs epetraFwdOpApplyAs;
68  EAdjointEpetraOp epetraFwdOpAdjointSupport;
69  eFwdOp.getEpetraOpView(outArg(epetraFwdOp), outArg(epetraFwdOpTransp),
70  outArg(epetraFwdOpApplyAs), outArg(epetraFwdOpAdjointSupport) );
71  const Epetra_RowMatrix &eRMOp =
72  Teuchos::dyn_cast<const Epetra_RowMatrix>(*epetraFwdOp);
73  const Epetra_Map &map = eRMOp.OperatorDomainMap();
75  e_diag = Teuchos::rcp(new Epetra_Vector(map));
76  eRMOp.ExtractDiagonalCopy(*e_diag);
78  space = create_VectorSpace(Teuchos::rcp(new Epetra_Map(map)));
80  diag = create_Vector(e_diag,space);
81  Teuchos::set_extra_data<RCP<const LinearOpSourceBase<double> > >(
82  fwdOpSrc, "Thyra::DiagonalEpetraLinearOpWithSolveFactory::fwdOpSrc",
83  Teuchos::inOutArg(diag)
84  );
86  Teuchos::rcp_implicit_cast<const VectorBase<double> >(diag)
87  );
88  // Above cast is questionable but should be okay based on use.
89 }
90 
91 
94  ,RCP<const LinearOpSourceBase<double> > *fwdOpSrc
95  ,RCP<const PreconditionerBase<double> > *prec
96  ,RCP<const LinearOpSourceBase<double> > *approxFwdOpSrc
97  ,ESupportSolveUse * /* supportSolveUse */
98  ) const
99 {
100  using Teuchos::get_extra_data;
101  TEUCHOS_TEST_FOR_EXCEPT(Op==NULL);
105  diag = diagOp.getDiag();
106  if( fwdOpSrc ) {
107  if(diag.get()) {
108  *fwdOpSrc =
109  get_extra_data<RCP<const LinearOpSourceBase<double> > >(
110  diag,"Thyra::DiagonalEpetraLinearOpWithSolveFactory::fwdOpSrc"
111  );
112  }
113  }
114  else {
115  *fwdOpSrc = Teuchos::null;
116  }
117  if(prec) *prec = Teuchos::null; // We never keep a preconditioner!
118  if(approxFwdOpSrc) *approxFwdOpSrc = Teuchos::null; // We never keep a preconditioner!
119 }
120 
121 
122 // Overridden from ParameterListAcceptor
123 
124 
126  RCP<Teuchos::ParameterList> const& /* paramList */
127  )
128 {}
129 
130 
133 {
134  return Teuchos::null;
135 }
136 
137 
140 {
141  return Teuchos::null;
142 }
143 
144 
147 {
148  return Teuchos::null;
149 }
150 
151 
154 {
155  return Teuchos::null;
156 }
157 
158 
159 } // namespace Thyra
void setParameterList(Teuchos::RCP< Teuchos::ParameterList > const &paramList)
Base class for all linear operators that can support a high-level solve operation.
EOpTransp
Enumeration for determining how a linear operator is applied. `*.
Teuchos::RCP< const Teuchos::ParameterList > getParameterList() const
RCP< const VectorSpaceBase< double > > create_VectorSpace(const RCP< const Epetra_Map > &epetra_map)
Create an VectorSpaceBase object given an Epetra_Map object.
RCP< const VectorBase< Scalar > > getDiag() const
Abstract base class for all LinearOpBase objects that can return an Epetra_Operator view of themselve...
bool isCompatible(const LinearOpSourceBase< double > &fwdOpSrc) const
T_To & dyn_cast(T_From &from)
T * get() const
virtual const Epetra_Map & OperatorDomainMap() const =0
virtual int ExtractDiagonalCopy(Epetra_Vector &Diagonal) const =0
Simple interface class to access a precreated preconditioner as one or more linear operators objects ...
void initialize(const RCP< const VectorSpaceBase< Scalar > > &space)
Initialize given a vector space which allocates a vector internally.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< LinearOpWithSolveBase< double > > createOp() const
RCP< VectorBase< double > > create_Vector(const RCP< Epetra_Vector > &epetra_v, const RCP< const VectorSpaceBase< double > > &space=Teuchos::null)
Create a non-const VectorBase object from a non-const Epetra_Vector object.
Concrete LinearOpWithSolveBase subclass for diagonal linear operators.
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
Uninitialize a LinearOpWithSolveBase object and return its remembered forward linear operator and pot...
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...
void initializeOp(const Teuchos::RCP< const LinearOpSourceBase< double > > &fwdOpSrc, LinearOpWithSolveBase< double > *Op, const ESupportSolveUse supportSolveUse) const
Abstract interface for finite-dimensional dense vectors.
Base interface for objects that can return a linear operator.
ESupportSolveUse
Enum that specifies how a LinearOpWithSolveBase object will be used for solves after it is constructe...
EApplyEpetraOpAs
Determine how the apply an Epetra_Operator as a linear operator.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
virtual Teuchos::RCP< const LinearOpBase< Scalar > > getOp() const =0
Return a const left preconditioner linear operator if one is designed or targeted to be applied on th...
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)