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 //
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 test instantiates the Belos classes using a std::complex scalar type
43 // and checks functionality.
44 //
45 
46 #include "BelosConfigDefs.hpp"
47 #include "BelosMVOPTester.hpp"
50 #include "BelosOutputManager.hpp"
51 
52 #ifdef HAVE_MPI
53 #include <mpi.h>
54 #endif
55 
56 // I/O for Harwell-Boeing files
57 #ifdef HAVE_BELOS_TRIUTILS
58 #include "Trilinos_Util_iohb.h"
59 #endif
60 
61 #include "MyMultiVec.hpp"
62 #include "MyOperator.hpp"
63 #include "MyBetterOperator.hpp"
64 
65 using namespace Teuchos;
66 
67 int main(int argc, char *argv[])
68 {
69  bool ierr, gerr;
70  gerr = true;
71 
72 #ifdef HAVE_MPI
73  // Initialize MPI and setup an Epetra communicator
74  MPI_Init(&argc,&argv);
75 #endif
76 
77  bool success = false;
78  bool verbose = false;
79  try {
80  int MyPID;
81 #ifdef HAVE_MPI
82  MPI_Comm_rank(MPI_COMM_WORLD, &MyPID);
83 #else
84  MyPID = 0;
85 #endif
86  (void) MyPID; // forestall "set but not used" warnings
87 
88  std::string filename("mhd1280b.cua");
89 
90  // number of global elements
91  int blockSize = 5;
92 
93  CommandLineProcessor cmdp(false,true);
94  cmdp.setOption("verbose","quiet",&verbose,"Print messages and results.");
95  cmdp.setOption("debug","quiet",&verbose,"Print messages and results.");
96  cmdp.setOption("filename",&filename,"Filename for Harwell-Boeing test matrix.");
97  if (cmdp.parse(argc,argv) != CommandLineProcessor::PARSE_SUCCESSFUL) {
98 #ifdef HAVE_MPI
99  MPI_Finalize();
100 #endif
101  return -1;
102  }
103 
104 #ifdef HAVE_COMPLEX
105  typedef std::complex<double> ST;
106 #elif HAVE_COMPLEX_H
107  typedef std::complex<double> ST;
108 #else
109  typedef double ST;
110  // no std::complex. quit with failure.
111  if (verbose && MyPID==0) {
112  std::cout << "Not compiled with std::complex support." << std::endl;
113  if (verbose && MyPID==0) {
114  std::cout << "End Result: TEST FAILED" << std::endl;
115  }
116 #ifdef HAVE_MPI
117  MPI_Finalize();
118 #endif
119  return -1;
120  }
121 #endif
122 
123  // Issue several useful typedefs;
124  typedef Belos::MultiVec<ST> MV;
125  typedef Belos::Operator<ST> OP;
126  typedef Belos::MultiVecTraits<ST,MV> MVT;
127  //typedef Belos::OperatorTraits<ST,MV,OP> OPT;
128 
129  // Create an output manager to handle the I/O from the solver
131  = rcp( new Belos::OutputManager<ST>() );
132  if (verbose) {
133  MyOM->setVerbosity( Belos::Warnings );
134  }
135 
136 #ifndef HAVE_BELOS_TRIUTILS
137  std::cout << "This test requires Triutils. Please configure with --enable-triutils." << std::endl;
138 #ifdef EPETRA_MPI
139  MPI_Finalize() ;
140 #endif
141  MyOM->print(Belos::Warnings,"End Result: TEST FAILED\n");
142  return -1;
143 #endif
144 
145  // Get the data from the HB file
146  int info;
147  int dim,dim2,nnz;
148  double *dvals;
149  int *colptr,*rowind;
150  nnz = -1;
151  info = readHB_newmat_double(filename.c_str(),&dim,&dim2,&nnz,&colptr,&rowind,&dvals);
152  if (info == 0 || nnz < 0) {
153  MyOM->stream(Belos::Warnings)
154  << "Warning reading '" << filename << "'" << std::endl
155  << "End Result: TEST FAILED" << std::endl;
156 #ifdef HAVE_MPI
157  MPI_Finalize();
158 #endif
159  return -1;
160  }
161  // Convert interleaved doubles to std::complex values
162  std::vector<ST> cvals(nnz);
163  for (int ii=0; ii<nnz; ii++) {
164  cvals[ii] = ST(dvals[ii*2],dvals[ii*2+1]);
165  }
166  // Build the problem matrix
168  = rcp( new MyBetterOperator<ST>(dim,colptr,nnz,rowind,&cvals[0]) );
169 
170 
171  // Create a MyMultiVec for cloning
172  std::vector<ScalarTraits<ST>::magnitudeType> v(blockSize);
173  RCP< MyMultiVec<ST> > ivec = rcp( new MyMultiVec<ST>(dim,blockSize) );
174  MVT::MvNorm(*ivec,v);
175 
176  // Create a MyOperator for testing against
177  RCP<MyOperator<ST> > A2 = rcp( new MyOperator<ST>(dim) );
178 
179  // test the multivector and its adapter
180  ierr = Belos::TestMultiVecTraits<ST,MV>(MyOM,ivec);
181  gerr &= ierr;
182  if (ierr) {
183  MyOM->print(Belos::Warnings, "*** MyMultiVec<std::complex> PASSED TestMultiVecTraits()\n");
184  }
185  else {
186  MyOM->print(Belos::Warnings, "*** MyMultiVec<std::complex> FAILED TestMultiVecTraits() ***\n\n");
187  }
188 
189  // test the operator and its adapter
190  ierr = Belos::TestOperatorTraits<ST,MV,OP>(MyOM,ivec,A2);
191  gerr &= ierr;
192  if (ierr) {
193  MyOM->print(Belos::Warnings,"*** MyOperator<std::complex> PASSED TestOperatorTraits()\n");
194  }
195  else {
196  MyOM->print(Belos::Warnings,"*** MyOperator<std::complex> FAILED TestOperatorTraits() ***\n\n");
197  }
198 
199  // test the operator and its adapter
200  ierr = Belos::TestOperatorTraits<ST,MV,OP>(MyOM,ivec,A1);
201  gerr &= ierr;
202  if (ierr) {
203  MyOM->print(Belos::Warnings,"*** MyBetterOperator<std::complex> PASSED TestOperatorTraits()\n");
204  }
205  else {
206  MyOM->print(Belos::Warnings,"*** MyBetterOperator<std::complex> FAILED TestOperatorTraits() ***\n\n");
207  }
208 
209  // Clean up.
210  free( dvals );
211  free( colptr );
212  free( rowind );
213 
214  success = gerr;
215  if (success) {
216  MyOM->print(Belos::Warnings,"End Result: TEST PASSED\n");
217  } else {
218  MyOM->print(Belos::Warnings,"End Result: TEST FAILED\n");
219  }
220  }
221  TEUCHOS_STANDARD_CATCH_STATEMENTS(verbose, std::cerr, success);
222 
223 #ifdef HAVE_MPI
224  MPI_Finalize();
225 #endif
226 
227  return ( success ? EXIT_SUCCESS : EXIT_FAILURE );
228 }
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:65
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:65
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.