Belos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
test_belos_projected_least_squares_solver.cpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Belos: Block Linear Solvers Package
4 //
5 // Copyright 2004-2016 NTESS and the Belos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
11 #include <Teuchos_DefaultComm.hpp>
16 
17 int
18 main (int argc, char *argv[])
19 {
25  using Teuchos::RCP;
26  using Teuchos::rcp;
27  using std::endl;
28 
29  typedef double scalar_type;
33 
34  Teuchos::oblackholestream blackHole;
35  // Initialize MPI using Teuchos wrappers, if Trilinos was built with
36  // MPI support. Otherwise, initialize a communicator with one
37  // process.
38  Teuchos::GlobalMPISession mpiSession (&argc, &argv, &blackHole);
39  const int myRank = mpiSession.getRank();
40  // Output stream only prints on MPI Proc 0.
41  std::ostream& out = (myRank == 0) ? std::cout : blackHole;
42 
43  // Command-line arguments
44  std::string robustnessLevel ("None");
45  bool testBlockGivens = false;
46  bool testGivensRotations = false;
47  bool verbose = false;
48  bool debug = false;
49  int testProblemSize = 10;
50 
51  // Parse command-line arguments
52  CommandLineProcessor cmdp (false,true);
53  cmdp.setOption ("robustness", &robustnessLevel,
54  "Robustness level: \"None\", \"Some\", or \"Lots\".");
55  cmdp.setOption ("testGivensRotations", "dontTestGivensRotations",
56  &testGivensRotations,
57  "Test the implementation of Givens rotations.");
58  cmdp.setOption ("testBlockGivens", "dontTestBlockGivens", &testBlockGivens,
59  "Test the panel version of the Givens rotations - based "
60  "update.");
61  cmdp.setOption ("verbose", "quiet", &verbose, "Print messages and results.");
62  cmdp.setOption ("debug", "release", &debug, "Print copious debug output.");
63  cmdp.setOption ("testProblemSize", &testProblemSize,
64  "Number of columns in the projected least-squares test "
65  "problem.");
66  if (cmdp.parse (argc,argv) != CommandLineProcessor::PARSE_SUCCESSFUL) {
67  out << "End Result: TEST FAILED" << endl;
68  return EXIT_FAILURE;
69  }
70  // Verbose output stream only prints on MPI Proc 0, and only in
71  // verbose mode.
72  std::ostream& verboseOut = verbose ? out : blackHole;
73 
74  // Robustness level for triangular solves.
75  const ERobustness robustness = robustnessStringToEnum (robustnessLevel);
76  verboseOut << "-- Robustness level: " << robustnessLevel << endl;
77 
78  bool success = true; // Innocent until proven guilty.
79 
80  // Seed the pseudorandom number generator with the same seed each
81  // time, to ensure repeatable results.
82  STS::seedrandom (0);
83  ProjectedLeastSquaresSolver<scalar_type> solver (out, robustness);
84 
85  if (testGivensRotations) {
86  solver.testGivensRotations (verboseOut);
87  }
88  if (testProblemSize > 0) {
89  const bool extraVerbose = debug;
90  success = success &&
91  solver.testUpdateColumn (verboseOut, testProblemSize,
92  testBlockGivens, extraVerbose);
93  success = success &&
94  solver.testTriangularSolves (verboseOut, testProblemSize,
95  robustness, extraVerbose);
96  }
97 
98  if (success) {
99  out << "End Result: TEST PASSED" << endl;
100  return EXIT_SUCCESS;
101  } else {
102  out << "End Result: TEST FAILED" << endl;
103  return EXIT_FAILURE;
104  }
105 }
int main(int argc, char *argv[])
std::string robustnessEnumToString(const ERobustness x)
Convert the given ERobustness enum value to a string.
ERobustness robustnessStringToEnum(const std::string &x)
Convert the given robustness string value to an ERobustness enum.
Methods for solving GMRES&#39; projected least-squares problem.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
ERobustness
Robustness level of projected least-squares solver operations.
Methods for solving GMRES&#39; projected least-squares problem.