Belos Package Browser (Single Doxygen Collection)  Development
test_tfqmr_complex_diag.cpp
Go to the documentation of this file.
1 //@HEADER
2 // ************************************************************************
3 //
4 // Belos: Block Linear Solvers Package
5 // Copyright 2004 Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ************************************************************************
40 //@HEADER
41 //
42 // This driver reads a problem from a Harwell-Boeing (HB) file.
43 // The right-hand-side from the HB file is used instead of random vectors.
44 // The initial guesses are all set to zero.
45 //
46 // NOTE: No preconditioner is used in this case.
47 //
48 #include "BelosConfigDefs.hpp"
49 #include "BelosLinearProblem.hpp"
50 #include "BelosTFQMRSolMgr.hpp"
52 #include "Teuchos_CommandLineProcessor.hpp"
53 #include "Teuchos_ParameterList.hpp"
54 #include "Teuchos_StandardCatchMacros.hpp"
55 
56 #ifdef HAVE_MPI
57 #include <mpi.h>
58 #endif
59 
60 #include "MyMultiVec.hpp"
61 #include "MyOperator.hpp"
62 
63 using namespace Teuchos;
64 
65 int main(int argc, char *argv[]) {
66  //
67 #ifdef HAVE_COMPLEX
68  typedef std::complex<double> ST;
69 #elif HAVE_COMPLEX_H
70  typedef std::complex<double> ST;
71 #else
72  std::cout << "Not compiled with std::complex support." << std::endl;
73  std::cout << "End Result: TEST FAILED" << std::endl;
74  return EXIT_FAILURE;
75 #endif
76 
77  typedef ScalarTraits<ST> SCT;
78  typedef SCT::magnitudeType MT;
79  typedef Belos::MultiVec<ST> MV;
80  typedef Belos::Operator<ST> OP;
81  typedef Belos::MultiVecTraits<ST,MV> MVT;
83  ST one = SCT::one();
84  ST zero = SCT::zero();
85 
86  Teuchos::GlobalMPISession session(&argc, &argv, NULL);
87  int MyPID = session.getRank();
88  //
89  using Teuchos::RCP;
90  using Teuchos::rcp;
91 
92  bool success = false;
93  bool verbose = false;
94  try {
95  bool norm_failure = false;
96  bool proc_verbose = false;
97  bool pseudo = false; // use pseudo block TFQMR to solve this linear system.
98  int frequency = -1; // how often residuals are printed by solver
99  int blocksize = 1;
100  int numrhs = 1;
101  int maxrestarts = 15;
102  int length = 50;
103  MT tol = 1.0e-5; // relative residual tolerance
104 
105  CommandLineProcessor cmdp(false,true);
106  cmdp.setOption("verbose","quiet",&verbose,"Print messages and results.");
107  cmdp.setOption("pseudo","regular",&pseudo,"Use pseudo-block TFQMR to solve the linear systems.");
108  cmdp.setOption("frequency",&frequency,"Solvers frequency for printing residuals (#iters).");
109  cmdp.setOption("tol",&tol,"Relative residual tolerance used by TFQMR solver.");
110  cmdp.setOption("num-rhs",&numrhs,"Number of right-hand sides to be solved for.");
111  cmdp.setOption("num-restarts",&maxrestarts,"Maximum number of restarts allowed for the TFQMR solver.");
112  cmdp.setOption("blocksize",&blocksize,"Block size used by TFQMR.");
113  cmdp.setOption("subspace-length",&length,"Maximum dimension of block-subspace used by TFQMR solver.");
114  if (cmdp.parse(argc,argv) != CommandLineProcessor::PARSE_SUCCESSFUL) {
115  return EXIT_FAILURE;
116  }
117 
118  proc_verbose = verbose && (MyPID==0); /* Only print on the zero processor */
119  if (proc_verbose) {
120  std::cout << Belos::Belos_Version() << std::endl << std::endl;
121  }
122  if (!verbose)
123  frequency = -1; // reset frequency if test is not verbose
124 
125 #ifndef HAVE_BELOS_TRIUTILS
126  std::cout << "This test requires Triutils. Please configure with --enable-triutils." << std::endl;
127  if (MyPID==0) {
128  std::cout << "End Result: TEST FAILED" << std::endl;
129  }
130  return EXIT_FAILURE;
131 #endif
132 
133  // Get the data from the HB file
134  int dim=100;
135 
136  // Build the problem matrix
137  std::vector<ST> diag( dim, (ST)4.0 );
138  RCP< MyOperator<ST> > A
139  = rcp( new MyOperator<ST>( diag ) );
140  //
141  // ********Other information used by block solver***********
142  // *****************(can be user specified)******************
143  //
144  int maxits = dim/blocksize; // maximum number of iterations to run
145  //
146  ParameterList belosList;
147  belosList.set( "Num Blocks", length ); // Maximum number of blocks in Krylov factorization
148  belosList.set( "Block Size", blocksize ); // Blocksize to be used by iterative solver
149  belosList.set( "Maximum Iterations", maxits ); // Maximum number of iterations allowed
150  belosList.set( "Maximum Restarts", maxrestarts ); // Maximum number of restarts allowed
151  belosList.set( "Convergence Tolerance", tol ); // Relative convergence tolerance requested
152  if (verbose) {
153  belosList.set( "Verbosity", Belos::Errors + Belos::Warnings +
155  if (frequency > 0)
156  belosList.set( "Output Frequency", frequency );
157  }
158  else
159  belosList.set( "Verbosity", Belos::Errors + Belos::Warnings );
160  //
161  // Construct the right-hand side and solution multivectors.
162  // NOTE: The right-hand side will be constructed such that the solution is
163  // a vectors of one.
164  //
165  RCP<MyMultiVec<ST> > soln = rcp( new MyMultiVec<ST>(dim,numrhs) );
166  RCP<MyMultiVec<ST> > rhs = rcp( new MyMultiVec<ST>(dim,numrhs) );
167  MVT::MvInit( *rhs, 1.0 );
168  MVT::MvInit( *soln, zero );
169  //
170  // Construct an unpreconditioned linear problem instance.
171  //
172  RCP<Belos::LinearProblem<ST,MV,OP> > problem =
173  rcp( new Belos::LinearProblem<ST,MV,OP>( A, soln, rhs ) );
174  bool set = problem->setProblem();
175  if (set == false) {
176  if (proc_verbose)
177  std::cout << std::endl << "ERROR: Belos::LinearProblem failed to set up correctly!" << std::endl;
178  return EXIT_FAILURE;
179  }
180 
181  //
182  // *******************************************************************
183  // *************Start the TFQMR iteration***********************
184  // *******************************************************************
185  //
186  Teuchos::RCP< Belos::SolverManager<ST,MV,OP> > solver;
187  if (pseudo)
188  solver = Teuchos::rcp( new Belos::PseudoBlockTFQMRSolMgr<ST,MV,OP>( problem, Teuchos::rcp(&belosList,false) ) );
189  else
190  solver = Teuchos::rcp( new Belos::TFQMRSolMgr<ST,MV,OP>( problem, Teuchos::rcp(&belosList,false) ) );
191 
192  //
193  // **********Print out information about problem*******************
194  //
195  if (proc_verbose) {
196  std::cout << std::endl << std::endl;
197  std::cout << "Dimension of matrix: " << dim << std::endl;
198  std::cout << "Number of right-hand sides: " << numrhs << std::endl;
199  std::cout << "Block size used by solver: " << blocksize << std::endl;
200  std::cout << "Max number of TFQMR iterations: " << maxits << std::endl;
201  std::cout << "Relative residual tolerance: " << tol << std::endl;
202  std::cout << std::endl;
203  }
204  //
205  // Perform solve
206  //
207  Belos::ReturnType ret = solver->solve();
208  //
209  // Compute actual residuals.
210  //
211  RCP<MyMultiVec<ST> > temp = rcp( new MyMultiVec<ST>(dim,numrhs) );
212  OPT::Apply( *A, *soln, *temp );
213  MVT::MvAddMv( one, *rhs, -one, *temp, *temp );
214  std::vector<MT> norm_num(numrhs), norm_denom(numrhs);
215  MVT::MvNorm( *temp, norm_num );
216  MVT::MvNorm( *rhs, norm_denom );
217  for (int i=0; i<numrhs; ++i) {
218  if (proc_verbose)
219  std::cout << "Relative residual "<<i<<" : " << norm_num[i] / norm_denom[i] << std::endl;
220  if ( norm_num[i] / norm_denom[i] > tol ) {
221  norm_failure = true;
222  }
223  }
224 
225  success = ret==Belos::Converged && !norm_failure;
226  if (success) {
227  if (proc_verbose)
228  std::cout << "End Result: TEST PASSED" << std::endl;
229  } else {
230  if (proc_verbose)
231  std::cout << "End Result: TEST FAILED" << std::endl;
232  }
233  }
234  TEUCHOS_STANDARD_CATCH_STATEMENTS(verbose, std::cerr, success);
235 
236  return ( success ? EXIT_SUCCESS : EXIT_FAILURE );
237 } // end test_bl_gmres_complex_hb.cpp
Belos header file which uses auto-configuration information to include necessary C++ headers.
Class which describes the linear problem to be solved by the iterative solver.
The Belos::PseudoBlockTFQMRSolMgr provides a solver manager for the pseudo-block TFQMR linear solver.
The Belos::TFQMRSolMgr provides a solver manager for the TFQMR linear solver.
A linear system to solve, and its associated information.
Traits class which defines basic operations on multivectors.
Interface for multivectors used by Belos' linear solvers.
Class which defines basic traits for the operator type.
Alternative run-time polymorphic interface for operators.
The Belos::PseudoBlockTFQMRSolMgr provides a powerful and fully-featured solver manager over the pseu...
The Belos::TFQMRSolMgr provides a powerful and fully-featured solver manager over the TFQMR linear so...
Simple example of a user's defined Belos::MultiVec class.
Definition: MyMultiVec.hpp:66
Simple example of a user's defined Belos::Operator class.
Definition: MyOperator.hpp:66
@ StatusTestDetails
Definition: BelosTypes.hpp:261
@ TimingDetails
Definition: BelosTypes.hpp:260
@ Warnings
Definition: BelosTypes.hpp:256
ReturnType
Whether the Belos solve converged for all linear systems.
Definition: BelosTypes.hpp:155
@ Converged
Definition: BelosTypes.hpp:156
std::string Belos_Version()
int main(int argc, char *argv[])