Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Tempus_ObserverTest.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 
9 #include "Teuchos_UnitTestHarness.hpp"
10 #include "Teuchos_XMLParameterListHelpers.hpp"
11 #include "Teuchos_TimeMonitor.hpp"
12 
13 #include "Thyra_VectorStdOps.hpp"
14 
15 #include "Tempus_IntegratorBasic.hpp"
16 #include "Tempus_IntegratorObserverLogging.hpp"
17 #include "Tempus_IntegratorObserverComposite.hpp"
18 
19 #include "Tempus_StepperRKObserverLogging.hpp"
20 #include "Tempus_StepperRKObserverComposite.hpp"
21 
22 #include "../TestModels/SinCosModel.hpp"
23 #include "../TestModels/VanDerPolModel.hpp"
24 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
25 
26 #include <vector>
27 
28 namespace Tempus_Test {
29 
30 using Teuchos::RCP;
31 using Teuchos::rcp;
32 using Teuchos::ParameterList;
33 using Teuchos::sublist;
34 using Teuchos::getParametersFromXmlFile;
35 
39 
40 #define TEST_INTEGRATOROBSERVERLOGGING
41 #define TEST_INTEGRATOROBSERVERCOMPOSITE
42 #define TEST_STEPPERRKOBSERVERLOGGING
43 
44 
45 // ************************************************************
46 // ************************************************************
47 #ifdef TEST_INTEGRATOROBSERVERLOGGING
48 TEUCHOS_UNIT_TEST(Observer, IntegratorObserverLogging)
49 {
50  // Read params from .xml file
51  RCP<ParameterList> pList =
52  getParametersFromXmlFile("Tempus_Observer_SinCos.xml");
53 
54  // Setup the SinCosModel
55  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
56  RCP<SinCosModel<double> > model =
57  Teuchos::rcp(new SinCosModel<double> (scm_pl));
58 
59  // Setup the Integrator and reset initial time step
60  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
61  RCP<Tempus::IntegratorBasic<double> > integrator =
62  Tempus::integratorBasic<double>(pl, model);
63 
64  RCP<Tempus::IntegratorObserverLogging<double> > loggingObs =
66  integrator->setObserver(loggingObs);
67 
68  // Integrate to timeMax
69  bool integratorStatus = integrator->advanceTime();
70  TEST_ASSERT(integratorStatus)
71 
72  // Test if at 'Final Time'
73  double time = integrator->getTime();
74  double timeFinal = pl->sublist("Demo Integrator")
75  .sublist("Time Step Control").get<double>("Final Time");
76  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
77 
78  // Construct the reference counter and order for comparison.
79  std::map<std::string,int> refCounters;
80  std::list<std::string> refOrder;
81 
82  refCounters[loggingObs->nameObserveStartIntegrator_ ] = 1;
83  refCounters[loggingObs->nameObserveStartTimeStep_ ] = 10;
84  refCounters[loggingObs->nameObserveNextTimeStep_ ] = 10;
85  refCounters[loggingObs->nameObserveBeforeTakeStep_ ] = 10;
86  refCounters[loggingObs->nameObserveAfterTakeStep_ ] = 10;
87  refCounters[loggingObs->nameObserveAfterCheckTimeStep_] = 10;
88  refCounters[loggingObs->nameObserveEndTimeStep_ ] = 10;
89  refCounters[loggingObs->nameObserveEndIntegrator_ ] = 1;
90 
91  refOrder.push_back(loggingObs->nameObserveStartIntegrator_ );
92  for (int i=0 ; i<10; ++i) {
93  refOrder.push_back(loggingObs->nameObserveStartTimeStep_ );
94  refOrder.push_back(loggingObs->nameObserveNextTimeStep_ );
95  refOrder.push_back(loggingObs->nameObserveBeforeTakeStep_ );
96  refOrder.push_back(loggingObs->nameObserveAfterTakeStep_ );
97  refOrder.push_back(loggingObs->nameObserveAfterCheckTimeStep_);
98  refOrder.push_back(loggingObs->nameObserveEndTimeStep_ );
99  }
100  refOrder.push_back(loggingObs->nameObserveEndIntegrator_ );
101 
102  const std::map<std::string,int>& counters = *(loggingObs->getCounters());
103  const std::list<std::string>& order = *(loggingObs->getOrder());
104 
105  // Compare against reference.
106  TEST_EQUALITY(
107  counters.find(loggingObs->nameObserveStartIntegrator_ )->second,
108  refCounters.find(loggingObs->nameObserveStartIntegrator_ )->second);
109  TEST_EQUALITY(
110  counters.find(loggingObs->nameObserveStartTimeStep_ )->second,
111  refCounters.find(loggingObs->nameObserveStartTimeStep_ )->second);
112  TEST_EQUALITY(
113  counters.find(loggingObs->nameObserveNextTimeStep_ )->second,
114  refCounters.find(loggingObs->nameObserveNextTimeStep_ )->second);
115  TEST_EQUALITY(
116  counters.find(loggingObs->nameObserveBeforeTakeStep_ )->second,
117  refCounters.find(loggingObs->nameObserveBeforeTakeStep_ )->second);
118  TEST_EQUALITY(
119  counters.find(loggingObs->nameObserveAfterTakeStep_ )->second,
120  refCounters.find(loggingObs->nameObserveAfterTakeStep_ )->second);
121  TEST_EQUALITY(
122  counters.find(loggingObs->nameObserveAfterCheckTimeStep_)->second,
123  refCounters.find(loggingObs->nameObserveAfterCheckTimeStep_)->second);
124  TEST_EQUALITY(
125  counters.find(loggingObs->nameObserveEndTimeStep_ )->second,
126  refCounters.find(loggingObs->nameObserveEndTimeStep_ )->second);
127  TEST_EQUALITY(
128  counters.find(loggingObs->nameObserveEndIntegrator_ )->second,
129  refCounters.find(loggingObs->nameObserveEndIntegrator_ )->second);
130 
131  TEUCHOS_ASSERT(order.size() == refOrder.size());
132  std::list<std::string>::const_iterator orderIter = order.begin();
133  std::list<std::string>::const_iterator refOrderIter = refOrder.begin();
134  for ( ; orderIter != order.end(); ++orderIter,++refOrderIter) {
135  //std::cout << *orderIter << std::endl;
136  TEST_EQUALITY(*orderIter, *refOrderIter);
137  }
138 
139  // Test the reset.
140  loggingObs->resetLogCounters();
141  TEST_EQUALITY(
142  counters.find(loggingObs->nameObserveStartIntegrator_ )->second, 0);
143  TEST_EQUALITY(
144  counters.find(loggingObs->nameObserveStartTimeStep_ )->second, 0);
145  TEST_EQUALITY(
146  counters.find(loggingObs->nameObserveNextTimeStep_ )->second, 0);
147  TEST_EQUALITY(
148  counters.find(loggingObs->nameObserveBeforeTakeStep_ )->second, 0);
149  TEST_EQUALITY(
150  counters.find(loggingObs->nameObserveAfterTakeStep_ )->second, 0);
151  TEST_EQUALITY(
152  counters.find(loggingObs->nameObserveAfterCheckTimeStep_)->second, 0);
153  TEST_EQUALITY(
154  counters.find(loggingObs->nameObserveEndTimeStep_ )->second, 0);
155  TEST_EQUALITY(
156  counters.find(loggingObs->nameObserveEndIntegrator_ )->second, 0);
157  TEST_EQUALITY(order.size(), 0);
158 
159  Teuchos::TimeMonitor::summarize();
160 }
161 #endif // TEST_INTEGRATOROBSERVERLOGGING
162 
163 #ifdef TEST_INTEGRATOROBSERVERCOMPOSITE
164 TEUCHOS_UNIT_TEST( Observer, IntegratorObserverComposite)
165 {
166 
167  // Read params from .xml file
168  RCP<ParameterList> pList =
169  getParametersFromXmlFile("Tempus_Observer_SinCos.xml");
170 
171  // Setup the SinCosModel
172  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
173  RCP<SinCosModel<double> > model =
174  Teuchos::rcp(new SinCosModel<double> (scm_pl));
175 
176  // Setup the Integrator and reset initial time step
177  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
178  RCP<Tempus::IntegratorBasic<double> > integrator =
179  Tempus::integratorBasic<double>(pl, model);
180 
181  RCP<Tempus::IntegratorObserverLogging<double> > loggingObs =
183 
184  // creating another logging observer
185  RCP<Tempus::IntegratorObserverLogging<double> > loggingObs2 =
187 
188  RCP<Tempus::IntegratorObserverComposite<double> > compObs =
190 
191  compObs->addObserver(loggingObs);
192  compObs->addObserver(loggingObs2);
193 
194  compObs->observeStartIntegrator(*integrator);
195  compObs->observeStartTimeStep(*integrator);
196  compObs->observeBeforeTakeStep(*integrator);
197  compObs->observeAfterTakeStep(*integrator);
198  compObs->observeAfterCheckTimeStep(*integrator);
199  compObs->observeEndTimeStep(*integrator);
200 
201 
202  for (int i=0 ; i<10; ++i) {
203  compObs->observeStartTimeStep(*integrator);
204  compObs->observeBeforeTakeStep(*integrator);
205  compObs->observeAfterTakeStep(*integrator);
206  compObs->observeAfterCheckTimeStep(*integrator);
207  compObs->observeEndTimeStep(*integrator);
208  }
209  compObs->observeEndIntegrator(*integrator);
210 
211  const std::map<std::string,int>& counters = *(loggingObs->getCounters());
212 
213  TEST_EQUALITY(
214  counters.find(loggingObs->nameObserveStartIntegrator_ )->second, 1);
215  TEST_EQUALITY(
216  counters.find(loggingObs->nameObserveStartTimeStep_ )->second,11);
217  TEST_EQUALITY(
218  counters.find(loggingObs->nameObserveBeforeTakeStep_ )->second,11);
219  TEST_EQUALITY(
220  counters.find(loggingObs->nameObserveAfterTakeStep_ )->second,11);
221  TEST_EQUALITY(
222  counters.find(loggingObs->nameObserveAfterCheckTimeStep_ )->second,11);
223  TEST_EQUALITY(
224  counters.find(loggingObs->nameObserveEndTimeStep_ )->second,11);
225  TEST_EQUALITY(
226  counters.find(loggingObs->nameObserveEndIntegrator_ )->second, 1);
227  TEST_EQUALITY(
228  counters.find(loggingObs2->nameObserveStartIntegrator_ )->second, 1);
229  TEST_EQUALITY(
230  counters.find(loggingObs2->nameObserveStartTimeStep_ )->second,11);
231  TEST_EQUALITY(
232  counters.find(loggingObs2->nameObserveBeforeTakeStep_ )->second,11);
233  TEST_EQUALITY(
234  counters.find(loggingObs2->nameObserveAfterTakeStep_ )->second,11);
235  TEST_EQUALITY(
236  counters.find(loggingObs2->nameObserveAfterCheckTimeStep_)->second,11);
237  TEST_EQUALITY(
238  counters.find(loggingObs2->nameObserveEndTimeStep_ )->second,11);
239  TEST_EQUALITY(
240  counters.find(loggingObs2->nameObserveEndIntegrator_ )->second, 1);
241 }
242 #endif
243 
244 
245 // ************************************************************
246 // ************************************************************
247 #ifdef TEST_STEPPERRKOBSERVERLOGGING
248 TEUCHOS_UNIT_TEST(Observer, StepperRKObserverLogging)
249 {
250 
251  // number of stages for the RK method
252  const int erkStageCount = 4;
253 
254  // Read params from .xml file
255  RCP<ParameterList> pList =
256  getParametersFromXmlFile("Tempus_Observer_SinCos.xml");
257 
258  // Setup the SinCosModel
259  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
260  RCP<SinCosModel<double> > model =
261  Teuchos::rcp(new SinCosModel<double> (scm_pl));
262 
263  // Setup the Integrator and reset initial time step
264  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
265  pl->sublist("Demo Integrator").set("Stepper Name", "Demo RK Stepper");
266  RCP<Tempus::IntegratorBasic<double> > integrator =
267  Tempus::integratorBasic<double>(pl, model);
268 
269  RCP<Tempus::StepperRKObserverLogging<double> > loggingObs =
271 
272  const auto stepper_observer =
273  Teuchos::rcp_dynamic_cast<Tempus::StepperRKObserverComposite<double> >
274  (integrator->getStepper()->getObserver(), true);
275 
276  stepper_observer->clearObservers();
277  stepper_observer->addObserver(loggingObs);
278 
279  // Integrate to timeMax
280  bool integratorStatus = integrator->advanceTime();
281  TEST_ASSERT(integratorStatus)
282 
283  // Test if at 'Final Time'
284  double time = integrator->getTime();
285  const int numberTimeStep = integrator->getIndex();
286  double timeFinal = pl->sublist("Demo Integrator")
287  .sublist("Time Step Control").get<double>("Final Time");
288  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
289 
290  // Construct the reference counter and order for comparison.
291  std::map<std::string,int> refCounters;
292  std::list<std::string> refOrder;
293 
294  refCounters[loggingObs->nameObserveBeginTakeStep_ ] = numberTimeStep;
295  refCounters[loggingObs->nameObserveBeginStage_ ] = erkStageCount * numberTimeStep;
296  refCounters[loggingObs->nameObserveBeforeImplicitExplicitly_] = erkStageCount * numberTimeStep;
297  refCounters[loggingObs->nameObserveBeforeSolve_ ] = erkStageCount * numberTimeStep;
298  refCounters[loggingObs->nameObserveAfterSolve_ ] = erkStageCount * numberTimeStep;
299  refCounters[loggingObs->nameObserveBeforeExplicit_ ] = erkStageCount * numberTimeStep;
300  refCounters[loggingObs->nameObserveEndStage_ ] = erkStageCount * numberTimeStep;
301  refCounters[loggingObs->nameObserveEndTakeStep_ ] = numberTimeStep;
302 
303  for (int i=0 ; i< numberTimeStep; ++i) {
304  refOrder.push_back(loggingObs->nameObserveBeginTakeStep_ );
305  for (int j=0; j<erkStageCount; ++j){
306  refOrder.push_back(loggingObs->nameObserveBeginStage_ );
307  refOrder.push_back(loggingObs->nameObserveBeforeImplicitExplicitly_ );
308  refOrder.push_back(loggingObs->nameObserveBeforeSolve_ );
309  refOrder.push_back(loggingObs->nameObserveAfterSolve_ );
310  refOrder.push_back(loggingObs->nameObserveBeforeExplicit_ );
311  refOrder.push_back(loggingObs->nameObserveEndStage_ );
312  }
313  refOrder.push_back(loggingObs->nameObserveEndTakeStep_ );
314  }
315 
316  const std::map<std::string,int>& counters = *(loggingObs->getCounters());
317  const std::list<std::string>& order = *(loggingObs->getOrder());
318 
319  // Compare against reference.
320  TEST_EQUALITY(
321  counters.find(loggingObs->nameObserveBeginTakeStep_ )->second,
322  refCounters.find(loggingObs->nameObserveBeginTakeStep_ )->second);
323  TEST_EQUALITY(
324  counters.find(loggingObs->nameObserveBeginStage_ )->second,
325  refCounters.find(loggingObs->nameObserveBeginStage_ )->second);
326  TEST_EQUALITY(
327  counters.find(loggingObs->nameObserveBeforeImplicitExplicitly_ )->second,
328  refCounters.find(loggingObs->nameObserveBeforeImplicitExplicitly_ )->second);
329  TEST_EQUALITY(
330  counters.find(loggingObs->nameObserveBeforeSolve_ )->second,
331  refCounters.find(loggingObs->nameObserveBeforeSolve_ )->second);
332  TEST_EQUALITY(
333  counters.find(loggingObs->nameObserveAfterSolve_ )->second,
334  refCounters.find(loggingObs->nameObserveAfterSolve_ )->second);
335  TEST_EQUALITY(
336  counters.find(loggingObs->nameObserveBeforeExplicit_ )->second,
337  refCounters.find(loggingObs->nameObserveBeforeExplicit_ )->second);
338  TEST_EQUALITY(
339  counters.find(loggingObs->nameObserveEndStage_ )->second,
340  refCounters.find(loggingObs->nameObserveEndStage_ )->second);
341  TEST_EQUALITY(
342  counters.find(loggingObs->nameObserveEndTakeStep_ )->second,
343  refCounters.find(loggingObs->nameObserveEndTakeStep_ )->second);
344 
345  TEUCHOS_ASSERT(order.size() == refOrder.size());
346 
347  std::list<std::string>::const_iterator orderIter = order.begin();
348  std::list<std::string>::const_iterator refOrderIter = refOrder.begin();
349 
350  for ( ; orderIter != order.end(); ++orderIter,++refOrderIter) {
351  TEST_EQUALITY(*orderIter, *refOrderIter);
352  }
353 
354  // Test the reset.
355  loggingObs->resetLogCounters();
356  TEST_EQUALITY(
357  counters.find(loggingObs->nameObserveBeginTakeStep_)->second, 0);
358  TEST_EQUALITY(
359  counters.find(loggingObs->nameObserveBeginStage_)->second, 0);
360  TEST_EQUALITY(
361  counters.find(loggingObs->nameObserveBeforeImplicitExplicitly_)->second, 0);
362  TEST_EQUALITY(
363  counters.find(loggingObs->nameObserveBeforeSolve_)->second, 0);
364  TEST_EQUALITY(
365  counters.find(loggingObs->nameObserveAfterSolve_)->second, 0);
366  TEST_EQUALITY(
367  counters.find(loggingObs->nameObserveBeforeExplicit_)->second, 0);
368  TEST_EQUALITY(
369  counters.find(loggingObs->nameObserveEndStage_)->second, 0);
370  TEST_EQUALITY(
371  counters.find(loggingObs->nameObserveEndTakeStep_)->second, 0);
372  TEST_EQUALITY(order.size(), 0);
373 
374 
375  Teuchos::TimeMonitor::summarize();
376 }
377 #endif
378 
379 
380 } // namespace Tempus_Test
This observer logs calls to observer functions. This observer simply logs and counts the calls to eac...
Sine-Cosine model problem from Rythmos. This is a canonical Sine-Cosine differential equation with a...
This observer logs calls to observer functions. This observer simply logs and counts the calls to eac...
TEUCHOS_UNIT_TEST(BackwardEuler, SinCos_ASA)
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Solution state for integrators and steppers. SolutionState contains the metadata for solutions and th...