Belos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosStatusTestUserOutput.hpp
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 
11 #ifndef BELOS_STATUS_TEST_USER_OUTPUT_HPP
12 #define BELOS_STATUS_TEST_USER_OUTPUT_HPP
13 
19 #include <vector>
20 #include "BelosConfigDefs.hpp"
21 #include "BelosTypes.hpp"
22 #include "BelosIteration.hpp"
23 
24 #include "BelosStatusTest.hpp"
25 #include "BelosStatusTestCombo.hpp"
29 
30 namespace Belos {
31 
37 template <class ScalarType, class MV, class OP>
38 class StatusTestUserOutput : public StatusTestOutput<ScalarType,MV,OP> {
39 
44 
45  public:
47 
48 
66  Teuchos::RCP<std::map<std::string,Teuchos::RCP<StatusTest<ScalarType,MV,OP> > > > taggedTests,
67  int mod = 1,
68  int printStates = Passed)
69  : printer_(printer),
70  taggedTests_(taggedTests),
71  state_(Undefined),
72  headerPrinted_(false),
73  stateTest_(printStates),
74  modTest_(mod),
75  lastNumIters_(-1),
76  comboType_(0),
77  blockSize_(1),
78  currNumRHS_(0),
79  currLSNum_(0),
80  numLSDgts_(1),
81  numIterDgts_(1)
82  {
83  // Set the input test.
84  setChild(test);
85  }
86 
88  virtual ~StatusTestUserOutput() {};
90 
92 
93 
111  {
112  TEUCHOS_TEST_FOR_EXCEPTION(iterTest_ == Teuchos::null,StatusTestError,"StatusTestUserOutput::checkStatus(): iteration test pointer is null.");
113  TEUCHOS_TEST_FOR_EXCEPTION(resTestVec_.size() == 0,StatusTestError,"StatusTestUserOutput::checkStatus(): residual test pointer is null.");
114  state_ = test_->checkStatus(solver);
115 
116  // Update some information for the header, if it has not printed or the linear system has changed.
117  LinearProblem<ScalarType,MV,OP> currProb = solver->getProblem();
118  //if (!headerPrinted_ || currLSNum_ != currProb.getLSNumber()) {
119  if (currLSNum_ != currProb.getLSNumber()) {
120  currLSNum_ = currProb.getLSNumber();
121  blockSize_ = solver->getBlockSize();
122  currIdx_ = currProb.getLSIndex();
123  currNumRHS_ = currIdx_.size();
124  numLSDgts_ = (int)std::floor((double)MVT::GetNumberVecs(*(currProb.getRHS())))+1;
125  numIterDgts_ = (int)std::floor(std::log10((double)iterTest_->getMaxIters()))+1;
126  }
127  // Print out current iteration information if it hasn't already been printed, or the status has changed
128  if (((iterTest_->getNumIters() % modTest_ == 0) && (iterTest_->getNumIters()!=lastNumIters_)) || (state_ == Passed)) {
129  lastNumIters_ = iterTest_->getNumIters();
130  if ( (state_ & stateTest_) == state_) {
131  if ( printer_->isVerbosity(StatusTestDetails) ) {
132  print( printer_->stream(StatusTestDetails) );
133  }
134  else if ( printer_->isVerbosity(Debug) ) {
135  print( printer_->stream(Debug) );
136  }
137  }
138  }
139 
140  return state_;
141  }
142 
145  return state_;
146  }
148 
149 
151 
152 
155  void setOutputManager(const Teuchos::RCP<OutputManager<ScalarType> > &printer) { printer_ = printer; }
156 
159  void setOutputFrequency(int mod) { modTest_ = mod; }
160 
166 
167  // First check to see if this test is a combination test
168  Teuchos::RCP<StatusTestCombo_t> comboTest = Teuchos::rcp_dynamic_cast<StatusTestCombo_t>(test);
169  TEUCHOS_TEST_FOR_EXCEPTION(comboTest == Teuchos::null,StatusTestError,"StatusTestUserOutput::setChild: The parameter \"test\" must be a Belos::StatusTestCombo.");
170  std::vector<Teuchos::RCP<StatusTest<ScalarType,MV,OP> > > tmpVec = comboTest->getStatusTests();
171 
172  // Get the number of tests.
173  int numTests = tmpVec.size();
174 
175  // Find the maximum iteration and residual tests
176  for (int i=0; i<numTests; ++i) {
177 
178  // Check if this is a maximum iteration test.
179  Teuchos::RCP<StatusTestMaxIters_t> tmpItrTest = Teuchos::rcp_dynamic_cast<StatusTestMaxIters_t>(tmpVec[i]);
180  if (tmpItrTest != Teuchos::null) {
181  iterTest_ = tmpItrTest;
182  continue;
183  }
184 
185  // Check if this is a single residual test
186  // this should only be the case if there are no user-specified tests when using the hard-coded solver-
187  // specific tests only
188  Teuchos::RCP<StatusTestResNorm_t> tmpResTest = Teuchos::rcp_dynamic_cast<StatusTestResNorm_t>(tmpVec[i]);
189  // If the residual status test is a single test, put in the vector
190  if (tmpResTest != Teuchos::null) {
191  resTestVec_.resize( 1 );
192  resTestVec_[0] = tmpResTest;
193  resTestNamesVec_.resize( 1 );
194  resTestNamesVec_[0] = "IMPLICIT RES";
195  continue;
196  }
197 
198  // Check if the residual test is a combination of several StatusTestResNorm objects.
199  // This should be the standard: we have a combo of the solver-specific tests and user-specific tests
200  // The user-specific tests are probably a combo again.
201  Teuchos::RCP<StatusTestCombo_t> tmpComboTest = Teuchos::rcp_dynamic_cast<StatusTestCombo_t>(tmpVec[i]);
202  TEUCHOS_TEST_FOR_EXCEPTION(tmpComboTest == Teuchos::null,StatusTestError,"StatusTestUserOutput(): test must be Belos::StatusTest[MaxIters|ResNorm|Combo].");
203  tmpVec = tmpComboTest->getStatusTests();
204  comboType_ = tmpComboTest->getComboType();
205 
206  // Add only status tests which are not in the user-specified list of tagged status tests
207  // More specifically: we want to add the implicit residual test
208  typename std::map<std::string,Teuchos::RCP<StatusTest<ScalarType,MV,OP> > >::iterator it;
209  for (size_t j=0; j<tmpVec.size(); ++j) {
210  tmpResTest = Teuchos::rcp_dynamic_cast<StatusTestResNorm_t>(tmpVec[j]);
211 
212  if(tmpResTest == Teuchos::null) continue;
213 
214  bool bFound = false;
215  for(it = taggedTests_->begin(); it != taggedTests_->end(); ++it) {
216  if(tmpVec[j] == it->second) { bFound = true; break; }
217  }
218  if(!bFound) {
219  resTestVec_.push_back(tmpResTest);
220  resTestNamesVec_.push_back("IMPLICIT RES");
221  }
222  }
223 
224  // add all tagged tests (the ordering is by the Tag names in alphabetical ordering)
225  for(it = taggedTests_->begin(); it != taggedTests_->end(); ++it) {
226  resTestVec_.push_back(it->second);
227  resTestNamesVec_.push_back(it->first);
228  }
229  }
230 
231  // Keep the pointer to the new test and reset the state to Undefined.
232  test_ = test;
233  state_ = Undefined;
234 
235  }
236 
239  return test_;
240  }
241 
244  void setSolverDesc(const std::string& solverDesc) { solverDesc_ = solverDesc; }
245 
248  void setPrecondDesc(const std::string& precondDesc) { precondDesc_ = precondDesc; }
250 
251 
253 
254 
259  void reset() {
260  state_ = Undefined;
261  test_->reset();
262  lastNumIters_ = -1;
263  headerPrinted_ = false;
264  }
265 
267 
270  void resetNumCalls() {}
271 
273 
275 
276 
278  void print(std::ostream& os, int indent = 0) const {
279  std::string ind(indent,' ');
280  std::string starLine(55,'*');
281  std::string starFront(5,'*');
282 
283  std::ios_base::fmtflags osFlags(os.flags());
284 
285  os.setf(std::ios::scientific, std::ios::floatfield);
286  os.precision(6);
287 
288  // Print header if this is the first call to this output status test.
289  if (!headerPrinted_) {
290  os << std::endl << ind << starLine << std::endl;
291  os << ind << starFront << " Belos Iterative Solver: " << solverDesc_ << std::endl;
292  if (precondDesc_ != "")
293  os << ind << starFront << " Preconditioner: " << precondDesc_ << std::endl;
294  os << ind << starFront << " Maximum Iterations: " << iterTest_->getMaxIters() << std::endl;
295  os << ind << starFront << " Block Size: " << blockSize_ << std::endl;
296  os << ind << starFront << " Status tests: " << std::endl;
297  test_->print(os,indent + 3);
298  os << ind << starLine << std::endl;
299  os.setf(std::ios_base::right, std::ios_base::adjustfield);
300  std::string indheader( 7 + numIterDgts_, ' ' );
301  os << indheader;
302  for (int i=0; i<currNumRHS_; ++i) {
303  if ( i > 0 && currIdx_[i]!=-1 ) {
304  // Put in space where 'Iter :' is in the previous lines
305  os << ind << indheader;
306  }
307  os << "[" << std::setw(numLSDgts_) << currIdx_[i]+1 << "] : ";
308  for (size_t j=0; j<resTestVec_.size(); ++j) {
309  os << std::setw(15) << resTestNamesVec_[j];
310  }
311  os << std::endl;
312  }
313  headerPrinted_ = true;
314  }
315 
316  // Print out residuals for each residual test.
317  os.setf(std::ios_base::right, std::ios_base::adjustfield);
318  std::string ind2( 7 + numIterDgts_, ' ' );
319  os << ind << "Iter " << std::setw(numIterDgts_) << iterTest_->getNumIters() << ", ";
320  for (int i=0; i<currNumRHS_; ++i) {
321  if ( i > 0 && currIdx_[i]!=-1 ) {
322  // Put in space where 'Iter :' is in the previous lines
323  os << ind << ind2;
324  }
325  os << "[" << std::setw(numLSDgts_) << currIdx_[i]+1 << "] : ";
326  for (size_t j=0; j<resTestVec_.size(); ++j) {
327  if ( resTestVec_[j]->getStatus() != Undefined && currIdx_[i]!=-1 ) {
328  // distinguish between ResNormTest derived and others
329  Teuchos::RCP<StatusTestResNorm_t> tempResTest = Teuchos::rcp_dynamic_cast<StatusTestResNorm_t>(resTestVec_[j]);
330  if(tempResTest != Teuchos::null)
331  os << std::setw(15) << (*tempResTest->getTestValue())[currIdx_[i]];
332  else {
333  if(resTestVec_[j]->getStatus() == Belos::Passed)
334  os << std::setw(15) << "Passed";
335  else if(resTestVec_[j]->getStatus() == Belos::Failed)
336  os << std::setw(15) << "Failed";
337  else os << std::setw(15) << "Undefined";
338  }
339  } else {
340  os << std::setw(15) << "---";
341  }
342  }
343  os << std::endl;
344  }
345  // reset os format
346  os.flags(osFlags);
347  }
348 
350 
351  private:
352  // Output manager.
354 
355  // Overall status test.
357 
358  // tagged tests
360 
361  // Iteration test (as passed in).
363 
365  std::vector<Teuchos::RCP<StatusTest<ScalarType,MV,OP> > > resTestVec_;
366 
368  std::vector<std::string> resTestNamesVec_;
369 
370  std::string solverDesc_;
371  std::string precondDesc_;
372  std::vector<int> currIdx_;
373  StatusType state_;
374  mutable bool headerPrinted_;
375  int stateTest_, modTest_;
376  int lastNumIters_, comboType_;
377  int blockSize_;
378  int currNumRHS_, currLSNum_;
379  int numLSDgts_, numIterDgts_;
380 };
381 
382 } // end of Belos namespace
383 
384 #endif /* BELOS_STATUS_TEST_USER_OUTPUT_HPP */
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 reset()
Informs the status test that it should reset its internal configuration to the uninitialized state...
virtual int getBlockSize() const =0
Get the blocksize to be used by the iterative solver in solving this linear problem.
Exception thrown to signal error in a status test during Belos::StatusTest::checkStatus().
void print(std::ostream &os, int indent=0) const
Output formatted description of stopping test to output stream.
A special StatusTest for printing other status tests in a simple format.
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)
An abstract class of StatusTest for stopping criteria using residual norms.
void setSolverDesc(const std::string &solverDesc)
Set a short solver description for output clarity.
Belos::StatusTest class for specifying a maximum number of iterations.
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.
StatusType
Whether the StatusTest wants iteration to stop.
Definition: BelosTypes.hpp:157
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 resetNumCalls()
Informs the outputting status test that it should reset the number of calls to zero.
A linear system to solve, and its associated information.
int getLSNumber() const
The number of linear systems that have been set.
virtual ~StatusTestUserOutput()
Destructor.
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > getChild() const
Get child test.
StatusTestUserOutput(const Teuchos::RCP< OutputManager< ScalarType > > &printer, Teuchos::RCP< StatusTest< ScalarType, MV, OP > > test, Teuchos::RCP< std::map< std::string, Teuchos::RCP< StatusTest< ScalarType, MV, OP > > > > taggedTests, int mod=1, int printStates=Passed)
Constructor.
virtual const LinearProblem< ScalarType, MV, OP > & getProblem() const =0
Get a constant reference to the linear problem.
virtual const std::vector< MagnitudeType > * getTestValue() const =0
Returns the test value, , computed in most recent call to CheckStatus.
Belos::StatusTest abstract class for specifying a residual norm stopping criteria.
void setPrecondDesc(const std::string &precondDesc)
Set a short preconditioner description for output clarity.
void setChild(Teuchos::RCP< StatusTest< ScalarType, MV, OP > > test)
Set child test, which must be a combination of a Belos::StatusTestMaxIters AND a combination of Belos...
void setOutputManager(const Teuchos::RCP< OutputManager< ScalarType > > &printer)
Set the output manager.
StatusType getStatus() const
Return the result of the most recent checkStatus call, or undefined if it has not been run...
StatusType checkStatus(Iteration< ScalarType, MV, OP > *solver)
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...
A virtual base class for StatusTest that print other status tests.
void setOutputFrequency(int mod)
Set how often the child test is printed.

Generated on Thu Oct 24 2024 09:25:34 for Belos by doxygen 1.8.5