Belos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
test_tfqmr_complex_hb.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"
22 
23 #ifdef HAVE_MPI
24 #include <mpi.h>
25 #endif
26 
27 // I/O for Harwell-Boeing files
28 #ifdef HAVE_BELOS_TRIUTILS
29 #include "Trilinos_Util_iohb.h"
30 #endif
31 
32 #include "MyMultiVec.hpp"
33 #include "MyBetterOperator.hpp"
34 #include "MyOperator.hpp"
35 
36 using namespace Teuchos;
37 
38 int main(int argc, char *argv[]) {
39  //
40  typedef std::complex<double> ST;
41  typedef ScalarTraits<ST> SCT;
42  typedef SCT::magnitudeType MT;
43  typedef Belos::MultiVec<ST> MV;
44  typedef Belos::Operator<ST> OP;
45  typedef Belos::MultiVecTraits<ST,MV> MVT;
47  ST one = SCT::one();
48  ST zero = SCT::zero();
49 
50  int info = 0;
51  bool norm_failure = false;
52 
53  Teuchos::GlobalMPISession session(&argc, &argv, NULL);
54  //
55  int MyPID = session.getRank();
56 
57  using Teuchos::RCP;
58  using Teuchos::rcp;
59 
60  bool success = false;
61  bool verbose = false;
62  try {
63  bool proc_verbose = false;
64  int frequency = -1; // how often residuals are printed by solver
65  int numrhs = 1;
66  std::string filename("mhd1280b.cua");
67  MT tol = 1.0e-5; // relative residual tolerance
68 
69  CommandLineProcessor cmdp(false, true);
70  cmdp.setOption("verbose","quiet",&verbose,"Print messages and results.");
71  cmdp.setOption("frequency",&frequency,"Solvers frequency for printing residuals (#iters).");
72  cmdp.setOption("filename",&filename,"Filename for Harwell-Boeing test matrix.");
73  cmdp.setOption("tol",&tol,"Relative residual tolerance used by TFQMR solver.");
74  cmdp.setOption("num-rhs",&numrhs,"Number of right-hand sides to be solved for.");
75  if (cmdp.parse(argc,argv) != CommandLineProcessor::PARSE_SUCCESSFUL) {
76  return EXIT_FAILURE;
77  }
78 
79  proc_verbose = verbose && (MyPID==0); /* Only print on the zero processor */
80  if (proc_verbose) {
81  std::cout << Belos::Belos_Version() << std::endl << std::endl;
82  }
83  if (!verbose)
84  frequency = -1; // reset frequency if test is not verbose
85 
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,dim2,nnz;
97  MT *dvals;
98  int *colptr,*rowind;
99  ST *cvals;
100  nnz = -1;
101  info = readHB_newmat_double(filename.c_str(),&dim,&dim2,&nnz,
102  &colptr,&rowind,&dvals);
103  if (info == 0 || nnz < 0) {
104  if (MyPID==0) {
105  std::cout << "Error reading '" << filename << "'" << std::endl;
106  std::cout << "End Result: TEST FAILED" << std::endl;
107  }
108  return EXIT_FAILURE;
109  }
110  // Convert interleaved doubles to std::complex values
111  cvals = new ST[nnz];
112  for (int ii=0; ii<nnz; ii++) {
113  cvals[ii] = ST(dvals[ii*2],dvals[ii*2+1]);
114  }
115  // Build the problem matrix
117  = rcp( new MyBetterOperator<ST>(dim,colptr,nnz,rowind,cvals) );
118  //
119  // ********Other information used by block solver***********
120  // *****************(can be user specified)******************
121  //
122  int maxits = dim; // maximum number of iterations to run
123  //
124  ParameterList belosList;
125  belosList.set( "Maximum Iterations", maxits ); // Maximum number of iterations allowed
126  belosList.set( "Convergence Tolerance", tol ); // Relative convergence tolerance requested
127  if (verbose) {
128  belosList.set( "Verbosity", Belos::Errors + Belos::Warnings +
130  if (frequency > 0)
131  belosList.set( "Output Frequency", frequency );
132  }
133  else
134  belosList.set( "Verbosity", Belos::Errors + Belos::Warnings );
135  //
136  // Construct the right-hand side and solution multivectors.
137  // NOTE: The right-hand side will be constructed such that the solution is
138  // a vectors of one.
139  //
140  RCP<MyMultiVec<ST> > soln = rcp( new MyMultiVec<ST>(dim,numrhs) );
141  RCP<MyMultiVec<ST> > rhs = rcp( new MyMultiVec<ST>(dim,numrhs) );
142  MVT::MvRandom( *soln );
143  OPT::Apply( *A, *soln, *rhs );
144  MVT::MvInit( *soln, zero );
145  //
146  // Construct an unpreconditioned linear problem instance.
147  //
149  rcp( new Belos::LinearProblem<ST,MV,OP>( A, soln, rhs ) );
150  bool set = problem->setProblem();
151  if (set == false) {
152  if (proc_verbose)
153  std::cout << std::endl << "ERROR: Belos::LinearProblem failed to set up correctly!" << std::endl;
154  return EXIT_FAILURE;
155  }
156  //
157  // *******************************************************************
158  // *************Start the TFQMR iteration***********************
159  // *******************************************************************
160  //
161  Belos::TFQMRSolMgr<ST,MV,OP> solver( problem, rcp(&belosList,false) );
162 
163  //
164  // **********Print out information about problem*******************
165  //
166  if (proc_verbose) {
167  std::cout << std::endl << std::endl;
168  std::cout << "Dimension of matrix: " << dim << std::endl;
169  std::cout << "Number of right-hand sides: " << numrhs << std::endl;
170  std::cout << "Max number of TFQMR iterations: " << maxits << std::endl;
171  std::cout << "Relative residual tolerance: " << tol << std::endl;
172  std::cout << std::endl;
173  }
174  //
175  // Perform solve
176  //
177  Belos::ReturnType ret = solver.solve();
178  //
179  // Compute actual residuals.
180  //
181  RCP<MyMultiVec<ST> > temp = rcp( new MyMultiVec<ST>(dim,numrhs) );
182  OPT::Apply( *A, *soln, *temp );
183  MVT::MvAddMv( one, *rhs, -one, *temp, *temp );
184  std::vector<MT> norm_num(numrhs), norm_denom(numrhs);
185  MVT::MvNorm( *temp, norm_num );
186  MVT::MvNorm( *rhs, norm_denom );
187  for (int i=0; i<numrhs; ++i) {
188  if (proc_verbose)
189  std::cout << "Relative residual "<<i<<" : " << norm_num[i] / norm_denom[i] << std::endl;
190  if ( norm_num[i] / norm_denom[i] > tol ) {
191  norm_failure = true;
192  }
193  }
194 
195  // Clean up.
196  delete [] dvals;
197  delete [] colptr;
198  delete [] rowind;
199  delete [] cvals;
200 
201  success = ret==Belos::Converged && !norm_failure;
202  if (success) {
203  if (proc_verbose)
204  std::cout << "End Result: TEST PASSED" << std::endl;
205  } else {
206  if (proc_verbose)
207  std::cout << "End Result: TEST FAILED" << std::endl;
208  }
209  }
210  TEUCHOS_STANDARD_CATCH_STATEMENTS(verbose, std::cerr, success);
211 
212  return ( success ? EXIT_SUCCESS : EXIT_FAILURE );
213 } // end test_tfqmr_complex_hb.cpp
std::string Belos_Version()
int main(int argc, char *argv[])
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
ReturnType solve() override
This method performs possibly repeated calls to the underlying linear solver&#39;s iterate() routine unti...
Traits class which defines basic operations on multivectors.
Simple example of a user&#39;s defined Belos::MultiVec class.
Definition: MyMultiVec.hpp:33
std::string filename
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
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...
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.
Simple example of a user&#39;s defined Belos::Operator class.