42 #include "sillyPowerMethod.hpp"
44 #include "Thyra_TestingTools.hpp"
52 #include "Epetra_CrsMatrix.h"
53 #include "Epetra_Map.h"
71 if(crsMatrix.MyGlobalRow(0)) {
72 const int numRowNz = crsMatrix.NumGlobalEntries(0);
74 int returndNumRowNz;
double rowValues[2];
int rowIndexes[2];
75 crsMatrix.ExtractGlobalRowCopy( 0, numRowNz, returndNumRowNz, &rowValues[0], &rowIndexes[0] );
77 for(
int k = 0; k < numRowNz; ++k)
if (rowIndexes[k] == 0) rowValues[k] *= diagScale;
78 crsMatrix.ReplaceGlobalValues( 0, numRowNz, rowValues, rowIndexes );
85 int main(
int argc,
char *argv[])
88 using Teuchos::outArg;
101 MPI_Comm mpiComm = MPI_COMM_WORLD;
109 out = Teuchos::VerboseObjectBase::getDefaultOStream();
118 CommandLineProcessor clp;
119 clp.throwExceptions(
false);
120 clp.addOutputSetupOptions(
true);
123 clp.setOption(
"global-dim", &globalDim,
"Global dimension of the linear system." );
126 clp.setOption(
"dump-all",
"no-dump", &dumpAll,
"Determines if quantities are dumped or not." );
128 CommandLineProcessor::EParseCommandLineReturn
129 parse_return = clp.parse(argc,argv);
130 if (parse_return != CommandLineProcessor::PARSE_SUCCESSFUL)
135 *out <<
"\n***\n*** Running power method example using Epetra implementation\n***\n" << std::scientific;
150 *out <<
"\n(1) Constructing tridagonal Epetra matrix A of global dimension = " << globalDim <<
" ...\n";
160 RCP<Thyra::LinearOpBase<double> >
161 A = Thyra::nonconstEpetraLinearOp(A_epetra);
163 if (dumpAll) *out <<
"\nA =\n" << *A;
168 *out <<
"\n(2) Running the power method on matrix A ...\n";
170 double tolerance = 1e-3;
171 int maxNumIters = 10*globalDim;
172 result = sillyPowerMethod<double>(*A, maxNumIters, tolerance, outArg(lambda), *out);
173 if(!result) success =
false;
174 *out <<
"\n Estimate of dominate eigenvalue lambda = " << lambda << std::endl;
179 *out <<
"\n(3) Scale the diagonal of A by a factor of 10 ...\n";
185 *out <<
"\n(4) Running the power method again on matrix A ...\n";
186 result = sillyPowerMethod<double>(*A, maxNumIters, tolerance, outArg(lambda), *out);
187 if(!result) success =
false;
188 *out <<
"\n Estimate of dominate eigenvalue lambda = " << lambda << std::endl;
194 *out <<
"\nCongratulations! All of the tests checked out!\n";
196 *out <<
"\nOh no! At least one of the tests failed!\n";
198 return success ? 0 : 1;
Teuchos::RCP< Epetra_Operator > createTridiagEpetraLinearOp(const int globalDim, const double diagScale, const bool verbose, std::ostream &out)
This function generates a tridiagonal linear operator using Epetra.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
T_To & dyn_cast(T_From &from)
#define TEUCHOS_STANDARD_CATCH_STATEMENTS(VERBOSE, ERR_STREAM, SUCCESS_FLAG)
int main(int argc, char *argv[])
void scaleFirstDiagElement(const double diagScale, Thyra::LinearOpBase< double > *A)
Teuchos::RCP< Epetra_Operator > get_Epetra_Operator(LinearOpBase< double > &op)
Full specialization for Scalar=double.
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)