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 
24 #ifdef HAVE_MPI
25 #include <mpi.h>
26 #endif
27 
28 // I/O for Harwell-Boeing files
29 #ifdef HAVE_BELOS_TRIUTILS
30 #include "Trilinos_Util_iohb.h"
31 #endif
32 
33 #include "MyMultiVec.hpp"
34 #include "MyBetterOperator.hpp"
35 #include "MyOperator.hpp"
36 
37 using namespace Teuchos;
38 
39 int main(int argc, char *argv[]) {
40  //
41  typedef std::complex<double> ST;
42  typedef ScalarTraits<ST> SCT;
43  typedef SCT::magnitudeType MT;
44  typedef Belos::MultiVec<ST> MV;
45  typedef Belos::Operator<ST> OP;
46  typedef Belos::MultiVecTraits<ST,MV> MVT;
48  ST one = SCT::one();
49  ST zero = SCT::zero();
50 
51  Teuchos::GlobalMPISession session(&argc, &argv, NULL);
52  //
53  using Teuchos::RCP;
54  using Teuchos::rcp;
55 
56  bool success = false;
57  bool verbose = false;
58  try {
59  int info = 0;
60  int MyPID = 0;
61  bool pseudo = false; // use pseudo block CG to solve this linear system.
62  bool norm_failure = false;
63  bool proc_verbose = false;
64  bool use_single_red = false;
65  int frequency = -1; // how often residuals are printed by solver
66  int blocksize = 1;
67  int numrhs = 1;
68  std::string filename("mhd1280b.cua");
69  MT tol = 1.0e-5; // relative residual tolerance
70 
71  CommandLineProcessor cmdp(false,true);
72  cmdp.setOption("verbose","quiet",&verbose,"Print messages and results.");
73  cmdp.setOption("pseudo","regular",&pseudo,"Use pseudo-block CG to solve the linear systems.");
74  cmdp.setOption("frequency",&frequency,"Solvers frequency for printing residuals (#iters).");
75  cmdp.setOption("filename",&filename,"Filename for Harwell-Boeing test matrix.");
76  cmdp.setOption("tol",&tol,"Relative residual tolerance used by CG solver.");
77  cmdp.setOption("num-rhs",&numrhs,"Number of right-hand sides to be solved for.");
78  cmdp.setOption("blocksize",&blocksize,"Block size used by CG .");
79  cmdp.setOption("use-single-red","use-standard-red",&use_single_red,"Use single-reduction variant of CG iteration.");
80  if (cmdp.parse(argc,argv) != CommandLineProcessor::PARSE_SUCCESSFUL) {
81  return -1;
82  }
83 
84  proc_verbose = verbose && (MyPID==0); /* Only print on the zero processor */
85  if (proc_verbose) {
86  std::cout << Belos::Belos_Version() << std::endl << std::endl;
87  }
88  if (!verbose)
89  frequency = -1; // reset frequency if test is not verbose
90 
91 
92 #ifndef HAVE_BELOS_TRIUTILS
93  std::cout << "This test requires Triutils. Please configure with --enable-triutils." << std::endl;
94  if (MyPID==0) {
95  std::cout << "End Result: TEST FAILED" << std::endl;
96  }
97  return -1;
98 #endif
99 
100  // Get the data from the HB file
101  int dim,dim2,nnz;
102  MT *dvals;
103  int *colptr,*rowind;
104  ST *cvals;
105  nnz = -1;
106  info = readHB_newmat_double(filename.c_str(),&dim,&dim2,&nnz,
107  &colptr,&rowind,&dvals);
108  if (info == 0 || nnz < 0) {
109  if (MyPID==0) {
110  std::cout << "Error reading '" << filename << "'" << std::endl;
111  std::cout << "End Result: TEST FAILED" << std::endl;
112  }
113  return -1;
114  }
115  // Convert interleaved doubles to std::complex values
116  cvals = new ST[nnz];
117  for (int ii=0; ii<nnz; ii++) {
118  cvals[ii] = ST(dvals[ii*2],dvals[ii*2+1]);
119  }
120  // Build the problem matrix
122  = rcp( new MyBetterOperator<ST>(dim,colptr,nnz,rowind,cvals) );
123  //
124  // ********Other information used by block solver***********
125  // *****************(can be user specified)******************
126  //
127  int maxits = dim/blocksize; // maximum number of iterations to run
128  //
129  ParameterList belosList;
130  belosList.set( "Block Size", blocksize ); // Blocksize to be used by iterative solver
131  belosList.set( "Maximum Iterations", maxits ); // Maximum number of iterations allowed
132  belosList.set( "Convergence Tolerance", tol ); // Relative convergence tolerance requested
133  if ((blocksize == 1) && use_single_red)
134  belosList.set( "Use Single Reduction", use_single_red ); // Use single reduction CG iteration
135  if (verbose) {
136  belosList.set( "Verbosity", Belos::Errors + Belos::Warnings +
138  if (frequency > 0)
139  belosList.set( "Output Frequency", frequency );
140  }
141  else
142  belosList.set( "Verbosity", Belos::Errors + Belos::Warnings );
143  //
144  // Construct the right-hand side and solution multivectors.
145  // NOTE: The right-hand side will be constructed such that the solution is
146  // a vectors of one.
147  //
148  RCP<MyMultiVec<ST> > soln = rcp( new MyMultiVec<ST>(dim,numrhs) );
149  RCP<MyMultiVec<ST> > rhs = rcp( new MyMultiVec<ST>(dim,numrhs) );
150  MVT::MvRandom( *soln );
151  OPT::Apply( *A, *soln, *rhs );
152  MVT::MvInit( *soln, zero );
153  //
154  // Construct an unpreconditioned linear problem instance.
155  //
157  rcp( new Belos::LinearProblem<ST,MV,OP>( A, soln, rhs ) );
158  bool set = problem->setProblem();
159  if (set == false) {
160  if (proc_verbose)
161  std::cout << std::endl << "ERROR: Belos::LinearProblem failed to set up correctly!" << std::endl;
162  return -1;
163  }
164 
165  //
166  // *******************************************************************
167  // *************Start the block CG iteration***********************
168  // *******************************************************************
169  //
171  if (pseudo)
172  solver = Teuchos::rcp( new Belos::PseudoBlockCGSolMgr<ST,MV,OP>( problem, Teuchos::rcp(&belosList,false) ) );
173  else
174  solver = Teuchos::rcp( new Belos::BlockCGSolMgr<ST,MV,OP>( problem, Teuchos::rcp(&belosList,false) ) );
175 
176  //
177  // **********Print out information about problem*******************
178  //
179  if (proc_verbose) {
180  std::cout << std::endl << std::endl;
181  std::cout << "Dimension of matrix: " << dim << std::endl;
182  std::cout << "Number of right-hand sides: " << numrhs << std::endl;
183  std::cout << "Block size used by solver: " << blocksize << std::endl;
184  std::cout << "Max number of CG iterations: " << maxits << std::endl;
185  std::cout << "Relative residual tolerance: " << tol << std::endl;
186  std::cout << std::endl;
187  }
188  //
189  // Perform solve
190  //
191  Belos::ReturnType ret = solver->solve();
192  //
193  // Compute actual residuals.
194  //
195  RCP<MyMultiVec<ST> > temp = rcp( new MyMultiVec<ST>(dim,numrhs) );
196  OPT::Apply( *A, *soln, *temp );
197  MVT::MvAddMv( one, *rhs, -one, *temp, *temp );
198  std::vector<MT> norm_num(numrhs), norm_denom(numrhs);
199  MVT::MvNorm( *temp, norm_num );
200  MVT::MvNorm( *rhs, norm_denom );
201  for (int i=0; i<numrhs; ++i) {
202  if (proc_verbose)
203  std::cout << "Relative residual "<<i<<" : " << norm_num[i] / norm_denom[i] << std::endl;
204  if ( norm_num[i] / norm_denom[i] > tol ) {
205  norm_failure = true;
206  }
207  }
208 
209  // Test achievedTol output
210  MT ach_tol = solver->achievedTol();
211  if (proc_verbose)
212  std::cout << "Achieved tol : "<<ach_tol<<std::endl;
213 
214  // Clean up.
215  delete [] dvals;
216  delete [] colptr;
217  delete [] rowind;
218  delete [] cvals;
219 
220  success = ret==Belos::Converged && !norm_failure;
221 
222  if (success) {
223  if (proc_verbose)
224  std::cout << "End Result: TEST PASSED" << std::endl;
225  } else {
226  if (proc_verbose)
227  std::cout << "End Result: TEST FAILED" << std::endl;
228  }
229  }
230  TEUCHOS_STANDARD_CATCH_STATEMENTS(verbose, std::cerr, success);
231 
232  return ( success ? EXIT_SUCCESS : EXIT_FAILURE );
233 } // 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.