00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039 #include "Thyra_BelosLinearOpWithSolveFactory.hpp"
00040 #include "Thyra_LinearOpWithSolveFactoryHelpers.hpp"
00041 #include "Thyra_EpetraLinearOp.hpp"
00042
00043
00044 #include "Epetra_SerialComm.h"
00045 #include "EpetraExt_readEpetraLinearSystem.h"
00046
00047
00048 #ifdef HAVE_BELOS_IFPACK
00049 #include "Thyra_IfpackPreconditionerFactory.hpp"
00050 #endif
00051
00052
00053 #include "Teuchos_ParameterList.hpp"
00054 #include "Teuchos_VerboseObject.hpp"
00055
00056 int main(int argc, char* argv[])
00057 {
00058
00059
00060
00061 Teuchos::RefCountPtr<Teuchos::FancyOStream>
00062 out = Teuchos::VerboseObjectBase::getDefaultOStream();
00063
00064
00065
00066 int blockSize = 2;
00067 int maxIterations = 400;
00068 int maxRestarts = 25;
00069 int gmresKrylovLength = 25;
00070 int outputFrequency = 1;
00071 bool outputMaxResOnly = true;
00072 double maxResid = 1e-6;
00073
00074 Teuchos::RefCountPtr<Teuchos::ParameterList>
00075 belosLOWSFPL = Teuchos::rcp( new Teuchos::ParameterList() );
00076
00077 belosLOWSFPL->set("Solver Type","GMRES");
00078 belosLOWSFPL->set("Max Iters",int(maxIterations));
00079 belosLOWSFPL->set("Default Rel Res Norm",double(maxResid));
00080 belosLOWSFPL->set("Max Restarts",int(maxRestarts));
00081 belosLOWSFPL->set("Block Size",int(blockSize));
00082 belosLOWSFPL->sublist("GMRES").set("Max Number of Krylov Vectors",int(gmresKrylovLength*blockSize));
00083 belosLOWSFPL->sublist("GMRES").set("Variant","Standard");
00084 Teuchos::ParameterList &outputterSL = belosLOWSFPL->sublist("Outputter");
00085 outputterSL.set("Output Frequency",int(outputFrequency));
00086 outputterSL.set("Output Max Res Only",bool(outputMaxResOnly));
00087
00088 #ifdef HAVE_BELOS_IFPACK
00089
00090
00091
00092 Teuchos::ParameterList &ifpackPFSL = belosLOWSFPL->sublist("IfpackPreconditionerFactory");
00093
00094 ifpackPFSL.set("Overlap",int(2));
00095 ifpackPFSL.set("Prec Type","ILUT");
00096 #endif
00097
00098
00099
00100 bool success = true;
00101
00102
00103 int numRhs = 5;
00104
00105
00106 std::string matrixFile = "orsirr1.hb";
00107
00108
00109 Epetra_SerialComm comm;
00110 Teuchos::RefCountPtr<Epetra_CrsMatrix> epetra_A;
00111 EpetraExt::readEpetraLinearSystem( matrixFile, comm, &epetra_A );
00112
00113
00114 Teuchos::RefCountPtr<Thyra::LinearOpBase<double> >
00115 A = Teuchos::rcp(new Thyra::EpetraLinearOp(epetra_A));
00116
00117
00118 Teuchos::RefCountPtr<const Thyra::VectorSpaceBase<double> > domain = A->domain();
00119
00120
00121 Teuchos::RefCountPtr<Thyra::LinearOpWithSolveFactoryBase<double> >
00122 belosLOWSFactory = Teuchos::rcp(new Thyra::BelosLinearOpWithSolveFactory<double>());
00123
00124 #ifdef HAVE_BELOS_IFPACK
00125
00126
00127 belosLOWSFactory->setPreconditionerFactory(
00128 Teuchos::rcp(new Thyra::IfpackPreconditionerFactory())
00129 ,"IfpackPreconditionerFactory"
00130 );
00131 #endif
00132
00133
00134 belosLOWSFactory->setParameterList( belosLOWSFPL );
00135
00136
00137 belosLOWSFactory->setVerbLevel(Teuchos::VERB_LOW);
00138
00139
00140 Teuchos::RefCountPtr<Thyra::LinearOpWithSolveBase<double> >
00141 nsA = belosLOWSFactory->createOp();
00142
00143
00144 Thyra::initializeOp<double>( *belosLOWSFactory, A, &*nsA );
00145
00146
00147 Teuchos::RefCountPtr< Thyra::MultiVectorBase<double> >
00148 b = Thyra::createMembers(domain, numRhs);
00149 Thyra::seed_randomize<double>(0);
00150 Thyra::randomize(-1.0, 1.0, &*b);
00151
00152
00153 Teuchos::RefCountPtr< Thyra::MultiVectorBase<double> >
00154 x = Thyra::createMembers(domain, numRhs);
00155 Thyra::assign(&*x, 0.0);
00156
00157
00158
00159 Thyra::SolveStatus<double> solveStatus;
00160 solveStatus = Thyra::solve( *nsA, Thyra::NONCONJ_ELE, *b, &*x );
00161
00162
00163 *out << "\nBelos LOWS Status: "<< solveStatus << endl;
00164
00165
00166
00167
00168 std::vector<double> norm_b(numRhs), norm_res(numRhs);
00169 Teuchos::RefCountPtr< Thyra::MultiVectorBase<double> >
00170 y = Thyra::createMembers(domain, numRhs);
00171
00172
00173 Thyra::norms_2( *b, &norm_b[0] );
00174
00175
00176 A->apply( Thyra::NONCONJ_ELE, *x, &*y );
00177
00178
00179 Thyra::update( -1.0, *b, &*y );
00180
00181
00182 Thyra::norms_2( *y, &norm_res[0] );
00183
00184
00185 double rel_res = 0.0;
00186 *out << "Final relative residual norms" << endl;
00187 for (int i=0; i<numRhs; ++i) {
00188 rel_res = norm_res[i]/norm_b[i];
00189 if (rel_res > maxResid)
00190 success = false;
00191 *out << "RHS " << i+1 << " : "
00192 << std::setw(16) << std::right << rel_res << endl;
00193 }
00194
00195 return ( success ? 0 : 1 );
00196 }