Thyra Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Thyra_EpetraLinearOp.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Thyra: Interfaces and Support for Abstract Numerical Algorithms
5 // Copyright (2004) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Roscoe A. Bartlett (bartlettra@ornl.gov)
38 //
39 // ***********************************************************************
40 // @HEADER
41 
42 #ifndef THYRA_EPETRA_LINEAR_OP_HPP
43 #define THYRA_EPETRA_LINEAR_OP_HPP
44 
45 #include "Thyra_LinearOpBase.hpp"
47 #include "Thyra_ScaledLinearOpBase.hpp"
48 #include "Thyra_RowStatLinearOpBase.hpp"
49 #include "Thyra_SpmdVectorSpaceBase.hpp"
50 
51 #include "Epetra_RowMatrix.h"
52 
53 
54 namespace Thyra {
55 
56 
79  : virtual public LinearOpBase<double>,
80  virtual public ScaledLinearOpBase<double>,
81  virtual public RowStatLinearOpBase<double>,
82  virtual public EpetraLinearOpBase
83 {
84 public:
85 
88 
94 
147  void initialize(
148  const RCP<Epetra_Operator> &op,
149  EOpTransp opTrans = NOTRANS,
152  const RCP<const VectorSpaceBase<double> > &range = Teuchos::null,
153  const RCP<const VectorSpaceBase<double> > &domain = Teuchos::null
154  );
155 
204  void partiallyInitialize(
205  const RCP<const VectorSpaceBase<double> > &range,
206  const RCP<const VectorSpaceBase<double> > &domain,
207  const RCP<Epetra_Operator> &op,
208  EOpTransp opTrans = NOTRANS,
211  );
212 
221  void setFullyInitialized(bool isFullyInitialized = true);
222 
230  void uninitialize(
231  RCP<Epetra_Operator> *op= NULL,
232  EOpTransp *opTrans = NULL,
233  EApplyEpetraOpAs *applyAs = NULL,
234  EAdjointEpetraOp *adjointSupport = NULL,
235  RCP<const VectorSpaceBase<double> > *range = NULL,
236  RCP<const VectorSpaceBase<double> > *domain = NULL
237  );
238 
248 
258 
261 
264 
266 
269 
272  const Ptr<RCP<Epetra_Operator> > &epetraOp,
273  const Ptr<EOpTransp> &epetraOpTransp,
274  const Ptr<EApplyEpetraOpAs> &epetraOpApplyAs,
275  const Ptr<EAdjointEpetraOp> &epetraOpAdjointSupport
276  );
278  void getEpetraOpView(
279  const Ptr<RCP<const Epetra_Operator> > &epetraOp,
280  const Ptr<EOpTransp> &epetraOpTransp,
281  const Ptr<EApplyEpetraOpAs> &epetraOpApplyAs,
282  const Ptr<EAdjointEpetraOp> &epetraOpAdjointSupport
283  ) const;
284 
286 
289 
292 
295 
298 
300 
303 
305  std::string description() const;
307  void describe(
308  FancyOStream &out,
309  const Teuchos::EVerbosityLevel verbLevel
310  ) const;
311 
313 
314 protected:
315 
318 
320  bool opSupportedImpl(EOpTransp M_trans) const;
321 
323  void applyImpl(
324  const EOpTransp M_trans,
325  const MultiVectorBase<double> &X,
326  const Ptr<MultiVectorBase<double> > &Y,
327  const double alpha,
328  const double beta
329  ) const;
330 
332 
335 
337  virtual bool supportsScaleLeftImpl() const;
338 
340  virtual bool supportsScaleRightImpl() const;
341 
343  virtual void scaleLeftImpl(const VectorBase<double> &row_scaling);
344 
346  virtual void scaleRightImpl(const VectorBase<double> &col_scaling);
347 
349 
352 
354  virtual bool rowStatIsSupportedImpl(
355  const RowStatLinearOpBaseUtils::ERowStat rowStat) const;
356 
358  virtual void getRowStatImpl(
359  const RowStatLinearOpBaseUtils::ERowStat rowStat,
360  const Ptr<VectorBase<double> > &rowStatVec) const;
361 
363 
366 
379  const RCP<Epetra_Operator> &op,
380  EOpTransp op_trans
381  ) const;
382 
394  allocateRange(
395  const RCP<Epetra_Operator> &op,
396  EOpTransp op_trans
397  ) const;
398 
400 
401 private:
402 
403  // ////////////////////////////////////
404  // Private data members
405 
409  EOpTransp opTrans_;
414 
415  // ////////////////////////////////////
416  // Private member functions
417 
418  const Epetra_Map& getRangeMap() const;
419  const Epetra_Map& getDomainMap() const;
420 
428  void computeAbsRowSum(Epetra_Vector & rowStatVec_in) const;
429 
430 }; // end class EpetraLinearOp
431 
432 
437 RCP<EpetraLinearOp> nonconstEpetraLinearOp();
438 
439 
445 partialNonconstEpetraLinearOp(
446  const RCP<const VectorSpaceBase<double> > &range,
447  const RCP<const VectorSpaceBase<double> > &domain,
448  const RCP<Epetra_Operator> &op,
449  EOpTransp opTrans = NOTRANS,
452  );
453 
454 
461 nonconstEpetraLinearOp(
462  const RCP<Epetra_Operator> &op,
463  EOpTransp opTrans = NOTRANS,
466  const RCP< const VectorSpaceBase<double> > &range = Teuchos::null,
467  const RCP< const VectorSpaceBase<double> > &domain = Teuchos::null
468  );
469 
470 
477 epetraLinearOp(
478  const RCP<const Epetra_Operator> &op,
479  EOpTransp opTrans = NOTRANS,
482  const RCP<const VectorSpaceBase<double> > &range = Teuchos::null,
483  const RCP<const VectorSpaceBase<double> > &domain = Teuchos::null
484  );
485 
486 
493 nonconstEpetraLinearOp(
494  const RCP<Epetra_Operator> &op,
495  const std::string &label,
496  EOpTransp opTrans = NOTRANS,
499  const RCP<const VectorSpaceBase<double> > &range = Teuchos::null,
500  const RCP<const VectorSpaceBase<double> > &domain = Teuchos::null
501  );
502 
503 
510 epetraLinearOp(
511  const RCP<const Epetra_Operator> &op,
512  const std::string &label,
513  EOpTransp opTrans = NOTRANS,
516  const RCP< const SpmdVectorSpaceBase<double> > &range = Teuchos::null,
517  const RCP< const SpmdVectorSpaceBase<double> > &domain = Teuchos::null
518  );
519 
520 
521 } // end namespace Thyra
522 
523 
524 #endif // THYRA_EPETRA_LINEAR_OP_HPP
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
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
Abstract base class for all LinearOpBase objects that can return an Epetra_Operator view of themselve...
RCP< const SpmdVectorSpaceBase< double > > domain_
std::string description() const
void getEpetraOpView(const Ptr< RCP< const Epetra_Operator > > &epetraOp, const Ptr< EOpTransp > &epetraOpTransp, const Ptr< EApplyEpetraOpAs > &epetraOpApplyAs, const Ptr< EAdjointEpetraOp > &epetraOpAdjointSupport) const
RCP< const SpmdVectorSpaceBase< double > > range_
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)
EAdjointEpetraOp adjointSupport_
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)
void computeAbsRowSum(Epetra_Vector &rowStatVec_in) const
Compute the absolute row sum for this matrix.
const Epetra_Map & getRangeMap() const
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
basic_FancyOStream< char > FancyOStream
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.
const Epetra_Map & getDomainMap() const
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(...)
virtual void scaleRightImpl(const VectorBase< double > &col_scaling)
RCP< Epetra_Operator > op_
RCP< const LinearOpBase< double > > clone() const
RCP< Epetra_RowMatrix > rowMatrix_
RCP< Epetra_Operator > epetra_op()