belos_epetra_thyra_lowsf_hb.cpp

Go to the documentation of this file.
00001 
00002 // @HEADER
00003 // ***********************************************************************
00004 //
00005 //                 Belos: Block Linear Solvers Package
00006 //                 Copyright (2004) Sandia Corporation
00007 //
00008 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
00009 // license for use of this work by or on behalf of the U.S. Government.
00010 //
00011 // This library is free software; you can redistribute it and/or modify
00012 // it under the terms of the GNU Lesser General Public License as
00013 // published by the Free Software Foundation; either version 2.1 of the
00014 // License, or (at your option) any later version.
00015 //
00016 // This library is distributed in the hope that it will be useful, but
00017 // WITHOUT ANY WARRANTY; without even the implied warranty of
00018 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00019 // Lesser General Public License for more details.
00020 //
00021 // You should have received a copy of the GNU Lesser General Public
00022 // License along with this library; if not, write to the Free Software
00023 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
00024 // USA
00025 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
00026 //
00027 // ***********************************************************************
00028 // @HEADER
00029 //
00030 // This driver reads a problem from a Harwell-Boeing (HB) file into an 
00031 // Epetra_CrsMatrix.  This matrix is then converted into a Thyra linear operator
00032 // through the Thyra-Epetra adapters.
00033 // The right-hand-side from the HB file is used instead of random vectors.
00034 // The initial guesses are all set to zero. 
00035 //
00036 //
00037 
00038 // Thyra includes
00039 #include "Thyra_BelosLinearOpWithSolveFactory.hpp"
00040 #include "Thyra_LinearOpWithSolveFactoryHelpers.hpp"
00041 #include "Thyra_EpetraLinearOp.hpp"
00042 
00043 // Epetra includes
00044 #include "Epetra_SerialComm.h"
00045 #include "EpetraExt_readEpetraLinearSystem.h"
00046 
00047 // Ifpack includes
00048 #ifdef HAVE_BELOS_IFPACK
00049   #include "Thyra_IfpackPreconditionerFactory.hpp"
00050 #endif
00051 
00052 // Teuchos includes
00053 #include "Teuchos_ParameterList.hpp"
00054 #include "Teuchos_VerboseObject.hpp"
00055 
00056 int main(int argc, char* argv[])
00057 {
00058   //
00059   // Get a default output stream from the Teuchos::VerboseObjectBase
00060   //
00061   Teuchos::RefCountPtr<Teuchos::FancyOStream>
00062     out = Teuchos::VerboseObjectBase::getDefaultOStream();
00063   //
00064   // Set the parameters for the Belos LOWS Factory and create a parameter list.
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   // Set the parameters for the Ifpack Preconditioner Factory and create parameter list
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   // Whether the linear solver succeeded.
00099   // (this will be set during the residual check at the end)
00100   bool success = true;
00101 
00102   // Number of random right-hand sides we will be solving for.
00103   int numRhs = 5;
00104 
00105   // Name of input matrix file
00106   std::string matrixFile = "orsirr1.hb";
00107 
00108   // Read in the matrix file (can be *.mtx, *.hb, etc.)
00109   Epetra_SerialComm comm;
00110   Teuchos::RefCountPtr<Epetra_CrsMatrix> epetra_A;
00111   EpetraExt::readEpetraLinearSystem( matrixFile, comm, &epetra_A );
00112   
00113   // Create a Thyra linear operator (A) using the Epetra_CrsMatrix (epetra_A).
00114   Teuchos::RefCountPtr<Thyra::LinearOpBase<double> > 
00115     A = Teuchos::rcp(new Thyra::EpetraLinearOp(epetra_A));
00116 
00117   // Get the domain space for the Thyra linear operator 
00118   Teuchos::RefCountPtr<const Thyra::VectorSpaceBase<double> > domain = A->domain();
00119 
00120   // Create the Belos LOWS factory.
00121   Teuchos::RefCountPtr<Thyra::LinearOpWithSolveFactoryBase<double> >
00122     belosLOWSFactory = Teuchos::rcp(new Thyra::BelosLinearOpWithSolveFactory<double>());
00123 
00124 #ifdef HAVE_BELOS_IFPACK
00125 
00126   // Set the preconditioner factory for the LOWS factory.
00127   belosLOWSFactory->setPreconditionerFactory(
00128                Teuchos::rcp(new Thyra::IfpackPreconditionerFactory())
00129                ,"IfpackPreconditionerFactory"
00130                );
00131 #endif
00132 
00133   // Set the parameter list to specify the behavior of the factory.
00134   belosLOWSFactory->setParameterList( belosLOWSFPL );
00135 
00136   // Set the output stream and the verbosity level (prints to std::cout by defualt)
00137   belosLOWSFactory->setVerbLevel(Teuchos::VERB_LOW);
00138 
00139   // Create a BelosLinearOpWithSolve object from the Belos LOWS factory.
00140   Teuchos::RefCountPtr<Thyra::LinearOpWithSolveBase<double> >
00141     nsA = belosLOWSFactory->createOp();
00142 
00143   // Initialize the BelosLinearOpWithSolve object with the Thyra linear operator.
00144   Thyra::initializeOp<double>( *belosLOWSFactory, A, &*nsA );
00145 
00146   // Create a right-hand side with numRhs vectors in it and randomize it.
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   // Create an initial vector with numRhs vectors in it and initialize it to zero.
00153   Teuchos::RefCountPtr< Thyra::MultiVectorBase<double> >
00154     x = Thyra::createMembers(domain, numRhs);
00155   Thyra::assign(&*x, 0.0);
00156 
00157   // Perform solve using the linear operator to get the approximate solution of Ax=b,
00158   // where b is the right-hand side and x is the left-hand side.
00159   Thyra::SolveStatus<double> solveStatus;
00160   solveStatus = Thyra::solve( *nsA, Thyra::NONCONJ_ELE, *b, &*x );
00161 
00162   // Print out status of solve.
00163   *out << "\nBelos LOWS Status: "<< solveStatus << endl;
00164 
00165   //
00166   // Compute residual and double check convergence.
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   // Compute the column norms of the right-hand side b.
00173   Thyra::norms_2( *b, &norm_b[0] );
00174 
00175   // Compute y=A*x, where x is the solution from the linear solver.
00176   A->apply( Thyra::NONCONJ_ELE, *x, &*y );
00177   
00178   // Compute A*x-b = y-b
00179   Thyra::update( -1.0, *b, &*y );
00180 
00181   // Compute the column norms of A*x-b.
00182   Thyra::norms_2( *y, &norm_res[0] );
00183 
00184   // Print out the final relative residual norms.
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 }

Generated on Thu Sep 18 12:30:20 2008 for Belos Package Browser (Single Doxygen Collection) by doxygen 1.3.9.1