Belos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
test_tfqmr_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"
18 #include "BelosTFQMRSolMgr.hpp"
23 
24 #ifdef HAVE_MPI
25 #include <mpi.h>
26 #endif
27 
28 #include "MyMultiVec.hpp"
29 #include "MyOperator.hpp"
30 
31 using namespace Teuchos;
32 
33 int main(int argc, char *argv[]) {
34  //
35  typedef std::complex<double> ST;
36  typedef ScalarTraits<ST> SCT;
37  typedef SCT::magnitudeType MT;
38  typedef Belos::MultiVec<ST> MV;
39  typedef Belos::Operator<ST> OP;
40  typedef Belos::MultiVecTraits<ST,MV> MVT;
42  ST one = SCT::one();
43  ST zero = SCT::zero();
44 
45  Teuchos::GlobalMPISession session(&argc, &argv, NULL);
46  int MyPID = session.getRank();
47  //
48  using Teuchos::RCP;
49  using Teuchos::rcp;
50 
51  bool success = false;
52  bool verbose = false;
53  try {
54  bool norm_failure = false;
55  bool proc_verbose = false;
56  bool pseudo = false; // use pseudo block TFQMR to solve this linear system.
57  int frequency = -1; // how often residuals are printed by solver
58  int blocksize = 1;
59  int numrhs = 1;
60  int maxrestarts = 15;
61  int length = 50;
62  MT tol = 1.0e-5; // relative residual tolerance
63 
64  CommandLineProcessor cmdp(false,true);
65  cmdp.setOption("verbose","quiet",&verbose,"Print messages and results.");
66  cmdp.setOption("pseudo","regular",&pseudo,"Use pseudo-block TFQMR to solve the linear systems.");
67  cmdp.setOption("frequency",&frequency,"Solvers frequency for printing residuals (#iters).");
68  cmdp.setOption("tol",&tol,"Relative residual tolerance used by TFQMR solver.");
69  cmdp.setOption("num-rhs",&numrhs,"Number of right-hand sides to be solved for.");
70  cmdp.setOption("num-restarts",&maxrestarts,"Maximum number of restarts allowed for the TFQMR solver.");
71  cmdp.setOption("blocksize",&blocksize,"Block size used by TFQMR.");
72  cmdp.setOption("subspace-length",&length,"Maximum dimension of block-subspace used by TFQMR solver.");
73  if (cmdp.parse(argc,argv) != CommandLineProcessor::PARSE_SUCCESSFUL) {
74  return EXIT_FAILURE;
75  }
76 
77  proc_verbose = verbose && (MyPID==0); /* Only print on the zero processor */
78  if (proc_verbose) {
79  std::cout << Belos::Belos_Version() << std::endl << std::endl;
80  }
81  if (!verbose)
82  frequency = -1; // reset frequency if test is not verbose
83 
84 #ifndef HAVE_BELOS_TRIUTILS
85  std::cout << "This test requires Triutils. Please configure with --enable-triutils." << std::endl;
86  if (MyPID==0) {
87  std::cout << "End Result: TEST FAILED" << std::endl;
88  }
89  return EXIT_FAILURE;
90 #endif
91 
92  // Get the data from the HB file
93  int dim=100;
94 
95  // Build the problem matrix
96  std::vector<ST> diag( dim, (ST)4.0 );
98  = rcp( new MyOperator<ST>( diag ) );
99  //
100  // ********Other information used by block solver***********
101  // *****************(can be user specified)******************
102  //
103  int maxits = dim/blocksize; // maximum number of iterations to run
104  //
105  ParameterList belosList;
106  belosList.set( "Num Blocks", length ); // Maximum number of blocks in Krylov factorization
107  belosList.set( "Block Size", blocksize ); // Blocksize to be used by iterative solver
108  belosList.set( "Maximum Iterations", maxits ); // Maximum number of iterations allowed
109  belosList.set( "Maximum Restarts", maxrestarts ); // Maximum number of restarts allowed
110  belosList.set( "Convergence Tolerance", tol ); // Relative convergence tolerance requested
111  if (verbose) {
112  belosList.set( "Verbosity", Belos::Errors + Belos::Warnings +
114  if (frequency > 0)
115  belosList.set( "Output Frequency", frequency );
116  }
117  else
118  belosList.set( "Verbosity", Belos::Errors + Belos::Warnings );
119  //
120  // Construct the right-hand side and solution multivectors.
121  // NOTE: The right-hand side will be constructed such that the solution is
122  // a vectors of one.
123  //
124  RCP<MyMultiVec<ST> > soln = rcp( new MyMultiVec<ST>(dim,numrhs) );
125  RCP<MyMultiVec<ST> > rhs = rcp( new MyMultiVec<ST>(dim,numrhs) );
126  MVT::MvInit( *rhs, 1.0 );
127  MVT::MvInit( *soln, zero );
128  //
129  // Construct an unpreconditioned linear problem instance.
130  //
132  rcp( new Belos::LinearProblem<ST,MV,OP>( A, soln, rhs ) );
133  bool set = problem->setProblem();
134  if (set == false) {
135  if (proc_verbose)
136  std::cout << std::endl << "ERROR: Belos::LinearProblem failed to set up correctly!" << std::endl;
137  return EXIT_FAILURE;
138  }
139 
140  //
141  // *******************************************************************
142  // *************Start the TFQMR iteration***********************
143  // *******************************************************************
144  //
146  if (pseudo)
147  solver = Teuchos::rcp( new Belos::PseudoBlockTFQMRSolMgr<ST,MV,OP>( problem, Teuchos::rcp(&belosList,false) ) );
148  else
149  solver = Teuchos::rcp( new Belos::TFQMRSolMgr<ST,MV,OP>( problem, Teuchos::rcp(&belosList,false) ) );
150 
151  //
152  // **********Print out information about problem*******************
153  //
154  if (proc_verbose) {
155  std::cout << std::endl << std::endl;
156  std::cout << "Dimension of matrix: " << dim << std::endl;
157  std::cout << "Number of right-hand sides: " << numrhs << std::endl;
158  std::cout << "Block size used by solver: " << blocksize << std::endl;
159  std::cout << "Max number of TFQMR iterations: " << maxits << std::endl;
160  std::cout << "Relative residual tolerance: " << tol << std::endl;
161  std::cout << std::endl;
162  }
163  //
164  // Perform solve
165  //
166  Belos::ReturnType ret = solver->solve();
167  //
168  // Compute actual residuals.
169  //
170  RCP<MyMultiVec<ST> > temp = rcp( new MyMultiVec<ST>(dim,numrhs) );
171  OPT::Apply( *A, *soln, *temp );
172  MVT::MvAddMv( one, *rhs, -one, *temp, *temp );
173  std::vector<MT> norm_num(numrhs), norm_denom(numrhs);
174  MVT::MvNorm( *temp, norm_num );
175  MVT::MvNorm( *rhs, norm_denom );
176  for (int i=0; i<numrhs; ++i) {
177  if (proc_verbose)
178  std::cout << "Relative residual "<<i<<" : " << norm_num[i] / norm_denom[i] << std::endl;
179  if ( norm_num[i] / norm_denom[i] > tol ) {
180  norm_failure = true;
181  }
182  }
183 
184  success = ret==Belos::Converged && !norm_failure;
185  if (success) {
186  if (proc_verbose)
187  std::cout << "End Result: TEST PASSED" << std::endl;
188  } else {
189  if (proc_verbose)
190  std::cout << "End Result: TEST FAILED" << std::endl;
191  }
192  }
193  TEUCHOS_STANDARD_CATCH_STATEMENTS(verbose, std::cerr, success);
194 
195  return ( success ? EXIT_SUCCESS : EXIT_FAILURE );
196 } // end test_bl_gmres_complex_hb.cpp
std::string Belos_Version()
int main(int argc, char *argv[])
The Belos::PseudoBlockTFQMRSolMgr provides a solver manager for the pseudo-block TFQMR linear solver...
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
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)
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
The Belos::PseudoBlockTFQMRSolMgr provides a powerful and fully-featured solver manager over the pseu...
ReturnType
Whether the Belos solve converged for all linear systems.
Definition: BelosTypes.hpp:123
The Belos::TFQMRSolMgr provides a powerful and fully-featured solver manager over the TFQMR linear so...
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.
Belos header file which uses auto-configuration information to include necessary C++ headers...
The Belos::TFQMRSolMgr provides a solver manager for the TFQMR linear solver.