Tempus Version of the Day
Time Integration
Loading...
Searching...
No Matches
Tempus_Test_BackwardEuler_CDR.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#include "Teuchos_DefaultComm.hpp"
14
15#include "Tempus_config.hpp"
16#include "Tempus_IntegratorBasic.hpp"
17#include "Tempus_StepperBackwardEuler.hpp"
18
19#ifdef TEMPUS_ENABLE_EPETRA_STACK
20#include "../TestModels/CDR_Model.hpp"
21#ifdef Tempus_ENABLE_MPI
22#include "Epetra_MpiComm.h"
23#else
24#include "Epetra_SerialComm.h"
25#endif
26#endif
27#ifdef TEMPUS_ENABLE_TPETRA_STACK
28#include "../TestModels/CDR_Model_Tpetra.hpp"
29#include "Tpetra_Core.hpp"
30#endif
31#include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
32
33#include "Stratimikos_DefaultLinearSolverBuilder.hpp"
34
35#include <vector>
36#include <fstream>
37#include <sstream>
38#include <limits>
39
40namespace Tempus_Test {
41
42using Teuchos::getParametersFromXmlFile;
43using Teuchos::ParameterList;
44using Teuchos::RCP;
45using Teuchos::rcp;
46using Teuchos::rcp_const_cast;
47using Teuchos::sublist;
48
52
53// ************************************************************
54// ************************************************************
55template <typename SC, typename Model, typename Comm>
56void CDR_Test(const Comm& comm, const int commSize, Teuchos::FancyOStream& out,
57 bool& success)
58{
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;
66 double dt = 0.2;
67 for (int n = 0; n < nTimeStepSizes; n++) {
68 // Read params from .xml file
69 RCP<ParameterList> pList =
70 getParametersFromXmlFile("Tempus_BackwardEuler_CDR.xml");
71
72 // Create CDR Model
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)");
79
80 auto model = rcp(new Model(comm, num_elements, left_end, right_end,
81 a_convection, k_source));
82
83 // Set the factory
84 ::Stratimikos::DefaultLinearSolverBuilder builder;
85
86 auto p = rcp(new ParameterList);
87 p->set("Linear Solver Type", "Belos");
88 p->set("Preconditioner Type", "None");
89 builder.setParameterList(p);
90
91 auto lowsFactory = builder.createLinearSolveStrategy("");
92
93 model->set_W_factory(lowsFactory);
94
95 // Set the step size
96 dt /= 2;
97
98 // Setup the Integrator and reset initial time step
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);
104
105 // Integrate to timeMax
106 bool integratorStatus = integrator->advanceTime();
107 TEST_ASSERT(integratorStatus)
108
109 // Test if at 'Final Time'
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);
116
117 // Store off the final solution and step size
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);
125
126 // Output finest temporal solution for plotting
127 // This only works for ONE MPI process
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";
132 const double dx =
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
142 << ", F=BLOCK\n";
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 << " ";
146 }
147 ftmp << std::endl;
148 for (int j = 0; j < num_elements + 1; j++)
149 ftmp << get_ele(*x, j) << " ";
150 ftmp << std::endl;
151 }
152 ftmp.close();
153 }
154 }
155
156 // Check the order and intercept
157 double xSlope = 0.0;
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,
162 xDotSlope, out);
163
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);
168 // At small dt, slopes should be equal to order.
169 // double order = stepper->getOrder();
170 // TEST_FLOATING_EQUALITY( xSlope, order, 0.01 );
171 // TEST_FLOATING_EQUALITY( xDotSlope, order, 0.01 );
172
173 // Write fine mesh solution at final time
174 // This only works for ONE MPI process
175 if (commSize == 1) {
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");
182
183 const Thyra::VectorBase<double>& x = *(solutions[solutions.size() - 1]);
184
185 std::ofstream ftmp("Tempus_BackwardEuler_CDR-Solution.dat");
186 for (int n = 0; n < num_elements + 1; n++) {
187 const double dx =
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;
191 }
192 ftmp.close();
193 }
194
195 Teuchos::TimeMonitor::summarize();
196}
197
198#ifdef TEMPUS_ENABLE_EPETRA_STACK
199// ************************************************************
200// ************************************************************
201TEUCHOS_UNIT_TEST(BackwardEuler, CDR)
202{
203 // Create a communicator for Epetra objects
204 RCP<Epetra_Comm> comm;
205#ifdef Tempus_ENABLE_MPI
206 comm = rcp(new Epetra_MpiComm(MPI_COMM_WORLD));
207#else
208 comm = rcp(new Epetra_SerialComm);
209#endif
210
211 CDR_Test<double, Tempus_Test::CDR_Model<double>>(comm, comm->NumProc(), out,
212 success);
213}
214#endif
215
216#ifdef TEMPUS_ENABLE_TPETRA_STACK
217// ************************************************************
218// ************************************************************
219TEUCHOS_UNIT_TEST(BackwardEuler, CDR_Tpetra)
220{
221 // Get default Tpetra template types
222 using SC = Tpetra::Vector<>::scalar_type;
223 using LO = Tpetra::Vector<>::local_ordinal_type;
224 using GO = Tpetra::Vector<>::global_ordinal_type;
225 using Node = Tpetra::Vector<>::node_type;
226
227 auto comm = Tpetra::getDefaultComm();
228
229 CDR_Test<SC, Tempus_Test::CDR_Model_Tpetra<SC, LO, GO, Node>>(
230 comm, comm->getSize(), out, success);
231}
232#endif
233
234} // namespace Tempus_Test
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Solution state for integrators and steppers.
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)
TEUCHOS_UNIT_TEST(BackwardEuler, SinCos_ASA)
void CDR_Test(const Comm &comm, const int commSize, Teuchos::FancyOStream &out, bool &success)