56void CDR_Test(
const Comm& comm,
const int commSize, Teuchos::FancyOStream& out,
59 RCP<Tempus::IntegratorBasic<double>> integrator;
60 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
61 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
62 std::vector<double> StepSize;
63 std::vector<double> xErrorNorm;
64 std::vector<double> xDotErrorNorm;
65 const int nTimeStepSizes = 5;
67 for (
int n = 0; n < nTimeStepSizes; n++) {
69 RCP<ParameterList> pList =
70 getParametersFromXmlFile(
"Tempus_BackwardEuler_CDR.xml");
73 RCP<ParameterList> model_pl = sublist(pList,
"CDR Model",
true);
74 const auto num_elements = model_pl->get<
int>(
"num elements");
75 const auto left_end = model_pl->get<SC>(
"left end");
76 const auto right_end = model_pl->get<SC>(
"right end");
77 const auto a_convection = model_pl->get<SC>(
"a (convection)");
78 const auto k_source = model_pl->get<SC>(
"k (source)");
80 auto model = rcp(
new Model(comm, num_elements, left_end, right_end,
81 a_convection, k_source));
84 ::Stratimikos::DefaultLinearSolverBuilder builder;
86 auto p = rcp(
new ParameterList);
87 p->set(
"Linear Solver Type",
"Belos");
88 p->set(
"Preconditioner Type",
"None");
89 builder.setParameterList(p);
91 auto lowsFactory = builder.createLinearSolveStrategy(
"");
93 model->set_W_factory(lowsFactory);
99 RCP<ParameterList> pl = sublist(pList,
"Tempus",
true);
100 pl->sublist(
"Demo Integrator")
101 .sublist(
"Time Step Control")
102 .set(
"Initial Time Step", dt);
103 integrator = Tempus::createIntegratorBasic<double>(pl, model);
106 bool integratorStatus = integrator->advanceTime();
107 TEST_ASSERT(integratorStatus)
110 double time = integrator->getTime();
111 double timeFinal = pl->sublist(
"Demo Integrator")
112 .sublist(
"Time Step Control")
113 .get<
double>(
"Final Time");
114 double tol = 100.0 * std::numeric_limits<double>::epsilon();
115 TEST_FLOATING_EQUALITY(time, timeFinal, tol);
118 StepSize.push_back(dt);
119 auto solution = Thyra::createMember(model->get_x_space());
120 Thyra::copy(*(integrator->getX()), solution.ptr());
121 solutions.push_back(solution);
122 auto solutionDot = Thyra::createMember(model->get_x_space());
123 Thyra::copy(*(integrator->getXDot()), solutionDot.ptr());
124 solutionsDot.push_back(solutionDot);
128 if ((n == nTimeStepSizes - 1) && (commSize == 1)) {
129 std::ofstream ftmp(
"Tempus_BackwardEuler_CDR.dat");
130 ftmp <<
"TITLE=\"Backward Euler Solution to CDR\"\n"
131 <<
"VARIABLES=\"z\",\"T\"\n";
133 std::fabs(left_end - right_end) /
static_cast<double>(num_elements);
134 RCP<const SolutionHistory<double>> solutionHistory =
135 integrator->getSolutionHistory();
136 int nStates = solutionHistory->getNumStates();
137 for (
int i = 0; i < nStates; i++) {
138 RCP<const SolutionState<double>> solutionState = (*solutionHistory)[i];
139 RCP<const Thyra::VectorBase<double>> x = solutionState->getX();
140 double ttime = solutionState->getTime();
141 ftmp <<
"ZONE T=\"Time=" << ttime <<
"\", I=" << num_elements + 1
143 for (
int j = 0; j < num_elements + 1; j++) {
144 const double x_coord = left_end +
static_cast<double>(j) * dx;
145 ftmp << x_coord <<
" ";
148 for (
int j = 0; j < num_elements + 1; j++)
149 ftmp << get_ele(*x, j) <<
" ";
158 double xDotSlope = 0.0;
159 RCP<Tempus::Stepper<double>> stepper = integrator->getStepper();
160 writeOrderError(
"Tempus_BackwardEuler_CDR-Error.dat", stepper, StepSize,
161 solutions, xErrorNorm, xSlope, solutionsDot, xDotErrorNorm,
164 TEST_FLOATING_EQUALITY(xSlope, 1.32213, 0.01);
165 TEST_FLOATING_EQUALITY(xErrorNorm[0], 0.116919, 1.0e-4);
166 TEST_FLOATING_EQUALITY(xDotSlope, 1.32052, 0.01);
167 TEST_FLOATING_EQUALITY(xDotErrorNorm[0], 0.449888, 1.0e-4);
176 RCP<ParameterList> pList =
177 getParametersFromXmlFile(
"Tempus_BackwardEuler_CDR.xml");
178 RCP<ParameterList> model_pl = sublist(pList,
"CDR Model",
true);
179 const int num_elements = model_pl->get<
int>(
"num elements");
180 const double left_end = model_pl->get<
double>(
"left end");
181 const double right_end = model_pl->get<
double>(
"right end");
185 std::ofstream ftmp(
"Tempus_BackwardEuler_CDR-Solution.dat");
186 for (
int n = 0; n < num_elements + 1; n++) {
188 std::fabs(left_end - right_end) /
static_cast<double>(num_elements);
189 const double x_coord = left_end +
static_cast<double>(n) * dx;
190 ftmp << x_coord <<
" " << Thyra::get_ele(x, n) << std::endl;
195 Teuchos::TimeMonitor::summarize();
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
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::FancyOStream &out)