Belos Package Browser (Single Doxygen Collection) Development
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 the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
00008 // the U.S. Government retains certain rights in this software.
00009 //
00010 // Redistribution and use in source and binary forms, with or without
00011 // modification, are permitted provided that the following conditions are
00012 // met:
00013 //
00014 // 1. Redistributions of source code must retain the above copyright
00015 // notice, this list of conditions and the following disclaimer.
00016 //
00017 // 2. Redistributions in binary form must reproduce the above copyright
00018 // notice, this list of conditions and the following disclaimer in the
00019 // documentation and/or other materials provided with the distribution.
00020 //
00021 // 3. Neither the name of the Corporation nor the names of the
00022 // contributors may be used to endorse or promote products derived from
00023 // this software without specific prior written permission.
00024 //
00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00036 //
00037 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
00038 //
00039 // ************************************************************************
00040 //@HEADER
00041 //
00042 // This driver reads a problem from a Harwell-Boeing (HB) file.
00043 // The right-hand-side from the HB file is used instead of random vectors.
00044 // The initial guesses are all set to zero. 
00045 //
00046 // NOTE: No preconditioner is used in this case. 
00047 //
00048 #include "BelosConfigDefs.hpp"
00049 #include "BelosLinearProblem.hpp"
00050 #include "BelosBlockGmresSolMgr.hpp"
00051 #include "BelosPseudoBlockGmresSolMgr.hpp"
00052 #include "Teuchos_CommandLineProcessor.hpp"
00053 #include "Teuchos_ParameterList.hpp"
00054 
00055 #ifdef HAVE_MPI
00056 #include <mpi.h>
00057 #endif
00058 
00059 // I/O for Harwell-Boeing files
00060 #ifdef HAVE_BELOS_TRIUTILS
00061 #include "iohb.h"
00062 #endif
00063 
00064 #include "MyMultiVec.hpp"
00065 #include "MyBetterOperator.hpp"
00066 #include "MyOperator.hpp"
00067 
00068 namespace Belos {
00069   class MPIFinalize {
00070   public:
00071     ~MPIFinalize() {
00072 #ifdef HAVE_MPI
00073       MPI_Finalize();
00074 #endif
00075     }
00076   };
00077 }
00078 
00079 using namespace Teuchos;
00080 
00081 int main(int argc, char *argv[]) {
00082   //
00083 #ifdef HAVE_COMPLEX
00084   typedef std::complex<double> ST;
00085 #elif HAVE_COMPLEX_H
00086   typedef std::complex<double> ST;
00087 #else
00088   std::cout << "Not compiled with std::complex support." << std::endl;
00089   std::cout << "End Result: TEST FAILED" << std::endl;
00090   return -1;
00091 #endif
00092 
00093   typedef ScalarTraits<ST>                 SCT;
00094   typedef SCT::magnitudeType                MT;
00095   typedef Belos::MultiVec<ST>               MV;
00096   typedef Belos::Operator<ST>               OP;
00097   typedef Belos::MultiVecTraits<ST,MV>     MVT;
00098   typedef Belos::OperatorTraits<ST,MV,OP>  OPT;
00099   ST one  = SCT::one();
00100   ST zero = SCT::zero();  
00101 
00102   int info = 0;
00103   int MyPID = 0;
00104   bool norm_failure = false;
00105 
00106 #ifdef HAVE_MPI 
00107   MPI_Init(&argc,&argv);
00108   MPI_Comm_rank(MPI_COMM_WORLD, &MyPID);
00109   Belos::MPIFinalize mpiFinalize; // Will call finalize with *any* return
00110   (void)mpiFinalize;
00111 #endif
00112   //
00113   using Teuchos::RCP;
00114   using Teuchos::rcp;
00115 
00116   bool verbose = false, proc_verbose = false;
00117   bool pseudo = false;   // use pseudo block GMRES to solve this linear system.
00118   int frequency = -1;  // how often residuals are printed by solver
00119   int blocksize = 1;
00120   int numrhs = 1;
00121   int maxrestarts = 15;
00122   int length = 50;
00123   std::string filename("mhd1280b.cua");
00124   MT tol = 1.0e-5;  // relative residual tolerance
00125 
00126   CommandLineProcessor cmdp(false,true);
00127   cmdp.setOption("verbose","quiet",&verbose,"Print messages and results.");
00128   cmdp.setOption("pseudo","regular",&pseudo,"Use pseudo-block GMRES to solve the linear systems.");
00129   cmdp.setOption("frequency",&frequency,"Solvers frequency for printing residuals (#iters).");
00130   cmdp.setOption("filename",&filename,"Filename for Harwell-Boeing test matrix.");
00131   cmdp.setOption("tol",&tol,"Relative residual tolerance used by GMRES solver.");
00132   cmdp.setOption("num-rhs",&numrhs,"Number of right-hand sides to be solved for.");
00133   cmdp.setOption("num-restarts",&maxrestarts,"Maximum number of restarts allowed for the GMRES solver.");
00134   cmdp.setOption("blocksize",&blocksize,"Block size used by GMRES.");
00135   cmdp.setOption("subspace-length",&length,"Maximum dimension of block-subspace used by GMRES solver.");
00136   if (cmdp.parse(argc,argv) != CommandLineProcessor::PARSE_SUCCESSFUL) {
00137     return -1;
00138   }
00139   
00140   proc_verbose = verbose && (MyPID==0);  /* Only print on the zero processor */
00141   if (proc_verbose) {
00142     std::cout << Belos::Belos_Version() << std::endl << std::endl;
00143   }
00144   if (!verbose)
00145     frequency = -1;  // reset frequency if test is not verbose
00146   
00147   
00148 #ifndef HAVE_BELOS_TRIUTILS
00149   std::cout << "This test requires Triutils. Please configure with --enable-triutils." << std::endl;
00150   if (MyPID==0) {
00151     std::cout << "End Result: TEST FAILED" << std::endl;  
00152   }
00153   return -1;
00154 #endif
00155   
00156   // Get the data from the HB file
00157   int dim,dim2,nnz;
00158   MT *dvals;
00159   int *colptr,*rowind;
00160   ST *cvals;
00161   nnz = -1;
00162   info = readHB_newmat_double(filename.c_str(),&dim,&dim2,&nnz,
00163                               &colptr,&rowind,&dvals);
00164   if (info == 0 || nnz < 0) {
00165     if (MyPID==0) {
00166       std::cout << "Error reading '" << filename << "'" << std::endl;
00167       std::cout << "End Result: TEST FAILED" << std::endl;
00168     }
00169     return -1;
00170   }
00171   // Convert interleaved doubles to std::complex values
00172   cvals = new ST[nnz];
00173   for (int ii=0; ii<nnz; ii++) {
00174     cvals[ii] = ST(dvals[ii*2],dvals[ii*2+1]);
00175   }
00176   // Build the problem matrix
00177   RCP< MyBetterOperator<ST> > A 
00178     = rcp( new MyBetterOperator<ST>(dim,colptr,nnz,rowind,cvals) );
00179   //
00180   // ********Other information used by block solver***********
00181   // *****************(can be user specified)******************
00182   //
00183   int maxits = dim/blocksize; // maximum number of iterations to run
00184   //
00185   ParameterList belosList;
00186   belosList.set( "Num Blocks", length );                 // Maximum number of blocks in Krylov factorization
00187   belosList.set( "Block Size", blocksize );              // Blocksize to be used by iterative solver
00188   belosList.set( "Maximum Iterations", maxits );         // Maximum number of iterations allowed
00189   belosList.set( "Maximum Restarts", maxrestarts );      // Maximum number of restarts allowed
00190   belosList.set( "Convergence Tolerance", tol );         // Relative convergence tolerance requested
00191   if (verbose) {
00192     belosList.set( "Verbosity", Belos::Errors + Belos::Warnings + 
00193        Belos::TimingDetails + Belos::StatusTestDetails );
00194     if (frequency > 0)
00195       belosList.set( "Output Frequency", frequency );
00196   }
00197   else
00198     belosList.set( "Verbosity", Belos::Errors + Belos::Warnings );
00199   //
00200   // Construct the right-hand side and solution multivectors.
00201   // NOTE:  The right-hand side will be constructed such that the solution is
00202   // a vectors of one.
00203   //
00204   RCP<MyMultiVec<ST> > soln = rcp( new MyMultiVec<ST>(dim,numrhs) );
00205   RCP<MyMultiVec<ST> > rhs = rcp( new MyMultiVec<ST>(dim,numrhs) );
00206   MVT::MvRandom( *soln );
00207   OPT::Apply( *A, *soln, *rhs );
00208   MVT::MvInit( *soln, zero );
00209   //
00210   //  Construct an unpreconditioned linear problem instance.
00211   //
00212   RCP<Belos::LinearProblem<ST,MV,OP> > problem = 
00213     rcp( new Belos::LinearProblem<ST,MV,OP>( A, soln, rhs ) );
00214   bool set = problem->setProblem();
00215   if (set == false) {
00216     if (proc_verbose)
00217       std::cout << std::endl << "ERROR:  Belos::LinearProblem failed to set up correctly!" << std::endl;
00218     return -1;
00219   }
00220   //
00221   // *******************************************************************
00222   // *************Start the block Gmres iteration***********************
00223   // *******************************************************************
00224   //
00225   Teuchos::RCP< Belos::SolverManager<ST,MV,OP> > solver;
00226   if (pseudo)
00227     solver = Teuchos::rcp( new Belos::PseudoBlockGmresSolMgr<ST,MV,OP>( problem, Teuchos::rcp(&belosList,false) ) );
00228   else
00229     solver = Teuchos::rcp( new Belos::BlockGmresSolMgr<ST,MV,OP>( problem, Teuchos::rcp(&belosList,false) ) );
00230   //
00231   // **********Print out information about problem*******************
00232   //
00233   if (proc_verbose) {
00234     std::cout << std::endl << std::endl;
00235     std::cout << "Dimension of matrix: " << dim << std::endl;
00236     std::cout << "Number of right-hand sides: " << numrhs << std::endl;
00237     std::cout << "Block size used by solver: " << blocksize << std::endl;
00238     std::cout << "Max number of Gmres iterations: " << maxits << std::endl; 
00239     std::cout << "Relative residual tolerance: " << tol << std::endl;
00240     std::cout << std::endl;
00241   }
00242   //
00243   // Perform solve
00244   //
00245   Belos::ReturnType ret = solver->solve();
00246   //
00247   // Compute actual residuals.
00248   //
00249   RCP<MyMultiVec<ST> > temp = rcp( new MyMultiVec<ST>(dim,numrhs) );
00250   OPT::Apply( *A, *soln, *temp );
00251   MVT::MvAddMv( one, *rhs, -one, *temp, *temp );
00252   std::vector<MT> norm_num(numrhs), norm_denom(numrhs);
00253   MVT::MvNorm( *temp, norm_num );
00254   MVT::MvNorm( *rhs, norm_denom );
00255   for (int i=0; i<numrhs; ++i) {
00256     if (proc_verbose) 
00257       std::cout << "Relative residual "<<i<<" : " << norm_num[i] / norm_denom[i] << std::endl;
00258     if ( norm_num[i] / norm_denom[i] > tol ) {
00259       norm_failure = true;
00260     }
00261   }
00262   
00263   // Clean up.
00264   delete [] dvals;
00265   delete [] colptr;
00266   delete [] rowind;
00267   delete [] cvals;
00268 
00269   if ( ret!=Belos::Converged || norm_failure ) {
00270     if (proc_verbose)
00271       std::cout << "End Result: TEST FAILED" << std::endl;  
00272     return -1;
00273   }
00274   //
00275   // Default return value
00276   //
00277   if (proc_verbose)
00278     std::cout << "End Result: TEST PASSED" << std::endl;
00279   return 0;
00280 
00281   //
00282 } // end test_bl_gmres_complex_hb.cpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines