Belos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
test_bl_gmres_complex_diag.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 //
10 // This driver reads a problem from a Harwell-Boeing (HB) file.
11 // The right-hand-side from the HB file is used instead of random vectors.
12 // The initial guesses are all set to zero.
13 //
14 // NOTE: No preconditioner is used in this case.
15 //
16 #include "BelosConfigDefs.hpp"
17 #include "BelosLinearProblem.hpp"
24 
25 #ifdef HAVE_MPI
26 #include <mpi.h>
27 #endif
28 
29 #include "MyMultiVec.hpp"
30 #include "MyOperator.hpp"
31 
32 using namespace Teuchos;
33 
34 int main(int argc, char *argv[]) {
35  //
36  typedef std::complex<double> ST;
37  typedef ScalarTraits<ST> SCT;
38  typedef SCT::magnitudeType MT;
39  typedef Belos::MultiVec<ST> MV;
40  typedef Belos::Operator<ST> OP;
41  typedef Belos::MultiVecTraits<ST,MV> MVT;
43  ST one = SCT::one();
44  ST zero = SCT::zero();
45 
46  Teuchos::GlobalMPISession session(&argc, &argv, NULL);
47  int MyPID = session.getRank();
48  //
49  using Teuchos::RCP;
50  using Teuchos::rcp;
51 
52  bool success = false;
53  bool verbose = false;
54  try {
55  bool norm_failure = false;
56  bool proc_verbose = false;
57  bool pseudo = false; // use pseudo block GMRES to solve this linear system.
58  int frequency = -1; // how often residuals are printed by solver
59  int blocksize = 1;
60  int numrhs = 1;
61  int maxrestarts = 15;
62  int length = 50;
63  std::string ortho("DGKS"); // The Belos documentation obscures the fact that
64  MT tol = 1.0e-5; // relative residual tolerance
65 
66  CommandLineProcessor cmdp(false,true);
67  cmdp.setOption("verbose","quiet",&verbose,"Print messages and results.");
68  cmdp.setOption("pseudo","regular",&pseudo,"Use pseudo-block GMRES to solve the linear systems.");
69  cmdp.setOption("frequency",&frequency,"Solvers frequency for printing residuals (#iters).");
70  cmdp.setOption("tol",&tol,"Relative residual tolerance used by GMRES solver.");
71  cmdp.setOption("num-rhs",&numrhs,"Number of right-hand sides to be solved for.");
72  cmdp.setOption("num-restarts",&maxrestarts,"Maximum number of restarts allowed for the GMRES solver.");
73  cmdp.setOption("blocksize",&blocksize,"Block size used by GMRES.");
74  cmdp.setOption("subspace-length",&length,"Maximum dimension of block-subspace used by GMRES solver.");
75  cmdp.setOption("ortho-type",&ortho,"Orthogonalization type, either DGKS, ICGS or IMGS");
76  if (cmdp.parse(argc,argv) != CommandLineProcessor::PARSE_SUCCESSFUL) {
77  return EXIT_FAILURE;
78  }
79 
80  proc_verbose = verbose && (MyPID==0); /* Only print on the zero processor */
81  if (proc_verbose) {
82  std::cout << Belos::Belos_Version() << std::endl << std::endl;
83  }
84  if (!verbose)
85  frequency = -1; // reset frequency if test is not verbose
86 
87 #ifndef HAVE_BELOS_TRIUTILS
88  std::cout << "This test requires Triutils. Please configure with --enable-triutils." << std::endl;
89  if (MyPID==0) {
90  std::cout << "End Result: TEST FAILED" << std::endl;
91  }
92  return EXIT_FAILURE;
93 #endif
94 
95  // Get the data from the HB file
96  int dim=100;
97 
98  // Build the problem matrix
99  std::vector<ST> diag( dim, (ST)4.0 );
101  = rcp( new MyOperator<ST>( diag ) );
102  //
103  // ********Other information used by block solver***********
104  // *****************(can be user specified)******************
105  //
106  int maxits = dim/blocksize; // maximum number of iterations to run
107  //
108  ParameterList belosList;
109  belosList.set( "Num Blocks", length ); // Maximum number of blocks in Krylov factorization
110  belosList.set( "Block Size", blocksize ); // Blocksize to be used by iterative solver
111  belosList.set( "Maximum Iterations", maxits ); // Maximum number of iterations allowed
112  belosList.set( "Maximum Restarts", maxrestarts ); // Maximum number of restarts allowed
113  belosList.set( "Convergence Tolerance", tol ); // Relative convergence tolerance requested
114  belosList.set( "Orthogonalization", ortho ); // Orthogonalization type
115  if (verbose) {
116  belosList.set( "Verbosity", Belos::Errors + Belos::Warnings +
118  if (frequency > 0)
119  belosList.set( "Output Frequency", frequency );
120  }
121  else
122  belosList.set( "Verbosity", Belos::Errors + Belos::Warnings );
123  //
124  // Construct the right-hand side and solution multivectors.
125  // NOTE: The right-hand side will be constructed such that the solution is
126  // a vectors of one.
127  //
128  RCP<MyMultiVec<ST> > soln = rcp( new MyMultiVec<ST>(dim,numrhs) );
129  RCP<MyMultiVec<ST> > rhs = rcp( new MyMultiVec<ST>(dim,numrhs) );
130  MVT::MvInit( *rhs, 1.0 );
131  MVT::MvInit( *soln, zero );
132  //
133  // Construct an unpreconditioned linear problem instance.
134  //
136  rcp( new Belos::LinearProblem<ST,MV,OP>( A, soln, rhs ) );
137  bool set = problem->setProblem();
138  if (set == false) {
139  if (proc_verbose)
140  std::cout << std::endl << "ERROR: Belos::LinearProblem failed to set up correctly!" << std::endl;
141  return EXIT_FAILURE;
142  }
143 
144  // Use a debugging status test to save absolute residual history.
145  // Debugging status tests are peer to the native status tests that are called whenever convergence is checked.
147 
148  //
149  // *******************************************************************
150  // *************Start the block Gmres iteration***********************
151  // *******************************************************************
152  //
154  if (pseudo)
155  solver = Teuchos::rcp( new Belos::PseudoBlockGmresSolMgr<ST,MV,OP>( problem, Teuchos::rcp(&belosList,false) ) );
156  else
157  solver = Teuchos::rcp( new Belos::BlockGmresSolMgr<ST,MV,OP>( problem, Teuchos::rcp(&belosList,false) ) );
158 
159  solver->setDebugStatusTest( Teuchos::rcp(&debugTest, false) );
160 
161  //
162  // **********Print out information about problem*******************
163  //
164  if (proc_verbose) {
165  std::cout << std::endl << std::endl;
166  std::cout << "Dimension of matrix: " << dim << std::endl;
167  std::cout << "Number of right-hand sides: " << numrhs << std::endl;
168  std::cout << "Block size used by solver: " << blocksize << std::endl;
169  std::cout << "Max number of Gmres iterations: " << maxits << std::endl;
170  std::cout << "Relative residual tolerance: " << tol << std::endl;
171  std::cout << std::endl;
172  }
173  //
174  // Perform solve
175  //
176  Belos::ReturnType ret = solver->solve();
177  //
178  // Compute actual residuals.
179  //
180  RCP<MyMultiVec<ST> > temp = rcp( new MyMultiVec<ST>(dim,numrhs) );
181  OPT::Apply( *A, *soln, *temp );
182  MVT::MvAddMv( one, *rhs, -one, *temp, *temp );
183  std::vector<MT> norm_num(numrhs), norm_denom(numrhs);
184  MVT::MvNorm( *temp, norm_num );
185  MVT::MvNorm( *rhs, norm_denom );
186  for (int i=0; i<numrhs; ++i) {
187  if (proc_verbose)
188  std::cout << "Relative residual "<<i<<" : " << norm_num[i] / norm_denom[i] << std::endl;
189  if ( norm_num[i] / norm_denom[i] > tol ) {
190  norm_failure = true;
191  }
192  }
193 
194  // Print absolute residual norm logging.
195  const std::vector<MT> residualLog = debugTest.getLogResNorm();
196  if (numrhs==1 && proc_verbose && residualLog.size())
197  {
198  std::cout << "Absolute residual 2-norm [ " << residualLog.size() << " ] : ";
199  for (unsigned int i=0; i<residualLog.size(); i++)
200  std::cout << residualLog[i] << " ";
201  std::cout << std::endl;
202  std::cout << "Final abs 2-norm / rhs 2-norm : " << residualLog[residualLog.size()-1] / norm_denom[0] << std::endl;
203  }
204 
205  success = ret==Belos::Converged && !norm_failure;
206  if (success) {
207  if (proc_verbose)
208  std::cout << "End Result: TEST PASSED" << std::endl;
209  } else {
210  if (proc_verbose)
211  std::cout << "End Result: TEST FAILED" << std::endl;
212  }
213  }
214  TEUCHOS_STANDARD_CATCH_STATEMENTS(verbose, std::cerr, success);
215 
216  return ( success ? EXIT_SUCCESS : EXIT_FAILURE );
217 } // end test_bl_gmres_complex_hb.cpp
std::string Belos_Version()
int main(int argc, char *argv[])
const std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > & getLogResNorm() const
Returns the log of the absolute residual norm from the iteration.
A Belos::StatusTest debugging class for storing the absolute residual norms generated during a solve...
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Interface to Block GMRES and Flexible GMRES.
Belos::StatusTest debugging class for storing the absolute residual norms generated during a solve...
Traits class which defines basic operations on multivectors.
Simple example of a user&#39;s defined Belos::MultiVec class.
Definition: MyMultiVec.hpp:33
Alternative run-time polymorphic interface for operators.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
The Belos::BlockGmresSolMgr provides a solver manager for the BlockGmres linear solver.
void setOption(const char option_true[], const char option_false[], bool *option_val, const char documentation[]=NULL)
#define TEUCHOS_STANDARD_CATCH_STATEMENTS(VERBOSE, ERR_STREAM, SUCCESS_FLAG)
A linear system to solve, and its associated information.
const double tol
Class which describes the linear problem to be solved by the iterative solver.
EParseCommandLineReturn parse(int argc, char *argv[], std::ostream *errout=&std::cerr) const
Interface to standard and &quot;pseudoblock&quot; GMRES.
ReturnType
Whether the Belos solve converged for all linear systems.
Definition: BelosTypes.hpp:123
Simple example of a user&#39;s defined Belos::Operator class.
Definition: MyOperator.hpp:33
Interface for multivectors used by Belos&#39; linear solvers.
Class which defines basic traits for the operator type.
The Belos::PseudoBlockGmresSolMgr provides a solver manager for the BlockGmres linear solver...
Belos header file which uses auto-configuration information to include necessary C++ headers...