00001 #include "test_single_belos_thyra_solver.hpp"
00002
00003 #ifndef __sun
00004
00005 #include "Thyra_BelosLinearOpWithSolveFactory.hpp"
00006 #include "Thyra_LinearOpWithSolveFactoryHelpers.hpp"
00007 #include "Thyra_EpetraLinearOp.hpp"
00008 #include "Thyra_LinearOpTester.hpp"
00009 #include "Thyra_LinearOpWithSolveBase.hpp"
00010 #include "Thyra_LinearOpWithSolveTester.hpp"
00011 #include "EpetraExt_readEpetraLinearSystem.h"
00012 #include "Epetra_SerialComm.h"
00013 #include "Teuchos_ParameterList.hpp"
00014 #ifdef HAVE_BELOS_IFPACK
00015 # include "Thyra_IfpackPreconditionerFactory.hpp"
00016 #endif
00017
00018 #endif // __sun
00019
00020 bool Thyra::test_single_belos_thyra_solver(
00021 const std::string matrixFile
00022 ,const bool testTranspose
00023 ,const bool usePreconditioner
00024 ,const int numRhs
00025 ,const int numRandomVectors
00026 ,const double maxFwdError
00027 ,const double maxResid
00028 ,const double maxSolutionError
00029 ,const bool showAllTests
00030 ,const bool dumpAll
00031 ,Teuchos::ParameterList *belosLOWSFPL
00032 ,Teuchos::FancyOStream *out_arg
00033 )
00034 {
00035 using Teuchos::rcp;
00036 using Teuchos::OSTab;
00037 bool result, success = true;
00038
00039 Teuchos::RefCountPtr<Teuchos::FancyOStream> out = Teuchos::rcp(out_arg,false);
00040
00041 try {
00042
00043 #ifndef __sun
00044
00045 if(out.get()) {
00046 *out << "\n***"
00047 << "\n*** Testing Thyra::BelosLinearOpWithSolveFactory (and Thyra::BelosLinearOpWithSolve)"
00048 << "\n***\n"
00049 << "\nEchoing input options:"
00050 << "\n matrixFile = " << matrixFile
00051 << "\n testTranspose = " << testTranspose
00052 << "\n usePreconditioner = " << usePreconditioner
00053 << "\n numRhs = " << numRhs
00054 << "\n numRandomVectors = " << numRandomVectors
00055 << "\n maxFwdError = " << maxFwdError
00056 << "\n maxResid = " << maxResid
00057 << "\n showAllTests = " << showAllTests
00058 << "\n dumpAll = " << dumpAll
00059 << std::endl;
00060 }
00061
00062 if(out.get()) *out << "\nA) Reading in an epetra matrix A from the file \'"<<matrixFile<<"\' ...\n";
00063
00064 Epetra_SerialComm comm;
00065 Teuchos::RefCountPtr<Epetra_CrsMatrix> epetra_A;
00066 EpetraExt::readEpetraLinearSystem( matrixFile, comm, &epetra_A );
00067
00068 Teuchos::RefCountPtr<LinearOpBase<double> > A = Teuchos::rcp(new EpetraLinearOp(epetra_A));
00069
00070 if(out.get() && dumpAll) *out << "\ndescribe(A) =\n" << describe(*A,Teuchos::VERB_EXTREME);
00071
00072 if(out.get()) *out << "\nB) Creating a BelosLinearOpWithSolveFactory object opFactory ...\n";
00073
00074 Teuchos::RefCountPtr<LinearOpWithSolveFactoryBase<double> >
00075 lowsFactory;
00076 if(1) {
00077 Teuchos::RefCountPtr<BelosLinearOpWithSolveFactory<double> >
00078 belosLowsFactory = Teuchos::rcp(new BelosLinearOpWithSolveFactory<double>());
00079 lowsFactory = belosLowsFactory;
00080 }
00081
00082 if(out.get()) {
00083 *out << "\nlowsFactory.getValidParameters() before setting preconditioner factory:\n";
00084 lowsFactory->getValidParameters()->print(*OSTab(out).getOStream(),0,true,false);
00085 }
00086
00087 if(usePreconditioner) {
00088 #ifdef HAVE_BELOS_IFPACK
00089 if(out.get()) {
00090 *out << "\nSetting an Ifpack preconditioner factory ...\n";
00091 }
00092 lowsFactory->setPreconditionerFactory(
00093 Teuchos::rcp(new IfpackPreconditionerFactory())
00094 ,"Ifpack"
00095 );
00096 #else
00097 TEST_FOR_EXCEPT(usePreconditioner);
00098 #endif
00099 }
00100
00101 if(out.get()) {
00102 *out << "\nlowsFactory.getValidParameters() after setting preconditioner factory:\n";
00103 lowsFactory->getValidParameters()->print(*OSTab(out).getOStream(),0,true,false);
00104 *out << "\nbelosLOWSFPL before setting parameters:\n";
00105 belosLOWSFPL->print(*OSTab(out).getOStream(),0,true);
00106 }
00107
00108 lowsFactory->setParameterList(Teuchos::rcp(belosLOWSFPL,false));
00109
00110 if(out.get()) {
00111 *out << "\nbelosLOWSFPL after setting parameters:\n";
00112 belosLOWSFPL->print(*OSTab(out).getOStream(),0,true);
00113 }
00114
00115 if(out.get()) *out << "\nC) Creating a BelosLinearOpWithSolve object nsA from A ...\n";
00116
00117 Teuchos::RefCountPtr<LinearOpWithSolveBase<double> >
00118 nsA = lowsFactory->createOp();
00119
00120 Thyra::initializeOp<double>( *lowsFactory, A, &*nsA );
00121
00122 if(out.get()) *out << "\nD) Testing the LinearOpBase interface of nsA ...\n";
00123
00124 LinearOpTester<double> linearOpTester;
00125 linearOpTester.check_adjoint(testTranspose);
00126 linearOpTester.num_rhs(numRhs);
00127 linearOpTester.num_random_vectors(numRandomVectors);
00128 linearOpTester.set_all_error_tol(maxFwdError);
00129 linearOpTester.set_all_warning_tol(1e-2*maxFwdError);
00130 linearOpTester.show_all_tests(showAllTests);
00131 linearOpTester.dump_all(dumpAll);
00132 Thyra::seed_randomize<double>(0);
00133 result = linearOpTester.check(*nsA,out.get());
00134 if(!result) success = false;
00135
00136 if(out.get()) *out << "\nE) Testing the LinearOpWithSolveBase interface of nsA ...\n";
00137
00138 LinearOpWithSolveTester<double> linearOpWithSolveTester;
00139 linearOpWithSolveTester.num_rhs(numRhs);
00140 linearOpWithSolveTester.turn_off_all_tests();
00141 linearOpWithSolveTester.check_forward_default(true);
00142 linearOpWithSolveTester.check_forward_residual(true);
00143 if(testTranspose) {
00144 linearOpWithSolveTester.check_adjoint_default(true);
00145 linearOpWithSolveTester.check_adjoint_residual(true);
00146 }
00147 else {
00148 linearOpWithSolveTester.check_adjoint_default(false);
00149 linearOpWithSolveTester.check_adjoint_residual(false);
00150 }
00151 linearOpWithSolveTester.set_all_solve_tol(maxResid);
00152 linearOpWithSolveTester.set_all_slack_error_tol(maxResid);
00153 linearOpWithSolveTester.set_all_slack_warning_tol(1e+1*maxResid);
00154 linearOpWithSolveTester.forward_default_residual_error_tol(2*maxResid);
00155 linearOpWithSolveTester.forward_default_solution_error_error_tol(maxSolutionError);
00156 linearOpWithSolveTester.adjoint_default_residual_error_tol(2*maxResid);
00157 linearOpWithSolveTester.adjoint_default_solution_error_error_tol(maxSolutionError);
00158 linearOpWithSolveTester.show_all_tests(showAllTests);
00159 linearOpWithSolveTester.dump_all(dumpAll);
00160 Thyra::seed_randomize<double>(0);
00161 result = linearOpWithSolveTester.check(*nsA,out.get());
00162 if(!result) success = false;
00163
00164 if(out.get()) {
00165 *out << "\nbelosLOWSFPL after solving:\n";
00166 belosLOWSFPL->print(*OSTab(out).getOStream(),0,true);
00167 }
00168
00169 #else // __sun
00170
00171 if(out.get()) *out << "\nTest failed since is was not even compiled since __sun was defined!\n";
00172 success = false;
00173
00174 #endif // __sun
00175
00176 }
00177 catch( const std::exception &excpt ) {
00178 if(out.get()) *out << std::flush;
00179 std::cerr << "*** Caught standard exception : " << excpt.what() << std::endl;
00180 success = false;
00181 }
00182
00183 return success;
00184
00185 }