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 "BelosBlockGmresSolMgr.hpp"
00038 #include "BelosPseudoBlockGmresSolMgr.hpp"
00039 #include "Teuchos_CommandLineProcessor.hpp"
00040 #include "Teuchos_ParameterList.hpp"
00041 
00042 #ifdef HAVE_MPI
00043 #include <mpi.h>
00044 #endif
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 namespace Belos {
00056   class MPIFinalize {
00057   public:
00058     ~MPIFinalize() {
00059 #ifdef HAVE_MPI
00060       MPI_Finalize();
00061 #endif
00062     }
00063   };
00064 }
00065 
00066 using namespace Teuchos;
00067 
00068 int main(int argc, char *argv[]) {
00069   //
00070 #ifdef HAVE_COMPLEX
00071   typedef std::complex<double> ST;
00072 #elif HAVE_COMPLEX_H
00073   typedef std::complex<double> ST;
00074 #else
00075   std::cout << "Not compiled with std::complex support." << std::endl;
00076   std::cout << "End Result: TEST FAILED" << std::endl;
00077   return -1;
00078 #endif
00079 
00080   typedef ScalarTraits<ST>                 SCT;
00081   typedef SCT::magnitudeType                MT;
00082   typedef Belos::MultiVec<ST>               MV;
00083   typedef Belos::Operator<ST>               OP;
00084   typedef Belos::MultiVecTraits<ST,MV>     MVT;
00085   typedef Belos::OperatorTraits<ST,MV,OP>  OPT;
00086   ST one  = SCT::one();
00087   ST zero = SCT::zero();  
00088 
00089   int info = 0;
00090   int MyPID = 0;
00091   bool norm_failure = false;
00092 
00093 #ifdef HAVE_MPI 
00094   MPI_Init(&argc,&argv);
00095   MPI_Comm_rank(MPI_COMM_WORLD, &MyPID);
00096   Belos::MPIFinalize mpiFinalize; // Will call finalize with *any* return
00097   (void)mpiFinalize;
00098 #endif
00099   //
00100   using Teuchos::RCP;
00101   using Teuchos::rcp;
00102 
00103   bool verbose = false, proc_verbose = false;
00104   bool pseudo = false;   // use pseudo block GMRES to solve this linear system.
00105   int frequency = -1;  // how often residuals are printed by solver
00106   int blocksize = 1;
00107   int numrhs = 1;
00108   int maxrestarts = 15;
00109   int length = 50;
00110   std::string filename("mhd1280b.cua");
00111   MT tol = 1.0e-5;  // relative residual tolerance
00112 
00113   CommandLineProcessor cmdp(false,true);
00114   cmdp.setOption("verbose","quiet",&verbose,"Print messages and results.");
00115   cmdp.setOption("pseudo","regular",&pseudo,"Use pseudo-block GMRES to solve the linear systems.");
00116   cmdp.setOption("frequency",&frequency,"Solvers frequency for printing residuals (#iters).");
00117   cmdp.setOption("filename",&filename,"Filename for Harwell-Boeing test matrix.");
00118   cmdp.setOption("tol",&tol,"Relative residual tolerance used by GMRES solver.");
00119   cmdp.setOption("num-rhs",&numrhs,"Number of right-hand sides to be solved for.");
00120   cmdp.setOption("num-restarts",&maxrestarts,"Maximum number of restarts allowed for the GMRES solver.");
00121   cmdp.setOption("blocksize",&blocksize,"Block size used by GMRES.");
00122   cmdp.setOption("subspace-length",&length,"Maximum dimension of block-subspace used by GMRES solver.");
00123   if (cmdp.parse(argc,argv) != CommandLineProcessor::PARSE_SUCCESSFUL) {
00124     return -1;
00125   }
00126   
00127   proc_verbose = verbose && (MyPID==0);  /* Only print on the zero processor */
00128   if (proc_verbose) {
00129     std::cout << Belos::Belos_Version() << std::endl << std::endl;
00130   }
00131   if (!verbose)
00132     frequency = -1;  // reset frequency if test is not verbose
00133   
00134   
00135 #ifndef HAVE_BELOS_TRIUTILS
00136   std::cout << "This test requires Triutils. Please configure with --enable-triutils." << std::endl;
00137   if (MyPID==0) {
00138     std::cout << "End Result: TEST FAILED" << std::endl;  
00139   }
00140   return -1;
00141 #endif
00142   
00143   // Get the data from the HB file
00144   int dim,dim2,nnz;
00145   MT *dvals;
00146   int *colptr,*rowind;
00147   ST *cvals;
00148   nnz = -1;
00149   info = readHB_newmat_double(filename.c_str(),&dim,&dim2,&nnz,
00150                               &colptr,&rowind,&dvals);
00151   if (info == 0 || nnz < 0) {
00152     if (MyPID==0) {
00153       std::cout << "Error reading '" << filename << "'" << std::endl;
00154       std::cout << "End Result: TEST FAILED" << std::endl;
00155     }
00156     return -1;
00157   }
00158   // Convert interleaved doubles to std::complex values
00159   cvals = new ST[nnz];
00160   for (int ii=0; ii<nnz; ii++) {
00161     cvals[ii] = ST(dvals[ii*2],dvals[ii*2+1]);
00162   }
00163   // Build the problem matrix
00164   RCP< MyBetterOperator<ST> > A 
00165     = rcp( new MyBetterOperator<ST>(dim,colptr,nnz,rowind,cvals) );
00166   //
00167   // ********Other information used by block solver***********
00168   // *****************(can be user specified)******************
00169   //
00170   int maxits = dim/blocksize; // maximum number of iterations to run
00171   //
00172   ParameterList belosList;
00173   belosList.set( "Num Blocks", length );                 // Maximum number of blocks in Krylov factorization
00174   belosList.set( "Block Size", blocksize );              // Blocksize to be used by iterative solver
00175   belosList.set( "Maximum Iterations", maxits );         // Maximum number of iterations allowed
00176   belosList.set( "Maximum Restarts", maxrestarts );      // Maximum number of restarts allowed
00177   belosList.set( "Convergence Tolerance", tol );         // Relative convergence tolerance requested
00178   if (verbose) {
00179     belosList.set( "Verbosity", Belos::Errors + Belos::Warnings + 
00180        Belos::TimingDetails + Belos::StatusTestDetails );
00181     if (frequency > 0)
00182       belosList.set( "Output Frequency", frequency );
00183   }
00184   else
00185     belosList.set( "Verbosity", Belos::Errors + Belos::Warnings );
00186   //
00187   // Construct the right-hand side and solution multivectors.
00188   // NOTE:  The right-hand side will be constructed such that the solution is
00189   // a vectors of one.
00190   //
00191   RCP<MyMultiVec<ST> > soln = rcp( new MyMultiVec<ST>(dim,numrhs) );
00192   RCP<MyMultiVec<ST> > rhs = rcp( new MyMultiVec<ST>(dim,numrhs) );
00193   MVT::MvRandom( *soln );
00194   OPT::Apply( *A, *soln, *rhs );
00195   MVT::MvInit( *soln, zero );
00196   //
00197   //  Construct an unpreconditioned linear problem instance.
00198   //
00199   RCP<Belos::LinearProblem<ST,MV,OP> > problem = 
00200     rcp( new Belos::LinearProblem<ST,MV,OP>( A, soln, rhs ) );
00201   bool set = problem->setProblem();
00202   if (set == false) {
00203     if (proc_verbose)
00204       std::cout << std::endl << "ERROR:  Belos::LinearProblem failed to set up correctly!" << std::endl;
00205     return -1;
00206   }
00207   //
00208   // *******************************************************************
00209   // *************Start the block Gmres iteration***********************
00210   // *******************************************************************
00211   //
00212   Teuchos::RCP< Belos::SolverManager<ST,MV,OP> > solver;
00213   if (pseudo)
00214     solver = Teuchos::rcp( new Belos::PseudoBlockGmresSolMgr<ST,MV,OP>( problem, Teuchos::rcp(&belosList,false) ) );
00215   else
00216     solver = Teuchos::rcp( new Belos::BlockGmresSolMgr<ST,MV,OP>( problem, Teuchos::rcp(&belosList,false) ) );
00217   //
00218   // **********Print out information about problem*******************
00219   //
00220   if (proc_verbose) {
00221     std::cout << std::endl << std::endl;
00222     std::cout << "Dimension of matrix: " << dim << std::endl;
00223     std::cout << "Number of right-hand sides: " << numrhs << std::endl;
00224     std::cout << "Block size used by solver: " << blocksize << std::endl;
00225     std::cout << "Max number of Gmres iterations: " << maxits << std::endl; 
00226     std::cout << "Relative residual tolerance: " << tol << std::endl;
00227     std::cout << std::endl;
00228   }
00229   //
00230   // Perform solve
00231   //
00232   Belos::ReturnType ret = solver->solve();
00233   //
00234   // Compute actual residuals.
00235   //
00236   RCP<MyMultiVec<ST> > temp = rcp( new MyMultiVec<ST>(dim,numrhs) );
00237   OPT::Apply( *A, *soln, *temp );
00238   MVT::MvAddMv( one, *rhs, -one, *temp, *temp );
00239   std::vector<MT> norm_num(numrhs), norm_denom(numrhs);
00240   MVT::MvNorm( *temp, norm_num );
00241   MVT::MvNorm( *rhs, norm_denom );
00242   for (int i=0; i<numrhs; ++i) {
00243     if (proc_verbose) 
00244       std::cout << "Relative residual "<<i<<" : " << norm_num[i] / norm_denom[i] << std::endl;
00245     if ( norm_num[i] / norm_denom[i] > tol ) {
00246       norm_failure = true;
00247     }
00248   }
00249   
00250   // Clean up.
00251   delete [] dvals;
00252   delete [] colptr;
00253   delete [] rowind;
00254   delete [] cvals;
00255 
00256   if ( ret!=Belos::Converged || norm_failure ) {
00257     if (proc_verbose)
00258       std::cout << "End Result: TEST FAILED" << std::endl;  
00259     return -1;
00260   }
00261   //
00262   // Default return value
00263   //
00264   if (proc_verbose)
00265     std::cout << "End Result: TEST PASSED" << std::endl;
00266   return 0;
00267 
00268   //
00269 } // end test_bl_gmres_complex_hb.cpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines
Generated on Wed Apr 13 10:05:15 2011 for Belos Package Browser (Single Doxygen Collection) by  doxygen 1.6.3