Thyra Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
EpetraExtDiagScalingTransformer_UnitTests.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_DefaultScaledAdjointLinearOp.hpp"
12 #include "Thyra_DefaultMultipliedLinearOp.hpp"
13 #include "Thyra_DefaultDiagonalLinearOp.hpp"
14 #include "Thyra_VectorStdOps.hpp"
15 #include "Thyra_TestingTools.hpp"
16 #include "Thyra_LinearOpTester.hpp"
18 #include "Thyra_EpetraLinearOp.hpp"
20 #include "EpetraExt_readEpetraLinearSystem.h"
26 #ifdef _MSC_VER
27 // This is required for operator not to be defined
28 #include <iso646.h>
29 #endif
30 
31 #ifdef HAVE_MPI
32 # include "Epetra_MpiComm.h"
33 #else
34 # include "Epetra_SerialComm.h"
35 #endif
36 
38 
39 
40 namespace {
41 
42 
43 using Teuchos::null;
44 using Teuchos::RCP;
46 using Thyra::epetraExtDiagScalingTransformer;
47 using Thyra::VectorBase;
48 using Thyra::LinearOpBase;
49 using Thyra::createMember;
50 using Thyra::LinearOpTester;
51 using Thyra::adjoint;
52 using Thyra::multiply;
53 using Thyra::diagonal;
54 
55 
56 std::string matrixFile = "";
57 
58 
60 {
62  "matrix-file", &matrixFile,
63  "Defines the Epetra_CrsMatrix to read in." );
64 }
65 
66 // helper function to excercise all different versions of B*D*G
68 buildADOperator(int scenario,const Teuchos::RCP<const Thyra::LinearOpBase<double> > & A,
69  const double vecScale, std::ostream & out)
70 {
71  // Create the implicit operator
72  double scalea=-7.0;
73  double scaled=10.0;
74 
75  RCP<const LinearOpBase<double> > M;
76  RCP<VectorBase<double> > d;
77  if(scenario<=2)
78  d = createMember(A->domain(), "d");
79  else
80  d = createMember(A->range(), "d");
81  V_S( d.ptr(), vecScale ); // ToDo: Set ton != 1.0 and generalize
82 
83  // create an operator based on requested scenario
84  // (these are the same as in EpetraExt)
85  switch(scenario) {
86  case 1:
87  M = multiply( scale(scalea,A), scale(scaled,diagonal(d)), "M" );
88  out << " Testing A*D" << std::endl;
89  break;
90  case 2:
91  M = multiply( scale(scaled,diagonal(d)), scale(scalea,A), "M" );
92  out << " Testing D*A" << std::endl;
93  break;
94  case 3:
95  M = multiply( A, scale(scaled,diagonal(d)), "M" );
96  out << " Testing adj(A)*D" << std::endl;
97  break;
98  case 4:
99  M = multiply( scale(scaled,diagonal(d)), A, "M" );
100  out << " Testing D*adj(A)" << std::endl;
101  break;
102  default:
103  TEUCHOS_ASSERT(false);
104  }
105  out << "\nM = " << *M;
106 
107  return M;
108 }
109 
110 TEUCHOS_UNIT_TEST( EpetraExtDiagScalingTransformer, isCompatible)
111 {
112  out << "\nReading linear system in Epetra format from the file \'"<<matrixFile<<"\'\n";
113 
114 #ifdef HAVE_MPI
115  Epetra_MpiComm comm(MPI_COMM_WORLD);
116 #else
117  Epetra_SerialComm comm;
118 #endif
119 
120  RCP<Epetra_CrsMatrix> epetra_A;
121  RCP<Epetra_CrsMatrix> epetra_B;
122  EpetraExt::readEpetraLinearSystem( matrixFile, comm, &epetra_A, NULL, NULL, NULL );
123  EpetraExt::readEpetraLinearSystem( matrixFile, comm, &epetra_B, NULL, NULL, NULL );
124  const RCP<const Thyra::LinearOpBase<double> > A = Thyra::epetraLinearOp(epetra_A);
125  const RCP<const Thyra::LinearOpBase<double> > B = Thyra::epetraLinearOp(epetra_B);
126 
127  RCP<VectorBase<double> > d = createMember(B->domain(), "d");
128  V_S( d.ptr(), 3.0 ); // ToDo: Set ton != 1.0 and generalize
129 
130  const RCP<EpetraExtDiagScalingTransformer> BD_transformer =
132 
133  TEST_ASSERT(BD_transformer != null);
134 
135  RCP<const LinearOpBase<double> > M;
136 
137  M = multiply(A,B);
138  TEST_ASSERT(not BD_transformer->isCompatible(*M));
139 
140  M = multiply(A,scale(3.9,diagonal(d)));
141  TEST_ASSERT(BD_transformer->isCompatible(*M));
142 
143  M = multiply(scale(3.0,diagonal(d)),scale(9.0,B));
144  TEST_ASSERT(BD_transformer->isCompatible(*M));
145 
146  M = multiply(B,scale(3.0,diagonal(d)),scale(9.0,B));
147  TEST_ASSERT(not BD_transformer->isCompatible(*M));
148 }
149 
150 TEUCHOS_UNIT_TEST( EpetraExtDiagScalingTransformer, reapply_scaling)
151 {
152  //
153  // A) Read in problem matrices
154  //
155 
156  out << "\nReading linear system in Epetra format from the file \'"<<matrixFile<<"\'\n";
157 
158 #ifdef HAVE_MPI
159  Epetra_MpiComm comm(MPI_COMM_WORLD);
160 #else
161  Epetra_SerialComm comm;
162 #endif
163  RCP<Epetra_CrsMatrix> epetra_A;
164  EpetraExt::readEpetraLinearSystem( matrixFile, comm, &epetra_A, NULL, NULL, NULL );
165 
166  //
167  // B) Create the Thyra wrapped version
168  //
169 
170  const RCP<const Thyra::LinearOpBase<double> > A = Thyra::epetraLinearOp(epetra_A);
171 
172  RCP<const Thyra::LinearOpBase<double> > M;
173 
174  // test scale on right
175  for(int scenario=1;scenario<=2;scenario++) {
176  out << "**********************************************************" << std::endl;
177  out << "**************** Starting Scenario = " << scenario << std::endl;
178  out << "**********************************************************" << std::endl;
179 
180  //
181  // D) Check compatibility
182  //
183  const RCP<EpetraExtDiagScalingTransformer> BD_transformer =
185 
186  TEST_ASSERT(BD_transformer != null);
187 
188  //
189  // E) Do the transformation (twice)
190  //
191 
192  const RCP<LinearOpBase<double> > M_explicit = BD_transformer->createOutputOp();
193  TEST_ASSERT(M_explicit != null);
194 
195  // do trans
196  M = buildADOperator(scenario,A,4.5,out);
197  BD_transformer->transform( *M, M_explicit.ptr() );
198  const RCP<const Epetra_Operator> eOp_explicit0 = Thyra::get_Epetra_Operator(*M_explicit);
199 
200  M = buildADOperator(scenario,A,7.5,out);
201  BD_transformer->transform( *M, M_explicit.ptr() );
202  const RCP<const Epetra_Operator> eOp_explicit1 = Thyra::get_Epetra_Operator(*M_explicit);
203 
204  // Epetra operator should not change!
205  TEST_EQUALITY(eOp_explicit0,eOp_explicit1);
206 
207  out << "\nM_explicit = " << *M_explicit;
208 
209  //
210  // F) Check the explicit operator
211  //
212 
213  LinearOpTester<double> M_explicit_tester;
214  M_explicit_tester.show_all_tests(true);;
215 
216  const bool result = M_explicit_tester.compare( *M, *M_explicit, Teuchos::inoutArg(out) );
217  if (!result) success = false;
218  }
219 }
220 
221 TEUCHOS_UNIT_TEST( EpetraExtDiagScalingTransformer, basic_BDG_GDB)
222 {
223 
224  //
225  // A) Read in problem matrices
226  //
227 
228  out << "\nReading linear system in Epetra format from the file \'"<<matrixFile<<"\'\n";
229 
230 #ifdef HAVE_MPI
231  Epetra_MpiComm comm(MPI_COMM_WORLD);
232 #else
233  Epetra_SerialComm comm;
234 #endif
235  RCP<Epetra_CrsMatrix> epetra_A;
236  EpetraExt::readEpetraLinearSystem( matrixFile, comm, &epetra_A, NULL, NULL, NULL );
237 
238  //
239  // B) Create the Thyra wrapped version
240  //
241 
242  const RCP<const Thyra::LinearOpBase<double> > A = Thyra::epetraLinearOp(epetra_A);
243 
244  //
245  // C) Create implicit B*D*B operator
246  //
247 
248  // build scenario=1 -> B*D*B, scenario=2-> B*D*B', scenario=3 -> B'*D*B
249  for(int scenario=1;scenario<=4;scenario++) {
250  out << "**********************************************************" << std::endl;
251  out << "**************** Starting Scenario = " << scenario << std::endl;
252  out << "**********************************************************" << std::endl;
253 
254  RCP<const Thyra::LinearOpBase<double> > M = buildADOperator(scenario,A,4.5,out);
255 
256  //
257  // D) Check compatibility
258  //
259  const RCP<EpetraExtDiagScalingTransformer> BD_transformer =
261 
262  TEST_ASSERT(BD_transformer != null);
263 
264  BD_transformer->isCompatible(*M);
265 
266  //
267  // E) Do the transformation
268  //
269 
270  const RCP<LinearOpBase<double> > M_explicit = BD_transformer->createOutputOp();
271  BD_transformer->transform( *M, M_explicit.ptr() );
272 
273  out << "\nM_explicit = " << *M_explicit;
274 
275  //
276  // F) Check the explicit operator
277  //
278 
279  LinearOpTester<double> M_explicit_tester;
280  M_explicit_tester.show_all_tests(true);
281 
282  const bool result = M_explicit_tester.compare( *M, *M_explicit, Teuchos::inoutArg(out) );
283  if (!result) success = false;
284  }
285 
286 }
287 
288 } // namespace
#define TEST_ASSERT(v1)
static CommandLineProcessor & getCLP()
TEUCHOS_UNIT_TEST(EpetraOperatorWrapper, basic)
RCP< EpetraExtDiagScalingTransformer > epetraExtDiagScalingTransformer()
Nonmember constructor.
void setOption(const char option_true[], const char option_false[], bool *option_val, const char documentation[]=NULL)
Teuchos::RCP< Epetra_Operator > get_Epetra_Operator(LinearOpBase< double > &op)
Full specialization for Scalar=double.
Transformer subclass for diagonally scaling a Epetra/Thyra operator.
#define TEST_EQUALITY(v1, v2)
#define TEUCHOS_ASSERT(assertion_test)