Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_RKButcherTableau.hpp
Go to the documentation of this file.
1 //@HEADER
2 // *****************************************************************************
3 // Tempus: Time Integration and Sensitivity Analysis Package
4 //
5 // Copyright 2017 NTESS and the Tempus contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 //@HEADER
9 
10 #ifndef Tempus_RKButcherTableau_hpp
11 #define Tempus_RKButcherTableau_hpp
12 
13 // disable clang warnings
14 #ifdef __clang__
15 #pragma clang system_header
16 #endif
17 
20 
21 #include "Tempus_config.hpp"
23 
24 namespace Tempus {
25 
43 template <class Scalar>
45  : virtual public Teuchos::Describable,
46  virtual public Teuchos::VerboseObject<RKButcherTableau<Scalar> > {
47  public:
48  RKButcherTableau(std::string stepperType,
52  const int order, const int orderMin, const int orderMax,
55  bool checkC = true)
56  : description_(stepperType)
57  {
58  const int numStages = A.numRows();
62  TEUCHOS_ASSERT(order > 0);
63  A_ = A;
64  b_ = b;
65  c_ = c;
66  order_ = order;
69  this->set_isImplicit();
70  this->set_isDIRK();
71 
72  // Consistency check on b
74  Scalar sumb = ST::zero();
75  for (size_t i = 0; i < this->numStages(); i++) sumb += b_(i);
77  std::abs(ST::one() - sumb) > 1.0e-08, std::runtime_error,
78  "Error - Butcher Tableau b fails to satisfy Sum(b_i) = 1.\n"
79  << " Sum(b_i) = " << sumb << "\n");
80 
81  // Consistency check on c
82  if (checkC) {
83  for (size_t i = 0; i < this->numStages(); i++) {
84  Scalar sumai = ST::zero();
85  for (size_t j = 0; j < this->numStages(); j++) sumai += A_(i, j);
86  bool failed = false;
87  if (std::abs(sumai) > 1.0e-08)
88  failed = (std::abs((sumai - c_(i)) / sumai) > 1.0e-08);
89  else
90  failed = (std::abs(c_(i)) > 1.0e-08);
91 
93  failed, std::runtime_error,
94  "Error - Butcher Tableau fails to satisfy c_i = Sum_j(a_ij).\n"
95  << " Stage i = " << i << "\n"
96  << " c_i = " << c_(i) << "\n"
97  << " Sum_j(a_ij) = " << sumai << "\n"
98  << " This may be OK as some valid tableaus do not satisfy\n"
99  << " this condition. If so, construct this RKButcherTableau\n"
100  << " with checkC = false.\n");
101  }
102  }
103 
104  if (bstar.length() > 0) {
106  isEmbedded_ = true;
107  }
108  else {
109  isEmbedded_ = false;
110  }
111  bstar_ = bstar;
112  }
113 
115  virtual std::size_t numStages() const { return A_.numRows(); }
118  {
119  return A_;
120  }
123  {
124  return b_;
125  }
128  {
129  return bstar_;
130  }
133  {
134  return c_;
135  }
137  virtual int order() const { return order_; }
139  virtual int orderMin() const { return orderMin_; }
141  virtual int orderMax() const { return orderMax_; }
143  virtual bool isImplicit() const { return isImplicit_; }
145  virtual bool isDIRK() const { return isDIRK_; }
147  virtual bool isEmbedded() const { return isEmbedded_; }
149  virtual bool isTVD() const { return isTVD_; }
151  virtual Scalar getTVDCoeff() const { return tvdCoeff_; }
153  virtual void setTVDCoeff(const Scalar a) { tvdCoeff_ = a; }
154  virtual void setTVD(bool a) { isTVD_ = a; }
155  virtual void setEmbedded(bool a) { isEmbedded_ = a; }
156 
157  /* \brief Redefined from Teuchos::Describable */
159  virtual std::string description() const { return description_; }
160 
161  virtual void describe(Teuchos::FancyOStream& out,
162  const Teuchos::EVerbosityLevel verbLevel) const
163  {
164  out.setOutputToRootOnly(0);
165 
166  if (verbLevel != Teuchos::VERB_NONE) {
167  out << this->description() << std::endl;
168  out << "number of Stages = " << this->numStages() << std::endl;
169  out << "A = " << printMat(this->A()) << std::endl;
170  out << "b = " << printMat(this->b()) << std::endl;
171  out << "c = " << printMat(this->c()) << std::endl;
172  out << "bstar = " << printMat(this->bstar()) << std::endl;
173  out << "order = " << this->order() << std::endl;
174  out << "orderMin = " << this->orderMin() << std::endl;
175  out << "orderMax = " << this->orderMax() << std::endl;
176  out << "isImplicit = " << this->isImplicit() << std::endl;
177  out << "isDIRK = " << this->isDIRK() << std::endl;
178  out << "isEmbedded = " << this->isEmbedded() << std::endl;
179  if (this->isTVD())
180  out << "TVD Coeff = " << this->getTVDCoeff() << std::endl;
181  }
182  }
184 
185  bool operator==(const RKButcherTableau& t) const
186  {
187  const Scalar relTol = 1.0e-15;
188  if (A_->numRows() != t.A_->numRows() || A_->numCols() != t.A_->numCols()) {
189  return false;
190  }
191  else {
192  int i, j;
193  for (i = 0; i < A_->numRows(); i++) {
194  for (j = 0; j < A_->numCols(); j++) {
195  if (std::abs((t.A_(i, j) - A_(i, j)) / A_(i, j)) > relTol)
196  return false;
197  }
198  }
199  }
200 
201  if (b_->length() != t.b_->length() || b_->length() != t.c_->length() ||
202  b_->length() != t.bstar_->length()) {
203  return false;
204  }
205  else {
206  int i;
207  for (i = 0; i < A_->numRows(); i++) {
208  if (std::abs((t.b_(i) - b_(i)) / b_(i)) > relTol) return false;
209  if (std::abs((t.c_(i) - c_(i)) / c_(i)) > relTol) return false;
210  if (std::abs((t.bstar_(i) - bstar_(i)) / bstar_(i)) > relTol)
211  return false;
212  }
213  }
214  return true;
215  }
216 
217  bool operator!=(const RKButcherTableau& t) const { return !((*this) == t); }
218 
219  private:
221 
222  protected:
224  {
225  isImplicit_ = false;
226  for (size_t i = 0; i < this->numStages(); i++)
227  for (size_t j = i; j < this->numStages(); j++)
228  if (A_(i, j) != 0.0) isImplicit_ = true;
229  }
230 
232  void set_isDIRK()
233  {
234  isDIRK_ = true;
235  bool nonZero = false;
236  for (size_t i = 0; i < this->numStages(); i++) {
237  if (A_(i, i) != 0.0) nonZero = true;
238  for (size_t j = i + 1; j < this->numStages(); j++)
239  if (A_(i, j) != 0.0) isDIRK_ = false;
240  }
241  if (nonZero == false) isDIRK_ = false;
242  }
243 
244  std::string description_;
245 
249  int order_;
253  bool isDIRK_;
255  bool isTVD_ = false;
257  Scalar tvdCoeff_ = 0;
258 };
259 
260 } // namespace Tempus
261 
262 #endif // Tempus_RKButcherTableau_hpp
RKButcherTableau(std::string stepperType, const Teuchos::SerialDenseMatrix< int, Scalar > &A, const Teuchos::SerialDenseVector< int, Scalar > &b, const Teuchos::SerialDenseVector< int, Scalar > &c, const int order, const int orderMin, const int orderMax, const Teuchos::SerialDenseVector< int, Scalar > &bstar=Teuchos::SerialDenseVector< int, Scalar >(), bool checkC=true)
bool operator!=(const RKButcherTableau &t) const
SerialBandDenseMatrixPrinter< OrdinalType, ScalarType > printMat(const SerialBandDenseMatrix< OrdinalType, ScalarType > &obj)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual void setTVDCoeff(const Scalar a)
set TVD coefficient of RK method
virtual bool isTVD() const
Return true if the RK method is TVD.
Teuchos::SerialDenseVector< int, Scalar > bstar_
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
virtual bool isImplicit() const
Return true if the RK method is implicit.
virtual bool isEmbedded() const
Return true if the RK method has embedded capabilities.
virtual bool isDIRK() const
Return true if the RK method is Diagonally Implicit.
virtual int orderMax() const
Return the maximum order.
virtual const Teuchos::SerialDenseVector< int, Scalar > & c() const
Return the vector of stage positions.
virtual Scalar getTVDCoeff() const
Return TVD coefficient of RK method.
void set_isDIRK()
DIRK is defined as if a_ij = 0 for j&gt;i and a_ii != 0 for at least one i.
virtual const Teuchos::SerialDenseVector< int, Scalar > & b() const
Return the vector of quadrature weights.
virtual std::string description() const
OrdinalType length() const
basic_FancyOStream & setOutputToRootOnly(const int rootRank)
OrdinalType numCols() const
virtual const Teuchos::SerialDenseMatrix< int, Scalar > & A() const
Return the matrix coefficients.
virtual std::size_t numStages() const
Return the number of stages.
Teuchos::SerialDenseVector< int, Scalar > c_
#define TEUCHOS_ASSERT(assertion_test)
Teuchos::SerialDenseMatrix< int, Scalar > A_
#define TEUCHOS_ASSERT_EQUALITY(val1, val2)
virtual const Teuchos::SerialDenseVector< int, Scalar > & bstar() const
Return the vector of quadrature weights for embedded methods.
virtual int order() const
Return the order.
Teuchos::SerialDenseVector< int, Scalar > b_
virtual int orderMin() const
Return the minimum order.
bool operator==(const RKButcherTableau &t) const
OrdinalType numRows() const