13 #include "Teuchos_DefaultComm.hpp"
15 #include "Tempus_config.hpp"
16 #include "Tempus_IntegratorBasic.hpp"
17 #include "Tempus_StepperBDF2.hpp"
19 #ifdef TEMPUS_ENABLE_EPETRA_STACK
20 #include "../TestModels/CDR_Model.hpp"
21 #ifdef Tempus_ENABLE_MPI
22 #include "Epetra_MpiComm.h"
24 #include "Epetra_SerialComm.h"
27 #ifdef TEMPUS_ENABLE_TPETRA_STACK
28 #include "../TestModels/CDR_Model_Tpetra.hpp"
29 #include "Tpetra_Core.hpp"
31 #include "../TestModels/VanDerPolModel.hpp"
32 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
34 #include "Stratimikos_DefaultLinearSolverBuilder.hpp"
41 namespace Tempus_Test {
43 using Teuchos::getParametersFromXmlFile;
47 using Teuchos::rcp_const_cast;
48 using Teuchos::sublist;
56 template <
typename SC,
typename Model,
typename Comm>
60 RCP<Tempus::IntegratorBasic<double>> integrator;
61 std::vector<RCP<Thyra::VectorBase<double>>> solutions;
62 std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
63 std::vector<double> StepSize;
66 auto pList = getParametersFromXmlFile(
"Tempus_BDF2_CDR.xml");
70 auto pl = sublist(pList,
"Tempus",
true);
71 auto dt = pl->sublist(
"Demo Integrator")
72 .sublist(
"Time Step Control")
73 .get<
double>(
"Initial Time Step");
75 auto model_pl = sublist(pList,
"CDR Model",
true);
77 const auto nTimeStepSizes =
78 model_pl->get<
int>(
"Number of Time Step Sizes", 5);
80 for (
int n = 0; n < nTimeStepSizes; n++) {
82 const int num_elements = model_pl->get<
int>(
"num elements");
83 const auto left_end = model_pl->get<SC>(
"left end");
84 const auto right_end = model_pl->get<SC>(
"right end");
85 const auto a_convection = model_pl->get<SC>(
"a (convection)");
86 const auto k_source = model_pl->get<SC>(
"k (source)");
88 auto model =
rcp(
new Model(comm, num_elements, left_end, right_end,
89 a_convection, k_source));
92 ::Stratimikos::DefaultLinearSolverBuilder builder;
94 auto p =
rcp(
new ParameterList);
95 p->set(
"Linear Solver Type",
"Belos");
96 p->set(
"Preconditioner Type",
"None");
97 builder.setParameterList(p);
99 auto lowsFactory = builder.createLinearSolveStrategy(
"");
101 model->set_W_factory(lowsFactory);
107 pl->sublist(
"Demo Integrator")
108 .sublist(
"Time Step Control")
109 .set(
"Initial Time Step", dt);
110 integrator = Tempus::createIntegratorBasic<double>(pl, model);
113 bool integratorStatus = integrator->advanceTime();
117 double time = integrator->getTime();
118 double timeFinal = pl->sublist("Demo Integrator")
119 .sublist("Time Step Control")
120 .get<
double>("Final Time");
121 double tol = 100.0 * std::numeric_limits<
double>::epsilon();
122 TEST_FLOATING_EQUALITY(time, timeFinal, tol);
125 StepSize.push_back(dt);
126 auto solution = Thyra::createMember(model->get_x_space());
127 Thyra::copy(*(integrator->getX()), solution.ptr());
128 solutions.push_back(solution);
129 auto solutionDot = Thyra::createMember(model->get_x_space());
130 Thyra::copy(*(integrator->getXDot()), solutionDot.ptr());
131 solutionsDot.push_back(solutionDot);
135 if ((n == nTimeStepSizes - 1) && (commSize == 1)) {
136 std::ofstream ftmp(
"Tempus_BDF2_CDR.dat");
137 ftmp <<
"TITLE=\"BDF2 Solution to CDR\"\n"
138 <<
"VARIABLES=\"z\",\"T\"\n";
140 std::fabs(left_end - right_end) /
static_cast<double>(num_elements);
141 auto solutionHistory = integrator->getSolutionHistory();
142 int nStates = solutionHistory->getNumStates();
143 for (
int i = 0; i < nStates; i++) {
144 auto solutionState = (*solutionHistory)[i];
145 auto x = solutionState->getX();
146 auto ttime = solutionState->getTime();
147 ftmp <<
"ZONE T=\"Time=" << ttime <<
"\", I=" << num_elements + 1
149 for (
int j = 0; j < num_elements + 1; j++) {
150 const auto x_coord = left_end +
static_cast<double>(j) * dx;
151 ftmp << x_coord <<
" ";
154 for (
int j = 0; j < num_elements + 1; j++)
155 ftmp << get_ele(*x, j) <<
" ";
163 if (nTimeStepSizes > 2) {
165 double xDotSlope = 0.0;
166 std::vector<double> xErrorNorm;
167 std::vector<double> xDotErrorNorm;
168 auto stepper = integrator->getStepper();
169 auto order = stepper->getOrder();
170 writeOrderError(
"Tempus_BDF2_CDR-Error.dat", stepper, StepSize, solutions,
171 xErrorNorm, xSlope, solutionsDot, xDotErrorNorm, xDotSlope,
185 auto pListCDR = getParametersFromXmlFile(
"Tempus_BDF2_CDR.xml");
186 auto model_pl_CDR = sublist(pListCDR,
"CDR Model",
true);
187 const auto num_elements = model_pl_CDR->get<
int>(
"num elements");
188 const auto left_end = model_pl_CDR->get<
double>(
"left end");
189 const auto right_end = model_pl_CDR->get<
double>(
"right end");
191 const auto& x = *(solutions[solutions.size() - 1]);
193 std::ofstream ftmp(
"Tempus_BDF2_CDR-Solution.dat");
194 for (
int n = 0; n < num_elements + 1; n++) {
196 std::fabs(left_end - right_end) /
static_cast<double>(num_elements);
197 const auto x_coord = left_end +
static_cast<double>(n) * dx;
198 ftmp << x_coord <<
" " << Thyra::get_ele(x, n) << std::endl;
206 #ifdef TEMPUS_ENABLE_EPETRA_STACK
212 RCP<Epetra_Comm> comm;
213 #ifdef Tempus_ENABLE_MPI
214 comm =
rcp(
new Epetra_MpiComm(MPI_COMM_WORLD));
216 comm =
rcp(
new Epetra_SerialComm);
219 CDR_Test<double, Tempus_Test::CDR_Model<double>>(comm, comm->NumProc(), out,
224 #ifdef TEMPUS_ENABLE_TPETRA_STACK
230 using SC = Tpetra::Vector<>::scalar_type;
231 using LO = Tpetra::Vector<>::local_ordinal_type;
232 using GO = Tpetra::Vector<>::global_ordinal_type;
233 using Node = Tpetra::Vector<>::node_type;
235 auto comm = Tpetra::getDefaultComm();
237 CDR_Test<SC, Tempus_Test::CDR_Model_Tpetra<SC, LO, GO, Node>>(
238 comm, comm->getSize(), out, success);
#define TEST_COMPARE(v1, comp, v2)
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)
#define TEST_FLOATING_EQUALITY(v1, v2, tol)
TEUCHOS_UNIT_TEST(BackwardEuler, SinCos_ASA)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
static void summarize(Ptr< const Comm< int > > comm, std::ostream &out=std::cout, const bool alwaysWriteLocal=false, const bool writeGlobalStats=true, const bool writeZeroTimers=true, const ECounterSetOp setOp=Intersection, const std::string &filter="", const bool ignoreZeroTimers=false)
void CDR_Test(const Comm &comm, const int commSize, Teuchos::FancyOStream &out, bool &success)
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Solution state for integrators and steppers.