Tempus Version of the Day
Time Integration
Loading...
Searching...
No Matches
Tempus_UnitTest_ERK_ForwardEuler.cpp
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
11
12#include "../TestModels/DahlquistTestModel.hpp"
13
14namespace Tempus_Unit_Test {
15
16using Teuchos::RCP;
17using Teuchos::rcp;
18using Teuchos::rcp_const_cast;
19using Teuchos::rcp_dynamic_cast;
20
21// ************************************************************
22// ************************************************************
23TEUCHOS_UNIT_TEST(ERK_ForwardEuler, Default_Construction)
24{
25 auto stepper = rcp(new Tempus::StepperERK_ForwardEuler<double>());
27
28 // Test stepper properties.
29 TEUCHOS_ASSERT(stepper->getOrder() == 1);
30 const auto rk_fe = stepper->getTableau();
31
32 TEUCHOS_ASSERT(rk_fe->isTVD());
33 TEUCHOS_ASSERT(rk_fe->getTVDCoeff() == 1);
34}
35
36// ************************************************************
37// ************************************************************
38TEUCHOS_UNIT_TEST(ERK_ForwardEuler, StepperFactory_Construction)
39{
40 auto model = rcp(new Tempus_Test::SinCosModel<double>());
41 testFactoryConstruction("RK Forward Euler", model);
42}
43
44// ************************************************************
45// ************************************************************
46TEUCHOS_UNIT_TEST(ERK_ForwardEuler, AppAction)
47{
48 auto stepper = rcp(new Tempus::StepperERK_ForwardEuler<double>());
49 auto model = rcp(new Tempus_Test::SinCosModel<double>());
50 testRKAppAction(stepper, model, out, success);
51}
52
53// ************************************************************
54// ************************************************************
55TEUCHOS_UNIT_TEST(ERK_ForwardEuler, FSAL)
56{
57 auto stepper = rcp(new Tempus::StepperERK_ForwardEuler<double>());
58 Teuchos::RCP<const Thyra::ModelEvaluator<double> > model =
59 rcp(new Tempus_Test::DahlquistTestModel<double>(-1.0, true));
60
61 stepper->setModel(model);
62 stepper->setICConsistency("Consistent");
63 stepper->setUseFSAL(true);
64 stepper->initialize();
65
66 // Create a SolutionHistory.
67 auto solutionHistory = Tempus::createSolutionHistoryME(model);
68
69 // Take one time step.
70 stepper->setInitialConditions(solutionHistory);
71 solutionHistory->initWorkingState();
72 double dt = 1.0;
73 solutionHistory->getWorkingState()->setTimeStep(dt);
74 solutionHistory->getWorkingState()->setTime(dt);
75 stepper->takeStep(solutionHistory);
76
77 // Test solution.
78 const double relTol = 1.0e-14;
79
80 // ICs
81 auto currentState = solutionHistory->getCurrentState();
82 const double x_0 = get_ele(*(currentState->getX()), 0);
83 const double xDot_0 = get_ele(*(currentState->getXDot()), 0);
84 TEST_FLOATING_EQUALITY(x_0, 1.0, relTol);
85 TEST_FLOATING_EQUALITY(xDot_0, -1.0, relTol);
86 TEST_ASSERT(std::abs(currentState->getTime()) < relTol);
87
88 // After one step.
89 auto workingState = solutionHistory->getWorkingState();
90 const double x_1 = get_ele(*(workingState->getX()), 0);
91 const double xDot_1 = get_ele(*(workingState->getXDot()), 0);
92 // out << "xDot_1 = " << xDot_1 << std::endl;
93 TEST_ASSERT(std::abs(x_1) < relTol);
94 TEST_ASSERT(std::abs(xDot_1) < relTol);
95 TEST_FLOATING_EQUALITY(workingState->getTime(), 1.0, relTol);
96}
97
98} // namespace Tempus_Unit_Test
Forward Euler Runge-Kutta Butcher Tableau.
The classic Dahlquist Test Problem.
Sine-Cosine model problem from Rythmos. This is a canonical Sine-Cosine differential equation.
void testExplicitRKAccessorsFullConstruction(const RCP< Tempus::StepperExplicitRK< double > > &stepper)
Unit test utility for ExplicitRK Stepper construction and accessors.
void testRKAppAction(const Teuchos::RCP< Tempus::StepperRKBase< double > > &stepper, const Teuchos::RCP< const Thyra::ModelEvaluator< double > > &model, Teuchos::FancyOStream &out, bool &success)
Unit test utility for Stepper RK AppAction.
TEUCHOS_UNIT_TEST(BackwardEuler, Default_Construction)
void testFactoryConstruction(std::string stepperType, const Teuchos::RCP< const Thyra::ModelEvaluator< double > > &model)
Unit test utility for Stepper construction through StepperFactory.
Teuchos::RCP< SolutionHistory< Scalar > > createSolutionHistoryME(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model)
Nonmember contructor from a Thyra ModelEvaluator.