Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_SolutionStateMetaData_decl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ****************************************************************************
3 // Tempus: Copyright (2017) Sandia Corporation
4 //
5 // Distributed under BSD 3-clause license (See accompanying file Copyright.txt)
6 // ****************************************************************************
7 // @HEADER
8 
9 #ifndef Tempus_SolutionStateMetaData_decl_hpp
10 #define Tempus_SolutionStateMetaData_decl_hpp
11 
12 // Teuchos
13 #include "Teuchos_VerboseObject.hpp"
14 #include "Teuchos_Describable.hpp"
15 // Tempus
16 #include "Tempus_config.hpp"
17 #include "Tempus_Types.hpp"
18 
19 
20 namespace Tempus {
21 
24 template<class Scalar>
26  public Teuchos::Describable,
27  public Teuchos::VerboseObject<Tempus::SolutionStateMetaData<Scalar> >
28 {
29 public:
30 
33 
36  const Scalar time,
37  const int iStep,
38  const Scalar dt,
39  const Scalar errorAbs,
40  const Scalar errorRel,
41  const int order,
42  const int nFailures,
43  const int nRunningFailures,
44  const int nConsecutiveFailures,
45  const Scalar tolRel,
46  const Scalar tolAbs,
47  const Scalar xNormL2,
48  const Scalar dxNormL2Rel,
49  const Scalar dxNormL2Abs,
50  const bool computeNorms,
51  const Status solutionStatus,
52  const bool output,
53  const bool outputScreen,
54  const bool isSynced,
55  const bool isInterpolated,
56  const Scalar accuracy);
57 
60 
63 
65  void copy(const Teuchos::RCP<const SolutionStateMetaData<Scalar> >& ssmd);
66 
69 
71 
72  Scalar getTime() const {return time_;}
73  int getIStep() const {return iStep_;}
74  Scalar getDt() const {return dt_;}
75  Scalar getErrorAbs() const {return errorAbs_;}
76  Scalar getErrorRel() const {return errorRel_;}
77  Scalar getOrder() const {return order_;}
78  int getNFailures() const {return nFailures_;}
79  int getNRunningFailures() const {return nRunningFailures_;}
81  Scalar getTolAbs() const {return tolAbs_;}
82  Scalar getTolRel() const {return tolRel_;}
83  Scalar getXNormL2() const {return xNormL2_;}
84  Scalar getDxNormL2Abs() const {return dxNormL2Abs_;}
85  Scalar getDxNormL2Rel() const {return dxNormL2Rel_;}
86  bool getComputeNorms() const {return computeNorms_;}
88  bool getOutput() const {return output_;}
89  bool getOutputScreen() const {return outputScreen_;}
90  bool getIsSynced() const {return isSynced_;}
91  bool getIsInterpolated() const {return isInterpolated_;}
92  Scalar getAccuracy() const {return accuracy_;}
93 
94  void setTime(Scalar time) {time_ = time;}
95  void setIStep(int iStep) {iStep_ = iStep;}
96  void setDt(Scalar dt) {dt_ = dt;}
97  void setErrorAbs (Scalar errorAbs) {errorAbs_ = errorAbs;}
98  void setErrorRel (Scalar errorRel) {errorRel_ = errorRel;}
99  void setOrder(Scalar order) {order_ = order;}
100  void setNFailures(int nFailures) {nFailures_ = nFailures;}
101  void setNRunningFailures(int nFailures) {nRunningFailures_ = nFailures;}
102  void setNConsecutiveFailures(int nConsecutiveFailures)
103  {nConsecutiveFailures_ = nConsecutiveFailures;}
104  void setTolRel (Scalar tolRel) {tolRel_ = tolRel;}
105  void setTolAbs (Scalar tolAbs) {tolAbs_ = tolAbs;}
106  void setXNormL2 (Scalar xNormL2) {xNormL2_ = xNormL2;}
107  void setDxNormL2Rel (Scalar dxNormL2Rel) {dxNormL2Rel_ = dxNormL2Rel;}
108  void setDxNormL2Abs (Scalar dxNormL2Abs) {dxNormL2Abs_ = dxNormL2Abs;}
109  void setComputeNorms (bool computeNorms) {computeNorms_ = computeNorms;}
110  void setSolutionStatus(Status solutionStatus)
111  {solutionStatus_ = solutionStatus;}
112  void setOutput(bool output) {output_ = output;}
113  void setOutputScreen(bool outputScreen) {outputScreen_ = outputScreen;}
114  void setIsSynced(bool isSynced) {isSynced_=isSynced;}
115  void setIsInterpolated(bool isInterpolated)
116  {isInterpolated_ = isInterpolated;}
117  void setAccuracy(Scalar accuracy) {accuracy_ = accuracy;}
119 
121 
122  virtual std::string description() const;
123  virtual void describe(Teuchos::FancyOStream &out,
124  const Teuchos::EVerbosityLevel verbLevel) const;
126 
127 protected:
128  Scalar time_;
129  int iStep_;
130  Scalar dt_;
131  Scalar errorAbs_;
132  Scalar errorRel_;
133  Scalar order_;
137  Scalar tolRel_;
138  Scalar tolAbs_;
139  Scalar xNormL2_;
140  Scalar dxNormL2Rel_;
141  Scalar dxNormL2Abs_;
143 
152  bool output_;
154 
158  bool isSynced_;
160  Scalar accuracy_;
161 
162 };
163 
164 
165 } // namespace Tempus
166 
167 #endif // Tempus_SolutionStateMetaData_decl_hpp
Teuchos::RCP< SolutionStateMetaData< Scalar > > clone() const
Clone constructor.
Status solutionStatus_
The solutionStatus is used to indicate.
int nFailures_
Total number of stepper failures.
Scalar dxNormL2Rel_
Relative L2-Norm of the change in solution, ||x_i-x_{i-1}||/||x_{i-1}||.
void copy(const Teuchos::RCP< const SolutionStateMetaData< Scalar > > &ssmd)
This is a deep copy.
Scalar dxNormL2Abs_
Absolute L2-Norm of the change in solution, ||x_i-x_{i-1}||.
int nRunningFailures_
Total number of running stepper failures.
void setNConsecutiveFailures(int nConsecutiveFailures)
int iStep_
Time step index for this solution.
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
bool isInterpolated_
F - soln is time integrated; T - soln is interpolated.
Status
Status for the Integrator, the Stepper and the SolutionState.
bool output_
SolutionState should be or has been outputted.
int nConsecutiveFailures_
Consecutive number of stepper failures.
Scalar errorRel_
Relative local truncation error.
Scalar dt_
Time step for this solution.
Scalar accuracy_
Interpolation accuracy of solution.
bool computeNorms_
flag to compute norms of solution
bool isSynced_
True - all of soln (x, xDot, xDotDot) is at the same time level. False - solution is at different tim...
Scalar errorAbs_
Absolute local truncation error.