Thyra Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Thyra_EpetraExtDiagScalingTransformer.cpp
Go to the documentation of this file.
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 
11 #include "Thyra_MultipliedLinearOpBase.hpp"
12 #include "Thyra_DiagonalLinearOpBase.hpp"
13 #include "Thyra_ScaledAdjointLinearOpBase.hpp"
14 #include "Thyra_EpetraLinearOp.hpp"
17 #include "Epetra_Map.h"
18 #include "Epetra_LocalMap.h"
19 #include "Epetra_SerialComm.h"
20 #include "Epetra_CrsMatrix.h"
21 #include "Teuchos_Assert.hpp"
22 #include "Teuchos_RCP.hpp"
23 
24 
25 namespace Thyra {
26 
27 
28 // Overridden from LinearOpTransformerBase
29 
30 
32  const LinearOpBase<double> &op_in) const
33 {
34  using Thyra::unwrap;
35  using Teuchos::rcp;
36  using Teuchos::rcp_dynamic_cast;
37  using Teuchos::dyn_cast;
38 
39  const MultipliedLinearOpBase<double> &multi_op =
40  dyn_cast<const MultipliedLinearOpBase<double> >(op_in);
41 
42  // this operation must only have two operands
43  if(multi_op.numOps()!=2)
44  return false;
45 
46  double scalar = 0.0;
47  EOpTransp transp = NOTRANS;
48 
49  // get properties of first operator: Transpose, scaler multiply...etc
51  unwrap( multi_op.getOp(0), &scalar, &transp, &A );
52  if(transp!=NOTRANS) return false;
53 
54  // get properties of second operator: Transpose, scaler multiply...etc
56  unwrap( multi_op.getOp(1), &scalar, &transp, &B );
57  if(transp!=NOTRANS) return false;
58 
59  // see if exactly one operator is a diagonal linear op
61  = rcp_dynamic_cast<const DiagonalLinearOpBase<double> >(A);
63  = rcp_dynamic_cast<const DiagonalLinearOpBase<double> >(B);
64 
65  if(dA==Teuchos::null && dB!=Teuchos::null)
66  return true;
67 
68  if(dA!=Teuchos::null && dB==Teuchos::null)
69  return true;
70 
71  return false;
72 }
73 
74 
77 {
78  return nonconstEpetraLinearOp();
79 }
80 
81 
83  const LinearOpBase<double> &op_in,
84  const Ptr<LinearOpBase<double> > &op_inout) const
85 {
86  using Thyra::unwrap;
87  using Teuchos::rcp;
88  using Teuchos::rcp_dynamic_cast;
89  using Teuchos::dyn_cast;
90 
91  //
92  // A) Get the component Thyra objects for M = op(A) * op(B)
93  //
94 
95  const MultipliedLinearOpBase<double> &multi_op =
96  dyn_cast<const MultipliedLinearOpBase<double> >(op_in);
97 
98  TEUCHOS_ASSERT(multi_op.numOps()==2);
99 
100  // get properties of first operator: Transpose, scaler multiply...etc
101  const RCP<const LinearOpBase<double> > op_A = multi_op.getOp(0);
102  double A_scalar = 0.0;
103  EOpTransp A_transp = NOTRANS;
105  unwrap( op_A, &A_scalar, &A_transp, &A );
106  TEUCHOS_ASSERT(A_transp==NOTRANS || A_transp==CONJTRANS); // sanity check
107 
108  // get properties of third operator: Transpose, scaler multiply...etc
109  const RCP<const LinearOpBase<double> > op_B = multi_op.getOp(1);
110  double B_scalar = 0.0;
111  EOpTransp B_transp = NOTRANS;
113  unwrap( op_B, &B_scalar, &B_transp, &B );
114  TEUCHOS_ASSERT(B_transp==NOTRANS || B_transp==CONJTRANS); // sanity check
115 
116  //
117  // B) Extract out the CrsMatrix and Diagonal objects
118  //
119 
120  // see if exactly one operator is a diagonal linear op
122  = rcp_dynamic_cast<const DiagonalLinearOpBase<double> >(A);
124  = rcp_dynamic_cast<const DiagonalLinearOpBase<double> >(B);
125 
128  if(dA==Teuchos::null) {
129  bool exactly_one_op_must_be_diagonal__dB_neq_null = dB!=Teuchos::null;
130  TEUCHOS_ASSERT(exactly_one_op_must_be_diagonal__dB_neq_null);
131 
132  // convert second operator to an Epetra_CrsMatrix
133  epetra_A = rcp_dynamic_cast<const Epetra_CrsMatrix>(get_Epetra_Operator(*A), true);
134  }
135  else if(dB==Teuchos::null) {
136  bool exactly_one_op_must_be_diagonal__dA_neq_null = dA!=Teuchos::null;
137  TEUCHOS_ASSERT(exactly_one_op_must_be_diagonal__dA_neq_null);
138 
139  // convert second operator to an Epetra_CrsMatrix
140  epetra_B = rcp_dynamic_cast<const Epetra_CrsMatrix>(get_Epetra_Operator(*B), true);
141  }
142  else {
143  bool exactly_one_op_must_be_diagonal=false;
144  TEUCHOS_ASSERT(exactly_one_op_must_be_diagonal);
145  }
146 
147  TEUCHOS_ASSERT( A_transp == NOTRANS ); // ToDo: Handle the transpose
148  TEUCHOS_ASSERT( B_transp == NOTRANS ); // ToDo: Handle the transpose
149 
150 
151  // get or allocate an object to use as the destination
152  EpetraLinearOp & thyra_epetra_op_inout = dyn_cast<EpetraLinearOp>(*op_inout);
153  RCP<Epetra_CrsMatrix> epetra_op =
154  rcp_dynamic_cast<Epetra_CrsMatrix>(thyra_epetra_op_inout.epetra_op());
155 
156  bool rightScale = dB!=Teuchos::null;
157 
158  if(rightScale) {
159  RCP<const Epetra_Vector> v = get_Epetra_Vector(epetra_A->OperatorDomainMap(),dB->getDiag());
160 
161  // if needed allocate a new operator, otherwise use old one assuming
162  // constant sparsity
163  if(epetra_op==Teuchos::null)
164  epetra_op = Teuchos::rcp(new Epetra_CrsMatrix(*epetra_A));
165  else
166  *epetra_op = *epetra_A;
167  epetra_op->RightScale(*v);
168  }
169  else {
170  RCP<const Epetra_Vector> v = get_Epetra_Vector(epetra_B->OperatorRangeMap(),dA->getDiag());
171 
172  // if needed allocate a new operator, otherwise use old one assuming
173  // constant sparsity
174  if(epetra_op==Teuchos::null)
175  epetra_op = Teuchos::rcp(new Epetra_CrsMatrix(*epetra_B));
176  else
177  *epetra_op = *epetra_B;
178  epetra_op->LeftScale(*v);
179  }
180 
181  epetra_op->Scale(A_scalar*B_scalar);
182 
183  // set output operator to use newly create epetra_op
184  thyra_epetra_op_inout.initialize(epetra_op);
185 }
186 
187 
188 } // namespace Thyra
Concrete LinearOpBase adapter subclass for Epetra_Operator object.
virtual bool isCompatible(const LinearOpBase< double > &op_in) 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.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
T_To & dyn_cast(T_From &from)
RCP< Epetra_Vector > get_Epetra_Vector(const Epetra_Map &map, const RCP< VectorBase< double > > &v)
Get a non-const Epetra_Vector view from a non-const VectorBase object if possible.
virtual RCP< LinearOpBase< double > > createOutputOp() const
Teuchos::RCP< Epetra_Operator > get_Epetra_Operator(LinearOpBase< double > &op)
Full specialization for Scalar=double.
virtual void transform(const LinearOpBase< double > &op_in, const Ptr< LinearOpBase< double > > &op_inout) const
#define TEUCHOS_ASSERT(assertion_test)
RCP< Epetra_Operator > epetra_op()