Belos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
test_minres_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 "BelosMinresSolMgr.hpp"
21 
22 #ifdef HAVE_MPI
23 #include <mpi.h>
24 #endif
25 
26 // I/O for Harwell-Boeing files
27 #ifdef HAVE_BELOS_TRIUTILS
28 #include "Trilinos_Util_iohb.h"
29 #endif
30 
31 #include "MyMultiVec.hpp"
32 #include "MyBetterOperator.hpp"
33 #include "MyOperator.hpp"
34 
35 using namespace Teuchos;
36 
37 int main(int argc, char *argv[]) {
38  //
39  typedef std::complex<double> ST;
40  typedef ScalarTraits<ST> SCT;
41  typedef SCT::magnitudeType MT;
42  typedef Belos::MultiVec<ST> MV;
43  typedef Belos::Operator<ST> OP;
44  typedef Belos::MultiVecTraits<ST,MV> MVT;
46  ST one = SCT::one();
47  ST zero = SCT::zero();
48 
49  int info = 0;
50  int MyPID = 0;
51  bool norm_failure = false;
52 
53  Teuchos::GlobalMPISession session(&argc, &argv, NULL);
54  //
55  using Teuchos::RCP;
56  using Teuchos::rcp;
57 
58 bool success = false;
59 bool verbose = false;
60 try {
61 bool proc_verbose = false;
62  int frequency = -1; // how often residuals are printed by solver
63  int blocksize = 1;
64  int numrhs = 1;
65  std::string filename("mhd1280b.cua");
66  MT tol = 1.0e-5; // relative residual tolerance
67 
68  CommandLineProcessor cmdp(false,true);
69  cmdp.setOption("verbose","quiet",&verbose,"Print messages and results.");
70  cmdp.setOption("frequency",&frequency,"Solvers frequency for printing residuals (#iters).");
71  cmdp.setOption("filename",&filename,"Filename for Harwell-Boeing test matrix.");
72  cmdp.setOption("tol",&tol,"Relative residual tolerance used by MINRES solver.");
73  cmdp.setOption("num-rhs",&numrhs,"Number of right-hand sides to be solved for.");
74  if (cmdp.parse(argc,argv) != CommandLineProcessor::PARSE_SUCCESSFUL) {
75  return -1;
76  }
77 
78  proc_verbose = verbose && (MyPID==0); /* Only print on the zero processor */
79  if (proc_verbose) {
80  std::cout << Belos::Belos_Version() << std::endl << std::endl;
81  }
82  if (!verbose)
83  frequency = -1; // reset frequency if test is not verbose
84 
85 
86 #ifndef HAVE_BELOS_TRIUTILS
87  std::cout << "This test requires Triutils. Please configure with --enable-triutils." << std::endl;
88  if (MyPID==0) {
89  std::cout << "End Result: TEST FAILED" << std::endl;
90  }
91  return -1;
92 #endif
93 
94  // Get the data from the HB file
95  int dim,dim2,nnz;
96  MT *dvals;
97  int *colptr,*rowind;
98  ST *cvals;
99  nnz = -1;
100  info = readHB_newmat_double(filename.c_str(),&dim,&dim2,&nnz,
101  &colptr,&rowind,&dvals);
102  if (info == 0 || nnz < 0) {
103  if (MyPID==0) {
104  std::cout << "Error reading '" << filename << "'" << std::endl;
105  std::cout << "End Result: TEST FAILED" << std::endl;
106  }
107  return -1;
108  }
109  // Convert interleaved doubles to std::complex values
110  cvals = new ST[nnz];
111  for (int ii=0; ii<nnz; ii++) {
112  cvals[ii] = ST(dvals[ii*2],dvals[ii*2+1]);
113  }
114  // Build the problem matrix
116  = rcp( new MyBetterOperator<ST>(dim,colptr,nnz,rowind,cvals) );
117  //
118  // ********Other information used by block solver***********
119  // *****************(can be user specified)******************
120  //
121  int maxits = dim/blocksize; // maximum number of iterations to run
122  //
123  ParameterList belosList;
124  belosList.set( "Maximum Iterations", maxits ); // Maximum number of iterations allowed
125  belosList.set( "Convergence Tolerance", tol ); // Relative convergence tolerance requested
126  if (verbose) {
127  belosList.set( "Verbosity", Belos::Errors + Belos::Warnings +
129  if (frequency > 0)
130  belosList.set( "Output Frequency", frequency );
131  }
132  else
133  belosList.set( "Verbosity", Belos::Errors + Belos::Warnings );
134  //
135  // Construct the right-hand side and solution multivectors.
136  // NOTE: The right-hand side will be constructed such that the solution is
137  // a vectors of one.
138  //
139  RCP<MyMultiVec<ST> > soln = rcp( new MyMultiVec<ST>(dim,numrhs) );
140  RCP<MyMultiVec<ST> > rhs = rcp( new MyMultiVec<ST>(dim,numrhs) );
141  MVT::MvRandom( *soln );
142  OPT::Apply( *A, *soln, *rhs );
143  MVT::MvInit( *soln, zero );
144  //
145  // Construct an unpreconditioned linear problem instance.
146  //
148  rcp( new Belos::LinearProblem<ST,MV,OP>( A, soln, rhs ) );
149  bool set = problem->setProblem();
150  if (set == false) {
151  if (proc_verbose)
152  std::cout << std::endl << "ERROR: Belos::LinearProblem failed to set up correctly!" << std::endl;
153  return -1;
154  }
155  //
156  // *******************************************************************
157  // *************Start the MINRES iteration***********************
158  // *******************************************************************
159  //
160  Belos::MinresSolMgr<ST,MV,OP> solver( problem, rcp(&belosList,false) );
161 
162  //
163  // **********Print out information about problem*******************
164  //
165  if (proc_verbose) {
166  std::cout << std::endl << std::endl;
167  std::cout << "Dimension of matrix: " << dim << std::endl;
168  std::cout << "Number of right-hand sides: " << numrhs << std::endl;
169  std::cout << "Block size used by solver: " << blocksize << std::endl;
170  std::cout << "Max number of MINRES 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 
203 if (success) {
204  if (proc_verbose)
205  std::cout << "End Result: TEST PASSED" << std::endl;
206 } else {
207 if (proc_verbose)
208  std::cout << "End Result: TEST FAILED" << std::endl;
209 }
210 }
211 TEUCHOS_STANDARD_CATCH_STATEMENTS(verbose, std::cerr, success);
212 
213 return ( success ? EXIT_SUCCESS : EXIT_FAILURE );
214 } // end test_minres_complex_hb.cpp
Solver manager for the MINRES linear solver.
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)
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)
MINRES linear solver solution manager.
#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
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...
ReturnType solve() override
Iterate until the status test tells us to stop.
Simple example of a user&#39;s defined Belos::Operator class.