test_bl_cg_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 "BelosStatusTestResNorm.hpp"
00040 #include "BelosStatusTestOutputter.hpp"
00041 #include "BelosStatusTestCombo.hpp"
00042 #include "BelosBlockCG.hpp"
00043 #include "Teuchos_CommandLineProcessor.hpp"
00044 
00045 // I/O for Harwell-Boeing files
00046 #ifdef HAVE_BELOS_TRIUTILS
00047 #include "iohb.h"
00048 #endif
00049 
00050 #include "MyMultiVec.hpp"
00051 #include "MyBetterOperator.hpp"
00052 
00053 using namespace Teuchos;
00054 
00055 int main(int argc, char *argv[]) {
00056   //
00057 #ifdef HAVE_COMPLEX
00058   typedef std::complex<double> ST;
00059 #elif HAVE_COMPLEX_H
00060   typedef ::complex<double> ST;
00061 #else
00062   typedef double ST;
00063   // no complex. quit with failure.
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 
00077   int info = 0;
00078   int MyPID = 0;
00079 
00080 #ifdef HAVE_MPI 
00081   // Initialize MPI 
00082   MPI_Init(&argc,&argv);  
00083   MPI_Comm_rank(MPI_COMM_WORLD, &MyPID);
00084 #endif
00085   //
00086   using Teuchos::RefCountPtr;
00087   using Teuchos::rcp;
00088 
00089   bool verbose = false, proc_verbose = false;
00090   int frequency = -1;  // how often residuals are printed by solver
00091   int numrhs = 1;  // total number of right-hand sides to solve for
00092   int blockSize = 1;  // blocksize used by solver
00093   std::string filename("mhd1280b.cua");
00094   MT tol = 1.0e-5;  // relative residual tolerance
00095 
00096   CommandLineProcessor cmdp(false,true);
00097   cmdp.setOption("verbose","quiet",&verbose,"Print messages and results.");
00098   cmdp.setOption("frequency",&frequency,"Solvers frequency for printing residuals (#iters).");
00099   cmdp.setOption("filename",&filename,"Filename for Harwell-Boeing test matrix.");
00100   cmdp.setOption("tol",&tol,"Relative residual tolerance used by CG solver.");
00101   cmdp.setOption("num-rhs",&numrhs,"Number of right-hand sides to be solved for.");
00102   cmdp.setOption("block-size",&blockSize,"Block size to be used by CG solver.");
00103   if (cmdp.parse(argc,argv) != CommandLineProcessor::PARSE_SUCCESSFUL) {
00104 #ifdef HAVE_MPI
00105     MPI_Finalize();
00106 #endif
00107     return -1;
00108   }
00109 
00110   proc_verbose = verbose && (MyPID==0);  /* Only print on the zero processor */
00111   if (proc_verbose) {
00112     cout << Belos::Belos_Version() << endl << endl;
00113   }
00114   if (!verbose)
00115     frequency = -1;  // reset frequency if test is not verbose
00116   //
00117   //
00118 
00119 #ifndef HAVE_BELOS_TRIUTILS
00120    cout << "This test requires Triutils. Please configure with --enable-triutils." << endl;
00121 #ifdef HAVE_MPI
00122   MPI_Finalize() ;
00123 #endif
00124   if (MyPID==0) {
00125     cout << "End Result: TEST FAILED" << endl;  
00126   }
00127   return -1;
00128 #endif
00129 
00130   // Get the data from the HB file
00131   int dim,dim2,nnz;
00132   double *dvals;
00133   int *colptr,*rowind;
00134   ST *cvals;
00135   nnz = -1;
00136   info = readHB_newmat_double(filename.c_str(),&dim,&dim2,&nnz,
00137                               &colptr,&rowind,&dvals);
00138   if (info == 0 || nnz < 0) {
00139     if (MyPID==0) {
00140       cout << "Error reading '" << filename << "'" << endl;
00141       cout << "End Result: TEST FAILED" << endl;
00142     }
00143 #ifdef HAVE_MPI
00144     MPI_Finalize();
00145 #endif
00146     return -1;
00147   }
00148   // Convert interleaved doubles to complex values
00149   cvals = new ST[nnz];
00150   for (int ii=0; ii<nnz; ii++) {
00151     cvals[ii] = ST(dvals[ii*2],dvals[ii*2+1]);
00152   }
00153   // Build the problem matrix
00154   RefCountPtr< MyBetterOperator<ST> > A 
00155     = rcp( new MyBetterOperator<ST>(dim,colptr,nnz,rowind,cvals) );
00156   //
00157   // ********Other information used by block solver***********
00158   // *****************(can be user specified)******************
00159   //
00160   int maxits = dim/blockSize; // maximum number of iterations to run
00161   //
00162   // Construct the right-hand side and solution multivectors.
00163   // NOTE:  The right-hand side will be constructed such that the solution is
00164   // a vectors of one.
00165   //
00166   RefCountPtr<MyMultiVec<ST> > soln = rcp( new MyMultiVec<ST>(dim,numrhs) );
00167   RefCountPtr<MyMultiVec<ST> > rhs = rcp( new MyMultiVec<ST>(dim,numrhs) );
00168 
00169   if (numrhs == 1) {
00170     MVT::MvInit( *soln, SCT::one() );
00171     OPT::Apply( *A, *soln, *rhs );
00172     MVT::MvInit( *soln, SCT::zero() );
00173   } else {
00174     MVT::MvRandom( *soln );
00175     OPT::Apply( *A, *soln, *rhs );
00176     MVT::MvInit( *soln, SCT::zero() );
00177   }
00178   //
00179   //  Construct an unpreconditioned linear problem instance.
00180   //
00181   RefCountPtr<Belos::LinearProblem<ST,MV,OP> > MyLP
00182     = rcp( new Belos::LinearProblem<ST,MV,OP>( A, soln, rhs ) );
00183   MyLP->SetBlockSize( blockSize );
00184   //
00185   // *******************************************************************
00186   // *************Start the block CG iteration*************************
00187   // *******************************************************************
00188   //
00189   // Create default output manager 
00190   RefCountPtr<Belos::OutputManager<ST> > MyOM 
00191     = rcp( new Belos::OutputManager<ST>( MyPID ) );
00192   // Set verbosity level
00193   if (verbose) {
00194     MyOM->SetVerbosity( Belos::Errors + Belos::Warnings 
00195       + Belos::TimingDetails + Belos::FinalSummary );
00196   }
00197 
00198   Belos::StatusTestMaxIters<ST,MV,OP> test1( maxits );
00199   Belos::StatusTestResNorm<ST,MV,OP> test2( tol );
00200   Belos::StatusTestOutputter<ST,MV,OP> test3( frequency, false );
00201   test3.set_resNormStatusTest( rcp(&test2,false) );
00202   test3.set_outputManager( MyOM );
00203       
00204   RefCountPtr<Belos::StatusTestCombo<ST,MV,OP> > MyTest
00205     = rcp( new Belos::StatusTestCombo<ST,MV,OP>( Belos::StatusTestCombo<ST,MV,OP>::OR, test1, test3 ) );
00206   //
00207   Belos::BlockCG<ST,MV,OP>
00208     MyBlockCG( MyLP, MyTest, MyOM );
00209   //
00210   // **********Print out information about problem*******************
00211   //
00212   if (proc_verbose) {
00213     cout << endl << endl;
00214     cout << "Dimension of matrix: " << dim << endl;
00215     cout << "Number of right-hand sides: " << numrhs << endl;
00216     cout << "Block size used by solver: " << blockSize << endl;
00217     cout << "Max number of CG iterations: " << maxits << endl; 
00218     cout << "Relative residual tolerance: " << tol << endl;
00219     cout << endl;
00220   }
00221   //
00222   //
00223   if (proc_verbose) {
00224     cout << endl << endl;
00225     cout << "Running Block CG -- please wait" << endl;
00226     cout << (numrhs+blockSize-1)/blockSize 
00227    << " pass(es) through the solver required to solve for " << endl; 
00228     cout << numrhs << " right-hand side(s) -- using a block size of " << blockSize
00229    << endl << endl;
00230   }
00231 
00232   MyBlockCG.Solve();
00233   
00234 #ifdef HAVE_MPI
00235   MPI_Finalize();
00236 #endif
00237   
00238   if (MyTest->GetStatus()!=Belos::Converged) {
00239     if (proc_verbose)
00240       cout << "End Result: TEST FAILED" << endl;  
00241     return -1;
00242   }
00243   //
00244   // Default return value
00245   //
00246   if (proc_verbose)
00247     cout << "End Result: TEST PASSED" << endl;
00248   return 0;
00249 
00250   //
00251 } // end test_bl_cg_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