42 std::vector<double> StepSize;
43 std::vector<double> ErrorNorm;
44 const int nTimeStepSizes = 7;
47 for (
int n = 0; n < nTimeStepSizes; n++) {
49 RCP<ParameterList> pList =
50 getParametersFromXmlFile(
"Tempus_PhysicsState_SinCos.xml");
57 RCP<ParameterList> scm_pl = sublist(pList,
"SinCosModel",
true);
59 RCP<SinCosModel<double> > model =
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);
74 Teuchos::RCP<Tempus::Stepper<double> > physicsStepper =
76 integrator->setStepper(physicsStepper);
77 order = integrator->getStepper()->getOrder();
83 RCP<Thyra::VectorBase<double> > x0 =
84 model->getNominalValues().get_x()->clone_v();
85 integrator->initializeSolutionHistory(0.0, x0);
89 RCP<PhysicsStateCounter<double> > pSC =
91 integrator->getSolutionHistory()->getCurrentState()->setPhysicsState(pSC);
94 bool integratorStatus = integrator->advanceTime();
95 TEST_ASSERT(integratorStatus)
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);
103 TEST_EQUALITY(pSC->getCounter(), (
int)(1.0 / dt));
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);
113 RCP<Thyra::VectorBase<double> > x = integrator->getX();
114 RCP<const Thyra::VectorBase<double> > x_exact =
115 model->getExactSolution(time).get_x();
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;
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);
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);
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;
162 Teuchos::TimeMonitor::summarize();
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...