Belos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
test_bl_cg_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 "BelosBlockCGSolMgr.hpp"
23 #include "Tpetra_Util_iohb.h"
24 
25 #ifdef HAVE_MPI
26 #include <mpi.h>
27 #endif
28 
29 #include "MyMultiVec.hpp"
30 #include "MyBetterOperator.hpp"
31 #include "MyOperator.hpp"
32 
33 using namespace Teuchos;
34 
35 int main(int argc, char *argv[]) {
36  //
37  typedef std::complex<double> ST;
38  typedef ScalarTraits<ST> SCT;
39  typedef SCT::magnitudeType MT;
40  typedef Belos::MultiVec<ST> MV;
41  typedef Belos::Operator<ST> OP;
42  typedef Belos::MultiVecTraits<ST,MV> MVT;
44  ST one = SCT::one();
45  ST zero = SCT::zero();
46 
47  Teuchos::GlobalMPISession session(&argc, &argv, NULL);
48  //
49  using Teuchos::RCP;
50  using Teuchos::rcp;
51 
52  bool success = false;
53  bool verbose = false;
54  try {
55  int info = 0;
56  int MyPID = 0;
57  bool pseudo = false; // use pseudo block CG to solve this linear system.
58  bool norm_failure = false;
59  bool proc_verbose = false;
60  bool use_single_red = false;
61  int frequency = -1; // how often residuals are printed by solver
62  int blocksize = 1;
63  int numrhs = 1;
64  std::string filename("mhd1280b.cua");
65  MT tol = 1.0e-5; // relative residual tolerance
66 
67  CommandLineProcessor cmdp(false,true);
68  cmdp.setOption("verbose","quiet",&verbose,"Print messages and results.");
69  cmdp.setOption("pseudo","regular",&pseudo,"Use pseudo-block CG to solve the linear systems.");
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 CG solver.");
73  cmdp.setOption("num-rhs",&numrhs,"Number of right-hand sides to be solved for.");
74  cmdp.setOption("blocksize",&blocksize,"Block size used by CG .");
75  cmdp.setOption("use-single-red","use-standard-red",&use_single_red,"Use single-reduction variant of CG iteration.");
76  if (cmdp.parse(argc,argv) != CommandLineProcessor::PARSE_SUCCESSFUL) {
77  return -1;
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  // Get the data from the HB file
88  int dim,dim2,nnz;
89  MT *dvals;
90  int *colptr,*rowind;
91  ST *cvals;
92  nnz = -1;
93  info = readHB_newmat_double(filename.c_str(),&dim,&dim2,&nnz,
94  &colptr,&rowind,&dvals);
95  if (info == 0 || nnz < 0) {
96  if (MyPID==0) {
97  std::cout << "Error reading '" << filename << "'" << std::endl;
98  std::cout << "End Result: TEST FAILED" << std::endl;
99  }
100  return -1;
101  }
102  // Convert interleaved doubles to std::complex values
103  cvals = new ST[nnz];
104  for (int ii=0; ii<nnz; ii++) {
105  cvals[ii] = ST(dvals[ii*2],dvals[ii*2+1]);
106  }
107  // Build the problem matrix
109  = rcp( new MyBetterOperator<ST>(dim,colptr,nnz,rowind,cvals) );
110  //
111  // ********Other information used by block solver***********
112  // *****************(can be user specified)******************
113  //
114  int maxits = dim/blocksize; // maximum number of iterations to run
115  //
116  ParameterList belosList;
117  belosList.set( "Block Size", blocksize ); // Blocksize to be used by iterative solver
118  belosList.set( "Maximum Iterations", maxits ); // Maximum number of iterations allowed
119  belosList.set( "Convergence Tolerance", tol ); // Relative convergence tolerance requested
120  if ((blocksize == 1) && use_single_red)
121  belosList.set( "Use Single Reduction", use_single_red ); // Use single reduction CG iteration
122  if (verbose) {
123  belosList.set( "Verbosity", Belos::Errors + Belos::Warnings +
125  if (frequency > 0)
126  belosList.set( "Output Frequency", frequency );
127  }
128  else
129  belosList.set( "Verbosity", Belos::Errors + Belos::Warnings );
130  //
131  // Construct the right-hand side and solution multivectors.
132  // NOTE: The right-hand side will be constructed such that the solution is
133  // a vectors of one.
134  //
135  RCP<MyMultiVec<ST> > soln = rcp( new MyMultiVec<ST>(dim,numrhs) );
136  RCP<MyMultiVec<ST> > rhs = rcp( new MyMultiVec<ST>(dim,numrhs) );
137  MVT::MvRandom( *soln );
138  OPT::Apply( *A, *soln, *rhs );
139  MVT::MvInit( *soln, zero );
140  //
141  // Construct an unpreconditioned linear problem instance.
142  //
144  rcp( new Belos::LinearProblem<ST,MV,OP>( A, soln, rhs ) );
145  bool set = problem->setProblem();
146  if (set == false) {
147  if (proc_verbose)
148  std::cout << std::endl << "ERROR: Belos::LinearProblem failed to set up correctly!" << std::endl;
149  return -1;
150  }
151 
152  //
153  // *******************************************************************
154  // *************Start the block CG iteration***********************
155  // *******************************************************************
156  //
158  if (pseudo)
159  solver = Teuchos::rcp( new Belos::PseudoBlockCGSolMgr<ST,MV,OP>( problem, Teuchos::rcp(&belosList,false) ) );
160  else
161  solver = Teuchos::rcp( new Belos::BlockCGSolMgr<ST,MV,OP>( problem, Teuchos::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 << "Block size used by solver: " << blocksize << std::endl;
171  std::cout << "Max number of CG iterations: " << maxits << std::endl;
172  std::cout << "Relative residual tolerance: " << tol << std::endl;
173  std::cout << std::endl;
174  }
175  //
176  // Perform solve
177  //
178  Belos::ReturnType ret = solver->solve();
179  //
180  // Compute actual residuals.
181  //
182  RCP<MyMultiVec<ST> > temp = rcp( new MyMultiVec<ST>(dim,numrhs) );
183  OPT::Apply( *A, *soln, *temp );
184  MVT::MvAddMv( one, *rhs, -one, *temp, *temp );
185  std::vector<MT> norm_num(numrhs), norm_denom(numrhs);
186  MVT::MvNorm( *temp, norm_num );
187  MVT::MvNorm( *rhs, norm_denom );
188  for (int i=0; i<numrhs; ++i) {
189  if (proc_verbose)
190  std::cout << "Relative residual "<<i<<" : " << norm_num[i] / norm_denom[i] << std::endl;
191  if ( norm_num[i] / norm_denom[i] > tol ) {
192  norm_failure = true;
193  }
194  }
195 
196  // Test achievedTol output
197  MT ach_tol = solver->achievedTol();
198  if (proc_verbose)
199  std::cout << "Achieved tol : "<<ach_tol<<std::endl;
200 
201  // Clean up.
202  delete [] dvals;
203  delete [] colptr;
204  delete [] rowind;
205  delete [] cvals;
206 
207  success = ret==Belos::Converged && !norm_failure;
208 
209  if (success) {
210  if (proc_verbose)
211  std::cout << "End Result: TEST PASSED" << std::endl;
212  } else {
213  if (proc_verbose)
214  std::cout << "End Result: TEST FAILED" << std::endl;
215  }
216  }
217  TEUCHOS_STANDARD_CATCH_STATEMENTS(verbose, std::cerr, success);
218 
219  return ( success ? EXIT_SUCCESS : EXIT_FAILURE );
220 } // end test_bl_cg_complex_hb.cpp
std::string Belos_Version()
int main(int argc, char *argv[])
The Belos::PseudoBlockCGSolMgr provides a solver manager for the BlockCG linear solver.
The Belos::BlockCGSolMgr provides a powerful and fully-featured solver manager over the CG and BlockC...
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.
The Belos::PseudoBlockCGSolMgr provides a powerful and fully-featured solver manager over the pseudo-...
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::BlockCGSolMgr provides a solver manager for the BlockCG linear solver. ...
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...
Simple example of a user&#39;s defined Belos::Operator class.