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_StepperExplicitRK.hpp"
20 
21 #include "Tempus_StepperRKObserverLogging.hpp"
22 #include "Tempus_StepperRKObserverComposite.hpp"
23 
24 #include "../TestModels/SinCosModel.hpp"
25 #include "../TestModels/VanDerPolModel.hpp"
26 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
27 
28 #include <vector>
29 
30 namespace Tempus_Test {
31 
32 using Teuchos::RCP;
33 using Teuchos::rcp;
34 using Teuchos::ParameterList;
35 using Teuchos::sublist;
36 using Teuchos::getParametersFromXmlFile;
37 
41 
42 
43 
44 // ************************************************************
45 // ************************************************************
46 TEUCHOS_UNIT_TEST(Observer, IntegratorObserverLogging)
47 {
48  // Read params from .xml file
49  RCP<ParameterList> pList =
50  getParametersFromXmlFile("Tempus_Observer_SinCos.xml");
51 
52  // Setup the SinCosModel
53  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
54  RCP<SinCosModel<double> > model =
55  Teuchos::rcp(new SinCosModel<double> (scm_pl));
56 
57  // Setup the Integrator and reset initial time step
58  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
59  RCP<Tempus::IntegratorBasic<double> > integrator =
60  Tempus::integratorBasic<double>(pl, model);
61 
62  RCP<Tempus::IntegratorObserverLogging<double> > loggingObs =
64  integrator->setObserver(loggingObs);
65 
66  // Integrate to timeMax
67  bool integratorStatus = integrator->advanceTime();
68  TEST_ASSERT(integratorStatus)
69 
70  // Test if at 'Final Time'
71  double time = integrator->getTime();
72  double timeFinal = pl->sublist("Demo Integrator")
73  .sublist("Time Step Control").get<double>("Final Time");
74  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
75 
76  // Construct the reference counter and order for comparison.
77  std::map<std::string,int> refCounters;
78  std::list<std::string> refOrder;
79 
80  refCounters[loggingObs->nameObserveStartIntegrator_ ] = 1;
81  refCounters[loggingObs->nameObserveStartTimeStep_ ] = 10;
82  refCounters[loggingObs->nameObserveNextTimeStep_ ] = 10;
83  refCounters[loggingObs->nameObserveBeforeTakeStep_ ] = 10;
84  refCounters[loggingObs->nameObserveAfterTakeStep_ ] = 10;
85  refCounters[loggingObs->nameObserveAfterCheckTimeStep_] = 10;
86  refCounters[loggingObs->nameObserveEndTimeStep_ ] = 10;
87  refCounters[loggingObs->nameObserveEndIntegrator_ ] = 1;
88 
89  refOrder.push_back(loggingObs->nameObserveStartIntegrator_ );
90  for (int i=0 ; i<10; ++i) {
91  refOrder.push_back(loggingObs->nameObserveStartTimeStep_ );
92  refOrder.push_back(loggingObs->nameObserveNextTimeStep_ );
93  refOrder.push_back(loggingObs->nameObserveBeforeTakeStep_ );
94  refOrder.push_back(loggingObs->nameObserveAfterTakeStep_ );
95  refOrder.push_back(loggingObs->nameObserveAfterCheckTimeStep_);
96  refOrder.push_back(loggingObs->nameObserveEndTimeStep_ );
97  }
98  refOrder.push_back(loggingObs->nameObserveEndIntegrator_ );
99 
100  const std::map<std::string,int>& counters = *(loggingObs->getCounters());
101  const std::list<std::string>& order = *(loggingObs->getOrder());
102 
103  // Compare against reference.
104  TEST_EQUALITY(
105  counters.find(loggingObs->nameObserveStartIntegrator_ )->second,
106  refCounters.find(loggingObs->nameObserveStartIntegrator_ )->second);
107  TEST_EQUALITY(
108  counters.find(loggingObs->nameObserveStartTimeStep_ )->second,
109  refCounters.find(loggingObs->nameObserveStartTimeStep_ )->second);
110  TEST_EQUALITY(
111  counters.find(loggingObs->nameObserveNextTimeStep_ )->second,
112  refCounters.find(loggingObs->nameObserveNextTimeStep_ )->second);
113  TEST_EQUALITY(
114  counters.find(loggingObs->nameObserveBeforeTakeStep_ )->second,
115  refCounters.find(loggingObs->nameObserveBeforeTakeStep_ )->second);
116  TEST_EQUALITY(
117  counters.find(loggingObs->nameObserveAfterTakeStep_ )->second,
118  refCounters.find(loggingObs->nameObserveAfterTakeStep_ )->second);
119  TEST_EQUALITY(
120  counters.find(loggingObs->nameObserveAfterCheckTimeStep_)->second,
121  refCounters.find(loggingObs->nameObserveAfterCheckTimeStep_)->second);
122  TEST_EQUALITY(
123  counters.find(loggingObs->nameObserveEndTimeStep_ )->second,
124  refCounters.find(loggingObs->nameObserveEndTimeStep_ )->second);
125  TEST_EQUALITY(
126  counters.find(loggingObs->nameObserveEndIntegrator_ )->second,
127  refCounters.find(loggingObs->nameObserveEndIntegrator_ )->second);
128 
129  TEUCHOS_ASSERT(order.size() == refOrder.size());
130  std::list<std::string>::const_iterator orderIter = order.begin();
131  std::list<std::string>::const_iterator refOrderIter = refOrder.begin();
132  for ( ; orderIter != order.end(); ++orderIter,++refOrderIter) {
133  //std::cout << *orderIter << std::endl;
134  TEST_EQUALITY(*orderIter, *refOrderIter);
135  }
136 
137  // Test the reset.
138  loggingObs->resetLogCounters();
139  TEST_EQUALITY(
140  counters.find(loggingObs->nameObserveStartIntegrator_ )->second, 0);
141  TEST_EQUALITY(
142  counters.find(loggingObs->nameObserveStartTimeStep_ )->second, 0);
143  TEST_EQUALITY(
144  counters.find(loggingObs->nameObserveNextTimeStep_ )->second, 0);
145  TEST_EQUALITY(
146  counters.find(loggingObs->nameObserveBeforeTakeStep_ )->second, 0);
147  TEST_EQUALITY(
148  counters.find(loggingObs->nameObserveAfterTakeStep_ )->second, 0);
149  TEST_EQUALITY(
150  counters.find(loggingObs->nameObserveAfterCheckTimeStep_)->second, 0);
151  TEST_EQUALITY(
152  counters.find(loggingObs->nameObserveEndTimeStep_ )->second, 0);
153  TEST_EQUALITY(
154  counters.find(loggingObs->nameObserveEndIntegrator_ )->second, 0);
155  TEST_EQUALITY(order.size(), 0);
156 
157  Teuchos::TimeMonitor::summarize();
158 }
159 
160 TEUCHOS_UNIT_TEST( Observer, IntegratorObserverComposite)
161 {
162 
163  // Read params from .xml file
164  RCP<ParameterList> pList =
165  getParametersFromXmlFile("Tempus_Observer_SinCos.xml");
166 
167  // Setup the SinCosModel
168  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
169  RCP<SinCosModel<double> > model =
170  Teuchos::rcp(new SinCosModel<double> (scm_pl));
171 
172  // Setup the Integrator and reset initial time step
173  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
174  RCP<Tempus::IntegratorBasic<double> > integrator =
175  Tempus::integratorBasic<double>(pl, model);
176 
177  RCP<Tempus::IntegratorObserverLogging<double> > loggingObs =
179 
180  // creating another logging observer
181  RCP<Tempus::IntegratorObserverLogging<double> > loggingObs2 =
183 
184  RCP<Tempus::IntegratorObserverComposite<double> > compObs =
186 
187  compObs->addObserver(loggingObs);
188  compObs->addObserver(loggingObs2);
189 
190  compObs->observeStartIntegrator(*integrator);
191  compObs->observeStartTimeStep(*integrator);
192  compObs->observeBeforeTakeStep(*integrator);
193  compObs->observeAfterTakeStep(*integrator);
194  compObs->observeAfterCheckTimeStep(*integrator);
195  compObs->observeEndTimeStep(*integrator);
196 
197 
198  for (int i=0 ; i<10; ++i) {
199  compObs->observeStartTimeStep(*integrator);
200  compObs->observeBeforeTakeStep(*integrator);
201  compObs->observeAfterTakeStep(*integrator);
202  compObs->observeAfterCheckTimeStep(*integrator);
203  compObs->observeEndTimeStep(*integrator);
204  }
205  compObs->observeEndIntegrator(*integrator);
206 
207  const std::map<std::string,int>& counters = *(loggingObs->getCounters());
208 
209  TEST_EQUALITY(
210  counters.find(loggingObs->nameObserveStartIntegrator_ )->second, 1);
211  TEST_EQUALITY(
212  counters.find(loggingObs->nameObserveStartTimeStep_ )->second,11);
213  TEST_EQUALITY(
214  counters.find(loggingObs->nameObserveBeforeTakeStep_ )->second,11);
215  TEST_EQUALITY(
216  counters.find(loggingObs->nameObserveAfterTakeStep_ )->second,11);
217  TEST_EQUALITY(
218  counters.find(loggingObs->nameObserveAfterCheckTimeStep_ )->second,11);
219  TEST_EQUALITY(
220  counters.find(loggingObs->nameObserveEndTimeStep_ )->second,11);
221  TEST_EQUALITY(
222  counters.find(loggingObs->nameObserveEndIntegrator_ )->second, 1);
223  TEST_EQUALITY(
224  counters.find(loggingObs2->nameObserveStartIntegrator_ )->second, 1);
225  TEST_EQUALITY(
226  counters.find(loggingObs2->nameObserveStartTimeStep_ )->second,11);
227  TEST_EQUALITY(
228  counters.find(loggingObs2->nameObserveBeforeTakeStep_ )->second,11);
229  TEST_EQUALITY(
230  counters.find(loggingObs2->nameObserveAfterTakeStep_ )->second,11);
231  TEST_EQUALITY(
232  counters.find(loggingObs2->nameObserveAfterCheckTimeStep_)->second,11);
233  TEST_EQUALITY(
234  counters.find(loggingObs2->nameObserveEndTimeStep_ )->second,11);
235  TEST_EQUALITY(
236  counters.find(loggingObs2->nameObserveEndIntegrator_ )->second, 1);
237 }
238 
239 
240 } // namespace Tempus_Test
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...