Tempus  Version of the Day
Time Integration
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 "../TestModels/SinCosModel.hpp"
20 #include "../TestModels/VanDerPolModel.hpp"
21 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
22 
23 #include <vector>
24 
25 namespace Tempus_Test {
26 
27 using Teuchos::RCP;
28 using Teuchos::rcp;
29 using Teuchos::ParameterList;
30 using Teuchos::sublist;
31 using Teuchos::getParametersFromXmlFile;
32 
36 
37 
38 // ************************************************************
39 // ************************************************************
40 TEUCHOS_UNIT_TEST(Observer, IntegratorObserverLogging)
41 {
42  // Read params from .xml file
43  RCP<ParameterList> pList =
44  getParametersFromXmlFile("Tempus_Observer_SinCos.xml");
45 
46  // Setup the SinCosModel
47  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
48  RCP<SinCosModel<double> > model =
49  Teuchos::rcp(new SinCosModel<double> (scm_pl));
50 
51  // Setup the Integrator and reset initial time step
52  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
53  RCP<Tempus::IntegratorBasic<double> > integrator =
54  Tempus::integratorBasic<double>(pl, model);
55 
56  RCP<Tempus::IntegratorObserverLogging<double> > loggingObs =
58  integrator->setObserver(loggingObs);
59 
60  // Integrate to timeMax
61  bool integratorStatus = integrator->advanceTime();
62  TEST_ASSERT(integratorStatus)
63 
64  // Test if at 'Final Time'
65  double time = integrator->getTime();
66  double timeFinal = pl->sublist("Demo Integrator")
67  .sublist("Time Step Control").get<double>("Final Time");
68  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
69 
70  // Construct the reference counter and order for comparison.
71  std::map<std::string,int> refCounters;
72  std::list<std::string> refOrder;
73 
74  refCounters[loggingObs->nameObserveStartIntegrator_ ] = 1;
75  refCounters[loggingObs->nameObserveStartTimeStep_ ] = 10;
76  refCounters[loggingObs->nameObserveNextTimeStep_ ] = 10;
77  refCounters[loggingObs->nameObserveBeforeTakeStep_ ] = 10;
78  refCounters[loggingObs->nameObserveAfterTakeStep_ ] = 10;
79  refCounters[loggingObs->nameObserveAcceptedTimeStep_] = 10;
80  refCounters[loggingObs->nameObserveEndIntegrator_ ] = 1;
81 
82  refOrder.push_back(loggingObs->nameObserveStartIntegrator_ );
83  for (int i=0 ; i<10; ++i) {
84  refOrder.push_back(loggingObs->nameObserveStartTimeStep_ );
85  refOrder.push_back(loggingObs->nameObserveNextTimeStep_ );
86  refOrder.push_back(loggingObs->nameObserveBeforeTakeStep_ );
87  refOrder.push_back(loggingObs->nameObserveAfterTakeStep_ );
88  refOrder.push_back(loggingObs->nameObserveAcceptedTimeStep_);
89  }
90  refOrder.push_back(loggingObs->nameObserveEndIntegrator_ );
91 
92  const std::map<std::string,int>& counters = *(loggingObs->getCounters());
93  const std::list<std::string>& order = *(loggingObs->getOrder());
94 
95  // Compare against reference.
96  TEST_EQUALITY(
97  counters.find(loggingObs->nameObserveStartIntegrator_ )->second,
98  refCounters.find(loggingObs->nameObserveStartIntegrator_ )->second);
99  TEST_EQUALITY(
100  counters.find(loggingObs->nameObserveStartTimeStep_ )->second,
101  refCounters.find(loggingObs->nameObserveStartTimeStep_ )->second);
102  TEST_EQUALITY(
103  counters.find(loggingObs->nameObserveNextTimeStep_ )->second,
104  refCounters.find(loggingObs->nameObserveNextTimeStep_ )->second);
105  TEST_EQUALITY(
106  counters.find(loggingObs->nameObserveBeforeTakeStep_ )->second,
107  refCounters.find(loggingObs->nameObserveBeforeTakeStep_ )->second);
108  TEST_EQUALITY(
109  counters.find(loggingObs->nameObserveAfterTakeStep_ )->second,
110  refCounters.find(loggingObs->nameObserveAfterTakeStep_ )->second);
111  TEST_EQUALITY(
112  counters.find(loggingObs->nameObserveAcceptedTimeStep_)->second,
113  refCounters.find(loggingObs->nameObserveAcceptedTimeStep_)->second);
114  TEST_EQUALITY(
115  counters.find(loggingObs->nameObserveEndIntegrator_ )->second,
116  refCounters.find(loggingObs->nameObserveEndIntegrator_ )->second);
117 
118  TEUCHOS_ASSERT(order.size() == refOrder.size());
119  std::list<std::string>::const_iterator orderIter = order.begin();
120  std::list<std::string>::const_iterator refOrderIter = refOrder.begin();
121  for ( ; orderIter != order.end(); ++orderIter,++refOrderIter) {
122  //std::cout << *orderIter << std::endl;
123  TEST_EQUALITY(*orderIter, *refOrderIter);
124  }
125 
126  // Test the reset.
127  loggingObs->resetLogCounters();
128  TEST_EQUALITY(
129  counters.find(loggingObs->nameObserveStartIntegrator_ )->second, 0);
130  TEST_EQUALITY(
131  counters.find(loggingObs->nameObserveStartTimeStep_ )->second, 0);
132  TEST_EQUALITY(
133  counters.find(loggingObs->nameObserveNextTimeStep_ )->second, 0);
134  TEST_EQUALITY(
135  counters.find(loggingObs->nameObserveBeforeTakeStep_ )->second, 0);
136  TEST_EQUALITY(
137  counters.find(loggingObs->nameObserveAfterTakeStep_ )->second, 0);
138  TEST_EQUALITY(
139  counters.find(loggingObs->nameObserveAcceptedTimeStep_)->second, 0);
140  TEST_EQUALITY(
141  counters.find(loggingObs->nameObserveEndIntegrator_ )->second, 0);
142  TEST_EQUALITY(order.size(), 0);
143 
144  Teuchos::TimeMonitor::summarize();
145 }
146 
147 TEUCHOS_UNIT_TEST( Observer, IntegratorObserverComposite) {
148 
149  // Read params from .xml file
150  RCP<ParameterList> pList =
151  getParametersFromXmlFile("Tempus_Observer_SinCos.xml");
152 
153  // Setup the SinCosModel
154  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
155  RCP<SinCosModel<double> > model =
156  Teuchos::rcp(new SinCosModel<double> (scm_pl));
157 
158  // Setup the Integrator and reset initial time step
159  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
160  RCP<Tempus::IntegratorBasic<double> > integrator =
161  Tempus::integratorBasic<double>(pl, model);
162 
163  RCP<Tempus::IntegratorObserverLogging<double> > loggingObs =
165 
166  // creating another logging observer
167  RCP<Tempus::IntegratorObserverLogging<double> > loggingObs2 =
169 
170  RCP<Tempus::IntegratorObserverComposite<double> > compObs =
172 
173  compObs->addObserver(loggingObs);
174  compObs->addObserver(loggingObs2);
175 
176  compObs->observeStartIntegrator(*integrator);
177  compObs->observeStartTimeStep(*integrator);
178  compObs->observeBeforeTakeStep(*integrator);
179  compObs->observeAfterTakeStep(*integrator);
180  compObs->observeAcceptedTimeStep(*integrator);
181 
182 
183  for (int i=0 ; i<10; ++i) {
184  compObs->observeStartTimeStep(*integrator);
185  compObs->observeBeforeTakeStep(*integrator);
186  compObs->observeAfterTakeStep(*integrator);
187  compObs->observeAcceptedTimeStep(*integrator);
188  }
189  compObs->observeEndIntegrator(*integrator);
190 
191  const std::map<std::string,int>& counters = *(loggingObs->getCounters());
192 
193  TEST_EQUALITY(
194  counters.find(loggingObs->nameObserveStartIntegrator_ )->second, 1);
195  TEST_EQUALITY(
196  counters.find(loggingObs->nameObserveStartTimeStep_ )->second,11);
197  TEST_EQUALITY(
198  counters.find(loggingObs->nameObserveBeforeTakeStep_ )->second,11);
199  TEST_EQUALITY(
200  counters.find(loggingObs->nameObserveAfterTakeStep_ )->second,11);
201  TEST_EQUALITY(
202  counters.find(loggingObs->nameObserveAcceptedTimeStep_ )->second,11);
203  TEST_EQUALITY(
204  counters.find(loggingObs->nameObserveEndIntegrator_ )->second, 1);
205  TEST_EQUALITY(
206  counters.find(loggingObs2->nameObserveStartIntegrator_ )->second, 1);
207  TEST_EQUALITY(
208  counters.find(loggingObs2->nameObserveStartTimeStep_ )->second,11);
209  TEST_EQUALITY(
210  counters.find(loggingObs2->nameObserveBeforeTakeStep_ )->second,11);
211  TEST_EQUALITY(
212  counters.find(loggingObs2->nameObserveAfterTakeStep_ )->second,11);
213  TEST_EQUALITY(
214  counters.find(loggingObs2->nameObserveAcceptedTimeStep_)->second,11);
215  TEST_EQUALITY(
216  counters.find(loggingObs2->nameObserveEndIntegrator_ )->second, 1);
217 }
218 
219 
220 } // namespace Tempus_Test
Tempus::IntegratorObserverComposite
This observer.
Definition: Tempus_IntegratorObserverComposite_decl.hpp:21
Tempus_Test
Definition: Tempus_BackwardEuler_ASA.cpp:34
Tempus::IntegratorObserverLogging
This observer logs calls to observer functions. This observer simply logs and counts the calls to eac...
Definition: Tempus_IntegratorObserverLogging_decl.hpp:23
Tempus::IntegratorBasic
Basic time integrator.
Definition: Tempus_IntegratorBasic_decl.hpp:35
Tempus_Test::SinCosModel
Sine-Cosine model problem from Rythmos. This is a canonical Sine-Cosine differential equation.
Definition: SinCosModel_decl.hpp:91
Tempus::SolutionState
Solution state for integrators and steppers. SolutionState contains the metadata for solutions and th...
Definition: Tempus_SolutionState_decl.hpp:56
Tempus_Test::TEUCHOS_UNIT_TEST
TEUCHOS_UNIT_TEST(BackwardEuler, SinCos_ASA)
Definition: Tempus_BackwardEuler_ASA.cpp:47
Tempus::SolutionHistory
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Definition: Tempus_Integrator.hpp:25