Tpetra Matrix/Vector Services Version of the Day
Tpetra_Power_Method_From_File.cpp
00001 #include <Teuchos_GlobalMPISession.hpp>
00002 #include <Teuchos_oblackholestream.hpp>
00003 #include <Teuchos_CommandLineProcessor.hpp>
00004 #include <Teuchos_Array.hpp>
00005 
00006 // I/O for Harwell-Boeing files
00007 #include <Tpetra_MatrixIO.hpp>
00008 
00009 #include "Tpetra_Power_Method.hpp"
00010 
00011 #include "Tpetra_DefaultPlatform.hpp"
00012 #include "Tpetra_Version.hpp"
00013 #include "Tpetra_Map.hpp"
00014 #include "Tpetra_MultiVector.hpp"
00015 #include "Tpetra_Vector.hpp"
00016 #include "Tpetra_CrsMatrix.hpp"
00017 
00018 int main(int argc, char *argv[]) {
00019   Teuchos::oblackholestream blackhole;
00020   Teuchos::GlobalMPISession mpiSession(&argc,&argv,&blackhole);
00021 
00022   //
00023   // Specify types used in this example
00024   // 
00025   typedef double                                                  Scalar;
00026   typedef Teuchos::ScalarTraits<Scalar>::magnitudeType            Magnitude;
00027   typedef int                                                     Ordinal;
00028   typedef Tpetra::DefaultPlatform::DefaultPlatformType            Platform;
00029   typedef Tpetra::DefaultPlatform::DefaultPlatformType::NodeType  Node;
00030   typedef Tpetra::CrsMatrix<Scalar,Ordinal,Ordinal,Node>          CrsMatrix;
00031   using Teuchos::RCP;
00032   using Teuchos::tuple;
00033 
00034   // 
00035   // Get the default communicator and node
00036   //
00037   Platform &platform = Tpetra::DefaultPlatform::getDefaultPlatform();
00038   RCP<const Teuchos::Comm<int> > comm = platform.getComm();
00039   RCP<Node>                      node = platform.getNode();
00040   const int myRank = comm->getRank();
00041 
00042   //
00043   // Get example parameters from command-line processor
00044   //  
00045   bool printMatrix = false;
00046   bool verbose = (myRank==0);
00047   int niters = 100;
00048   Magnitude tolerance = 1.0e-2;
00049   std::string filename("bcsstk14.hb");
00050   Teuchos::CommandLineProcessor cmdp(false,true);
00051   cmdp.setOption("verbose","quiet",&verbose,"Print messages and results.");
00052   cmdp.setOption("filename",&filename,"Filename for Harwell-Boeing test matrix.");
00053   cmdp.setOption("tolerance",&tolerance,"Relative residual tolerance used for solver.");
00054   cmdp.setOption("iterations",&niters,"Maximum number of iterations.");
00055   cmdp.setOption("printMatrix","noPrintMatrix",&printMatrix,"Print the full matrix after reading it.");
00056   if (cmdp.parse(argc,argv) != Teuchos::CommandLineProcessor::PARSE_SUCCESSFUL) {
00057     return -1;
00058   }
00059 
00060   // 
00061   // Say hello, print some communicator info
00062   //
00063   if (verbose) {
00064     std::cout << "\n" << Tpetra::version() << std::endl << std::endl;
00065   }
00066   std::cout << "Comm info: " << *comm;
00067 
00068   //
00069   // Read Tpetra::CrsMatrix from file
00070   //
00071   RCP<CrsMatrix> A;
00072   Tpetra::Utils::readHBMatrix(filename,comm,node,A);
00073   if (printMatrix) {
00074     RCP<Teuchos::FancyOStream> fos = Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
00075     A->describe(*fos, Teuchos::VERB_EXTREME);
00076   }
00077   else if (verbose) {
00078     std::cout << std::endl << A->description() << std::endl << std::endl;
00079   }
00080 
00081   //
00082   // Iterate
00083   //
00084   TpetraExamples::powerMethod<Scalar,Ordinal>(A, niters, tolerance, verbose);
00085 
00086   if (verbose) {
00087     std::cout << "\nEnd Result: TEST PASSED" << std::endl;
00088   }
00089   return 0;
00090 }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines