Thyra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Thyra_EpetraLinearOp.hpp
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 #ifndef THYRA_EPETRA_LINEAR_OP_HPP
11 #define THYRA_EPETRA_LINEAR_OP_HPP
12 
13 #include "Thyra_LinearOpBase.hpp"
14 #include "Thyra_EpetraLinearOpBase.hpp"
15 #include "Thyra_ScaledLinearOpBase.hpp"
16 #include "Thyra_RowStatLinearOpBase.hpp"
17 #include "Thyra_SpmdVectorSpaceBase.hpp"
18 
19 #include "Epetra_RowMatrix.h"
20 
21 
22 namespace Thyra {
23 
24 
47  : virtual public LinearOpBase<double>,
48  virtual public ScaledLinearOpBase<double>,
49  virtual public RowStatLinearOpBase<double>,
50  virtual public EpetraLinearOpBase
51 {
52 public:
53 
56 
62 
115  void initialize(
116  const RCP<Epetra_Operator> &op,
117  EOpTransp opTrans = NOTRANS,
120  const RCP<const VectorSpaceBase<double> > &range = Teuchos::null,
121  const RCP<const VectorSpaceBase<double> > &domain = Teuchos::null
122  );
123 
172  void partiallyInitialize(
173  const RCP<const VectorSpaceBase<double> > &range,
174  const RCP<const VectorSpaceBase<double> > &domain,
175  const RCP<Epetra_Operator> &op,
176  EOpTransp opTrans = NOTRANS,
179  );
180 
189  void setFullyInitialized(bool isFullyInitialized = true);
190 
198  void uninitialize(
199  RCP<Epetra_Operator> *op= NULL,
200  EOpTransp *opTrans = NULL,
201  EApplyEpetraOpAs *applyAs = NULL,
202  EAdjointEpetraOp *adjointSupport = NULL,
203  RCP<const VectorSpaceBase<double> > *range = NULL,
204  RCP<const VectorSpaceBase<double> > *domain = NULL
205  );
206 
216 
226 
229 
232 
234 
237 
240  const Ptr<RCP<Epetra_Operator> > &epetraOp,
241  const Ptr<EOpTransp> &epetraOpTransp,
242  const Ptr<EApplyEpetraOpAs> &epetraOpApplyAs,
243  const Ptr<EAdjointEpetraOp> &epetraOpAdjointSupport
244  );
246  void getEpetraOpView(
247  const Ptr<RCP<const Epetra_Operator> > &epetraOp,
248  const Ptr<EOpTransp> &epetraOpTransp,
249  const Ptr<EApplyEpetraOpAs> &epetraOpApplyAs,
250  const Ptr<EAdjointEpetraOp> &epetraOpAdjointSupport
251  ) const;
252 
254 
257 
260 
263 
266 
268 
271 
273  std::string description() const;
275  void describe(
276  FancyOStream &out,
277  const Teuchos::EVerbosityLevel verbLevel
278  ) const;
279 
281 
282 protected:
283 
286 
288  bool opSupportedImpl(EOpTransp M_trans) const;
289 
291  void applyImpl(
292  const EOpTransp M_trans,
293  const MultiVectorBase<double> &X,
294  const Ptr<MultiVectorBase<double> > &Y,
295  const double alpha,
296  const double beta
297  ) const;
298 
300 
303 
305  virtual bool supportsScaleLeftImpl() const;
306 
308  virtual bool supportsScaleRightImpl() const;
309 
311  virtual void scaleLeftImpl(const VectorBase<double> &row_scaling);
312 
314  virtual void scaleRightImpl(const VectorBase<double> &col_scaling);
315 
317 
320 
322  virtual bool rowStatIsSupportedImpl(
323  const RowStatLinearOpBaseUtils::ERowStat rowStat) const;
324 
326  virtual void getRowStatImpl(
327  const RowStatLinearOpBaseUtils::ERowStat rowStat,
328  const Ptr<VectorBase<double> > &rowStatVec) const;
329 
331 
334 
347  const RCP<Epetra_Operator> &op,
348  EOpTransp op_trans
349  ) const;
350 
362  allocateRange(
363  const RCP<Epetra_Operator> &op,
364  EOpTransp op_trans
365  ) const;
366 
368 
369 private:
370 
371  // ////////////////////////////////////
372  // Private data members
373 
374  bool isFullyInitialized_;
376  RCP<Epetra_RowMatrix> rowMatrix_;
377  EOpTransp opTrans_;
378  EApplyEpetraOpAs applyAs_;
379  EAdjointEpetraOp adjointSupport_;
382 
383  // ////////////////////////////////////
384  // Private member functions
385 
386  const Epetra_Map& getRangeMap() const;
387  const Epetra_Map& getDomainMap() const;
388 
396  void computeAbsRowSum(Epetra_Vector & rowStatVec_in) const;
397 
398 }; // end class EpetraLinearOp
399 
400 
405 RCP<EpetraLinearOp> nonconstEpetraLinearOp();
406 
407 
413 partialNonconstEpetraLinearOp(
414  const RCP<const VectorSpaceBase<double> > &range,
415  const RCP<const VectorSpaceBase<double> > &domain,
416  const RCP<Epetra_Operator> &op,
417  EOpTransp opTrans = NOTRANS,
420  );
421 
422 
429 nonconstEpetraLinearOp(
430  const RCP<Epetra_Operator> &op,
431  EOpTransp opTrans = NOTRANS,
434  const RCP< const VectorSpaceBase<double> > &range = Teuchos::null,
435  const RCP< const VectorSpaceBase<double> > &domain = Teuchos::null
436  );
437 
438 
445 epetraLinearOp(
446  const RCP<const Epetra_Operator> &op,
447  EOpTransp opTrans = NOTRANS,
450  const RCP<const VectorSpaceBase<double> > &range = Teuchos::null,
451  const RCP<const VectorSpaceBase<double> > &domain = Teuchos::null
452  );
453 
454 
461 nonconstEpetraLinearOp(
462  const RCP<Epetra_Operator> &op,
463  const std::string &label,
464  EOpTransp opTrans = NOTRANS,
467  const RCP<const VectorSpaceBase<double> > &range = Teuchos::null,
468  const RCP<const VectorSpaceBase<double> > &domain = Teuchos::null
469  );
470 
471 
478 epetraLinearOp(
479  const RCP<const Epetra_Operator> &op,
480  const std::string &label,
481  EOpTransp opTrans = NOTRANS,
484  const RCP< const SpmdVectorSpaceBase<double> > &range = Teuchos::null,
485  const RCP< const SpmdVectorSpaceBase<double> > &domain = Teuchos::null
486  );
487 
488 
489 } // end namespace Thyra
490 
491 
492 #endif // THYRA_EPETRA_LINEAR_OP_HPP
493 
494 #if defined(Thyra_SHOW_DEPRECATED_WARNINGS)
495 #ifdef __GNUC__
496 #warning "The ThyraEpetraAdapters package is deprecated"
497 #endif
498 #endif
499 
RCP< const SpmdVectorSpaceBase< double > > spmdDomain() const
Return a smart pointer to the SpmdVectorSpaceBase object for the domain.
virtual RCP< const SpmdVectorSpaceBase< double > > allocateDomain(const RCP< Epetra_Operator > &op, EOpTransp op_trans) const
Allocate the domain space of the operator.
Concrete LinearOpBase adapter subclass for Epetra_Operator object.
virtual bool supportsScaleRightImpl() const
EOpTransp
Enumeration for determining how a linear operator is applied. `*.
void partiallyInitialize(const RCP< const VectorSpaceBase< double > > &range, const RCP< const VectorSpaceBase< double > > &domain, const RCP< Epetra_Operator > &op, EOpTransp opTrans=NOTRANS, EApplyEpetraOpAs applyAs=EPETRA_OP_APPLY_APPLY, EAdjointEpetraOp adjointSupport=EPETRA_OP_ADJOINT_SUPPORTED)
Partially initialize.
RCP< const SpmdVectorSpaceBase< double > > spmdRange() const
Return a smart pointer to the SpmdVectorSpaceBase object for the range.
virtual bool rowStatIsSupportedImpl(const RowStatLinearOpBaseUtils::ERowStat rowStat) const
virtual bool supportsScaleLeftImpl() const
basic_FancyOStream< char > FancyOStream
Use the non-transposed operator.
Abstract base class for all LinearOpBase objects that can return an Epetra_Operator view of themselve...
std::string description() const
Abstract interface for objects that represent a space for vectors.
void getEpetraOpView(const Ptr< RCP< const Epetra_Operator > > &epetraOp, const Ptr< EOpTransp > &epetraOpTransp, const Ptr< EApplyEpetraOpAs > &epetraOpApplyAs, const Ptr< EAdjointEpetraOp > &epetraOpAdjointSupport) const
void initialize(const RCP< Epetra_Operator > &op, EOpTransp opTrans=NOTRANS, EApplyEpetraOpAs applyAs=EPETRA_OP_APPLY_APPLY, EAdjointEpetraOp adjointSupport=EPETRA_OP_ADJOINT_SUPPORTED, const RCP< const VectorSpaceBase< double > > &range=Teuchos::null, const RCP< const VectorSpaceBase< double > > &domain=Teuchos::null)
Fully initialize.
virtual void scaleLeftImpl(const VectorBase< double > &row_scaling)
Interface for a collection of column vectors called a multi-vector.
EAdjointEpetraOp
Determine if adjoints are supported on Epetra_Opeator or not.
EpetraLinearOp()
Construct to uninitialized.
void setFullyInitialized(bool isFullyInitialized=true)
Set to fully initialized.
void getNonconstEpetraOpView(const Ptr< RCP< Epetra_Operator > > &epetraOp, const Ptr< EOpTransp > &epetraOpTransp, const Ptr< EApplyEpetraOpAs > &epetraOpApplyAs, const Ptr< EAdjointEpetraOp > &epetraOpAdjointSupport)
Abstract interface for finite-dimensional dense vectors.
RCP< const VectorSpaceBase< double > > range() const
RCP< const VectorSpaceBase< double > > domain() const
void applyImpl(const EOpTransp M_trans, const MultiVectorBase< double > &X, const Ptr< MultiVectorBase< double > > &Y, const double alpha, const double beta) const
void describe(FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
Base class for all linear operators.
Applies left or right sclaing to the linear operator.
Interface for exxtracting row statistics as a VectorBase from a supporting LinearOpBase object...
EApplyEpetraOpAs
Determine how the apply an Epetra_Operator as a linear operator.
void uninitialize(RCP< Epetra_Operator > *op=NULL, EOpTransp *opTrans=NULL, EApplyEpetraOpAs *applyAs=NULL, EAdjointEpetraOp *adjointSupport=NULL, RCP< const VectorSpaceBase< double > > *range=NULL, RCP< const VectorSpaceBase< double > > *domain=NULL)
Set to uninitialized and optionally return the current state.
virtual void getRowStatImpl(const RowStatLinearOpBaseUtils::ERowStat rowStat, const Ptr< VectorBase< double > > &rowStatVec) const
virtual RCP< const SpmdVectorSpaceBase< double > > allocateRange(const RCP< Epetra_Operator > &op, EOpTransp op_trans) const
Allocate the range space of the operator.
bool opSupportedImpl(EOpTransp M_trans) const
Apply using Epetra_Operator::Apply(...)
Base abstract VectorSpaceBase class for all SPMD-based vector spaces.
virtual void scaleRightImpl(const VectorBase< double > &col_scaling)
RCP< const LinearOpBase< double > > clone() const
RCP< Epetra_Operator > epetra_op()