44 #include "Thyra_SpmdMultiVectorBase.hpp"
45 #include "Thyra_MultiVectorStdOps.hpp"
46 #include "Thyra_AssertOp.hpp"
53 #include "Epetra_Map.h"
54 #include "Epetra_Vector.h"
55 #include "Epetra_Operator.h"
56 #include "Epetra_CrsMatrix.h"
66 :isFullyInitialized_(false),
78 const RCP<
const VectorSpaceBase<double> > &range_in,
79 const RCP<
const VectorSpaceBase<double> > &domain_in
83 using Teuchos::rcp_dynamic_cast;
84 typedef SpmdVectorSpaceBase<double> SPMDVSB;
89 "Thyra::EpetraLinearOp::initialize(...): Error!" );
95 l_spmdRange = rcp_dynamic_cast<const SPMDVSB>(range_in,
true);
102 l_spmdDomain = rcp_dynamic_cast<const SPMDVSB>(domain_in,
true);
110 rowMatrix_ = Teuchos::rcp_dynamic_cast<Epetra_RowMatrix>(
op_);
121 const RCP<
const VectorSpaceBase<double> > &range_in,
122 const RCP<
const VectorSpaceBase<double> > &domain_in,
130 using Teuchos::rcp_dynamic_cast;
131 typedef SpmdVectorSpaceBase<double> SPMDVSB;
136 "Thyra::EpetraLinearOp::partiallyInitialize(...): Error!" );
138 "Thyra::EpetraLinearOp::partiallyInitialize(...): Error!" );
140 "Thyra::EpetraLinearOp::partiallyInitialize(...): Error!" );
144 l_spmdRange = rcp_dynamic_cast<const SPMDVSB>(range_in,
true);
146 l_spmdDomain = rcp_dynamic_cast<const SPMDVSB>(domain_in,
true);
151 rowMatrix_ = Teuchos::rcp_dynamic_cast<Epetra_RowMatrix>(
op_);
173 RCP<
const VectorSpaceBase<double> > *range_out,
174 RCP<
const VectorSpaceBase<double> > *domain_out
182 if(range_out) *range_out =
range_;
183 if(domain_out) *domain_out =
domain_;
230 const Ptr<EOpTransp> &epetraOpTransp,
231 const Ptr<EApplyEpetraOpAs> &epetraOpApplyAs,
232 const Ptr<EAdjointEpetraOp> &epetraOpAdjointSupport
244 const Ptr<EOpTransp> &epetraOpTransp,
245 const Ptr<EApplyEpetraOpAs> &epetraOpApplyAs,
246 const Ptr<EAdjointEpetraOp> &epetraOpAdjointSupport
277 return Teuchos::null;
286 std::ostringstream oss;
290 oss <<
",rangeDim="<<this->
range()->dim();
291 oss <<
",domainDim="<<this->
domain()->dim();
308 using Teuchos::rcp_dynamic_cast;
310 using Teuchos::describe;
319 <<
"rangeDim=" << this->
range()->dim()
320 <<
",domainDim=" << this->
domain()->dim()
333 csr_op = rcp_dynamic_cast<const Epetra_CrsMatrix>(
op_);
340 out <<
"op=NULL"<<
"\n";
362 const EOpTransp M_trans,
363 const MultiVectorBase<double> &X_in,
364 const Ptr<MultiVectorBase<double> > &Y_inout,
370 THYRA_FUNC_TIME_MONITOR(
"Thyra::EpetraLinearOp::euclideanApply");
372 const EOpTransp real_M_trans = real_trans(M_trans);
376 THYRA_ASSERT_LINEAR_OP_MULTIVEC_APPLY_SPACES(
377 "EpetraLinearOp::euclideanApply(...)", *
this, M_trans, X_in, Y_inout
381 Exceptions::OpNotSupported,
382 "EpetraLinearOp::apply(...): *this was informed that adjoints "
383 "are not supported when initialized."
388 const int numCols = XY_domain->dim();
401 THYRA_FUNC_TIME_MONITOR_DIFF(
402 "Thyra::EpetraLinearOp::euclideanApply: Convert MultiVectors", MultiVectors);
421 bool oldState =
op_->UseTranspose();
422 op_->SetUseTranspose(
429 THYRA_FUNC_TIME_MONITOR_DIFF(
430 "Thyra::EpetraLinearOp::euclideanApply: Apply", Apply);
434 THYRA_FUNC_TIME_MONITOR_DIFF(
435 "Thyra::EpetraLinearOp::euclideanApply: Apply(beta==0): Apply",
437 op_->Apply( *X, *Y );
440 THYRA_FUNC_TIME_MONITOR_DIFF(
441 "Thyra::EpetraLinearOp::euclideanApply: Apply(beta==0): ApplyInverse",
443 op_->ApplyInverse( *X, *Y );
452 THYRA_FUNC_TIME_MONITOR_DIFF(
453 "Thyra::EpetraLinearOp::euclideanApply: Apply(beta==0): Scale Y",
461 THYRA_FUNC_TIME_MONITOR_DIFF(
462 "Thyra::EpetraLinearOp::euclideanApply: Apply(beta!=0): Scale Y",
464 scale( beta, Y_inout );
467 THYRA_FUNC_TIME_MONITOR_DIFF(
468 "Thyra::EpetraLinearOp::euclideanApply: Apply(beta!=0): Y=0",
470 assign( Y_inout, 0.0 );
473 Epetra_MultiVector T(
op_->OperatorRangeMap(), numCols,
false);
478 THYRA_FUNC_TIME_MONITOR_DIFF(
479 "Thyra::EpetraLinearOp::euclideanApply: Apply(beta!=0): Apply",
484 THYRA_FUNC_TIME_MONITOR_DIFF(
485 "Thyra::EpetraLinearOp::euclideanApply: Apply(beta!=0): ApplyInverse",
487 op_->ApplyInverse( *X, T );
496 THYRA_FUNC_TIME_MONITOR_DIFF(
497 "Thyra::EpetraLinearOp::euclideanApply: Apply(beta!=0): Update Y",
513 op_->SetUseTranspose(oldState);
538 using Teuchos::rcpFromRef;
547 using Teuchos::rcpFromRef;
558 const RowStatLinearOpBaseUtils::ERowStat rowStat)
const
564 case RowStatLinearOpBaseUtils::ROW_STAT_INV_ROW_SUM:
565 case RowStatLinearOpBaseUtils::ROW_STAT_ROW_SUM:
575 const RowStatLinearOpBaseUtils::ERowStat rowStat,
576 const Ptr<VectorBase<double> > &rowStatVec_in
579 using Teuchos::rcpFromPtr;
583 case RowStatLinearOpBaseUtils::ROW_STAT_INV_ROW_SUM:
586 case RowStatLinearOpBaseUtils::ROW_STAT_ROW_SUM:
605 return Teuchos::rcp_dynamic_cast<
const SpmdVectorSpaceBase<double> >(
618 return Teuchos::rcp_dynamic_cast<
const SpmdVectorSpaceBase<double> >(
630 ?
op_->OperatorRangeMap() :
op_->OperatorDomainMap() );
638 ?
op_->OperatorDomainMap() :
op_->OperatorRangeMap() );
647 = Teuchos::rcp_dynamic_cast<Epetra_CrsMatrix>(
rowMatrix_);
650 Exceptions::OpNotSupported,
651 "EpetraLinearOp::computeAbsRowSum(...): wrapped matrix must be of type "
652 "Epetra_CrsMatrix for this method. Other operator types are not supported."
660 if (crsMatrix->Filled()) {
662 std::invalid_argument,
663 "EpetraLinearOp::computeAbsRowSum(...): Epetra_CrsMatrix must be filled"
668 double * xp = (
double*)x.Values();
669 if (crsMatrix->Graph().RangeMap().SameAs(x.Map()) && crsMatrix->Exporter() != 0) {
670 Epetra_Vector x_tmp(crsMatrix->RowMap());
671 x_tmp.PutScalar(0.0);
672 double * x_tmp_p = (
double*)x_tmp.Values();
673 for (i=0; i < crsMatrix->NumMyRows(); i++) {
675 double * RowValues = 0;
676 crsMatrix->ExtractMyRowView(i,NumEntries,RowValues);
677 for (j=0; j < NumEntries; j++) x_tmp_p[i] += std::abs(RowValues[j]);
681 else if (crsMatrix->Graph().RowMap().SameAs(x.Map())) {
682 for (i=0; i < crsMatrix->NumMyRows(); i++) {
684 double * RowValues = 0;
685 crsMatrix->ExtractMyRowView(i,NumEntries,RowValues);
687 for (j=0; j < NumEntries; j++) scale += std::abs(RowValues[j]);
703 Thyra::nonconstEpetraLinearOp()
710 Thyra::partialNonconstEpetraLinearOp(
711 const RCP<
const VectorSpaceBase<double> > &range,
712 const RCP<
const VectorSpaceBase<double> > &domain,
713 const RCP<Epetra_Operator> &op,
719 RCP<EpetraLinearOp> thyraEpetraOp =
Teuchos::rcp(
new EpetraLinearOp());
720 thyraEpetraOp->partiallyInitialize(
721 range, domain,op,opTrans, applyAs, adjointSupport
723 return thyraEpetraOp;
728 Thyra::nonconstEpetraLinearOp(
729 const RCP<Epetra_Operator> &op,
733 const RCP<
const VectorSpaceBase<double> > &range,
734 const RCP<
const VectorSpaceBase<double> > &domain
737 RCP<EpetraLinearOp> thyraEpetraOp =
Teuchos::rcp(
new EpetraLinearOp());
738 thyraEpetraOp->initialize(
739 op,opTrans, applyAs, adjointSupport, range, domain
741 return thyraEpetraOp;
746 Thyra::epetraLinearOp(
747 const RCP<const Epetra_Operator> &op,
751 const RCP<
const VectorSpaceBase<double> > &range,
752 const RCP<
const VectorSpaceBase<double> > &domain
755 RCP<EpetraLinearOp> thyraEpetraOp =
Teuchos::rcp(
new EpetraLinearOp());
756 thyraEpetraOp->initialize(
757 Teuchos::rcp_const_cast<Epetra_Operator>(op),
758 opTrans, applyAs, adjointSupport, range, domain
760 return thyraEpetraOp;
765 Thyra::nonconstEpetraLinearOp(
766 const RCP<Epetra_Operator> &op,
767 const std::string &label,
771 const RCP<
const VectorSpaceBase<double> > &range,
772 const RCP<
const VectorSpaceBase<double> > &domain
775 RCP<EpetraLinearOp> thyraEpetraOp =
Teuchos::rcp(
new EpetraLinearOp());
776 thyraEpetraOp->initialize(
777 op,opTrans, applyAs, adjointSupport, range, domain
779 thyraEpetraOp->setObjectLabel(label);
780 return thyraEpetraOp;
785 Thyra::epetraLinearOp(
786 const RCP<const Epetra_Operator> &op,
787 const std::string &label,
791 const RCP<
const SpmdVectorSpaceBase<double> > &range,
792 const RCP<
const SpmdVectorSpaceBase<double> > &domain
795 RCP<EpetraLinearOp> thyraEpetraOp =
Teuchos::rcp(
new EpetraLinearOp());
796 thyraEpetraOp->initialize(
797 Teuchos::rcp_const_cast<Epetra_Operator>(op),
798 opTrans, applyAs, adjointSupport, range, domain
800 thyraEpetraOp->setObjectLabel(label);
801 return thyraEpetraOp;