60 #ifdef HAVE_BELOS_TRIUTILS 
   61 #include "Trilinos_Util_iohb.h" 
   68 using namespace Teuchos;
 
   70 int main(
int argc, 
char *argv[]) {
 
   73   typedef std::complex<double> ST;
 
   75   typedef std::complex<double> ST;
 
   77   std::cout << 
"Not compiled with std::complex support." << std::endl;
 
   78   std::cout << 
"End Result: TEST FAILED" << std::endl;
 
   83   typedef SCT::magnitudeType                MT;
 
   89   ST zero = SCT::zero();
 
   92   bool norm_failure = 
false;
 
  101   bool success = 
false;
 
  102   bool verbose = 
false;
 
  104     bool proc_verbose = 
false;
 
  107     std::string 
filename(
"mhd1280b.cua");
 
  111     cmdp.
setOption(
"verbose",
"quiet",&verbose,
"Print messages and results.");
 
  112     cmdp.
setOption(
"frequency",&frequency,
"Solvers frequency for printing residuals (#iters).");
 
  113     cmdp.
setOption(
"filename",&filename,
"Filename for Harwell-Boeing test matrix.");
 
  114     cmdp.
setOption(
"tol",&tol,
"Relative residual tolerance used by TFQMR solver.");
 
  115     cmdp.
setOption(
"num-rhs",&numrhs,
"Number of right-hand sides to be solved for.");
 
  120     proc_verbose = verbose && (MyPID==0);  
 
  128 #ifndef HAVE_BELOS_TRIUTILS 
  129     std::cout << 
"This test requires Triutils. Please configure with --enable-triutils." << std::endl;
 
  131       std::cout << 
"End Result: TEST FAILED" << std::endl;
 
  142     info = readHB_newmat_double(filename.c_str(),&dim,&dim2,&nnz,
 
  143         &colptr,&rowind,&dvals);
 
  144     if (info == 0 || nnz < 0) {
 
  146         std::cout << 
"Error reading '" << filename << 
"'" << std::endl;
 
  147         std::cout << 
"End Result: TEST FAILED" << std::endl;
 
  153     for (
int ii=0; ii<nnz; ii++) {
 
  154       cvals[ii] = ST(dvals[ii*2],dvals[ii*2+1]);
 
  166     belosList.
set( 
"Maximum Iterations", maxits );         
 
  167     belosList.
set( 
"Convergence Tolerance", tol );         
 
  172         belosList.
set( 
"Output Frequency", frequency );
 
  183     MVT::MvRandom( *soln );
 
  184     OPT::Apply( *A, *soln, *rhs );
 
  185     MVT::MvInit( *soln, zero );
 
  191     bool set = problem->setProblem();
 
  194         std::cout << std::endl << 
"ERROR:  Belos::LinearProblem failed to set up correctly!" << std::endl;
 
  208       std::cout << std::endl << std::endl;
 
  209       std::cout << 
"Dimension of matrix: " << dim << std::endl;
 
  210       std::cout << 
"Number of right-hand sides: " << numrhs << std::endl;
 
  211       std::cout << 
"Max number of TFQMR iterations: " << maxits << std::endl;
 
  212       std::cout << 
"Relative residual tolerance: " << tol << std::endl;
 
  213       std::cout << std::endl;
 
  223     OPT::Apply( *A, *soln, *temp );
 
  224     MVT::MvAddMv( one, *rhs, -one, *temp, *temp );
 
  225     std::vector<MT> norm_num(numrhs), norm_denom(numrhs);
 
  226     MVT::MvNorm( *temp, norm_num );
 
  227     MVT::MvNorm( *rhs, norm_denom );
 
  228     for (
int i=0; i<numrhs; ++i) {
 
  230         std::cout << 
"Relative residual "<<i<<
" : " << norm_num[i] / norm_denom[i] << std::endl;
 
  231       if ( norm_num[i] / norm_denom[i] > tol ) {
 
  245         std::cout << 
"End Result: TEST PASSED" << std::endl;
 
  248         std::cout << 
"End Result: TEST FAILED" << std::endl;
 
  253   return ( success ? EXIT_SUCCESS : EXIT_FAILURE );
 
std::string Belos_Version()
 
int main(int argc, char *argv[])
 
ParameterList & set(std::string const &name, T const &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's iterate() routine unti...
 
Traits class which defines basic operations on multivectors. 
 
Simple example of a user's defined Belos::MultiVec class. 
 
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. 
 
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. 
 
The Belos::TFQMRSolMgr provides a powerful and fully-featured solver manager over the TFQMR linear so...
 
Interface for multivectors used by Belos' 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's defined Belos::Operator class.