Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_SolutionStateMetaData_impl.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_impl_hpp
10 #define Tempus_SolutionStateMetaData_impl_hpp
11 
12 
13 namespace Tempus {
14 
15 
16 template<class Scalar>
18  :time_ (0.0),
19  iStep_ (0),
20  dt_ (0.0),
21  errorAbs_ (0.0),
22  errorRel_ (0.0),
23  order_ (1),
24  nFailures_ (0),
25  nRunningFailures_(0),
26  nConsecutiveFailures_(0),
27  tolRel_ (1.0e-02),
28  tolAbs_ (0.0),
29  xNormL2_ (0.0),
30  dxNormL2Rel_ (0.0),
31  dxNormL2Abs_ (0.0),
32  computeNorms_ (false),
33  solutionStatus_(WORKING),
34  output_ (false),
35  outputScreen_ (false),
36  isSynced_ (true),
37  isInterpolated_(false),
38  accuracy_ (0.0)
39 {}
40 
41 template<class Scalar>
43  const Scalar time,
44  const int iStep,
45  const Scalar dt,
46  const Scalar errorAbs,
47  const Scalar errorRel,
48  const int order,
49  const int nFailures,
50  const int nRunningFailures,
51  const int nConsecutiveFailures,
52  const Scalar tolRel,
53  const Scalar tolAbs,
54  const Scalar xNormL2,
55  const Scalar dxNormL2Rel,
56  const Scalar dxNormL2Abs,
57  const bool computeNorms,
58  const Status solutionStatus,
59  const bool output,
60  const bool outputScreen,
61  const bool isSynced,
62  const bool isInterpolated,
63  const Scalar accuracy)
64  :time_ (time),
65  iStep_ (iStep),
66  dt_ (dt),
67  errorAbs_ (errorAbs),
68  errorRel_ (errorRel),
69  order_ (order),
70  nFailures_ (nFailures),
71  nRunningFailures_(nRunningFailures),
72  nConsecutiveFailures_(nConsecutiveFailures),
73  tolRel_ (tolRel),
74  tolAbs_ (tolAbs),
75  xNormL2_ (xNormL2),
76  dxNormL2Rel_ (dxNormL2Rel),
77  dxNormL2Abs_ (dxNormL2Abs),
78  computeNorms_ (computeNorms),
79  solutionStatus_(solutionStatus),
80  output_ (output),
81  outputScreen_ (outputScreen),
82  isSynced_ (isSynced),
83  isInterpolated_(isInterpolated),
84  accuracy_ (accuracy)
85 {}
86 
87 template<class Scalar>
89  :time_ (ssmd.time_),
90  iStep_ (ssmd.iStep_),
91  dt_ (ssmd.dt_),
92  errorAbs_ (ssmd.errorAbs_),
93  errorRel_ (ssmd.errorRel_),
94  order_ (ssmd.order_),
95  nFailures_ (ssmd.nFailures_),
96  nRunningFailures_(ssmd.nRunningFailures_),
97  nConsecutiveFailures_(ssmd.nConsecutiveFailures_),
98  tolRel_ (ssmd.tolRel_),
99  tolAbs_ (ssmd.tolAbs_),
100  dxNormL2Rel_ (ssmd.dxNormL2Rel_),
101  dxNormL2Abs_ (ssmd.dxNormL2Abs_),
102  computeNorms_ (ssmd.computeNorms_),
103  solutionStatus_(ssmd.solutionStatus_),
104  output_ (ssmd.output_),
105  outputScreen_ (ssmd.outputScreen_),
106  isSynced_ (ssmd.isSynced_),
107  isInterpolated_(ssmd.isInterpolated_),
108  accuracy_ (ssmd.accuracy_)
109 {}
110 
111 
112 template<class Scalar>
114 {
117  time_,
118  iStep_,
119  dt_,
120  errorAbs_,
121  errorRel_,
122  order_,
123  nFailures_,
124  nRunningFailures_,
125  nConsecutiveFailures_,
126  tolRel_,
127  tolAbs_,
128  xNormL2_,
129  dxNormL2Rel_,
130  dxNormL2Abs_,
131  computeNorms_,
132  solutionStatus_,
133  output_,
134  outputScreen_,
135  isSynced_,
136  isInterpolated_,
137  accuracy_));
138 
139  return md;
140 }
141 
142 
143 template<class Scalar>
146 {
147  time_ = ssmd->time_;
148  iStep_ = ssmd->iStep_;
149  dt_ = ssmd->dt_;
150  errorAbs_ = ssmd->errorAbs_;
151  errorRel_ = ssmd->errorRel_;
152  order_ = ssmd->order_;
153  nFailures_ = ssmd->nFailures_;
154  nRunningFailures_= ssmd->nRunningFailures_;
155  nConsecutiveFailures_ = ssmd->nConsecutiveFailures_;
156  tolRel_ = ssmd->tolRel_,
157  tolAbs_ = ssmd->tolAbs_,
158  xNormL2_ = ssmd->xNormL2_,
159  dxNormL2Rel_ = ssmd->dxNormL2Rel_,
160  dxNormL2Abs_ = ssmd->dxNormL2Abs_,
161  computeNorms_ = ssmd->computeNorms_,
162  solutionStatus_ = ssmd->solutionStatus_;
163  output_ = ssmd->output_;
164  outputScreen_ = ssmd->outputScreen_;
165  isSynced_ = ssmd->isSynced_;
166  isInterpolated_ = ssmd->isInterpolated_;
167  accuracy_ = ssmd->accuracy_;
168 }
169 
170 
171 template<class Scalar>
173 {
174  std::string name = "Tempus::SolutionStateMetaData";
175  return(name);
176 }
177 
178 
179 template<class Scalar>
182  const Teuchos::EVerbosityLevel verbLevel) const
183 {
184  if (verbLevel == Teuchos::VERB_EXTREME) {
185  auto l_out = Teuchos::fancyOStream( out.getOStream() );
186  l_out->setOutputToRootOnly(0);
187  *l_out << description() << "::describe:" << std::endl
188  << "time = " << time_ << std::endl
189  << "iStep = " << iStep_ << std::endl
190  << "dt = " << dt_ << std::endl
191  << "errorAbs = " << errorAbs_ << std::endl
192  << "errorRel = " << errorRel_ << std::endl
193  << "order = " << order_ << std::endl
194  << "nFailures = " << nFailures_ << std::endl
195  << "nRunningFailures = " << nRunningFailures_<< std::endl
196  << "nConsecutiveFailures = " << nConsecutiveFailures_ << std::endl
197  << "tolRel = " << tolRel_ << std::endl
198  << "tolAbs = " << tolAbs_ << std::endl
199  << "xNormL2 = " << xNormL2_ << std::endl
200  << "dxNormL2Rel = " << dxNormL2Rel_ << std::endl
201  << "dxNormL2Abs = " << dxNormL2Abs_ << std::endl
202  << "computeNorms = " << computeNorms_ << std::endl
203  << "solutionStatus = " << toString(solutionStatus_) << std::endl
204  << "output = " << output_ << std::endl
205  << "outputScreen = " << outputScreen_ << std::endl
206  << "isSynced = " << isSynced_ << std::endl
207  << "isInterpolated = " << isInterpolated_ << std::endl
208  << "accuracy = " << accuracy_ << std::endl;
209  }
210 }
211 
212 } // namespace Tempus
213 #endif // Tempus_SolutionStateMetaData_impl_hpp
Teuchos::RCP< SolutionStateMetaData< Scalar > > clone() const
Clone constructor.
const std::string toString(const Status status)
Convert Status to string.
void copy(const Teuchos::RCP< const SolutionStateMetaData< Scalar > > &ssmd)
This is a deep copy.
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Status
Status for the Integrator, the Stepper and the SolutionState.
RCP< std::basic_ostream< char_type, traits_type > > getOStream()