Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_UnitTest_EDIRK_TrapezoidalRule.cpp
Go to the documentation of this file.
1 // @HEADER
2 // ****************************************************************************
3 // Tempus: Copyright (2017) Sandia Corporation
4 //
5 // Distributed under BSD 3-clause license (See accompanying file Copyright.txt)
6 // ****************************************************************************
7 // @HEADER
8 
10 
12 
13 #include "../TestModels/DahlquistTestModel.hpp"
14 
15 namespace Tempus_Unit_Test {
16 
18 using Teuchos::RCP;
19 using Teuchos::rcp;
20 using Teuchos::rcp_const_cast;
21 using Teuchos::rcp_dynamic_cast;
22 using Teuchos::sublist;
23 
24 // ************************************************************
25 // ************************************************************
26 TEUCHOS_UNIT_TEST(EDIRK_TrapezoidalRule, Default_Construction)
27 {
30 
31  // Test stepper properties.
32  TEUCHOS_ASSERT(stepper->getOrder() == 2);
33 }
34 
35 // ************************************************************
36 // ************************************************************
37 TEUCHOS_UNIT_TEST(EDIRK_TrapezoidalRule, StepperFactory_Construction)
38 {
39  auto model = rcp(new Tempus_Test::SinCosModel<double>());
40  testFactoryConstruction("RK Trapezoidal Rule", model);
41 }
42 
43 // ************************************************************
44 // ************************************************************
45 TEUCHOS_UNIT_TEST(EDIRK_TrapezoidalRule, AppAction)
46 {
48  auto model = rcp(new Tempus_Test::SinCosModel<double>());
49  testRKAppAction(stepper, model, out, success);
50 }
51 
52 // ************************************************************
53 // ************************************************************
54 
55 class StepperRKModifierEDIRK_TrapezoidaTest
56  : virtual public Tempus::StepperRKModifierBase<double> {
57  public:
59  StepperRKModifierEDIRK_TrapezoidaTest(Teuchos::FancyOStream &Out,
60  bool &Success)
61  : out(Out), success(Success)
62  {
63  }
64 
65  // FSAL
67  virtual ~StepperRKModifierEDIRK_TrapezoidaTest() {}
68 
70  virtual void modify(
74  {
75  const double relTol = 1.0e-14;
76  auto stageNumber = stepper->getStageNumber();
77  Teuchos::SerialDenseVector<int, double> c = stepper->getTableau()->c();
78 
79  auto currentState = sh->getCurrentState();
80  auto workingState = sh->getWorkingState();
81  const double dt = workingState->getTimeStep();
82  double time = currentState->getTime();
83  if (stageNumber >= 0) time += c(stageNumber) * dt;
84 
85  auto x = workingState->getX();
86  auto xDot = workingState->getXDot();
87  if (xDot == Teuchos::null) xDot = stepper->getStepperXDot();
88 
89  switch (actLoc) {
91  {
92  auto DME = Teuchos::rcp_dynamic_cast<
94  stepper->getModel());
95  TEST_FLOATING_EQUALITY(DME->getLambda(), -1.0, relTol);
96  }
97  TEST_FLOATING_EQUALITY(dt, 1.0, relTol);
98 
99  const double x_0 = get_ele(*(x), 0);
100  const double xDot_0 = get_ele(*(xDot), 0);
101  TEST_FLOATING_EQUALITY(x_0, 1.0, relTol); // Should be x_0
102  TEST_FLOATING_EQUALITY(xDot_0, -1.0, relTol); // Should be xDot_0
103  TEST_ASSERT(std::abs(time) < relTol);
104  TEST_FLOATING_EQUALITY(dt, 1.0, relTol);
105  TEST_COMPARE(stageNumber, ==, -1);
106  break;
107  }
110  const double X_i = get_ele(*(x), 0);
111  const double f_i = get_ele(*(xDot), 0);
112 
113  if (stageNumber == 0) {
114  TEST_FLOATING_EQUALITY(X_i, 1.0, relTol);
115  TEST_ASSERT(std::abs(f_i) < relTol);
116  TEST_ASSERT(std::abs(time) < relTol);
117  }
118  else if (stageNumber == 1) {
119  TEST_FLOATING_EQUALITY(X_i, 1.0, relTol);
120  TEST_FLOATING_EQUALITY(f_i, -1.0, relTol);
121  TEST_FLOATING_EQUALITY(time, 1.0, relTol);
122  }
123  else {
124  TEUCHOS_TEST_FOR_EXCEPT(!(-1 < stageNumber && stageNumber < 2));
125  }
126 
127  break;
128  }
132  const double X_i = get_ele(*(x), 0);
133  const double f_i = get_ele(*(xDot), 0);
134 
135  if (stageNumber == 0) {
136  // X_i = 1, f_1 = 0
137  TEST_FLOATING_EQUALITY(X_i, 1.0, relTol);
138  TEST_FLOATING_EQUALITY(f_i, -1.0, relTol);
139  TEST_ASSERT(std::abs(time) < relTol);
140  }
141  else if (stageNumber == 1) {
142  // X_i = , f_i =
143  TEST_FLOATING_EQUALITY(X_i, 1.0 / 3.0, relTol);
144  TEST_FLOATING_EQUALITY(f_i, -1.0 / 3.0, relTol);
145  TEST_FLOATING_EQUALITY(time, 1.0, relTol);
146  }
147  else {
148  TEUCHOS_TEST_FOR_EXCEPT(!(-1 < stageNumber && stageNumber < 2));
149  }
150 
151  break;
152  }
154  const double x_1 = get_ele(*(x), 0);
155  time = workingState->getTime();
156  TEST_FLOATING_EQUALITY(x_1, 1.0 / 3.0, relTol); // Should be x_1
157  TEST_FLOATING_EQUALITY(time, 1.0, relTol);
158  TEST_FLOATING_EQUALITY(dt, 1.0, relTol);
159  TEST_COMPARE(stageNumber, ==, -1);
160 
161  if (stepper->getUseEmbedded() == true) {
162  TEST_FLOATING_EQUALITY(workingState->getTolRel(), 1.0, relTol);
163  TEST_ASSERT(std::abs(workingState->getTolAbs()) < relTol);
164  // e = 0 from doxygen above.
165  TEST_ASSERT(std::abs(workingState->getErrorRel()) < relTol);
166  }
167  }
168  }
169  }
170 
171  private:
173  bool &success;
174 };
175 
176 // ************************************************************
177 // ************************************************************
178 TEUCHOS_UNIT_TEST(DIRK, EDIRK_Trapezoida_Modifier)
179 {
183 
184  auto modifier = rcp(new StepperRKModifierEDIRK_TrapezoidaTest(out, success));
185 
186  stepper->setModel(model);
187  stepper->setAppAction(modifier);
188  stepper->setICConsistency("Consistent");
189  stepper->setUseFSAL(false);
190  stepper->initialize();
191 
192  // Create a SolutionHistory.
193  auto solutionHistory = Tempus::createSolutionHistoryME(model);
194 
195  // Take one time step.
196  stepper->setInitialConditions(solutionHistory);
197  solutionHistory->initWorkingState();
198  double dt = 1.0;
199  solutionHistory->getWorkingState()->setTimeStep(dt);
200  solutionHistory->getWorkingState()->setTime(dt);
201  stepper->takeStep(solutionHistory); // Primary testing occurs in
202  // modifierX during takeStep().
203  // Test stepper properties.
204  TEUCHOS_ASSERT(stepper->getOrder() == 2);
205 }
206 } // 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...
RK Trapezoidal Rule (A.K.A. RK Crank-Nicolson)
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)