9 #include "Teuchos_UnitTestHarness.hpp"
10 #include "Teuchos_XMLParameterListHelpers.hpp"
11 #include "Teuchos_TimeMonitor.hpp"
13 #include "Thyra_VectorStdOps.hpp"
15 #include "Tempus_config.hpp"
16 #include "Tempus_IntegratorBasic.hpp"
17 #include "Tempus_StepperLeapfrog.hpp"
19 #include "../TestModels/HarmonicOscillatorModel.hpp"
20 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
23 #ifdef Tempus_ENABLE_MPI
24 #include "Epetra_MpiComm.h"
26 #include "Epetra_SerialComm.h"
34 namespace Tempus_Test {
38 using Teuchos::rcp_const_cast;
39 using Teuchos::ParameterList;
40 using Teuchos::sublist;
41 using Teuchos::getParametersFromXmlFile;
56 RCP<ParameterList> pList =
57 getParametersFromXmlFile(
"Tempus_Leapfrog_SinCos.xml");
58 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
61 RCP<ParameterList> hom_pl = sublist(pList,
"HarmonicOscillatorModel",
true);
66 stepper->setModel(model);
67 stepper->initialize();
71 ParameterList tscPL = pl->sublist(
"Default Integrator")
72 .sublist(
"Time Step Control");
73 timeStepControl->setStepType (tscPL.get<std::string>(
"Integrator Step Type"));
74 timeStepControl->setInitIndex(tscPL.get<
int> (
"Initial Time Index"));
75 timeStepControl->setInitTime (tscPL.get<
double>(
"Initial Time"));
76 timeStepControl->setFinalTime(tscPL.get<
double>(
"Final Time"));
77 timeStepControl->setInitTimeStep(dt);
78 timeStepControl->initialize();
81 Thyra::ModelEvaluatorBase::InArgs<double> inArgsIC =
82 stepper->getModel()->getNominalValues();
83 auto icX = rcp_const_cast<Thyra::VectorBase<double> > (inArgsIC.get_x());
84 auto icXDot = rcp_const_cast<Thyra::VectorBase<double> > (inArgsIC.get_x_dot());
85 auto icXDotDot = rcp_const_cast<Thyra::VectorBase<double> > (inArgsIC.get_x_dot_dot());
86 auto icState = Tempus::createSolutionStateX<double>(icX, icXDot, icXDotDot);
87 icState->setTime (timeStepControl->getInitTime());
88 icState->setIndex (timeStepControl->getInitIndex());
89 icState->setTimeStep(0.0);
90 icState->setOrder (stepper->getOrder());
95 solutionHistory->setName(
"Forward States");
97 solutionHistory->setStorageLimit(2);
98 solutionHistory->addState(icState);
101 RCP<Tempus::IntegratorBasic<double> > integrator =
102 Tempus::integratorBasic<double>();
103 integrator->setStepperWStepper(stepper);
104 integrator->setTimeStepControl(timeStepControl);
105 integrator->setSolutionHistory(solutionHistory);
107 integrator->initialize();
111 bool integratorStatus = integrator->advanceTime();
112 TEST_ASSERT(integratorStatus)
116 double time = integrator->getTime();
117 double timeFinal =pl->sublist(
"Default Integrator")
118 .sublist(
"Time Step Control").get<
double>(
"Final Time");
119 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
122 RCP<Thyra::VectorBase<double> > x = integrator->getX();
123 RCP<const Thyra::VectorBase<double> > x_exact =
124 model->getExactSolution(time).get_x();
127 RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
128 Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
131 std::cout <<
" Stepper = " << stepper->description() << std::endl;
132 std::cout <<
" =========================" << std::endl;
133 std::cout <<
" Exact solution : " << get_ele(*(x_exact), 0) << std::endl;
134 std::cout <<
" Computed solution: " << get_ele(*(x ), 0) << std::endl;
135 std::cout <<
" Difference : " << get_ele(*(xdiff ), 0) << std::endl;
136 std::cout <<
" =========================" << std::endl;
137 TEST_FLOATING_EQUALITY(get_ele(*(x), 0), 0.167158, 1.0e-4 );
145 RCP<Tempus::IntegratorBasic<double> > integrator;
146 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
147 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
148 std::vector<double> StepSize;
149 std::vector<double> xErrorNorm;
150 std::vector<double> xDotErrorNorm;
151 const int nTimeStepSizes = 9;
155 RCP<ParameterList> pList =
156 getParametersFromXmlFile(
"Tempus_Leapfrog_SinCos.xml");
159 RCP<ParameterList> hom_pl = sublist(pList,
"HarmonicOscillatorModel",
true);
164 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
167 double dt =pl->sublist(
"Default Integrator")
168 .sublist(
"Time Step Control").get<
double>(
"Initial Time Step");
171 for (
int n=0; n<nTimeStepSizes; n++) {
175 std::cout <<
"\n \n time step #" << n
176 <<
" (out of " << nTimeStepSizes-1 <<
"), dt = " << dt <<
"\n";
177 pl->sublist(
"Default Integrator")
178 .sublist(
"Time Step Control").set(
"Initial Time Step", dt);
179 integrator = Tempus::integratorBasic<double>(pl, model);
182 bool integratorStatus = integrator->advanceTime();
183 TEST_ASSERT(integratorStatus)
186 time = integrator->getTime();
187 double timeFinal =pl->sublist(
"Default Integrator")
188 .sublist(
"Time Step Control").get<
double>(
"Final Time");
189 TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
192 if (n == nTimeStepSizes-1) {
193 RCP<const SolutionHistory<double> > solutionHistory =
194 integrator->getSolutionHistory();
195 writeSolution(
"Tempus_Leapfrog_SinCos.dat", solutionHistory);
198 for (
int i=0; i<solutionHistory->getNumStates(); i++) {
199 double time_i = (*solutionHistory)[i]->getTime();
201 rcp_const_cast<Thyra::VectorBase<double> > (
202 model->getExactSolution(time_i).get_x()),
203 rcp_const_cast<Thyra::VectorBase<double> > (
204 model->getExactSolution(time_i).get_x_dot()));
205 state->setTime((*solutionHistory)[i]->getTime());
206 solnHistExact->addState(state);
208 writeSolution(
"Tempus_Leapfrog_SinCos-Ref.dat", solnHistExact);
214 StepSize.push_back(dt);
215 auto solution = Thyra::createMember(model->get_x_space());
216 Thyra::copy(*(integrator->getX()),solution.ptr());
217 solutions.push_back(solution);
218 auto solutionDot = Thyra::createMember(model->get_x_space());
219 Thyra::copy(*(integrator->getXdot()),solutionDot.ptr());
220 solutionsDot.push_back(solutionDot);
221 if (n == nTimeStepSizes-1) {
222 StepSize.push_back(0.0);
223 auto solutionExact = Thyra::createMember(model->get_x_space());
224 Thyra::copy(*(model->getExactSolution(time).get_x()),solutionExact.ptr());
225 solutions.push_back(solutionExact);
226 auto solutionDotExact = Thyra::createMember(model->get_x_space());
227 Thyra::copy(*(model->getExactSolution(time).get_x_dot()),
228 solutionDotExact.ptr());
229 solutionsDot.push_back(solutionDotExact);
235 double xDotSlope = 0.0;
236 RCP<Tempus::Stepper<double> > stepper = integrator->getStepper();
237 double order = stepper->getOrder();
240 solutions, xErrorNorm, xSlope,
241 solutionsDot, xDotErrorNorm, xDotSlope);
243 TEST_FLOATING_EQUALITY( xSlope, order, 0.02 );
244 TEST_FLOATING_EQUALITY( xErrorNorm[0], 0.0157928, 1.0e-4 );
245 TEST_FLOATING_EQUALITY( xDotSlope, 1.09387, 0.01 );
246 TEST_FLOATING_EQUALITY( xDotErrorNorm[0], 0.563002, 1.0e-4 );
248 Teuchos::TimeMonitor::summarize();
Teuchos::RCP< SolutionState< Scalar > > createSolutionStateX(const Teuchos::RCP< Thyra::VectorBase< Scalar > > &x, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &xdot=Teuchos::null, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &xdotdot=Teuchos::null)
Nonmember constructor from non-const solution vectors, x.
void writeSolution(const std::string filename, Teuchos::RCP< const Tempus::SolutionHistory< Scalar > > solutionHistory)
Consider the ODE: where is a constant, is a constant damping parameter, is a constant forcing par...
void writeOrderError(const std::string filename, Teuchos::RCP< Tempus::Stepper< Scalar > > stepper, std::vector< Scalar > &StepSize, std::vector< Teuchos::RCP< Thyra::VectorBase< Scalar >>> &solutions, std::vector< Scalar > &xErrorNorm, Scalar &xSlope, std::vector< Teuchos::RCP< Thyra::VectorBase< Scalar >>> &solutionsDot, std::vector< Scalar > &xDotErrorNorm, Scalar &xDotSlope, std::vector< Teuchos::RCP< Thyra::VectorBase< Scalar >>> &solutionsDotDot, std::vector< Scalar > &xDotDotErrorNorm, Scalar &xDotDotSlope)
TEUCHOS_UNIT_TEST(BackwardEuler, SinCos_ASA)
TimeStepControl manages the time step size. There several mechanicisms that effect the time step size...
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Keep a fix number of states.
Solution state for integrators and steppers. SolutionState contains the metadata for solutions and th...