Tempus Version of the Day
Time Integration
Loading...
Searching...
No Matches
Tempus_PhysicsStateTest.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
10#include "Teuchos_UnitTestHarness.hpp"
11#include "Teuchos_XMLParameterListHelpers.hpp"
12#include "Teuchos_TimeMonitor.hpp"
13
14#include "Thyra_VectorStdOps.hpp"
15
16#include "Tempus_IntegratorBasic.hpp"
19
20#include "../TestModels/SinCosModel.hpp"
21#include "../TestModels/VanDerPolModel.hpp"
22#include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
23
24#include <fstream>
25#include <vector>
26
27namespace Tempus_Test {
28
29using Teuchos::getParametersFromXmlFile;
30using Teuchos::ParameterList;
31using Teuchos::RCP;
32using Teuchos::sublist;
33
37
38// ************************************************************
39// ************************************************************
40TEUCHOS_UNIT_TEST(PhysicsState, SinCos)
41{
42 std::vector<double> StepSize;
43 std::vector<double> ErrorNorm;
44 const int nTimeStepSizes = 7;
45 double dt = 0.2;
46 double order = 0.0;
47 for (int n = 0; n < nTimeStepSizes; n++) {
48 // Read params from .xml file
49 RCP<ParameterList> pList =
50 getParametersFromXmlFile("Tempus_PhysicsState_SinCos.xml");
51
52 // std::ofstream ftmp("PL.txt");
53 // pList->print(ftmp);
54 // ftmp.close();
55
56 // Setup the SinCosModel
57 RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
58 // RCP<SinCosModel<double> > model = sineCosineModel(scm_pl);
59 RCP<SinCosModel<double> > model =
60 Teuchos::rcp(new SinCosModel<double>(scm_pl));
61
62 dt /= 2;
63
64 // Setup the Integrator and reset initial time step
65 RCP<ParameterList> pl = sublist(pList, "Tempus", true);
66 pl->sublist("Demo Integrator")
67 .sublist("Time Step Control")
68 .set("Initial Time Step", dt);
69 RCP<Tempus::IntegratorBasic<double> > integrator =
70 Tempus::createIntegratorBasic<double>(pl, model);
71
72 // Replace Tempus::StepperForwardEuler with
73 // Tempus_Test::StepperPhysicsStateTest
74 Teuchos::RCP<Tempus::Stepper<double> > physicsStepper =
75 Teuchos::rcp(new StepperPhysicsStateTest<double>(model));
76 integrator->setStepper(physicsStepper);
77 order = integrator->getStepper()->getOrder();
78
79 // Initial Conditions
80 // During the Integrator construction, the initial SolutionState
81 // is set by default to model->getNominalVales().get_x(). However,
82 // the application can set it also by integrator->initializeSolutionHistory.
83 RCP<Thyra::VectorBase<double> > x0 =
84 model->getNominalValues().get_x()->clone_v();
85 integrator->initializeSolutionHistory(0.0, x0);
86
87 // Replace Tempus::PhysicsState with
88 // Tempus_Test::PhysicsStateCounter
89 RCP<PhysicsStateCounter<double> > pSC =
90 Teuchos::rcp(new PhysicsStateCounter<double>("PhysicsStateTest", 0));
91 integrator->getSolutionHistory()->getCurrentState()->setPhysicsState(pSC);
92
93 // Integrate to timeMax
94 bool integratorStatus = integrator->advanceTime();
95 TEST_ASSERT(integratorStatus)
96
97 // Test PhysicsState
98 Teuchos::RCP<Tempus::PhysicsState<double> > pS =
99 integrator->getSolutionHistory()->getCurrentState()->getPhysicsState();
100 TEST_EQUALITY(pS->getName(), "PhysicsStateTest");
101 pSC = Teuchos::rcp_dynamic_cast<PhysicsStateCounter<double> >(pS);
102 // out << " Counter = " << pSC->getCounter() << std::endl;
103 TEST_EQUALITY(pSC->getCounter(), (int)(1.0 / dt));
104
105 // Test if at 'Final Time'
106 double time = integrator->getTime();
107 double timeFinal = pl->sublist("Demo Integrator")
108 .sublist("Time Step Control")
109 .get<double>("Final Time");
110 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
111
112 // Time-integrated solution and the exact solution
113 RCP<Thyra::VectorBase<double> > x = integrator->getX();
114 RCP<const Thyra::VectorBase<double> > x_exact =
115 model->getExactSolution(time).get_x();
116
117 // Plot sample solution and exact solution
118 if (n == 0) {
119 std::ofstream ftmp("Tempus_ForwardEuler_SinCos.dat");
120 RCP<const SolutionHistory<double> > solutionHistory =
121 integrator->getSolutionHistory();
122 RCP<const Thyra::VectorBase<double> > x_exact_plot;
123 for (int i = 0; i < solutionHistory->getNumStates(); i++) {
124 RCP<const SolutionState<double> > solutionState = (*solutionHistory)[i];
125 double time_i = solutionState->getTime();
126 RCP<const Thyra::VectorBase<double> > x_plot = solutionState->getX();
127 x_exact_plot = model->getExactSolution(time_i).get_x();
128 ftmp << time_i << " " << Thyra::get_ele(*(x_plot), 0) << " "
129 << Thyra::get_ele(*(x_plot), 1) << " "
130 << Thyra::get_ele(*(x_exact_plot), 0) << " "
131 << Thyra::get_ele(*(x_exact_plot), 1) << std::endl;
132 }
133 ftmp.close();
134 }
135
136 // Calculate the error
137 RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
138 Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
139 StepSize.push_back(dt);
140 const double L2norm = Thyra::norm_2(*xdiff);
141 ErrorNorm.push_back(L2norm);
142 }
143
144 // Check the order and intercept
145 double slope = computeLinearRegressionLogLog<double>(StepSize, ErrorNorm);
146 out << " Stepper = ForwardEuler" << std::endl;
147 out << " =========================" << std::endl;
148 out << " Expected order: " << order << std::endl;
149 out << " Observed order: " << slope << std::endl;
150 out << " =========================" << std::endl;
151 TEST_FLOATING_EQUALITY(slope, order, 0.01);
152 TEST_FLOATING_EQUALITY(ErrorNorm[0], 0.051123, 1.0e-4);
153
154 std::ofstream ftmp("Tempus_ForwardEuler_SinCos-Error.dat");
155 double error0 = 0.8 * ErrorNorm[0];
156 for (int n = 0; n < nTimeStepSizes; n++) {
157 ftmp << StepSize[n] << " " << ErrorNorm[n] << " "
158 << error0 * (pow(StepSize[n] / StepSize[0], order)) << std::endl;
159 }
160 ftmp.close();
161
162 Teuchos::TimeMonitor::summarize();
163}
164
165} // namespace Tempus_Test
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Solution state for integrators and steppers.
PhysicsStateCounter is a simple PhysicsState that counts steps.
Sine-Cosine model problem from Rythmos. This is a canonical Sine-Cosine differential equation.
This is a Forward Euler time stepper to test the PhysicsState.
TEUCHOS_UNIT_TEST(BackwardEuler, SinCos_ASA)