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_prec_lap.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 creates a 2D Laplacian operator as a Tpetra::CisMatrix
00044 // and a preconditioner.  These matrices are then converted into Thyra 
00045 // linear operators through the Thyra-Tpetra adapters.
00046 //
00047 // The right-hand sides are all 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_PreconditionerFactoryHelpers.hpp"
00056 #include "Thyra_TpetraLinearOp.hpp"
00057 #include "Thyra_DefaultPreconditioner.hpp"
00058 
00059 // Tpetra includes
00060 #include "Tpetra_ElementSpace.hpp"
00061 #include "Tpetra_VectorSpace.hpp"
00062 #include "Tpetra_CisMatrix.hpp"
00063 
00064 // Teuchos includes
00065 #include "Teuchos_ParameterList.hpp"
00066 #include "Teuchos_VerboseObject.hpp"
00067 #include "Teuchos_GlobalMPISession.hpp"
00068 
00069 #ifdef HAVE_MPI
00070 #  include "Tpetra_MpiPlatform.hpp"
00071 #else
00072 #  include "Tpetra_SerialPlatform.hpp"
00073 #endif
00074 
00075 int main(int argc, char* argv[])
00076 {
00077   //
00078   // Get a default output stream from the Teuchos::VerboseObjectBase
00079   //
00080   Teuchos::RCP<Teuchos::FancyOStream>
00081     out = Teuchos::VerboseObjectBase::getDefaultOStream();
00082   
00083   Teuchos::GlobalMPISession mpiSession(&argc,&argv);
00084 
00085   //
00086   // Create the scalar-type typedefs
00087   //
00088 
00089 #ifdef HAVE_COMPLEX
00090   typedef std::complex<double> ST;  // Scalar-type typedef
00091 #elif HAVE_COMPLEX_H
00092   typedef std::complex<double> ST;     // Scalar-type typedef
00093 #else
00094   typedef double ST;                // Scalar-type typedef
00095 #endif
00096   typedef Teuchos::ScalarTraits<ST>::magnitudeType MT;  // Magnitude-type typdef
00097   typedef int OT;                   // Ordinal-type typdef
00098 
00099   ST one = Teuchos::ScalarTraits<ST>::one();
00100   ST zero = Teuchos::ScalarTraits<ST>::zero();
00101   
00102   //
00103   // Create 2D Laplacian using Tpetra::CisMatrix
00104   //
00105 
00106 #ifdef HAVE_MPI
00107   MPI_Comm mpiComm = MPI_COMM_WORLD;
00108   const Tpetra::MpiPlatform<OT,OT>  ordinalPlatform(mpiComm);
00109   const Tpetra::MpiPlatform<OT,ST>   scalarPlatform(mpiComm);
00110 #else
00111   const Tpetra::SerialPlatform<OT,OT>  ordinalPlatform;
00112   const Tpetra::SerialPlatform<OT,ST>   scalarPlatform;
00113 #endif
00114 
00115   // Declare global dimension of the linear operator
00116   OT globalDim = 500;
00117 
00118   // Create the element space and std::vector space
00119   const Tpetra::ElementSpace<OT> elementSpace(globalDim,0,ordinalPlatform);
00120   const Tpetra::VectorSpace<OT,ST> vectorSpace(elementSpace,scalarPlatform);
00121   
00122   // Allocate the Tpetra::CisMatrix object.
00123   RCP<Tpetra::CisMatrix<OT,ST> >
00124     tpetra_A = rcp(new Tpetra::CisMatrix<OT,ST>(vectorSpace));
00125 
00126   // Allocate the Tpetra::CisMatrix preconditioning object.
00127   RCP<Tpetra::CisMatrix<OT,ST> >
00128     tpetra_Prec = rcp(new Tpetra::CisMatrix<OT,ST>(vectorSpace));
00129 
00130   // Get the indexes of the rows on this processor
00131   const int numMyElements = vectorSpace.getNumMyEntries();
00132   const std::vector<int> &myGlobalElements = vectorSpace.elementSpace().getMyGlobalElements();
00133 
00134   // Fill the local matrix entries one row at a time.
00135   const ST offDiag = -1.0, diag = 2.0;
00136   int numEntries; ST values[3]; int indexes[3];
00137   ST prec_value[1]; int prec_index[1];
00138   prec_value[0] = one/diag;   // Diagonal preconditioning
00139   for( int k = 0; k < numMyElements; ++k ) {
00140     const int rowIndex = myGlobalElements[k];
00141     prec_index[0] = rowIndex;
00142     if( rowIndex == 0 ) {                     // First row
00143       numEntries = 2;
00144       values[0]  = diag;             values[1]  = offDiag;
00145       indexes[0] = 0;                indexes[1] = 1; 
00146     }
00147     else if( rowIndex == globalDim - 1 ) {    // Last row
00148       numEntries = 2;
00149       values[0]  = offDiag;         values[1]  = diag;
00150       indexes[0] = globalDim-2;     indexes[1] = globalDim-1; 
00151     }
00152     else {                                    // Middle rows
00153       numEntries = 3;
00154       values[0]  = offDiag;         values[1]  = diag;          values[2]  = offDiag;
00155       indexes[0] = rowIndex-1;      indexes[1] = rowIndex;      indexes[2] = rowIndex+1; 
00156     }
00157     tpetra_A->submitEntries(Tpetra::Insert,rowIndex,numEntries,values,indexes);
00158     tpetra_Prec->submitEntries(Tpetra::Insert,rowIndex,1,prec_value,prec_index);
00159   }
00160 
00161   // Finish the construction of the Tpetra::CisMatrix
00162   tpetra_A->fillComplete();
00163   tpetra_Prec->fillComplete();
00164 
00165   // Create a Thyra linear operator (A) using the Tpetra::CisMatrix (tpetra_A).
00166   RCP<Thyra::LinearOpBase<ST> >
00167     A = Teuchos::rcp( new Thyra::TpetraLinearOp<OT,ST>(
00168                       Teuchos::rcp_implicit_cast<Tpetra::Operator<OT,ST> >(tpetra_A) ) );
00169 
00170   // Create a Thyra linear operator (Prec) using the Tpetra::CisMatrix (tpetra_Prec)
00171   RCP<Thyra::LinearOpBase<ST> >
00172     Prec = Teuchos::rcp( new Thyra::TpetraLinearOp<OT,ST>(
00173                     Teuchos::rcp_implicit_cast<Tpetra::Operator<OT,ST> >(tpetra_Prec) ) );
00174 
00175   // Create a Thyra default preconditioner (DefPrec) using the Thyra linear operator (Prec)
00176   RCP<const Thyra::DefaultPreconditioner<ST> >
00177     DefPrec = Thyra::unspecifiedPrec<ST>(Prec);
00178 
00179   //
00180   // Set the parameters for the Belos LOWS Factory and create a parameter list.
00181   // NOTE:  All the code below only uses Thyra and is independent of Tpetra.
00182   //
00183   int             blockSize              = 1;  // can only be 1 right now.
00184   int             maxIterations          = globalDim;
00185   int             maxRestarts            = 0;
00186   int             gmresKrylovLength      = globalDim;
00187   int             outputFrequency        = 100;
00188   bool            outputMaxResOnly       = true;
00189   MT              maxResid               = 1e-3;
00190 
00191   Teuchos::RCP<Teuchos::ParameterList>
00192     belosLOWSFPL = Teuchos::rcp( new Teuchos::ParameterList() );
00193  
00194   belosLOWSFPL->set("Solver Type","Block GMRES");
00195 
00196   Teuchos::ParameterList& belosLOWSFPL_solver =
00197     belosLOWSFPL->sublist("Solver Types");
00198 
00199   Teuchos::ParameterList& belosLOWSFPL_gmres =
00200     belosLOWSFPL_solver.sublist("Block GMRES");
00201 
00202   belosLOWSFPL_gmres.set("Maximum Iterations",int(maxIterations));
00203   belosLOWSFPL_gmres.set("Convergence Tolerance",double(maxResid));
00204   belosLOWSFPL_gmres.set("Maximum Restarts",int(maxRestarts));
00205   belosLOWSFPL_gmres.set("Block Size",int(blockSize));
00206   belosLOWSFPL_gmres.set("Num Blocks",int(gmresKrylovLength));
00207   belosLOWSFPL_gmres.set("Output Frequency",int(outputFrequency));
00208   belosLOWSFPL_gmres.set("Show Maximum Residual Norm Only",bool(outputMaxResOnly));
00209  
00210   // Whether the linear solver succeeded.
00211   // (this will be set during the residual check at the end)
00212   bool success = true;
00213 
00214   // Number of random right-hand sides we will be solving for.
00215   int numRhs = 5;
00216 
00217   // Get the domain space for the Thyra linear operator 
00218   Teuchos::RCP<const Thyra::VectorSpaceBase<ST> > domain = A->domain();
00219 
00220   // Create the Belos LOWS factory.
00221   Teuchos::RCP<Thyra::LinearOpWithSolveFactoryBase<ST> >
00222     belosLOWSFactory = Teuchos::rcp(new Thyra::BelosLinearOpWithSolveFactory<ST>());
00223 
00224   // Set the parameter list to specify the behavior of the factory.
00225   belosLOWSFactory->setParameterList( belosLOWSFPL );
00226 
00227   // Set the output stream and the verbosity level (prints to std::cout by defualt)
00228   // NOTE:  Set to VERB_NONE for no output from the solver.
00229   belosLOWSFactory->setVerbLevel(Teuchos::VERB_LOW);
00230 
00231   // Create a BelosLinearOpWithSolve object from the Belos LOWS factory.
00232   Teuchos::RCP<Thyra::LinearOpWithSolveBase<ST> >
00233     nsA = belosLOWSFactory->createOp();
00234 
00235   // Initialize the BelosLinearOpWithSolve object with the Thyra linear operator (A)
00236   // and preconditioner (DefPrec).
00237   Thyra::initializeOp<ST>( *belosLOWSFactory, A, &*nsA );
00238 
00239   // Create a right-hand side with numRhs vectors in it and randomize it.
00240   Teuchos::RCP< Thyra::MultiVectorBase<ST> > 
00241     b = Thyra::createMembers(domain, numRhs);
00242   Thyra::seed_randomize<ST>(0);
00243   Thyra::randomize(-one, one, &*b);
00244 
00245   // Create an initial std::vector with numRhs vectors in it and initialize it to zero.
00246   Teuchos::RCP< Thyra::MultiVectorBase<ST> >
00247     x = Thyra::createMembers(domain, numRhs);
00248   Thyra::assign(&*x, zero);
00249 
00250   // Perform solve using the linear operator to get the approximate solution of Ax=b,
00251   // where b is the right-hand side and x is the left-hand side.
00252   Thyra::SolveStatus<ST> solveStatus;
00253   solveStatus = Thyra::solve( *nsA, Thyra::NONCONJ_ELE, *b, &*x );
00254 
00255   // Print out status of solve.
00256   *out << "\nBelos LOWS Status: "<< solveStatus << std::endl;
00257 
00258   //
00259   // Compute residual and ST check convergence.
00260   //
00261   std::vector<MT> norm_b(numRhs), norm_res(numRhs);
00262   Teuchos::RCP< Thyra::MultiVectorBase<ST> >
00263     y = Thyra::createMembers(domain, numRhs);
00264 
00265   // Compute the column norms of the right-hand side b.
00266   Thyra::norms_2( *b, &norm_b[0] );
00267 
00268   // Compute y=A*x, where x is the solution from the linear solver.
00269   A->apply( Thyra::NONCONJ_ELE, *x, &*y );
00270   
00271   // Compute A*x-b = y-b
00272   Thyra::update( -one, *b, &*y );
00273 
00274   // Compute the column norms of A*x-b.
00275   Thyra::norms_2( *y, &norm_res[0] );
00276 
00277   // Print out the final relative residual norms.
00278   MT rel_res = 0.0;
00279   *out << "Final relative residual norms" << std::endl;  
00280   for (int i=0; i<numRhs; ++i) {
00281     rel_res = norm_res[i]/norm_b[i];
00282     if (rel_res > maxResid)
00283       success = false;
00284     *out << "RHS " << i+1 << " : " 
00285          << std::setw(16) << std::right << rel_res << std::endl;
00286   }
00287 
00288   return ( success ? 0 : 1 );
00289 }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines