Belos  Version of the Day
 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 // 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_RESNORM_OUTPUT_HPP
12 #define BELOS_STATUS_TEST_RESNORM_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 
41 template <class ScalarType, class MV, class OP>
42 class StatusTestResNormOutput : public StatusTestOutput<ScalarType,MV,OP> {
43 
48 
49  public:
51 
52 
70  int mod = 1,
71  int printStates = Passed)
72  : printer_(printer),
73  state_(Undefined),
74  headerPrinted_(false),
75  stateTest_(printStates),
76  modTest_(mod),
77  lastNumIters_(-1),
78  comboType_(0),
79  numResTests_(0),
80  blockSize_(1),
81  currNumRHS_(0),
82  currLSNum_(0),
83  numLSDgts_(1),
84  numIterDgts_(1)
85  {
86  // Set the input test.
87  setChild(test);
88  }
89 
93 
95 
96 
114  {
115  TEUCHOS_TEST_FOR_EXCEPTION(iterTest_ == Teuchos::null,StatusTestError,"StatusTestResNormOutput::checkStatus(): iteration test pointer is null.");
116  TEUCHOS_TEST_FOR_EXCEPTION(resTestVec_.size() == 0,StatusTestError,"StatusTestResNormOutput::checkStatus(): residual test pointer is null.");
117  state_ = test_->checkStatus(solver);
118 
119  // Update some information for the header, if it has not printed or the linear system has changed.
120  LinearProblem<ScalarType,MV,OP> currProb = solver->getProblem();
121  //if (!headerPrinted_ || currLSNum_ != currProb.getLSNumber()) {
122  if (currLSNum_ != currProb.getLSNumber()) {
123  currLSNum_ = currProb.getLSNumber();
124  blockSize_ = solver->getBlockSize();
125  currIdx_ = currProb.getLSIndex();
126  currNumRHS_ = currIdx_.size();
127  numLSDgts_ = (int)std::floor((double)MVT::GetNumberVecs(*(currProb.getRHS())))+1;
128  numIterDgts_ = (int)std::floor(std::log10((double)iterTest_->getMaxIters()))+1;
129  }
130  // Print out current iteration information if it hasn't already been printed, or the status has changed
131  if (((iterTest_->getNumIters() % modTest_ == 0) && (iterTest_->getNumIters()!=lastNumIters_)) || (state_ == Passed)) {
132  lastNumIters_ = iterTest_->getNumIters();
133  if ( (state_ & stateTest_) == state_) {
134  if ( printer_->isVerbosity(StatusTestDetails) ) {
135  print( printer_->stream(StatusTestDetails) );
136  }
137  else if ( printer_->isVerbosity(Debug) ) {
138  print( printer_->stream(Debug) );
139  }
140  }
141  }
142 
143  return state_;
144  }
145 
148  return state_;
149  }
151 
152 
154 
155 
158  void setOutputManager(const Teuchos::RCP<OutputManager<ScalarType> > &printer) { printer_ = printer; }
159 
162  void setOutputFrequency(int mod) { modTest_ = mod; }
163 
169 
170  // First check to see if this test is a combination test
171  Teuchos::RCP<StatusTestCombo_t> comboTest = Teuchos::rcp_dynamic_cast<StatusTestCombo_t>(test);
172  TEUCHOS_TEST_FOR_EXCEPTION(comboTest == Teuchos::null,StatusTestError,"StatusTestResNormOutput(): test must be a Belos::StatusTestCombo.");
173  std::vector<Teuchos::RCP<StatusTest<ScalarType,MV,OP> > > tmpVec = comboTest->getStatusTests();
174 
175  // Get the number of tests.
176  int numTests = tmpVec.size();
177 
178  // Find the maximum iteration and residual tests
179  for (int i=0; i<numTests; ++i) {
180 
181  // Check if this is a maximum iteration test.
182  Teuchos::RCP<StatusTestMaxIters_t> tmpItrTest = Teuchos::rcp_dynamic_cast<StatusTestMaxIters_t>(tmpVec[i]);
183  if (tmpItrTest != Teuchos::null) {
184  iterTest_ = tmpItrTest;
185  continue;
186  }
187 
188  // Check if this is a single residual test
189  Teuchos::RCP<StatusTestResNorm_t> tmpResTest = Teuchos::rcp_dynamic_cast<StatusTestResNorm_t>(tmpVec[i]);
190  // If the residual status test is a single test, put in the vector
191  if (tmpResTest != Teuchos::null) {
192  numResTests_ = 1;
193  resTestVec_.resize( numResTests_ );
194  resTestVec_[0] = tmpResTest;
195  continue;
196  }
197 
198  // Check if the residual test is a combination of several StatusTestResNorm objects.
199  Teuchos::RCP<StatusTestCombo_t> tmpComboTest = Teuchos::rcp_dynamic_cast<StatusTestCombo_t>(tmpVec[i]);
200  TEUCHOS_TEST_FOR_EXCEPTION(tmpComboTest == Teuchos::null,StatusTestError,"StatusTestResNormOutput(): test must be Belos::StatusTest[MaxIters|ResNorm|Combo].");
201  tmpVec = tmpComboTest->getStatusTests();
202  comboType_ = tmpComboTest->getComboType();
203  numResTests_ = tmpVec.size();
204  resTestVec_.resize( numResTests_ );
205  for (int j=0; j<numResTests_; ++j) {
206  tmpResTest = Teuchos::rcp_dynamic_cast<StatusTestResNorm_t>(tmpVec[j]);
207  TEUCHOS_TEST_FOR_EXCEPTION(tmpResTest == Teuchos::null,StatusTestError,"StatusTestResNormOutput(): test must be a vector of Belos::StatusTestResNorm.");
208  resTestVec_[j] = tmpResTest;
209  }
210  }
211 
212  // Keep the pointer to the new test and reset the state to Undefined.
213  test_ = test;
214  state_ = Undefined;
215 
216  }
217 
220  return test_;
221  }
222 
225  void setSolverDesc(const std::string& solverDesc) { solverDesc_ = solverDesc; }
226 
229  void setPrecondDesc(const std::string& precondDesc) { precondDesc_ = precondDesc; }
231 
232 
234 
235 
240  void reset() {
241  state_ = Undefined;
242  test_->reset();
243  lastNumIters_ = -1;
244  headerPrinted_ = false;
245  }
246 
248 
251  void resetNumCalls() {}
252 
254 
256 
257 
259  void print(std::ostream& os, int indent = 0) const {
260  std::string ind(indent,' ');
261  std::string starLine(55,'*');
262  std::string starFront(5,'*');
263 
264  std::ios_base::fmtflags osFlags(os.flags());
265 
266  os.setf(std::ios::scientific, std::ios::floatfield);
267  os.precision(6);
268 
269  // Print header if this is the first call to this output status test.
270  if (!headerPrinted_) {
271  os << std::endl << ind << starLine << std::endl;
272  os << ind << starFront << " Belos Iterative Solver: " << solverDesc_ << std::endl;
273  if (precondDesc_ != "")
274  os << ind << starFront << " Preconditioner: " << precondDesc_ << std::endl;
275  os << ind << starFront << " Maximum Iterations: " << iterTest_->getMaxIters() << std::endl;
276  os << ind << starFront << " Block Size: " << blockSize_ << std::endl;
277  if (numResTests_ > 1) {
278  os << ind << starFront << " Residual Tests ("
279  << ((comboType_ == StatusTestCombo_t::OR) ? "OR" : (comboType_ == StatusTestCombo_t::AND) ? "AND" :"SEQ")
280  << "): " << std::endl;
281  } else {
282  os << ind << starFront << " Residual Test: " << std::endl;
283  }
284  for (int i=0; i<numResTests_; ++i) {
285  os << ind << starFront << " Test " << i+1 << " : " << resTestVec_[i]->description() << std::endl;
286  }
287  os << ind << starLine << std::endl;
288  headerPrinted_ = true;
289  }
290 
291  // Print out residuals for each residual test.
292  os.setf(std::ios_base::right, std::ios_base::adjustfield);
293  std::string ind2( 7 + numIterDgts_, ' ' );
294  os << ind << "Iter " << std::setw(numIterDgts_) << iterTest_->getNumIters() << ", ";
295  for (int i=0; i<currNumRHS_; ++i) {
296  if ( i > 0 && currIdx_[i]!=-1 ) {
297  // Put in space where 'Iter :' is in the previous lines
298  os << ind << ind2;
299  }
300  os << "[" << std::setw(numLSDgts_) << currIdx_[i]+1 << "] : ";
301  for (int j=0; j<numResTests_; ++j) {
302  if ( resTestVec_[j]->getStatus() != Undefined && currIdx_[i]!=-1 ) {
303  os << std::setw(15) << (*resTestVec_[j]->getTestValue())[currIdx_[i]];
304  } else {
305  os << std::setw(15) << "---";
306  }
307  }
308  os << std::endl;
309  }
310  // reset os format
311  os.flags(osFlags);
312  }
313 
315 
316  private:
317  // Output manager.
319 
320  // Overall status test.
322 
323  // Iteration test (as passed in).
325 
327  std::vector<Teuchos::RCP<StatusTestResNorm<ScalarType,MV,OP> > > resTestVec_;
328 
329  std::string solverDesc_;
330  std::string precondDesc_;
331  std::vector<int> currIdx_;
332  StatusType state_;
333  mutable bool headerPrinted_;
334  int stateTest_, modTest_;
335  int lastNumIters_, comboType_;
336  int numResTests_, blockSize_;
337  int currNumRHS_, currLSNum_;
338  int numLSDgts_, numIterDgts_;
339 };
340 
341 } // end of Belos namespace
342 
343 #endif /* BELOS_STATUS_TEST_RESNORM_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 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...
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)
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.
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.
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: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 setPrecondDesc(const std::string &precondDesc)
Set a short preconditioner description for output clarity.
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...
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.

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