Tempus Version of the Day
Time Integration
Loading...
Searching...
No Matches
Tempus_StepperExplicitRK_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_StepperExplicitRK_impl_hpp
11#define Tempus_StepperExplicitRK_impl_hpp
12
13#include "Thyra_VectorStdOps.hpp"
14
16
17namespace Tempus {
18
19template <class Scalar>
21{
22 this->setUseEmbedded(false);
23 this->setStageNumber(-1);
24 this->setAppAction(Teuchos::null);
25}
26
27template <class Scalar>
29 const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel,
30 bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck,
31 bool useEmbedded,
32 const Teuchos::RCP<StepperRKAppAction<Scalar> >& stepperRKAppAction)
33{
34 this->setUseFSAL(useFSAL);
35 this->setICConsistency(ICConsistency);
36 this->setICConsistencyCheck(ICConsistencyCheck);
37 this->setUseEmbedded(useEmbedded);
38 this->setStageNumber(-1);
39 this->setErrorNorm();
40
41 this->setAppAction(stepperRKAppAction);
42
43 if (appModel != Teuchos::null) {
44 this->setModel(appModel);
45 this->initialize();
46 }
47}
48
49template <class Scalar>
51 const Teuchos::RCP<SolutionHistory<Scalar> >& sh) const
52{
53 Scalar dt = Scalar(1.0e+99);
54 if (!this->getUseEmbedded()) return dt;
55
56 Teuchos::RCP<SolutionState<Scalar> > currentState = sh->getCurrentState();
57 const int order = currentState->getOrder();
58 const Scalar time = currentState->getTime();
59 const Scalar errorRel = currentState->getTolRel();
60 const Scalar errorAbs = currentState->getTolAbs();
61
62 Teuchos::RCP<Thyra::VectorBase<Scalar> > stageX, scratchX;
63 stageX = Thyra::createMember(this->appModel_->get_f_space());
64 scratchX = Thyra::createMember(this->appModel_->get_f_space());
65 Thyra::assign(stageX.ptr(), *(currentState->getX()));
66
67 std::vector<Teuchos::RCP<Thyra::VectorBase<Scalar> > > stageXDot(2);
68 for (int i = 0; i < 2; ++i) {
69 stageXDot[i] = Thyra::createMember(this->appModel_->get_f_space());
70 assign(stageXDot[i].ptr(), Teuchos::ScalarTraits<Scalar>::zero());
71 }
72
73 // A: one function evaluation at F(t_0, X_0)
74 typedef Thyra::ModelEvaluatorBase MEB;
75 MEB::InArgs<Scalar> inArgs = this->appModel_->getNominalValues();
76 MEB::OutArgs<Scalar> outArgs = this->appModel_->createOutArgs();
77 inArgs.set_x(stageX);
78 if (inArgs.supports(MEB::IN_ARG_t)) inArgs.set_t(time);
79 if (inArgs.supports(MEB::IN_ARG_x_dot)) inArgs.set_x_dot(Teuchos::null);
80 outArgs.set_f(stageXDot[0]); // K1
81 this->appModel_->evalModel(inArgs, outArgs);
82
83 this->stepperErrorNormCalculator_->setRelativeTolerance(errorRel);
84 this->stepperErrorNormCalculator_->setAbsoluteTolerance(errorAbs);
85
86 Scalar d0 = this->stepperErrorNormCalculator_->errorNorm(stageX);
87 Scalar d1 = this->stepperErrorNormCalculator_->errorNorm(stageXDot[0]);
88
89 // b) first guess for the step size
90 dt = Teuchos::as<Scalar>(0.01) * (d0 / d1);
91
92 // c) perform one explicit Euler step (X_1)
93 Thyra::Vp_StV(stageX.ptr(), dt, *(stageXDot[0]));
94
95 // compute F(t_0 + dt, X_1)
96 inArgs.set_x(stageX);
97 if (inArgs.supports(MEB::IN_ARG_t)) inArgs.set_t(time + dt);
98 if (inArgs.supports(MEB::IN_ARG_x_dot)) inArgs.set_x_dot(Teuchos::null);
99 outArgs.set_f(stageXDot[1]); // K2
100 this->appModel_->evalModel(inArgs, outArgs);
101
102 // d) compute estimate of the second derivative of the solution
103 // d2 = || f(t_0 + dt, X_1) - f(t_0, X_0) || / dt
104 Teuchos::RCP<Thyra::VectorBase<Scalar> > errX;
105 errX = Thyra::createMember(this->appModel_->get_f_space());
106 assign(errX.ptr(), Teuchos::ScalarTraits<Scalar>::zero());
107 Thyra::V_VmV(errX.ptr(), *(stageXDot[1]), *(stageXDot[0]));
108 Scalar d2 = this->stepperErrorNormCalculator_->errorNorm(errX) / dt;
109
110 // e) compute step size h_1 (from m = 0 order Taylor series)
111 Scalar max_d1_d2 = std::max(d1, d2);
112 Scalar h1 = std::pow((0.01 / max_d1_d2), (1.0 / (order + 1)));
113
114 // f) propose starting step size
115 dt = std::min(100 * dt, h1);
116 return dt;
117}
118
119template <class Scalar>
120Teuchos::RCP<const Teuchos::ParameterList>
122{
123 return this->getValidParametersBasicERK();
124}
125
126template <class Scalar>
127Teuchos::RCP<Teuchos::ParameterList>
129{
130 auto pl = this->getValidParametersBasic();
131 pl->template set<bool>(
132 "Use Embedded", this->getUseEmbedded(),
133 "'Whether to use Embedded Stepper (if available) or not\n"
134 " 'true' - Stepper will compute embedded solution and is adaptive.\n"
135 " 'false' - Stepper is not embedded(adaptive).\n");
136 pl->template set<std::string>("Description", this->getDescription());
137
138 return pl;
139}
140
141template <class Scalar>
143{
144 TEUCHOS_TEST_FOR_EXCEPTION(this->tableau_ == Teuchos::null, std::logic_error,
145 "Error - Need to set the tableau, before calling "
146 "StepperExplicitRK::initialize()\n");
147
148 TEUCHOS_TEST_FOR_EXCEPTION(
149 this->appModel_ == Teuchos::null, std::logic_error,
150 "Error - Need to set the model, setModel(), before calling "
151 "StepperExplicitRK::initialize()\n");
152
154}
155
156template <class Scalar>
158 const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& appModel)
159{
161
162 // Set the stage vectors
163 int numStages = this->tableau_->numStages();
164 stageXDot_.resize(numStages);
165 for (int i = 0; i < numStages; ++i) {
166 stageXDot_[i] = Thyra::createMember(this->appModel_->get_f_space());
167 assign(stageXDot_[i].ptr(), Teuchos::ScalarTraits<Scalar>::zero());
168 }
169
170 this->setEmbeddedMemory();
171 this->setErrorNorm();
172
173 this->isInitialized_ = false;
174}
175
176template <class Scalar>
178{
179 if (this->getModel() == Teuchos::null)
180 return; // Embedded memory will be set when setModel() is called.
181
182 if (this->tableau_->isEmbedded() && this->getUseEmbedded()) {
183 this->ee_ = Thyra::createMember(this->appModel_->get_f_space());
184 this->abs_u0 = Thyra::createMember(this->appModel_->get_f_space());
185 this->abs_u = Thyra::createMember(this->appModel_->get_f_space());
186 this->sc = Thyra::createMember(this->appModel_->get_f_space());
187 }
188 else {
189 this->ee_ = Teuchos::null;
190 this->abs_u0 = Teuchos::null;
191 this->abs_u = Teuchos::null;
192 this->sc = Teuchos::null;
193 }
194}
195
196template <class Scalar>
198 const Teuchos::RCP<SolutionHistory<Scalar> >& solutionHistory)
199{
200 if (this->getUseFSAL())
201 this->setStepperXDot(stageXDot_.back());
202 else
203 this->setStepperXDot(stageXDot_[0]);
204
206
207 auto xDot = solutionHistory->getCurrentState()->getXDot();
208 if (xDot != Teuchos::null && this->getUseFSAL())
209 Thyra::assign(this->getStepperXDot().ptr(), *(xDot));
210}
211
212template <class Scalar>
214 const Teuchos::RCP<SolutionHistory<Scalar> >& solutionHistory)
215{
216 this->checkInitialized();
217
218 using Teuchos::RCP;
219
220 TEMPUS_FUNC_TIME_MONITOR("Tempus::StepperExplicitRK::takeStep()");
221 {
222 TEUCHOS_TEST_FOR_EXCEPTION(
223 solutionHistory->getNumStates() < 2, std::logic_error,
224 "Error - StepperExplicitRK<Scalar>::takeStep(...)\n"
225 "Need at least two SolutionStates for ExplicitRK.\n"
226 " Number of States = "
227 << solutionHistory->getNumStates()
228 << "\n"
229 "Try setting in \"Solution History\" \"Storage Type\" = "
230 "\"Undo\"\n"
231 " or \"Storage Type\" = \"Static\" and \"Storage Limit\" = "
232 "\"2\"\n");
233
234 RCP<SolutionState<Scalar> > currentState =
235 solutionHistory->getCurrentState();
236 RCP<SolutionState<Scalar> > workingState =
237 solutionHistory->getWorkingState();
238 const Scalar dt = workingState->getTimeStep();
239 const Scalar time = currentState->getTime();
240
241 const int numStages = this->tableau_->numStages();
242 Teuchos::SerialDenseMatrix<int, Scalar> A = this->tableau_->A();
243 Teuchos::SerialDenseVector<int, Scalar> b = this->tableau_->b();
244 Teuchos::SerialDenseVector<int, Scalar> c = this->tableau_->c();
245
246 Thyra::assign(workingState->getX().ptr(), *(currentState->getX()));
247
248 RCP<StepperExplicitRK<Scalar> > thisStepper = Teuchos::rcpFromRef(*this);
249 this->stepperRKAppAction_->execute(
250 solutionHistory, thisStepper,
252
253 // Compute stage solutions
254 for (int i = 0; i < numStages; ++i) {
255 this->setStageNumber(i);
256 Thyra::assign(workingState->getX().ptr(), *(currentState->getX()));
257 for (int j = 0; j < i; ++j) {
258 if (A(i, j) != Teuchos::ScalarTraits<Scalar>::zero()) {
259 Thyra::Vp_StV(workingState->getX().ptr(), dt * A(i, j),
260 *stageXDot_[j]);
261 }
262 }
263 this->setStepperXDot(stageXDot_[i]);
264
265 this->stepperRKAppAction_->execute(
266 solutionHistory, thisStepper,
268 this->stepperRKAppAction_->execute(
269 solutionHistory, thisStepper,
271 this->stepperRKAppAction_->execute(
272 solutionHistory, thisStepper,
274 this->stepperRKAppAction_->execute(
275 solutionHistory, thisStepper,
277
278 if (i == 0 && this->getUseFSAL() &&
279 workingState->getNConsecutiveFailures() == 0) {
280 RCP<Thyra::VectorBase<Scalar> > tmp = stageXDot_[0];
281 stageXDot_[0] = stageXDot_.back();
282 stageXDot_.back() = tmp;
283 this->setStepperXDot(stageXDot_[0]);
284 }
285 else {
286 const Scalar ts = time + c(i) * dt;
287 auto p = Teuchos::rcp(new ExplicitODEParameters<Scalar>(dt));
288
289 // Evaluate xDot = f(x,t).
290 this->evaluateExplicitODE(stageXDot_[i], workingState->getX(), ts, p);
291 }
292
293 this->stepperRKAppAction_->execute(
294 solutionHistory, thisStepper,
296 }
297 // reset the stage number
298 this->setStageNumber(-1);
299
300 // Sum for solution: x_n = x_n-1 + Sum{ b(i) * dt*f(i) }
301 Thyra::assign((workingState->getX()).ptr(), *(currentState->getX()));
302 for (int i = 0; i < numStages; ++i) {
303 if (b(i) != Teuchos::ScalarTraits<Scalar>::zero()) {
304 Thyra::Vp_StV((workingState->getX()).ptr(), dt * b(i),
305 *(stageXDot_[i]));
306 }
307 }
308
309 if (this->getUseFSAL()) {
310 if (numStages == 1) {
311 const Scalar ts = time + dt;
312 auto p = Teuchos::rcp(new ExplicitODEParameters<Scalar>(dt));
313 // Evaluate xDot = f(x,t).
314 this->evaluateExplicitODE(stageXDot_[0], workingState->getX(), ts, p);
315 }
316 if (workingState->getXDot() != Teuchos::null)
317 Thyra::assign((workingState->getXDot()).ptr(), *(stageXDot_.back()));
318 }
319
320 // At this point, the stepper has passed.
321 // But when using adaptive time stepping, the embedded method
322 // can change the step status
323 workingState->setSolutionStatus(Status::PASSED);
324
325 if (this->tableau_->isEmbedded() && this->getUseEmbedded()) {
326 const Scalar tolRel = workingState->getTolRel();
327 const Scalar tolAbs = workingState->getTolAbs();
328
329 // update the tolerance
330 this->stepperErrorNormCalculator_->setRelativeTolerance(tolRel);
331 this->stepperErrorNormCalculator_->setAbsoluteTolerance(tolAbs);
332
333 // just compute the error weight vector
334 // (all that is needed is the error, and not the embedded solution)
335 Teuchos::SerialDenseVector<int, Scalar> errWght = b;
336 errWght -= this->tableau_->bstar();
337
338 // Compute local truncation error estimate: | u^{n+1} - \hat{u}^{n+1} |
339 // Sum for solution: ee_n = Sum{ (b(i) - bstar(i)) * dt*f(i) }
340 assign(this->ee_.ptr(), Teuchos::ScalarTraits<Scalar>::zero());
341 for (int i = 0; i < numStages; ++i) {
342 if (errWght(i) != Teuchos::ScalarTraits<Scalar>::zero()) {
343 Thyra::Vp_StV(this->ee_.ptr(), dt * errWght(i), *(stageXDot_[i]));
344 }
345 }
346
347 Scalar err = this->stepperErrorNormCalculator_->computeWRMSNorm(
348 currentState->getX(), workingState->getX(), this->ee_);
349 workingState->setErrorRel(err);
350
351 // Test if step should be rejected
352 if (std::isinf(err) || std::isnan(err) || err > Teuchos::as<Scalar>(1.0))
353 workingState->setSolutionStatus(Status::FAILED);
354 }
355
356 workingState->setOrder(this->getOrder());
357 workingState->computeNorms(currentState);
358 this->stepperRKAppAction_->execute(
359 solutionHistory, thisStepper,
361 }
362 return;
363}
364
371template <class Scalar>
372Teuchos::RCP<Tempus::StepperState<Scalar> >
374{
375 Teuchos::RCP<Tempus::StepperState<Scalar> > stepperState =
376 rcp(new StepperState<Scalar>(this->getStepperType()));
377 return stepperState;
378}
379
380template <class Scalar>
382 Teuchos::FancyOStream& out, const Teuchos::EVerbosityLevel verbLevel) const
383{
384 out.setOutputToRootOnly(0);
385 out << std::endl;
386 Stepper<Scalar>::describe(out, verbLevel);
387 StepperExplicit<Scalar>::describe(out, verbLevel);
388
389 out << "--- StepperExplicitRK ---\n";
390 if (this->tableau_ != Teuchos::null) this->tableau_->describe(out, verbLevel);
391 out << " tableau_ = " << this->tableau_ << std::endl;
392 out << " stepperRKAppAction_= " << this->stepperRKAppAction_ << std::endl;
393 out << " stageXDot_.size() = " << stageXDot_.size() << std::endl;
394 const int numStages = stageXDot_.size();
395 for (int i = 0; i < numStages; ++i)
396 out << " stageXDot_[" << i << "] = " << stageXDot_[i] << std::endl;
397 out << " useEmbedded_ = " << Teuchos::toString(this->useEmbedded_)
398 << std::endl;
399 out << " ee_ = " << this->ee_ << std::endl;
400 out << " abs_u0 = " << this->abs_u0 << std::endl;
401 out << " abs_u = " << this->abs_u << std::endl;
402 out << " sc = " << this->sc << std::endl;
403 out << "-------------------------" << std::endl;
404}
405
406template <class Scalar>
407bool StepperExplicitRK<Scalar>::isValidSetup(Teuchos::FancyOStream& out) const
408{
409 out.setOutputToRootOnly(0);
410 bool isValidSetup = true;
411
412 if (!Stepper<Scalar>::isValidSetup(out)) isValidSetup = false;
413 if (!StepperExplicit<Scalar>::isValidSetup(out)) isValidSetup = false;
414
415 if (this->tableau_ == Teuchos::null) {
416 isValidSetup = false;
417 out << "The tableau is not set!\n";
418 }
419
420 if (this->stepperRKAppAction_ == Teuchos::null) {
421 isValidSetup = false;
422 out << "The AppAction is not set!\n";
423 }
424
425 return isValidSetup;
426}
427
428} // namespace Tempus
429#endif // Tempus_StepperExplicitRK_impl_hpp
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
virtual void setup(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel, bool useFSAL, std::string ICConsistency, bool ICConsistencyCheck, bool useEmbedded, const Teuchos::RCP< StepperRKAppAction< Scalar > > &stepperRKAppAction)
Setup for constructor.
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
virtual void setInitialConditions(const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory)
Set the initial conditions and make them consistent.
virtual Scalar getInitTimeStep(const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory) const
virtual void setupDefault()
Default setup for constructor.
virtual void setModel(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel)
Set model.
virtual Teuchos::RCP< Tempus::StepperState< Scalar > > getDefaultStepperState()
Get a default (initial) StepperState.
virtual bool isValidSetup(Teuchos::FancyOStream &out) const
virtual void takeStep(const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory)
Take the specified timestep, dt, and return true if successful.
virtual Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Teuchos::RCP< Teuchos::ParameterList > getValidParametersBasicERK() const
virtual void initialize()
Initialize during construction and after changing input parameters.
Thyra Base interface for implicit time steppers.
virtual void setInitialConditions(const Teuchos::RCP< SolutionHistory< Scalar > > &solutionHistory)
Set the initial conditions, make them consistent, and set needed memory.
virtual void setModel(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &appModel)
Set model.
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
Application Action for StepperRKBase.
StepperState is a simple class to hold state information about the stepper.
Thyra Base interface for time steppers.
virtual void initialize()
Initialize after construction and changing input parameters.
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const