Tempus Version of the Day
Time Integration
Loading...
Searching...
No Matches
Tempus_UnitTest_NewmarkExplicitAForm.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 "Teuchos_XMLParameterListHelpers.hpp"
13
14#include "Tempus_StepperNewmarkExplicitAForm.hpp"
19
20#include "../TestModels/HarmonicOscillatorModel.hpp"
21
22namespace Tempus_Unit_Test {
23
24using Teuchos::ParameterList;
25using Teuchos::RCP;
26using Teuchos::rcp;
27using Teuchos::rcp_const_cast;
28using Teuchos::rcp_dynamic_cast;
29using Teuchos::sublist;
30
31// ************************************************************
32// ************************************************************
33class StepperNewmarkExplicitAFormModifierTest
35 public:
37 StepperNewmarkExplicitAFormModifierTest()
38 : testBEGIN_STEP(false),
39 testBEFORE_EXPLICIT_EVAL(false),
40 testAFTER_EXPLICIT_EVAL(false),
41 testEND_STEP(false),
42 testCurrentValue(-0.99),
43 testDt(-1.5),
44 testName("")
45 {
46 }
47
49 virtual ~StepperNewmarkExplicitAFormModifierTest() {}
50
52 virtual void modify(
53 Teuchos::RCP<Tempus::SolutionHistory<double>> sh,
56 double>::ACTION_LOCATION actLoc)
57 {
58 switch (actLoc) {
59 case StepperNewmarkExplicitAFormAppAction<double>::BEGIN_STEP: {
60 testBEGIN_STEP = true;
61 break;
62 }
63 case StepperNewmarkExplicitAFormAppAction<double>::BEFORE_EXPLICIT_EVAL: {
64 testBEFORE_EXPLICIT_EVAL = true;
65 testName = "Newmark Explicit A Form - Modifier";
66 stepper->setStepperName(testName);
67 break;
68 }
69 case StepperNewmarkExplicitAFormAppAction<double>::AFTER_EXPLICIT_EVAL: {
70 testAFTER_EXPLICIT_EVAL = true;
71 testDt = sh->getWorkingState()->getTimeStep();
72 // sh->getWorkingState()->setTimeStep(testDt);
73 break;
74 }
75 case StepperNewmarkExplicitAFormAppAction<double>::END_STEP: {
76 testEND_STEP = true;
77 auto x = sh->getWorkingState()->getX();
78 testCurrentValue = get_ele(*(x), 0);
79 break;
80 }
81 default:
82 TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
83 "Error - unknown action location.\n");
84 }
85 }
86
87 bool testBEGIN_STEP;
88 bool testBEFORE_EXPLICIT_EVAL;
89 bool testAFTER_EXPLICIT_EVAL;
90 bool testEND_STEP;
91 double testCurrentValue;
92 double testDt;
93 std::string testName;
94};
95
96// ************************************************************
97// ************************************************************
98class StepperNewmarkExplicitAFormModifierXTest
100 public:
102 StepperNewmarkExplicitAFormModifierXTest()
103 : testX_BEGIN_STEP(false),
104 testX_BEFORE_EXPLICIT_EVAL(false),
105 testX_AFTER_EXPLICIT_EVAL(false),
106 testX_END_STEP(false),
107 testX(-0.99),
108 testXDot(-0.99),
109 testDt(-1.5),
110 testTime(-1.5)
111 {
112 }
113
115 virtual ~StepperNewmarkExplicitAFormModifierXTest() {}
116
118 virtual void modify(
119 Teuchos::RCP<Thyra::VectorBase<double>> x, const double time,
120 const double dt,
122 double>::MODIFIER_TYPE modType)
123 {
124 switch (modType) {
125 case StepperNewmarkExplicitAFormModifierXBase<double>::X_BEGIN_STEP: {
126 testX_BEGIN_STEP = true;
127 testDt = dt;
128 testX = get_ele(*(x), 0);
129 break;
130 }
131 case StepperNewmarkExplicitAFormModifierXBase<
133 testX_BEFORE_EXPLICIT_EVAL = true;
134 testTime = time;
135 testX = get_ele(*(x), 0);
136 break;
137 }
138 case StepperNewmarkExplicitAFormModifierXBase<
140 testX_AFTER_EXPLICIT_EVAL = true;
141 testX = get_ele(*(x), 0);
142 break;
143 }
144 case StepperNewmarkExplicitAFormModifierXBase<double>::X_END_STEP: {
145 testX_END_STEP = true;
146 testTime = time;
147 testX = get_ele(*(x), 0);
148 break;
149 }
150 default:
151 TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
152 "Error - unknown action location.\n");
153 }
154 }
155
156 bool testX_BEGIN_STEP;
157 bool testX_BEFORE_EXPLICIT_EVAL;
158 bool testX_AFTER_EXPLICIT_EVAL;
159 bool testX_END_STEP;
160 double testX;
161 double testXDot;
162 double testDt;
163 double testTime;
164};
165
166TEUCHOS_UNIT_TEST(NewmarkExplicitAForm, AppAction_Modifier)
167{
168 using Teuchos::ParameterList;
169 using Teuchos::RCP;
170 using Teuchos::sublist;
171
172 RCP<Tempus::IntegratorBasic<double>> integrator;
173 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
174 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
175
176 // Read params from .xml file
177 RCP<ParameterList> pList = Teuchos::getParametersFromXmlFile(
178 "Tempus_NewmarkExplicitAForm_HarmonicOscillator_Damped.xml");
179
180 // Setup the HarmonicOscillatorModel
181 RCP<ParameterList> hom_pl = sublist(pList, "HarmonicOscillatorModel", true);
182 RCP<Tempus_Test::HarmonicOscillatorModel<double>> model =
183 Teuchos::rcp(new Tempus_Test::HarmonicOscillatorModel<double>(hom_pl));
184
185 // Setup the Integrator and reset initial time step
186 RCP<ParameterList> pl = sublist(pList, "Tempus", true);
187 RCP<ParameterList> stepperPL = sublist(pl, "Default Stepper", true);
188 stepperPL->remove("Zero Initial Guess");
189
190 double dt = pl->sublist("Default Integrator")
191 .sublist("Time Step Control")
192 .get<double>("Initial Time Step");
193 dt *= 2.0;
194
195 pl->sublist("Default Integrator")
196 .sublist("Time Step Control")
197 .set("Initial Time Step", dt);
198 integrator = Tempus::createIntegratorBasic<double>(pl, model);
199
200 RCP<Tempus::StepperNewmarkExplicitAForm<double>> stepper =
201 Teuchos::rcp_dynamic_cast<Tempus::StepperNewmarkExplicitAForm<double>>(
202 integrator->getStepper(), true);
203
204 auto modifier = rcp(new StepperNewmarkExplicitAFormModifierTest());
205 stepper->setAppAction(modifier);
206 stepper->initialize();
207 integrator->initialize();
208
209 // Integrate to timeMax
210 bool integratorStatus = integrator->advanceTime();
211 TEST_ASSERT(integratorStatus)
212
213 // Testing that each ACTION_LOCATION has been called.
214 TEST_COMPARE(modifier->testBEGIN_STEP, ==, true);
215 TEST_COMPARE(modifier->testBEFORE_EXPLICIT_EVAL, ==, true);
216 TEST_COMPARE(modifier->testAFTER_EXPLICIT_EVAL, ==, true);
217 TEST_COMPARE(modifier->testEND_STEP, ==, true);
218
219 // Testing that values can be set through the Modifier.
220 auto x = integrator->getX();
221 auto Dt = integrator->getTime();
222 TEST_FLOATING_EQUALITY(modifier->testDt, Dt, 1.0e-14);
223 TEST_FLOATING_EQUALITY(modifier->testCurrentValue, get_ele(*(x), 0), 1.0e-14);
224 TEST_COMPARE(modifier->testName, ==, stepper->getStepperName());
225}
226
227TEUCHOS_UNIT_TEST(NewmarkExplicitAForm, AppAction_ModifierX)
228{
229 using Teuchos::ParameterList;
230 using Teuchos::RCP;
231 using Teuchos::sublist;
232
233 RCP<Tempus::IntegratorBasic<double>> integrator;
234 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
235 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
236
237 // Read params from .xml file
238 RCP<ParameterList> pList = Teuchos::getParametersFromXmlFile(
239 "Tempus_NewmarkExplicitAForm_HarmonicOscillator_Damped.xml");
240
241 // Setup the HarmonicOscillatorModel
242 RCP<ParameterList> hom_pl = sublist(pList, "HarmonicOscillatorModel", true);
243 RCP<Tempus_Test::HarmonicOscillatorModel<double>> model =
244 Teuchos::rcp(new Tempus_Test::HarmonicOscillatorModel<double>(hom_pl));
245
246 // Setup the Integrator and reset initial time step
247 RCP<ParameterList> pl = sublist(pList, "Tempus", true);
248 RCP<ParameterList> stepperPL = sublist(pl, "Default Stepper", true);
249 stepperPL->remove("Zero Initial Guess");
250
251 double dt = pl->sublist("Default Integrator")
252 .sublist("Time Step Control")
253 .get<double>("Initial Time Step");
254 dt *= 2.0;
255
256 pl->sublist("Default Integrator")
257 .sublist("Time Step Control")
258 .set("Initial Time Step", dt);
259 integrator = Tempus::createIntegratorBasic<double>(pl, model);
260
261 RCP<Tempus::StepperNewmarkExplicitAForm<double>> stepper =
262 Teuchos::rcp_dynamic_cast<Tempus::StepperNewmarkExplicitAForm<double>>(
263 integrator->getStepper(), true);
264
265 auto modifierX = rcp(new StepperNewmarkExplicitAFormModifierXTest());
266 stepper->setAppAction(modifierX);
267 stepper->initialize();
268 integrator->initialize();
269
270 // Integrate to timeMax
271 bool integratorStatus = integrator->advanceTime();
272 TEST_ASSERT(integratorStatus)
273
274 // Testing that each ACTION_LOCATION has been called.
275 TEST_COMPARE(modifierX->testX_BEGIN_STEP, ==, true);
276 TEST_COMPARE(modifierX->testX_BEFORE_EXPLICIT_EVAL, ==, true);
277 TEST_COMPARE(modifierX->testX_AFTER_EXPLICIT_EVAL, ==, true);
278 TEST_COMPARE(modifierX->testX_END_STEP, ==, true);
279
280 // Testing that values can be set through the Modifier.
281 auto Dt = integrator->getTime();
282 TEST_FLOATING_EQUALITY(modifierX->testDt, Dt, 1.0e-14);
283
284 const auto x = integrator->getX();
285 TEST_FLOATING_EQUALITY(modifierX->testX, get_ele(*(x), 0), 1.0e-14);
286}
287
288} // namespace Tempus_Unit_Test
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
ACTION_LOCATION
Indicates the location of application action (see algorithm).
virtual void modify(Teuchos::RCP< Thyra::VectorBase< double > >, const double, const double, const MODIFIER_TYPE modType)=0
Modify solution based on the MODIFIER_TYPE.
MODIFIER_TYPE
Indicates the location of application action (see algorithm).
TEUCHOS_UNIT_TEST(BackwardEuler, Default_Construction)