Belos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
cxx_main_complex.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 test instantiates the Belos classes using a std::complex scalar type
11 // and checks functionality.
12 //
13 
14 #include "BelosConfigDefs.hpp"
15 #include "BelosMVOPTester.hpp"
18 #include "BelosOutputManager.hpp"
19 
20 #ifdef HAVE_MPI
21 #include <mpi.h>
22 #endif
23 
24 // I/O for Harwell-Boeing files
25 #ifdef HAVE_BELOS_TRIUTILS
26 #include "Trilinos_Util_iohb.h"
27 #endif
28 
29 #include "MyMultiVec.hpp"
30 #include "MyOperator.hpp"
31 #include "MyBetterOperator.hpp"
32 
33 using namespace Teuchos;
34 
35 int main(int argc, char *argv[])
36 {
37  bool ierr, gerr;
38  gerr = true;
39 
40 #ifdef HAVE_MPI
41  // Initialize MPI and setup an Epetra communicator
42  MPI_Init(&argc,&argv);
43 #endif
44 
45  bool success = false;
46  bool verbose = false;
47  try {
48  int MyPID;
49 #ifdef HAVE_MPI
50  MPI_Comm_rank(MPI_COMM_WORLD, &MyPID);
51 #else
52  MyPID = 0;
53 #endif
54  (void) MyPID; // forestall "set but not used" warnings
55 
56  std::string filename("mhd1280b.cua");
57 
58  // number of global elements
59  int blockSize = 5;
60 
61  CommandLineProcessor cmdp(false,true);
62  cmdp.setOption("verbose","quiet",&verbose,"Print messages and results.");
63  cmdp.setOption("debug","quiet",&verbose,"Print messages and results.");
64  cmdp.setOption("filename",&filename,"Filename for Harwell-Boeing test matrix.");
65  if (cmdp.parse(argc,argv) != CommandLineProcessor::PARSE_SUCCESSFUL) {
66 #ifdef HAVE_MPI
67  MPI_Finalize();
68 #endif
69  return -1;
70  }
71 
72  typedef std::complex<double> ST;
73 
74  // Issue several useful typedefs;
75  typedef Belos::MultiVec<ST> MV;
76  typedef Belos::Operator<ST> OP;
77  typedef Belos::MultiVecTraits<ST,MV> MVT;
78  //typedef Belos::OperatorTraits<ST,MV,OP> OPT;
79 
80  // Create an output manager to handle the I/O from the solver
82  = rcp( new Belos::OutputManager<ST>() );
83  if (verbose) {
84  MyOM->setVerbosity( Belos::Warnings );
85  }
86 
87 #ifndef HAVE_BELOS_TRIUTILS
88  std::cout << "This test requires Triutils. Please configure with --enable-triutils." << std::endl;
89 #ifdef EPETRA_MPI
90  MPI_Finalize() ;
91 #endif
92  MyOM->print(Belos::Warnings,"End Result: TEST FAILED\n");
93  return -1;
94 #endif
95 
96  // Get the data from the HB file
97  int info;
98  int dim,dim2,nnz;
99  double *dvals;
100  int *colptr,*rowind;
101  nnz = -1;
102  info = readHB_newmat_double(filename.c_str(),&dim,&dim2,&nnz,&colptr,&rowind,&dvals);
103  if (info == 0 || nnz < 0) {
104  MyOM->stream(Belos::Warnings)
105  << "Warning reading '" << filename << "'" << std::endl
106  << "End Result: TEST FAILED" << std::endl;
107 #ifdef HAVE_MPI
108  MPI_Finalize();
109 #endif
110  return -1;
111  }
112  // Convert interleaved doubles to std::complex values
113  std::vector<ST> cvals(nnz);
114  for (int ii=0; ii<nnz; ii++) {
115  cvals[ii] = ST(dvals[ii*2],dvals[ii*2+1]);
116  }
117  // Build the problem matrix
119  = rcp( new MyBetterOperator<ST>(dim,colptr,nnz,rowind,&cvals[0]) );
120 
121 
122  // Create a MyMultiVec for cloning
123  std::vector<ScalarTraits<ST>::magnitudeType> v(blockSize);
124  RCP< MyMultiVec<ST> > ivec = rcp( new MyMultiVec<ST>(dim,blockSize) );
125  MVT::MvNorm(*ivec,v);
126 
127  // Create a MyOperator for testing against
128  RCP<MyOperator<ST> > A2 = rcp( new MyOperator<ST>(dim) );
129 
130  // test the multivector and its adapter
131  ierr = Belos::TestMultiVecTraits<ST,MV>(MyOM,ivec);
132  gerr &= ierr;
133  if (ierr) {
134  MyOM->print(Belos::Warnings, "*** MyMultiVec<std::complex> PASSED TestMultiVecTraits()\n");
135  }
136  else {
137  MyOM->print(Belos::Warnings, "*** MyMultiVec<std::complex> FAILED TestMultiVecTraits() ***\n\n");
138  }
139 
140  // test the operator and its adapter
141  ierr = Belos::TestOperatorTraits<ST,MV,OP>(MyOM,ivec,A2);
142  gerr &= ierr;
143  if (ierr) {
144  MyOM->print(Belos::Warnings,"*** MyOperator<std::complex> PASSED TestOperatorTraits()\n");
145  }
146  else {
147  MyOM->print(Belos::Warnings,"*** MyOperator<std::complex> FAILED TestOperatorTraits() ***\n\n");
148  }
149 
150  // test the operator and its adapter
151  ierr = Belos::TestOperatorTraits<ST,MV,OP>(MyOM,ivec,A1);
152  gerr &= ierr;
153  if (ierr) {
154  MyOM->print(Belos::Warnings,"*** MyBetterOperator<std::complex> PASSED TestOperatorTraits()\n");
155  }
156  else {
157  MyOM->print(Belos::Warnings,"*** MyBetterOperator<std::complex> FAILED TestOperatorTraits() ***\n\n");
158  }
159 
160  // Clean up.
161  free( dvals );
162  free( colptr );
163  free( rowind );
164 
165  success = gerr;
166  if (success) {
167  MyOM->print(Belos::Warnings,"End Result: TEST PASSED\n");
168  } else {
169  MyOM->print(Belos::Warnings,"End Result: TEST FAILED\n");
170  }
171  }
172  TEUCHOS_STANDARD_CATCH_STATEMENTS(verbose, std::cerr, success);
173 
174 #ifdef HAVE_MPI
175  MPI_Finalize();
176 #endif
177 
178  return ( success ? EXIT_SUCCESS : EXIT_FAILURE );
179 }
Belos&#39;s basic output manager for sending information of select verbosity levels to the appropriate ou...
Class which manages the output and verbosity of the Belos solvers.
int main(int argc, char *argv[])
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)
#define TEUCHOS_STANDARD_CATCH_STATEMENTS(VERBOSE, ERR_STREAM, SUCCESS_FLAG)
Test routines for MultiVecTraits and OperatorTraits conformity.
EParseCommandLineReturn parse(int argc, char *argv[], std::ostream *errout=&std::cerr) const
Simple example of a user&#39;s defined Belos::Operator class.
Definition: MyOperator.hpp:33
Interface for multivectors used by Belos&#39; linear solvers.
Belos header file which uses auto-configuration information to include necessary C++ headers...
Simple example of a user&#39;s defined Belos::Operator class.