test_bl_gmres_complex_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 HB file is used instead of random vectors.
00031 // 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 "BelosBlockGmres.hpp"
00044 #include "Teuchos_CommandLineProcessor.hpp"
00045 
00046 // I/O for Harwell-Boeing files
00047 #ifdef HAVE_BELOS_TRIUTILS
00048 #include "iohb.h"
00049 #endif
00050 
00051 #include "MyMultiVec.hpp"
00052 #include "MyBetterOperator.hpp"
00053 #include "MyOperator.hpp"
00054 
00055 using namespace Teuchos;
00056 
00057 int main(int argc, char *argv[]) {
00058   //
00059 #ifdef HAVE_COMPLEX
00060   typedef std::complex<double> ST;
00061 #elif HAVE_COMPLEX_H
00062   typedef ::complex<double> ST;
00063 #else
00064   cout << "Not compiled with complex support." << endl;
00065   cout << "End Result: TEST FAILED" << endl;
00066   return -1;
00067   }
00068 #endif
00069   typedef ScalarTraits<ST>                 SCT;
00070   typedef SCT::magnitudeType                MT;
00071   typedef Belos::MultiVec<ST>               MV;
00072   typedef Belos::Operator<ST>               OP;
00073   typedef Belos::MultiVecTraits<ST,MV>     MVT;
00074   typedef Belos::OperatorTraits<ST,MV,OP>  OPT;
00075   ST one  = SCT::one();
00076   ST zero = SCT::zero();  
00077 
00078   int info = 0;
00079   int MyPID = 0;
00080   bool norm_failure = false;
00081 
00082 #ifdef HAVE_MPI 
00083   // Initialize MPI 
00084   MPI_Init(&argc,&argv);  
00085   MPI_Comm_rank(MPI_COMM_WORLD, &MyPID);
00086 #endif
00087   //
00088   using Teuchos::RefCountPtr;
00089   using Teuchos::rcp;
00090 
00091   bool verbose = false, proc_verbose = false;
00092   int frequency = -1;  // how often residuals are printed by solver
00093   int blocksize = 1;
00094   int numrhs = 1;
00095   int numrestarts = 15;
00096   int length = 50;
00097   std::string filename("mhd1280b.cua");
00098   MT tol = 1.0e-5;  // relative residual tolerance
00099 
00100   CommandLineProcessor cmdp(false,true);
00101   cmdp.setOption("verbose","quiet",&verbose,"Print messages and results.");
00102   cmdp.setOption("frequency",&frequency,"Solvers frequency for printing residuals (#iters).");
00103   cmdp.setOption("filename",&filename,"Filename for Harwell-Boeing test matrix.");
00104   cmdp.setOption("tol",&tol,"Relative residual tolerance used by GMRES solver.");
00105   cmdp.setOption("num-rhs",&numrhs,"Number of right-hand sides to be solved for.");
00106   cmdp.setOption("num-restarts",&numrestarts,"Number of restarts allowed for the GMRES solver.");
00107   cmdp.setOption("block-size",&blocksize,"Block size used by GMRES.");
00108   cmdp.setOption("subspace-length",&length,"Maximum dimension of block-subspace used by GMRES solver.");
00109   if (cmdp.parse(argc,argv) != CommandLineProcessor::PARSE_SUCCESSFUL) {
00110 #ifdef HAVE_MPI
00111     MPI_Finalize();
00112 #endif
00113     return -1;
00114   }
00115 
00116   proc_verbose = verbose && (MyPID==0);  /* Only print on the zero processor */
00117   if (proc_verbose) {
00118     cout << Belos::Belos_Version() << endl << endl;
00119   }
00120   if (!verbose)
00121     frequency = -1;  // reset frequency if test is not verbose
00122       
00123 
00124 #ifndef HAVE_BELOS_TRIUTILS
00125    cout << "This test requires Triutils. Please configure with --enable-triutils." << endl;
00126 #ifdef HAVE_MPI
00127   MPI_Finalize() ;
00128 #endif
00129   if (MyPID==0) {
00130     cout << "End Result: TEST FAILED" << endl;  
00131   }
00132   return -1;
00133 #endif
00134 
00135   // Get the data from the HB file
00136   int dim,dim2,nnz;
00137   MT *dvals;
00138   int *colptr,*rowind;
00139   ST *cvals;
00140   nnz = -1;
00141   info = readHB_newmat_double(filename.c_str(),&dim,&dim2,&nnz,
00142                               &colptr,&rowind,&dvals);
00143   if (info == 0 || nnz < 0) {
00144     if (MyPID==0) {
00145       cout << "Error reading '" << filename << "'" << endl;
00146       cout << "End Result: TEST FAILED" << endl;
00147     }
00148 #ifdef HAVE_MPI
00149     MPI_Finalize();
00150 #endif
00151     return -1;
00152   }
00153   // Convert interleaved doubles to complex values
00154   cvals = new ST[nnz];
00155   for (int ii=0; ii<nnz; ii++) {
00156     cvals[ii] = ST(dvals[ii*2],dvals[ii*2+1]);
00157   }
00158   // Build the problem matrix
00159   RefCountPtr< MyBetterOperator<ST> > A 
00160     = rcp( new MyBetterOperator<ST>(dim,colptr,nnz,rowind,cvals) );
00161   //
00162   // ********Other information used by block solver***********
00163   // *****************(can be user specified)******************
00164   //
00165   int maxits = dim/blocksize; // maximum number of iterations to run
00166   //
00167   RefCountPtr<ParameterList> My_PL = rcp( new ParameterList() );
00168   My_PL->set( "Length", length );  // Maximum number of blocks in Krylov factorization
00169   //
00170   // Construct the right-hand side and solution multivectors.
00171   // NOTE:  The right-hand side will be constructed such that the solution is
00172   // a vectors of one.
00173   //
00174   RefCountPtr<MyMultiVec<ST> > soln = rcp( new MyMultiVec<ST>(dim,numrhs) );
00175   RefCountPtr<MyMultiVec<ST> > rhs = rcp( new MyMultiVec<ST>(dim,numrhs) );
00176   MVT::MvRandom( *soln );
00177   OPT::Apply( *A, *soln, *rhs );
00178   MVT::MvInit( *soln, zero );
00179   //
00180   //  Construct an unpreconditioned linear problem instance.
00181   //
00182   RefCountPtr<Belos::LinearProblem<ST,MV,OP> > My_LP = 
00183     rcp( new Belos::LinearProblem<ST,MV,OP>( A, soln, rhs ) );
00184   My_LP->SetBlockSize( blocksize );
00185   //
00186   // *******************************************************************
00187   // *************Start the block Gmres iteration***********************
00188   // *******************************************************************
00189   //
00190   RefCountPtr<Belos::OutputManager<ST> > My_OM = 
00191     rcp( new Belos::OutputManager<ST>( MyPID ) );
00192   if (verbose)
00193     My_OM->SetVerbosity( Belos::Errors + Belos::Warnings 
00194        + Belos::TimingDetails + Belos::IterationDetails 
00195        + Belos::FinalSummary ); 
00196   //
00197   typedef Belos::StatusTestCombo<ST,MV,OP> StatusTestCombo_t;
00198   Belos::StatusTestMaxIters<ST,MV,OP> test1( maxits );
00199   Belos::StatusTestMaxRestarts<ST,MV,OP> test2( numrestarts );
00200   Belos::StatusTestResNorm<ST,MV,OP> test3( tol );
00201   Belos::StatusTestOutputter<ST,MV,OP> test4( frequency, false );
00202   test4.set_resNormStatusTest( rcp(&test3, false) );
00203   test4.set_outputManager( My_OM );
00204   
00205   RefCountPtr<StatusTestCombo_t > My_Test =
00206     rcp( new StatusTestCombo_t( StatusTestCombo_t::OR, test1, test2 ) );
00207   My_Test->AddStatusTest( test4 );
00208   //
00209   Belos::BlockGmres<ST,MV,OP>
00210     MyBlockGmres( My_LP, My_Test, My_OM, My_PL );
00211   //
00212   // **********Print out information about problem*******************
00213   //
00214   if (proc_verbose) {
00215     cout << endl << endl;
00216     cout << "Dimension of matrix: " << dim << endl;
00217     cout << "Number of right-hand sides: " << numrhs << endl;
00218     cout << "Block size used by solver: " << blocksize << endl;
00219     cout << "Max number of Gmres iterations: " << maxits << endl; 
00220     cout << "Relative residual tolerance: " << tol << endl;
00221     cout << endl;
00222   }
00223   //
00224   //
00225   if (proc_verbose) {
00226     cout << endl << endl;
00227     cout << "Running Block Gmres -- please wait" << endl;
00228     cout << (numrhs+blocksize-1)/blocksize 
00229    << " pass(es) through the solver required to solve for " << endl; 
00230     cout << numrhs << " right-hand side(s) -- using a block size of " << blocksize
00231    << endl << endl;
00232   }
00233 
00234   // 
00235   // Perform solve.
00236   //
00237   MyBlockGmres.Solve();
00238 
00239   RefCountPtr<MyMultiVec<ST> > temp = rcp( new MyMultiVec<ST>(dim,numrhs) );
00240   OPT::Apply( *A, *soln, *temp );
00241   MVT::MvAddMv( one, *rhs, -one, *temp, *temp );
00242   std::vector<MT> norm_num(numrhs), norm_denom(numrhs);
00243   MVT::MvNorm( *temp, &norm_num );
00244   MVT::MvNorm( *rhs, &norm_denom );
00245   for (int i=0; i<numrhs; ++i) {
00246     if (proc_verbose) 
00247       cout << "Relative residual "<<i<<" : " << norm_num[i] / norm_denom[i] << endl;
00248     if ( norm_num[i] / norm_denom[i] > tol ) {
00249       norm_failure = true;
00250     }
00251   }
00252   
00253   // Clean up.
00254   delete [] dvals;
00255   delete [] colptr;
00256   delete [] rowind;
00257   delete [] cvals;
00258 
00259 #ifdef HAVE_MPI
00260     MPI_Finalize();
00261 #endif
00262 
00263   if ( My_Test->GetStatus()!=Belos::Converged || norm_failure ) {
00264     if (proc_verbose)
00265       cout << "End Result: TEST FAILED" << endl;  
00266     return -1;
00267   }
00268   //
00269   // Default return value
00270   //
00271   if (proc_verbose)
00272     cout << "End Result: TEST PASSED" << endl;
00273   return 0;
00274   //
00275 } // end test_bl_gmres_complex_hb.cpp

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