17 #include "../TestModels/DahlquistTestModel.hpp"
19 namespace Tempus_Unit_Test {
23 using Teuchos::rcp_const_cast;
24 using Teuchos::rcp_dynamic_cast;
118 class StepperRKModifierBogackiShampineTest
124 : out(Out), success(Success)
129 virtual ~StepperRKModifierBogackiShampineTest() {}
137 const double relTol = 1.0e-14;
138 auto stageNumber = stepper->getStageNumber();
141 auto currentState = sh->getCurrentState();
142 auto workingState = sh->getWorkingState();
143 const double dt = workingState->getTimeStep();
144 double time = currentState->getTime();
145 if (stageNumber >= 0) time += c(stageNumber) * dt;
147 auto x = workingState->getX();
148 auto xDot = workingState->getXDot();
149 if (xDot == Teuchos::null) xDot = stepper->getStepperXDot();
153 auto DME = Teuchos::rcp_dynamic_cast<
155 stepper->getModel());
160 const double x_0 = get_ele(*(x), 0);
161 const double xDot_0 = get_ele(*(xDot), 0);
173 const double X_i = get_ele(*(x), 0);
174 const double f_i = get_ele(*(xDot), 0);
175 if (stageNumber == 0) {
179 !stepper->getUseFSAL()) {
186 else if (stageNumber == 1) {
196 else if (stageNumber == 2) {
207 else if (stageNumber == 3) {
214 else if (workingState->getNConsecutiveFailures() > 0) {
227 const double x_1 = get_ele(*(x), 0);
228 time = workingState->getTime();
234 if (stepper->getUseEmbedded() ==
true) {
236 TEST_ASSERT(std::abs(workingState->getTolAbs()) < relTol);
238 TEST_ASSERT(std::abs(workingState->getErrorRel()) < relTol);
243 "Error - unknown action location.\n");
259 auto modifier =
rcp(
new StepperRKModifierBogackiShampineTest(out, success));
261 stepper->setModel(model);
262 stepper->setAppAction(modifier);
263 stepper->setICConsistency(
"Consistent");
264 stepper->setUseFSAL(
false);
265 stepper->initialize();
271 stepper->setInitialConditions(solutionHistory);
272 solutionHistory->initWorkingState();
274 solutionHistory->getWorkingState()->setTimeStep(dt);
275 solutionHistory->getWorkingState()->setTime(dt);
276 stepper->takeStep(solutionHistory);
289 auto modifier =
rcp(
new StepperRKModifierBogackiShampineTest(out, success));
291 stepper->setModel(model);
292 stepper->setAppAction(modifier);
293 stepper->setICConsistency(
"Consistent");
294 stepper->setUseFSAL(
true);
295 stepper->initialize();
301 stepper->setInitialConditions(solutionHistory);
302 solutionHistory->initWorkingState();
304 solutionHistory->getWorkingState()->setTimeStep(dt);
305 solutionHistory->getWorkingState()->setTime(dt);
306 stepper->takeStep(solutionHistory);
319 auto modifier =
rcp(
new StepperRKModifierBogackiShampineTest(out, success));
321 stepper->setModel(model);
322 stepper->setAppAction(modifier);
323 stepper->setICConsistency(
"Consistent");
324 stepper->setUseFSAL(
true);
325 stepper->initialize();
331 stepper->setInitialConditions(solutionHistory);
332 solutionHistory->initWorkingState();
334 solutionHistory->getWorkingState()->setTimeStep(dt);
335 solutionHistory->getWorkingState()->setTime(dt);
336 solutionHistory->getWorkingState()->setNConsecutiveFailures(1);
337 stepper->takeStep(solutionHistory);
350 auto modifier =
rcp(
new StepperRKModifierBogackiShampineTest(out, success));
352 stepper->setModel(model);
353 stepper->setAppAction(modifier);
354 stepper->setICConsistency(
"Consistent");
355 stepper->setUseFSAL(
false);
356 stepper->setUseEmbedded(
true);
357 stepper->initialize();
363 stepper->setInitialConditions(solutionHistory);
364 solutionHistory->initWorkingState();
366 solutionHistory->getWorkingState()->setTimeStep(dt);
367 solutionHistory->getWorkingState()->setTime(dt);
368 solutionHistory->getWorkingState()->setTolRel(1.0);
369 solutionHistory->getWorkingState()->setTolAbs(0.0);
370 stepper->takeStep(solutionHistory);
383 class StepperRKModifierXBogackiShampineTest
389 : out(Out), success(Success)
394 virtual ~StepperRKModifierXBogackiShampineTest() {}
399 const double dt,
const int stageNumber,
403 const double relTol = 1.0e-14;
405 const double x_0 = get_ele(*(x), 0);
417 const double X_i = get_ele(*(x), 0);
418 if (stageNumber == 0) {
422 else if (stageNumber == 1) {
426 else if (stageNumber == 2) {
430 else if (stageNumber == 3) {
440 const double x_1 = get_ele(*(x), 0);
448 "Error - unknown action location.\n");
464 auto modifierX =
rcp(
new StepperRKModifierXBogackiShampineTest(out, success));
466 stepper->setModel(model);
467 stepper->setAppAction(modifierX);
468 stepper->setICConsistency(
"Consistent");
469 stepper->setUseFSAL(
false);
470 stepper->initialize();
476 stepper->setInitialConditions(solutionHistory);
477 solutionHistory->initWorkingState();
479 solutionHistory->getWorkingState()->setTimeStep(dt);
480 solutionHistory->getWorkingState()->setTime(dt);
481 stepper->takeStep(solutionHistory);
494 auto modifierX =
rcp(
new StepperRKModifierXBogackiShampineTest(out, success));
496 stepper->setModel(model);
497 stepper->setAppAction(modifierX);
498 stepper->setICConsistency(
"Consistent");
499 stepper->setUseFSAL(
true);
500 stepper->initialize();
506 stepper->setInitialConditions(solutionHistory);
507 solutionHistory->initWorkingState();
509 solutionHistory->getWorkingState()->setTimeStep(dt);
510 solutionHistory->getWorkingState()->setTime(dt);
511 stepper->takeStep(solutionHistory);
524 auto modifierX =
rcp(
new StepperRKModifierXBogackiShampineTest(out, success));
526 stepper->setModel(model);
527 stepper->setAppAction(modifierX);
528 stepper->setICConsistency(
"Consistent");
529 stepper->setUseFSAL(
true);
530 stepper->initialize();
536 stepper->setInitialConditions(solutionHistory);
537 solutionHistory->initWorkingState();
539 solutionHistory->getWorkingState()->setTimeStep(dt);
540 solutionHistory->getWorkingState()->setTime(dt);
541 solutionHistory->getWorkingState()->setNConsecutiveFailures(1);
542 stepper->takeStep(solutionHistory);
At the beginning of the stage.
At the beginning of the step.
Modify at the end of the step.
Modify after the implicit solve.
Base ModifierX for StepperRK.
Modify at the beginning of the stage.
The classic Dahlquist Test Problem.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
#define TEST_COMPARE(v1, comp, v2)
#define TEST_FLOATING_EQUALITY(v1, v2, tol)
Before the implicit solve.
Teuchos::RCP< SolutionHistory< Scalar > > createSolutionHistoryME(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model)
Nonmember contructor from a Thyra ModelEvaluator.
Modify before the implicit solve.
Base class for Runge-Kutta methods, ExplicitRK, DIRK and IMEX.
Modify at the beginning of the step.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Explicit RK Bogacki-Shampine Butcher Tableau.
TEUCHOS_UNIT_TEST(BackwardEuler, Default_Construction)
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Modify at the end of the stage.
MODIFIER_TYPE
Indicates the location of application action (see algorithm).
Modify before the explicit evaluation.
ACTION_LOCATION
Indicates the location of application action (see algorithm).
After the implicit solve.
Before the explicit evaluation.
#define TEUCHOS_ASSERT(assertion_test)
Base modifier for StepperRK.
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
virtual void modify(Teuchos::RCP< Thyra::VectorBase< double > >, const double, const double, const int, const MODIFIER_TYPE modType)=0
Modify solution based on the MODIFIER_TYPE.