test_tfqmr_hb.cpp

Go to the documentation of this file.
00001 // @HEADER
00002 // ***********************************************************************
00003 //
00004 //                 Belos: Block Linear Solvers Package
00005 //                 Copyright (2004) Sandia Corporation
00006 //
00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
00008 // license for use of this work by or on behalf of the U.S. Government.
00009 //
00010 // This library is free software; you can redistribute it and/or modify
00011 // it under the terms of the GNU Lesser General Public License as
00012 // published by the Free Software Foundation; either version 2.1 of the
00013 // License, or (at your option) any later version.
00014 //
00015 // This library is distributed in the hope that it will be useful, but
00016 // WITHOUT ANY WARRANTY; without even the implied warranty of
00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00018 // Lesser General Public License for more details.
00019 //
00020 // You should have received a copy of the GNU Lesser General Public
00021 // License along with this library; if not, write to the Free Software
00022 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
00023 // USA
00024 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
00025 //
00026 // ***********************************************************************
00027 // @HEADER
00028 //
00029 // This driver reads a problem from a Harwell-Boeing (HB) file.
00030 // The right-hand-side from the problem is being used instead of multiple
00031 // random right-hand-sides.  The initial guesses are all set to zero. 
00032 //
00033 // NOTE: No preconditioner is used in this case. 
00034 //
00035 #include "BelosConfigDefs.hpp"
00036 #include "BelosLinearProblem.hpp"
00037 #include "BelosOutputManager.hpp"
00038 #include "BelosStatusTestMaxIters.hpp"
00039 #include "BelosStatusTestMaxRestarts.hpp"
00040 #include "BelosStatusTestResNorm.hpp"
00041 #include "BelosStatusTestOutputter.hpp"
00042 #include "BelosStatusTestCombo.hpp"
00043 #include "BelosEpetraAdapter.hpp"
00044 #include "BelosTFQMR.hpp"
00045 #include "createEpetraProblem.hpp"
00046 #include "Epetra_CrsMatrix.h"
00047 #include "Teuchos_CommandLineProcessor.hpp"
00048 #include "Teuchos_ParameterList.hpp"
00049 #include "Epetra_Map.h"
00050 
00051 int main(int argc, char *argv[]) {
00052   //
00053 #ifdef EPETRA_MPI 
00054   MPI_Init(&argc,&argv);
00055   Belos::MPIFinalize mpiFinalize; // Will call finalize with *any* return
00056 #endif  
00057   //
00058   typedef double                            ST;
00059   typedef Teuchos::ScalarTraits<ST>        SCT;
00060   typedef SCT::magnitudeType                MT;
00061   typedef Epetra_MultiVector                MV;
00062   typedef Epetra_Operator                   OP;
00063   typedef Belos::MultiVecTraits<ST,MV>     MVT;
00064   typedef Belos::OperatorTraits<ST,MV,OP>  OPT;
00065 
00066   using Teuchos::ParameterList;
00067   using Teuchos::RefCountPtr;
00068   using Teuchos::rcp;
00069 
00070   bool verbose = false, proc_verbose = false;
00071   int frequency = -1;
00072   int numrhs = 1;
00073   int numiters = 10; 
00074   std::string filename("orsirr1.hb");
00075   MT tol = 1.0e-5;  // relative residual tolerance
00076 
00077   Teuchos::CommandLineProcessor cmdp(false,true);
00078   cmdp.setOption("verbose","quiet",&verbose,"Print messages and results.");
00079   cmdp.setOption("frequency",&frequency,"Solvers frequency for printing residuals (#iters).");
00080   cmdp.setOption("filename",&filename,"Filename for Harwell-Boeing test matrix.");
00081   cmdp.setOption("tol",&tol,"Relative residual tolerance used by TFQMR solver.");
00082   cmdp.setOption("num-iters",&numiters,"Number of iterations allowed for TFQMR solver.");
00083   if (cmdp.parse(argc,argv) != Teuchos::CommandLineProcessor::PARSE_SUCCESSFUL) {
00084     return -1;
00085   }
00086   if (!verbose)
00087     frequency = -1;  // reset frequency if test is not verbose
00088   //
00089   // Get the problem
00090   //
00091   int MyPID;
00092   RefCountPtr<Epetra_CrsMatrix> A;
00093   RefCountPtr<Epetra_MultiVector> B, X;
00094   int return_val =Belos::createEpetraProblem(filename,NULL,&A,&B,&X,&MyPID);
00095   if(return_val != 0) return return_val;
00096   const Epetra_Map &Map = A->RowMap();
00097   proc_verbose = verbose && (MyPID==0);  /* Only print on the zero processor */
00098   //
00099   // ********Other information used by block solver***********
00100   // *****************(can be user specified)******************
00101   //
00102   const int NumGlobalElements = B->GlobalLength();
00103   //
00104   // Construct an unpreconditioned linear problem instance.
00105   //
00106   Belos::LinearProblem<double,MV,OP>
00107   My_LP( A, X, B );
00108   My_LP.SetBlockSize( 1 );
00109   //
00110   // *******************************************************************
00111   // *************Start the block Tfqmr iteration*************************
00112   // *******************************************************************
00113   //
00114   Belos::OutputManager<double> My_OM( MyPID );
00115   if (verbose)
00116     My_OM.SetVerbosity( Belos::Errors + Belos::Warnings
00117       + Belos::TimingDetails + Belos::FinalSummary );
00118   
00119   typedef Belos::StatusTestCombo<double,MV,OP> StatusTestCombo_t;
00120   Belos::StatusTestMaxIters<double,MV,OP> test1( numiters );
00121   Belos::StatusTestResNorm<double,MV,OP> test2( tol );
00122   Belos::StatusTestOutputter<ST,MV,OP> test3( frequency, false );
00123   test3.set_resNormStatusTest( rcp(&test2,false) );
00124   test3.set_outputManager( rcp(&My_OM,false) );    
00125   StatusTestCombo_t My_Test( StatusTestCombo_t::OR, test1, test3 );
00126   
00127   Belos::TFQMR<double,MV,OP>
00128     MyTFQMR( rcp(&My_LP,false), rcp(&My_Test,false), rcp(&My_OM,false) );
00129   
00130   //
00131   // **********Print out information about problem*******************
00132   //
00133   if (proc_verbose) {
00134     cout << endl << endl;
00135     cout << "Dimension of matrix: " << NumGlobalElements << endl;
00136     cout << "Number of right-hand sides: " << numrhs << endl;
00137     cout << "Max number of TFQMR iterations: " << numiters << endl; 
00138     cout << "Relative residual tolerance: " << tol << endl;
00139     cout << endl;
00140   }
00141   //
00142   // Perform solve
00143   //
00144   MyTFQMR.Solve();
00145 
00146   //
00147   // Compute actual residuals.
00148   //
00149   std::vector<double> actual_resids( numrhs );
00150   std::vector<double> rhs_norm( numrhs );
00151   Epetra_MultiVector resid(Map, numrhs);
00152   OPT::Apply( *A, *X, resid );
00153   MVT::MvAddMv( -1.0, resid, 1.0, *B, resid ); 
00154   MVT::MvNorm( resid, &actual_resids );
00155   MVT::MvNorm( *B, &rhs_norm );
00156   if (proc_verbose) {
00157     cout<< "---------- Actual Residuals (normalized) ----------"<<endl<<endl;
00158     for ( int i=0; i<numrhs; i++) {
00159       cout<<"Problem "<<i<<" : \t"<< actual_resids[i]/rhs_norm[i] <<endl;
00160     }
00161   }
00162 
00163   if (My_Test.GetStatus()!=Belos::Converged) {
00164   if (proc_verbose)
00165           cout << "End Result: TEST FAILED" << endl;  
00166   return -1;
00167   }
00168   //
00169   // Default return value
00170   //
00171   if (proc_verbose)
00172     cout << "End Result: TEST PASSED" << endl;
00173   return 0;
00174   //
00175 } // end test_bl_tfqmr_hb.cpp

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