Stratimikos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines
belos_tpetra_thyra_lowsf_hb.cpp
Go to the documentation of this file.
00001 
00002 // @HEADER
00003 // ***********************************************************************
00004 // 
00005 //         Stratimikos: Thyra-based strategies for linear solvers
00006 //                Copyright (2006) Sandia Corporation
00007 // 
00008 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
00009 // license for use of this work by or on behalf of the U.S. Government.
00010 // 
00011 // Redistribution and use in source and binary forms, with or without
00012 // modification, are permitted provided that the following conditions are
00013 // met:
00014 //
00015 // 1. Redistributions of source code must retain the above copyright
00016 // notice, this list of conditions and the following disclaimer.
00017 //
00018 // 2. Redistributions in binary form must reproduce the above copyright
00019 // notice, this list of conditions and the following disclaimer in the
00020 // documentation and/or other materials provided with the distribution.
00021 //
00022 // 3. Neither the name of the Corporation nor the names of the
00023 // contributors may be used to endorse or promote products derived from
00024 // this software without specific prior written permission.
00025 //
00026 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00027 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00028 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00029 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00030 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00031 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00032 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00033 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00034 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00035 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00036 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00037 //
00038 // Questions? Contact Roscoe A. Bartlett (rabartl@sandia.gov) 
00039 // 
00040 // ***********************************************************************
00041 // @HEADER
00042 //
00043 // This driver reads a problem from a Harwell-Boeing (HB) file into an 
00044 // Tpetra::CisMatrix.  This matrix is then converted into a Thyra linear operator
00045 // through the Thyra-Tpetra adapters.
00046 //
00047 // The right-hand-side from the HB file is used instead of random vectors.
00048 // The initial guesses are all set to zero. 
00049 //
00050 //
00051 
00052 // Thyra includes
00053 #include "Thyra_BelosLinearOpWithSolveFactory.hpp"
00054 #include "Thyra_LinearOpWithSolveFactoryHelpers.hpp"
00055 #include "Thyra_TpetraLinearOp.hpp"
00056 
00057 // Tpetra includes
00058 #include "Tpetra_ElementSpace.hpp"
00059 #include "Tpetra_VectorSpace.hpp"
00060 #include "Tpetra_CisMatrix.hpp"
00061 
00062 // I/O for Harwell-Boeing files
00063 #ifdef HAVE_BELOS_TRIUTILS
00064 #include "iohb.h"
00065 #endif
00066 
00067 // Teuchos includes
00068 #include "Teuchos_ParameterList.hpp"
00069 #include "Teuchos_VerboseObject.hpp"
00070 #include "Teuchos_GlobalMPISession.hpp"
00071 
00072 #ifdef HAVE_MPI
00073 #  include "Tpetra_MpiPlatform.hpp"
00074 #else
00075 #  include "Tpetra_SerialPlatform.hpp"
00076 #endif
00077 
00078 // My Tpetra::Operator include
00079 #include "MyOperator.hpp"
00080 
00081 int main(int argc, char* argv[])
00082 {
00083   //
00084   // Get a default output stream from the Teuchos::VerboseObjectBase
00085   //
00086   Teuchos::RCP<Teuchos::FancyOStream>
00087     out = Teuchos::VerboseObjectBase::getDefaultOStream();
00088   
00089   Teuchos::GlobalMPISession mpiSession(&argc,&argv);
00090 
00091 #ifdef HAVE_COMPLEX
00092   typedef std::complex<double> ST;  // Scalar-type typedef
00093 #elif HAVE_COMPLEX_H
00094   typedef std::complex<double> ST;     // Scalar-type typedef
00095 #else
00096   typedef double ST;                // Scalar-type typedef
00097 #endif
00098   
00099   typedef Teuchos::ScalarTraits<ST>::magnitudeType MT;  // Magnitude-type typedef
00100   typedef int OT;                   // Ordinal-type typedef
00101   ST one = Teuchos::ScalarTraits<ST>::one(); 
00102   ST zero = Teuchos::ScalarTraits<ST>::zero(); 
00103   
00104 #ifdef HAVE_MPI
00105   MPI_Comm mpiComm = MPI_COMM_WORLD;
00106   const Tpetra::MpiPlatform<OT,OT>  ordinalPlatform(mpiComm);
00107   const Tpetra::MpiPlatform<OT,ST>   scalarPlatform(mpiComm);
00108 #else
00109   const Tpetra::SerialPlatform<OT,OT>  ordinalPlatform;
00110   const Tpetra::SerialPlatform<OT,ST>   scalarPlatform;
00111 #endif
00112   
00113   //
00114   // Get the data from the HB file
00115   //
00116   
00117   // Name of input matrix file
00118   std::string matrixFile = "mhd1280b.cua";
00119   
00120   int info=0;
00121   int dim,dim2,nnz;
00122   MT *dvals;
00123   int *colptr,*rowind;
00124   ST *cvals;
00125   nnz = -1;
00126   info = readHB_newmat_double(matrixFile.c_str(),&dim,&dim2,&nnz,
00127                               &colptr,&rowind,&dvals);
00128 
00129   if (info == 0 || nnz < 0) {
00130     *out << "Error reading '" << matrixFile << "'" << std::endl;
00131   }
00132 #ifdef HAVE_MPI
00133   MPI_Finalize();
00134 #endif
00135 
00136   // Convert interleaved doubles to std::complex values
00137   cvals = new ST[nnz];
00138   for (int ii=0; ii<nnz; ii++) {
00139     cvals[ii] = ST(dvals[ii*2],dvals[ii*2+1]);
00140   }
00141   
00142   // Declare global dimension of the linear operator
00143   OT globalDim = dim;
00144   
00145   // Create the element space and std::vector space
00146   const Tpetra::ElementSpace<OT> elementSpace(globalDim,0,ordinalPlatform);
00147   const Tpetra::VectorSpace<OT,ST> vectorSpace(elementSpace,scalarPlatform);
00148   
00149   // Create my implementation of a Tpetra::Operator
00150   RCP<Tpetra::Operator<OT,ST> >
00151     tpetra_A = rcp( new MyOperator<OT,ST>(vectorSpace,dim,colptr,nnz,rowind,cvals) );
00152 
00153   // Create a Thyra linear operator (A) using the Tpetra::CisMatrix (tpetra_A).
00154   RCP<Thyra::LinearOpBase<ST> >
00155     A = Teuchos::rcp( new Thyra::TpetraLinearOp<OT,ST>(tpetra_A) );
00156 
00157   //
00158   // Set the parameters for the Belos LOWS Factory and create a parameter list.
00159   //
00160   int             blockSize              = 1;
00161   int             maxIterations          = globalDim;
00162   int             maxRestarts            = 15;
00163   int             gmresKrylovLength      = 50;
00164   int             outputFrequency        = 100;
00165   bool            outputMaxResOnly       = true;
00166   MT              maxResid               = 1e-5;
00167 
00168   Teuchos::RCP<Teuchos::ParameterList>
00169     belosLOWSFPL = Teuchos::rcp( new Teuchos::ParameterList() );
00170  
00171   belosLOWSFPL->set("Solver Type","Block GMRES");
00172 
00173   Teuchos::ParameterList& belosLOWSFPL_solver =
00174     belosLOWSFPL->sublist("Solver Types");
00175 
00176   Teuchos::ParameterList& belosLOWSFPL_gmres =
00177     belosLOWSFPL_solver.sublist("Block GMRES");
00178 
00179   belosLOWSFPL_gmres.set("Maximum Iterations",int(maxIterations));
00180   belosLOWSFPL_gmres.set("Convergence Tolerance",MT(maxResid));
00181   belosLOWSFPL_gmres.set("Maximum Restarts",int(maxRestarts));
00182   belosLOWSFPL_gmres.set("Block Size",int(blockSize));
00183   belosLOWSFPL_gmres.set("Num Blocks",int(gmresKrylovLength));
00184   belosLOWSFPL_gmres.set("Output Frequency",int(outputFrequency));
00185   belosLOWSFPL_gmres.set("Show Maximum Residual Norm Only",bool(outputMaxResOnly));
00186  
00187   // Whether the linear solver succeeded.
00188   // (this will be set during the residual check at the end)
00189   bool success = true;
00190 
00191   // Number of random right-hand sides we will be solving for.
00192   int numRhs = 1;
00193 
00194   // Get the domain space for the Thyra linear operator 
00195   Teuchos::RCP<const Thyra::VectorSpaceBase<ST> > domain = A->domain();
00196 
00197   // Create the Belos LOWS factory.
00198   Teuchos::RCP<Thyra::LinearOpWithSolveFactoryBase<ST> >
00199     belosLOWSFactory = Teuchos::rcp(new Thyra::BelosLinearOpWithSolveFactory<ST>());
00200 
00201   // Set the parameter list to specify the behavior of the factory.
00202   belosLOWSFactory->setParameterList( belosLOWSFPL );
00203 
00204   // Set the output stream and the verbosity level (prints to std::cout by defualt)
00205   // NOTE:  Set to VERB_NONE for no output from the solver.
00206   belosLOWSFactory->setVerbLevel(Teuchos::VERB_LOW);
00207 
00208   // Create a BelosLinearOpWithSolve object from the Belos LOWS factory.
00209   Teuchos::RCP<Thyra::LinearOpWithSolveBase<ST> >
00210     nsA = belosLOWSFactory->createOp();
00211 
00212   // Initialize the BelosLinearOpWithSolve object with the Thyra linear operator.
00213   Thyra::initializeOp<ST>( *belosLOWSFactory, A, &*nsA );
00214 
00215   // Create a right-hand side with numRhs vectors in it.
00216   Teuchos::RCP< Thyra::MultiVectorBase<ST> > 
00217     b = Thyra::createMembers(domain, numRhs);
00218 
00219   // Create an initial std::vector with numRhs vectors in it and initialize it to one.
00220   Teuchos::RCP< Thyra::MultiVectorBase<ST> >
00221     x = Thyra::createMembers(domain, numRhs);
00222   Thyra::assign(&*x, one);
00223 
00224   // Initialize the right-hand side so that the solution is a std::vector of ones.
00225   A->apply( Thyra::NONCONJ_ELE, *x, &*b );
00226   Thyra::assign(&*x, zero);
00227 
00228   // Perform solve using the linear operator to get the approximate solution of Ax=b,
00229   // where b is the right-hand side and x is the left-hand side.
00230   Thyra::SolveStatus<ST> solveStatus;
00231   solveStatus = Thyra::solve( *nsA, Thyra::NONCONJ_ELE, *b, &*x );
00232 
00233   // Print out status of solve.
00234   *out << "\nBelos LOWS Status: "<< solveStatus << std::endl;
00235 
00236   //
00237   // Compute residual and ST check convergence.
00238   //
00239   std::vector<MT> norm_b(numRhs), norm_res(numRhs);
00240   Teuchos::RCP< Thyra::MultiVectorBase<ST> >
00241     y = Thyra::createMembers(domain, numRhs);
00242 
00243   // Compute the column norms of the right-hand side b.
00244   Thyra::norms_2( *b, &norm_b[0] );
00245 
00246   // Compute y=A*x, where x is the solution from the linear solver.
00247   A->apply( Thyra::NONCONJ_ELE, *x, &*y );
00248   
00249   // Compute A*x-b = y-b
00250   Thyra::update( -one, *b, &*y );
00251 
00252   // Compute the column norms of A*x-b.
00253   Thyra::norms_2( *y, &norm_res[0] );
00254 
00255   // Print out the final relative residual norms.
00256   MT rel_res = 0.0;
00257   *out << "Final relative residual norms" << std::endl;  
00258   for (int i=0; i<numRhs; ++i) {
00259     rel_res = norm_res[i]/norm_b[i];
00260     if (rel_res > maxResid)
00261       success = false;
00262     *out << "RHS " << i+1 << " : " 
00263    << std::setw(16) << std::right << rel_res << std::endl;
00264   }
00265 
00266   return ( success ? 0 : 1 );
00267 }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines