Tempus Version of the Day
Time Integration
Loading...
Searching...
No Matches
Tempus_SolutionStateMetaData_impl.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_SolutionStateMetaData_impl_hpp
11#define Tempus_SolutionStateMetaData_impl_hpp
12
13namespace Tempus {
14
15template <class Scalar>
17 : time_(0.0),
18 iStep_(0),
19 dt_(0.0),
20 errorAbs_(0.0),
21 errorRel_(0.0),
22 errorRelNm1_(0.0),
23 errorRelNm2_(0.0),
24 order_(1),
25 nFailures_(0),
26 nRunningFailures_(0),
27 nConsecutiveFailures_(0),
28 tolRel_(1.0e-02),
29 tolAbs_(0.0),
30 xNormL2_(0.0),
31 dxNormL2Rel_(0.0),
32 dxNormL2Abs_(0.0),
33 computeNorms_(false),
34 solutionStatus_(WORKING),
35 output_(false),
36 outputScreen_(false),
37 isSynced_(true),
38 isInterpolated_(false),
39 accuracy_(0.0)
40{
41}
42
43template <class Scalar>
45 const Scalar time, const int iStep, const Scalar dt, const Scalar errorAbs,
46 const Scalar errorRel, const Scalar errorRelNm1, const Scalar errorRelNm2,
47 const int order, const int nFailures, const int nRunningFailures,
48 const int nConsecutiveFailures, const Scalar tolRel, const Scalar tolAbs,
49 const Scalar xNormL2, const Scalar dxNormL2Rel, const Scalar dxNormL2Abs,
50 const bool computeNorms, const Status solutionStatus, const bool output,
51 const bool outputScreen, const bool isSynced, const bool isInterpolated,
52 const Scalar accuracy)
53 : time_(time),
54 iStep_(iStep),
55 dt_(dt),
56 errorAbs_(errorAbs),
57 errorRel_(errorRel),
58 errorRelNm1_(errorRelNm1),
59 errorRelNm2_(errorRelNm2),
60 order_(order),
61 nFailures_(nFailures),
62 nRunningFailures_(nRunningFailures),
63 nConsecutiveFailures_(nConsecutiveFailures),
64 tolRel_(tolRel),
65 tolAbs_(tolAbs),
66 xNormL2_(xNormL2),
67 dxNormL2Rel_(dxNormL2Rel),
68 dxNormL2Abs_(dxNormL2Abs),
69 computeNorms_(computeNorms),
70 solutionStatus_(solutionStatus),
71 output_(output),
72 outputScreen_(outputScreen),
73 isSynced_(isSynced),
74 isInterpolated_(isInterpolated),
75 accuracy_(accuracy)
76{
77}
78
79template <class Scalar>
82 : time_(ssmd.time_),
83 iStep_(ssmd.iStep_),
84 dt_(ssmd.dt_),
85 errorAbs_(ssmd.errorAbs_),
86 errorRel_(ssmd.errorRel_),
87 errorRelNm1_(ssmd.errorRelNm1_),
88 errorRelNm2_(ssmd.errorRelNm2_),
89 order_(ssmd.order_),
90 nFailures_(ssmd.nFailures_),
91 nRunningFailures_(ssmd.nRunningFailures_),
92 nConsecutiveFailures_(ssmd.nConsecutiveFailures_),
93 tolRel_(ssmd.tolRel_),
94 tolAbs_(ssmd.tolAbs_),
95 dxNormL2Rel_(ssmd.dxNormL2Rel_),
96 dxNormL2Abs_(ssmd.dxNormL2Abs_),
97 computeNorms_(ssmd.computeNorms_),
98 solutionStatus_(ssmd.solutionStatus_),
99 output_(ssmd.output_),
100 outputScreen_(ssmd.outputScreen_),
101 isSynced_(ssmd.isSynced_),
102 isInterpolated_(ssmd.isInterpolated_),
103 accuracy_(ssmd.accuracy_)
104{
105}
106
107template <class Scalar>
108Teuchos::RCP<SolutionStateMetaData<Scalar> >
110{
111 Teuchos::RCP<SolutionStateMetaData<Scalar> > md =
113 time_, iStep_, dt_, errorAbs_, errorRel_, errorRelNm1_, errorRelNm2_,
114 order_, nFailures_, nRunningFailures_, nConsecutiveFailures_, tolRel_,
115 tolAbs_, xNormL2_, dxNormL2Rel_, dxNormL2Abs_, computeNorms_,
116 solutionStatus_, output_, outputScreen_, isSynced_, isInterpolated_,
117 accuracy_));
118
119 return md;
120}
121
122template <class Scalar>
124 const Teuchos::RCP<const SolutionStateMetaData<Scalar> >& ssmd)
125{
126 time_ = ssmd->time_;
127 iStep_ = ssmd->iStep_;
128 dt_ = ssmd->dt_;
129 errorAbs_ = ssmd->errorAbs_;
130 errorRel_ = ssmd->errorRel_;
131 errorRelNm1_ = ssmd->errorRelNm1_;
132 errorRelNm2_ = ssmd->errorRelNm2_;
133 order_ = ssmd->order_;
134 nFailures_ = ssmd->nFailures_;
135 nRunningFailures_ = ssmd->nRunningFailures_;
136 nConsecutiveFailures_ = ssmd->nConsecutiveFailures_;
137 tolRel_ = ssmd->tolRel_, tolAbs_ = ssmd->tolAbs_, xNormL2_ = ssmd->xNormL2_,
138 dxNormL2Rel_ = ssmd->dxNormL2Rel_, dxNormL2Abs_ = ssmd->dxNormL2Abs_,
139 computeNorms_ = ssmd->computeNorms_, solutionStatus_ = ssmd->solutionStatus_;
140 output_ = ssmd->output_;
141 outputScreen_ = ssmd->outputScreen_;
142 isSynced_ = ssmd->isSynced_;
143 isInterpolated_ = ssmd->isInterpolated_;
144 accuracy_ = ssmd->accuracy_;
145}
146
147template <class Scalar>
149{
150 std::string name = "Tempus::SolutionStateMetaData";
151 return (name);
152}
153
154template <class Scalar>
156 Teuchos::FancyOStream& out, const Teuchos::EVerbosityLevel verbLevel) const
157{
158 auto l_out = Teuchos::fancyOStream(out.getOStream());
159 Teuchos::OSTab ostab(*l_out, 2, this->description());
160 l_out->setOutputToRootOnly(0);
161
162 *l_out << "\n--- " << this->description() << " ---" << std::endl;
163
164 if (verbLevel >= Teuchos::VERB_MEDIUM) {
165 *l_out << " time = " << time_ << std::endl
166 << " iStep = " << iStep_ << std::endl
167 << " dt = " << dt_ << std::endl
168 << " errorAbs = " << errorAbs_ << std::endl
169 << " errorRel = " << errorRel_ << std::endl
170 << " errorRelNm1 = " << errorRelNm1_ << std::endl
171 << " errorRelNm2 = " << errorRelNm2_ << std::endl
172 << " order = " << order_ << std::endl
173 << " nFailures = " << nFailures_ << std::endl
174 << " nRunningFailures = " << nRunningFailures_ << std::endl
175 << " nConsecutiveFailures = " << nConsecutiveFailures_ << std::endl
176 << " tolRel = " << tolRel_ << std::endl
177 << " tolAbs = " << tolAbs_ << std::endl
178 << " xNormL2 = " << xNormL2_ << std::endl
179 << " dxNormL2Rel = " << dxNormL2Rel_ << std::endl
180 << " dxNormL2Abs = " << dxNormL2Abs_ << std::endl
181 << " computeNorms = " << computeNorms_ << std::endl
182 << " solutionStatus = " << toString(solutionStatus_) << std::endl
183 << " output = " << output_ << std::endl
184 << " outputScreen = " << outputScreen_ << std::endl
185 << " isSynced = " << isSynced_ << std::endl
186 << " isInterpolated = " << isInterpolated_ << std::endl
187 << " accuracy = " << accuracy_ << std::endl;
188 }
189 *l_out << std::string(this->description().length() + 8, '-') << std::endl;
190}
191
192} // namespace Tempus
193#endif // Tempus_SolutionStateMetaData_impl_hpp
Teuchos::RCP< SolutionStateMetaData< Scalar > > clone() const
Clone constructor.
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
Status
Status for the Integrator, the Stepper and the SolutionState.
const std::string toString(const Status status)
Convert Status to string.