Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_UnitTest_SDIRK_21Pair.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 
13 
14 #include "../TestModels/DahlquistTestModel.hpp"
15 
16 namespace Tempus_Unit_Test {
17 
19 using Teuchos::RCP;
20 using Teuchos::rcp;
21 using Teuchos::rcp_const_cast;
22 using Teuchos::rcp_dynamic_cast;
23 using Teuchos::sublist;
24 
25 // ************************************************************
26 // ************************************************************
27 TEUCHOS_UNIT_TEST(SDIRK_21Pair, Default_Construction)
28 {
29  auto stepper = rcp(new Tempus::StepperSDIRK_21Pair<double>());
31 
32  // Test stepper properties.
33  TEUCHOS_ASSERT(stepper->getOrder() == 2);
34 }
35 
36 // ************************************************************
37 // ************************************************************
38 TEUCHOS_UNIT_TEST(SDIRK_21Pair, StepperFactory_Construction)
39 {
40  auto model = rcp(new Tempus_Test::SinCosModel<double>());
41  testFactoryConstruction("SDIRK 2(1) Pair", model);
42 }
43 
44 // ************************************************************
45 // ************************************************************
46 TEUCHOS_UNIT_TEST(SDIRK_21Pair, AppAction)
47 {
48  auto stepper = rcp(new Tempus::StepperSDIRK_21Pair<double>());
49  auto model = rcp(new Tempus_Test::SinCosModel<double>());
50  testRKAppAction(stepper, model, out, success);
51 }
52 
53 // ************************************************************
54 // ************************************************************
55 class StepperRKModifierSDIRK21Test
56  : virtual public Tempus::StepperRKModifierBase<double> {
57  public:
59  StepperRKModifierSDIRK21Test(Teuchos::FancyOStream &Out, bool &Success)
60  : out(Out), success(Success)
61  {
62  }
63 
65  virtual ~StepperRKModifierSDIRK21Test() {}
66 
68  virtual void modify(
72  {
73  const double relTol = 1.0e-14;
74  auto stageNumber = stepper->getStageNumber();
75  Teuchos::SerialDenseVector<int, double> c = stepper->getTableau()->c();
76 
77  auto currentState = sh->getCurrentState();
78  auto workingState = sh->getWorkingState();
79  const double dt = workingState->getTimeStep();
80  double time = currentState->getTime();
81  if (stageNumber >= 0) time += c(stageNumber) * dt;
82 
83  auto x = workingState->getX();
84  auto xDot = workingState->getXDot();
85  if (xDot == Teuchos::null) xDot = stepper->getStepperXDot();
86 
87  switch (actLoc) {
89  {
90  auto DME = Teuchos::rcp_dynamic_cast<
92  stepper->getModel());
93  TEST_FLOATING_EQUALITY(DME->getLambda(), -1.0, relTol);
94  }
95  TEST_FLOATING_EQUALITY(dt, 1.0, relTol);
96 
97  const double x_0 = get_ele(*(x), 0);
98  const double xDot_0 = get_ele(*(xDot), 0);
99  TEST_FLOATING_EQUALITY(x_0, 1.0, relTol); // Should be x_0
100  TEST_FLOATING_EQUALITY(xDot_0, -1.0, relTol); // Should be xDot_0
101  TEST_ASSERT(std::abs(time) < relTol);
102  TEST_FLOATING_EQUALITY(dt, 1.0, relTol);
103  TEST_COMPARE(stageNumber, ==, -1);
104  break;
105  }
108  const double X_i = get_ele(*(x), 0);
109  const double f_i = get_ele(*(xDot), 0);
110 
111  if (stageNumber == 0) {
112  TEST_FLOATING_EQUALITY(X_i, 1.0, relTol);
113  TEST_FLOATING_EQUALITY(time, 1.0, relTol);
114  TEST_ASSERT(std::abs(f_i) < relTol);
115  }
116  else if (stageNumber == 1) {
117  TEST_FLOATING_EQUALITY(X_i, 0.5, relTol);
118  TEST_FLOATING_EQUALITY(f_i, -1.0, relTol);
119  TEST_ASSERT(std::abs(time) < relTol);
120  }
121  else {
122  TEUCHOS_TEST_FOR_EXCEPT(!(-1 < stageNumber && stageNumber < 2));
123  }
124 
125  break;
126  }
130  const double X_i = get_ele(*(x), 0); // 0.5
131  const double f_i = get_ele(*(xDot), 0); // -0.5
132 
133  if (stageNumber == 0) {
134  // X_i = 1, f_1 = 0
135  TEST_FLOATING_EQUALITY(X_i, 0.5, relTol);
136  TEST_FLOATING_EQUALITY(time, 1.0, relTol);
137  TEST_FLOATING_EQUALITY(f_i, -0.5, relTol);
138  }
139  else if (stageNumber == 1) {
140  // X_i = , f_i =
141  TEST_FLOATING_EQUALITY(X_i, 0.75, relTol);
142  TEST_FLOATING_EQUALITY(f_i, -0.75, relTol);
143  TEST_ASSERT(std::abs(time) < relTol);
144  }
145  else {
146  TEUCHOS_TEST_FOR_EXCEPT(!(-1 < stageNumber && stageNumber < 2));
147  }
148 
149  break;
150  }
152  const double x_1 = get_ele(*(x), 0);
153  time = workingState->getTime();
154  TEST_FLOATING_EQUALITY(x_1, 3.0 / 8.0, relTol); // Should be x_1
155  TEST_FLOATING_EQUALITY(time, 1.0, relTol);
156  TEST_FLOATING_EQUALITY(dt, 1.0, relTol);
157  TEST_COMPARE(stageNumber, ==, -1);
158 
159  if (stepper->getUseEmbedded() == true) {
160  TEST_FLOATING_EQUALITY(workingState->getTolRel(), 1.0, relTol);
161  TEST_ASSERT(std::abs(workingState->getTolAbs()) < relTol);
162  // e = 0 from doxygen above.
163  TEST_ASSERT(std::abs(workingState->getErrorRel()) < relTol);
164  }
165  }
166  }
167  }
168 
169  private:
171  bool &success;
172 };
173 
174 // ************************************************************
175 // ************************************************************
176 TEUCHOS_UNIT_TEST(SDIRK_21Pair, Modifier)
177 {
178  auto stepper = rcp(new Tempus::StepperSDIRK_21Pair<double>());
181 
182  auto modifier = rcp(new StepperRKModifierSDIRK21Test(out, success));
183 
184  stepper->setModel(model);
185  stepper->setAppAction(modifier);
186  stepper->setICConsistency("Consistent");
187  stepper->setUseFSAL(false);
188  stepper->initialize();
189 
190  // Create a SolutionHistory.
191  auto solutionHistory = Tempus::createSolutionHistoryME(model);
192 
193  // Take one time step.
194  stepper->setInitialConditions(solutionHistory);
195  solutionHistory->initWorkingState();
196  double dt = 1.0;
197  solutionHistory->getWorkingState()->setTimeStep(dt);
198  solutionHistory->getWorkingState()->setTime(dt);
199  stepper->takeStep(solutionHistory); // Primary testing occurs in
200  // modifierX during takeStep().
201  // Test stepper properties.
202  TEUCHOS_ASSERT(stepper->getOrder() == 2);
203 }
204 
205 } // namespace Tempus_Unit_Test
void testDIRKAccessorsFullConstruction(const RCP< Tempus::StepperDIRK< double >> &stepper)
Unit test utility for ExplicitRK Stepper construction and accessors.
void testFactoryConstruction(std::string stepperType, const Teuchos::RCP< const Thyra::ModelEvaluator< double > > &model)
Unit test utility for Stepper construction through StepperFactory.
The classic Dahlquist Test Problem.
#define TEST_COMPARE(v1, comp, v2)
#define TEST_FLOATING_EQUALITY(v1, v2, tol)
#define TEST_ASSERT(v1)
Sine-Cosine model problem from Rythmos. This is a canonical Sine-Cosine differential equation with a...
Teuchos::RCP< SolutionHistory< Scalar > > createSolutionHistoryME(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model)
Nonmember contructor from a Thyra ModelEvaluator.
Base class for Runge-Kutta methods, ExplicitRK, DIRK and IMEX.
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_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
TEUCHOS_UNIT_TEST(BackwardEuler, Default_Construction)
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
ACTION_LOCATION
Indicates the location of application action (see algorithm).
#define TEUCHOS_ASSERT(assertion_test)
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)