Belos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosStatusTestResNormOutput.hpp
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 
43 #ifndef BELOS_STATUS_TEST_RESNORM_OUTPUT_HPP
44 #define BELOS_STATUS_TEST_RESNORM_OUTPUT_HPP
45 
51 #include <vector>
52 #include "BelosConfigDefs.hpp"
53 #include "BelosTypes.hpp"
54 #include "BelosIteration.hpp"
55 
56 #include "BelosStatusTest.hpp"
57 #include "BelosStatusTestCombo.hpp"
61 
62 namespace Belos {
63 
73 template <class ScalarType, class MV, class OP>
74 class StatusTestResNormOutput : public StatusTestOutput<ScalarType,MV,OP> {
75 
80 
81  public:
83 
84 
102  int mod = 1,
103  int printStates = Passed)
104  : printer_(printer),
105  state_(Undefined),
106  headerPrinted_(false),
107  stateTest_(printStates),
108  modTest_(mod),
109  lastNumIters_(-1),
110  comboType_(0),
111  numResTests_(0),
112  blockSize_(1),
113  currNumRHS_(0),
114  currLSNum_(0),
115  numLSDgts_(1),
116  numIterDgts_(1)
117  {
118  // Set the input test.
119  setChild(test);
120  }
121 
125 
127 
128 
146  {
147  TEUCHOS_TEST_FOR_EXCEPTION(iterTest_ == Teuchos::null,StatusTestError,"StatusTestResNormOutput::checkStatus(): iteration test pointer is null.");
148  TEUCHOS_TEST_FOR_EXCEPTION(resTestVec_.size() == 0,StatusTestError,"StatusTestResNormOutput::checkStatus(): residual test pointer is null.");
149  state_ = test_->checkStatus(solver);
150 
151  // Update some information for the header, if it has not printed or the linear system has changed.
152  LinearProblem<ScalarType,MV,OP> currProb = solver->getProblem();
153  //if (!headerPrinted_ || currLSNum_ != currProb.getLSNumber()) {
154  if (currLSNum_ != currProb.getLSNumber()) {
155  currLSNum_ = currProb.getLSNumber();
156  blockSize_ = solver->getBlockSize();
157  currIdx_ = currProb.getLSIndex();
158  currNumRHS_ = currIdx_.size();
159  numLSDgts_ = (int)std::floor((double)MVT::GetNumberVecs(*(currProb.getRHS())))+1;
160  numIterDgts_ = (int)std::floor(std::log10((double)iterTest_->getMaxIters()))+1;
161  }
162  // Print out current iteration information if it hasn't already been printed, or the status has changed
163  if (((iterTest_->getNumIters() % modTest_ == 0) && (iterTest_->getNumIters()!=lastNumIters_)) || (state_ == Passed)) {
164  lastNumIters_ = iterTest_->getNumIters();
165  if ( (state_ & stateTest_) == state_) {
166  if ( printer_->isVerbosity(StatusTestDetails) ) {
167  print( printer_->stream(StatusTestDetails) );
168  }
169  else if ( printer_->isVerbosity(Debug) ) {
170  print( printer_->stream(Debug) );
171  }
172  }
173  }
174 
175  return state_;
176  }
177 
180  return state_;
181  }
183 
184 
186 
187 
190  void setOutputManager(const Teuchos::RCP<OutputManager<ScalarType> > &printer) { printer_ = printer; }
191 
194  void setOutputFrequency(int mod) { modTest_ = mod; }
195 
201 
202  // First check to see if this test is a combination test
203  Teuchos::RCP<StatusTestCombo_t> comboTest = Teuchos::rcp_dynamic_cast<StatusTestCombo_t>(test);
204  TEUCHOS_TEST_FOR_EXCEPTION(comboTest == Teuchos::null,StatusTestError,"StatusTestResNormOutput(): test must be a Belos::StatusTestCombo.");
205  std::vector<Teuchos::RCP<StatusTest<ScalarType,MV,OP> > > tmpVec = comboTest->getStatusTests();
206 
207  // Get the number of tests.
208  int numTests = tmpVec.size();
209 
210  // Find the maximum iteration and residual tests
211  for (int i=0; i<numTests; ++i) {
212 
213  // Check if this is a maximum iteration test.
214  Teuchos::RCP<StatusTestMaxIters_t> tmpItrTest = Teuchos::rcp_dynamic_cast<StatusTestMaxIters_t>(tmpVec[i]);
215  if (tmpItrTest != Teuchos::null) {
216  iterTest_ = tmpItrTest;
217  continue;
218  }
219 
220  // Check if this is a single residual test
221  Teuchos::RCP<StatusTestResNorm_t> tmpResTest = Teuchos::rcp_dynamic_cast<StatusTestResNorm_t>(tmpVec[i]);
222  // If the residual status test is a single test, put in the vector
223  if (tmpResTest != Teuchos::null) {
224  numResTests_ = 1;
225  resTestVec_.resize( numResTests_ );
226  resTestVec_[0] = tmpResTest;
227  continue;
228  }
229 
230  // Check if the residual test is a combination of several StatusTestResNorm objects.
231  Teuchos::RCP<StatusTestCombo_t> tmpComboTest = Teuchos::rcp_dynamic_cast<StatusTestCombo_t>(tmpVec[i]);
232  TEUCHOS_TEST_FOR_EXCEPTION(tmpComboTest == Teuchos::null,StatusTestError,"StatusTestResNormOutput(): test must be Belos::StatusTest[MaxIters|ResNorm|Combo].");
233  tmpVec = tmpComboTest->getStatusTests();
234  comboType_ = tmpComboTest->getComboType();
235  numResTests_ = tmpVec.size();
236  resTestVec_.resize( numResTests_ );
237  for (int j=0; j<numResTests_; ++j) {
238  tmpResTest = Teuchos::rcp_dynamic_cast<StatusTestResNorm_t>(tmpVec[j]);
239  TEUCHOS_TEST_FOR_EXCEPTION(tmpResTest == Teuchos::null,StatusTestError,"StatusTestResNormOutput(): test must be a vector of Belos::StatusTestResNorm.");
240  resTestVec_[j] = tmpResTest;
241  }
242  }
243 
244  // Keep the pointer to the new test and reset the state to Undefined.
245  test_ = test;
246  state_ = Undefined;
247 
248  }
249 
252  return test_;
253  }
254 
257  void setSolverDesc(const std::string& solverDesc) { solverDesc_ = solverDesc; }
258 
261  void setPrecondDesc(const std::string& precondDesc) { precondDesc_ = precondDesc; }
263 
264 
266 
267 
272  void reset() {
273  state_ = Undefined;
274  test_->reset();
275  lastNumIters_ = -1;
276  headerPrinted_ = false;
277  }
278 
280 
283  void resetNumCalls() {}
284 
286 
288 
289 
291  void print(std::ostream& os, int indent = 0) const {
292  std::string ind(indent,' ');
293  std::string starLine(55,'*');
294  std::string starFront(5,'*');
295 
296  std::ios_base::fmtflags osFlags(os.flags());
297 
298  os.setf(std::ios::scientific, std::ios::floatfield);
299  os.precision(6);
300 
301  // Print header if this is the first call to this output status test.
302  if (!headerPrinted_) {
303  os << std::endl << ind << starLine << std::endl;
304  os << ind << starFront << " Belos Iterative Solver: " << solverDesc_ << std::endl;
305  if (precondDesc_ != "")
306  os << ind << starFront << " Preconditioner: " << precondDesc_ << std::endl;
307  os << ind << starFront << " Maximum Iterations: " << iterTest_->getMaxIters() << std::endl;
308  os << ind << starFront << " Block Size: " << blockSize_ << std::endl;
309  if (numResTests_ > 1) {
310  os << ind << starFront << " Residual Tests ("
311  << ((comboType_ == StatusTestCombo_t::OR) ? "OR" : (comboType_ == StatusTestCombo_t::AND) ? "AND" :"SEQ")
312  << "): " << std::endl;
313  } else {
314  os << ind << starFront << " Residual Test: " << std::endl;
315  }
316  for (int i=0; i<numResTests_; ++i) {
317  os << ind << starFront << " Test " << i+1 << " : " << resTestVec_[i]->description() << std::endl;
318  }
319  os << ind << starLine << std::endl;
320  headerPrinted_ = true;
321  }
322 
323  // Print out residuals for each residual test.
324  os.setf(std::ios_base::right, std::ios_base::adjustfield);
325  std::string ind2( 7 + numIterDgts_, ' ' );
326  os << ind << "Iter " << std::setw(numIterDgts_) << iterTest_->getNumIters() << ", ";
327  for (int i=0; i<currNumRHS_; ++i) {
328  if ( i > 0 && currIdx_[i]!=-1 ) {
329  // Put in space where 'Iter :' is in the previous lines
330  os << ind << ind2;
331  }
332  os << "[" << std::setw(numLSDgts_) << currIdx_[i]+1 << "] : ";
333  for (int j=0; j<numResTests_; ++j) {
334  if ( resTestVec_[j]->getStatus() != Undefined && currIdx_[i]!=-1 ) {
335  os << std::setw(15) << (*resTestVec_[j]->getTestValue())[currIdx_[i]];
336  } else {
337  os << std::setw(15) << "---";
338  }
339  }
340  os << std::endl;
341  }
342  // reset os format
343  os.flags(osFlags);
344  }
345 
347 
348  private:
349  // Output manager.
351 
352  // Overall status test.
354 
355  // Iteration test (as passed in).
357 
359  std::vector<Teuchos::RCP<StatusTestResNorm<ScalarType,MV,OP> > > resTestVec_;
360 
361  std::string solverDesc_;
362  std::string precondDesc_;
363  std::vector<int> currIdx_;
365  mutable bool headerPrinted_;
371 };
372 
373 } // end of Belos namespace
374 
375 #endif /* BELOS_STATUS_TEST_RESNORM_OUTPUT_HPP */
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > test_
Collection of types and exceptions used within the Belos solvers.
Belos&#39;s basic output manager for sending information of select verbosity levels to the appropriate ou...
void setChild(Teuchos::RCP< StatusTest< ScalarType, MV, OP > > test)
Set child test, which must be a combination of a Belos::StatusTestMaxIters AND a single or combinatio...
Teuchos::RCP< OutputManager< ScalarType > > printer_
void setOutputManager(const Teuchos::RCP< OutputManager< ScalarType > > &printer)
Set the output manager.
virtual int getBlockSize() const =0
Get the blocksize to be used by the iterative solver in solving this linear problem.
StatusType checkStatus(Iteration< ScalarType, MV, OP > *solver)
Teuchos::RCP< StatusTestMaxIters< ScalarType, MV, OP > > iterTest_
Exception thrown to signal error in a status test during Belos::StatusTest::checkStatus().
StatusTestResNormOutput(const Teuchos::RCP< OutputManager< ScalarType > > &printer, Teuchos::RCP< StatusTest< ScalarType, MV, OP > > test, int mod=1, int printStates=Passed)
Constructor.
std::vector< Teuchos::RCP< StatusTestResNorm< ScalarType, MV, OP > > > resTestVec_
Vector of residual status tests.
Virtual base class for StatusTest that printing status tests.
Pure virtual base class for defining the status testing capabilities of Belos.
ComboType getComboType() const
Return the type of combination (OR, AND, or SEQ).
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Belos::StatusTestMaxIters< ScalarType, MV, OP > StatusTestMaxIters_t
An abstract class of StatusTest for stopping criteria using residual norms.
Belos::StatusTest class for specifying a maximum number of iterations.
void reset()
Informs the status test that it should reset its internal configuration to the uninitialized state...
Pure virtual base class which describes the basic interface to the linear solver iteration.
Teuchos::RCP< const MV > getRHS() const
A pointer to the right-hand side B.
A pure virtual class for defining the status tests for the Belos iterative solvers.
void resetNumCalls()
Informs the outputting status test that it should reset the number of calls to zero.
void setSolverDesc(const std::string &solverDesc)
Set a short solver description for output clarity.
StatusType
Whether the StatusTest wants iteration to stop.
Definition: BelosTypes.hpp:189
static int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv.
Traits class which defines basic operations on multivectors.
Belos::StatusTest for logically combining several status tests.
A Belos::StatusTest class for specifying a maximum number of iterations.
const std::vector< int > & getLSIndex() const
(Zero-based) indices of the linear system(s) currently being solved.
void setPrecondDesc(const std::string &precondDesc)
Set a short preconditioner description for output clarity.
MultiVecTraits< ScalarType, MV > MVT
Belos::StatusTestCombo< ScalarType, MV, OP > StatusTestCombo_t
A linear system to solve, and its associated information.
int getLSNumber() const
The number of linear systems that have been set.
virtual const LinearProblem< ScalarType, MV, OP > & getProblem() const =0
Get a constant reference to the linear problem.
Belos::StatusTest abstract class for specifying a residual norm stopping criteria.
A special StatusTest for printing other status tests in a simple format.
StatusType getStatus() const
Return the result of the most recent checkStatus call, or undefined if it has not been run...
Belos::StatusTestResNorm< ScalarType, MV, OP > StatusTestResNorm_t
void test()
A class for extending the status testing capabilities of Belos via logical combinations.
st_vector getStatusTests()
Return the vector of status tests.
Belos header file which uses auto-configuration information to include necessary C++ headers...
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > getChild() const
Get child test.
A virtual base class for StatusTest that print other status tests.
void print(std::ostream &os, int indent=0) const
Output formatted description of stopping test to output stream.
void setOutputFrequency(int mod)
Set how often the child test is printed.