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 #include "Teuchos_StandardCatchMacros.hpp"
00055 #include "Teuchos_StandardCatchMacros.hpp"
00056 
00057 #ifdef HAVE_MPI
00058 #include <mpi.h>
00059 #endif
00060 
00061 // I/O for Harwell-Boeing files
00062 #ifdef HAVE_BELOS_TRIUTILS
00063 #include "iohb.h"
00064 #endif
00065 
00066 #include "MyMultiVec.hpp"
00067 #include "MyBetterOperator.hpp"
00068 #include "MyOperator.hpp"
00069 
00070 using namespace Teuchos;
00071 
00072 int main(int argc, char *argv[]) {
00073   //
00074 #ifdef HAVE_COMPLEX
00075   typedef std::complex<double> ST;
00076 #elif HAVE_COMPLEX_H
00077   typedef std::complex<double> ST;
00078 #else
00079   std::cout << "Not compiled with std::complex support." << std::endl;
00080   std::cout << "End Result: TEST FAILED" << std::endl;
00081   return EXIT_FAILURE;
00082 #endif
00083 
00084   typedef ScalarTraits<ST>                 SCT;
00085   typedef SCT::magnitudeType                MT;
00086   typedef Belos::MultiVec<ST>               MV;
00087   typedef Belos::Operator<ST>               OP;
00088   typedef Belos::MultiVecTraits<ST,MV>     MVT;
00089   typedef Belos::OperatorTraits<ST,MV,OP>  OPT;
00090   ST one  = SCT::one();
00091   ST zero = SCT::zero();
00092 
00093   Teuchos::GlobalMPISession session(&argc, &argv, NULL);
00094   int MyPID = session.getRank();
00095   //
00096   using Teuchos::RCP;
00097   using Teuchos::rcp;
00098 
00099   bool success = false;
00100   bool verbose = false;
00101   try {
00102     int info = 0;
00103     bool norm_failure = false;
00104     bool proc_verbose = false;
00105     bool pseudo = false;   // use pseudo block GMRES to solve this linear system.
00106     int frequency = -1;  // how often residuals are printed by solver
00107     int blocksize = 1;
00108     int numrhs = 1;
00109     int maxrestarts = 15;
00110     int length = 50;
00111     std::string filename("mhd1280b.cua");
00112     MT tol = 1.0e-5;  // relative residual tolerance
00113 
00114     CommandLineProcessor cmdp(false,true);
00115     cmdp.setOption("verbose","quiet",&verbose,"Print messages and results.");
00116     cmdp.setOption("pseudo","regular",&pseudo,"Use pseudo-block GMRES to solve the linear systems.");
00117     cmdp.setOption("frequency",&frequency,"Solvers frequency for printing residuals (#iters).");
00118     cmdp.setOption("filename",&filename,"Filename for Harwell-Boeing test matrix.");
00119     cmdp.setOption("tol",&tol,"Relative residual tolerance used by GMRES solver.");
00120     cmdp.setOption("num-rhs",&numrhs,"Number of right-hand sides to be solved for.");
00121     cmdp.setOption("num-restarts",&maxrestarts,"Maximum number of restarts allowed for the GMRES solver.");
00122     cmdp.setOption("blocksize",&blocksize,"Block size used by GMRES.");
00123     cmdp.setOption("subspace-length",&length,"Maximum dimension of block-subspace used by GMRES solver.");
00124     if (cmdp.parse(argc,argv) != CommandLineProcessor::PARSE_SUCCESSFUL) {
00125       return EXIT_FAILURE;
00126     }
00127 
00128     proc_verbose = verbose && (MyPID==0);  /* Only print on the zero processor */
00129     if (proc_verbose) {
00130       std::cout << Belos::Belos_Version() << std::endl << std::endl;
00131     }
00132     if (!verbose)
00133       frequency = -1;  // reset frequency if test is not verbose
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 EXIT_FAILURE;
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 EXIT_FAILURE;
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 EXIT_FAILURE;
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     success = ret==Belos::Converged && !norm_failure;
00257     if (success) {
00258       if (proc_verbose)
00259         std::cout << "End Result: TEST PASSED" << std::endl;
00260     } else {
00261       if (proc_verbose)
00262         std::cout << "End Result: TEST FAILED" << std::endl;
00263     }
00264   }
00265   TEUCHOS_STANDARD_CATCH_STATEMENTS(verbose, std::cerr, success);
00266 
00267   return ( success ? EXIT_SUCCESS : EXIT_FAILURE );
00268 } // end test_bl_gmres_complex_hb.cpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines