13 #include "Teuchos_DefaultComm.hpp"
15 #include "Tempus_config.hpp"
16 #include "Tempus_IntegratorBasic.hpp"
17 #include "Tempus_StepperBackwardEuler.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 "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
33 #include "Stratimikos_DefaultLinearSolverBuilder.hpp"
40 namespace Tempus_Test {
42 using Teuchos::getParametersFromXmlFile;
46 using Teuchos::rcp_const_cast;
47 using Teuchos::sublist;
55 template <
typename SC,
typename Model,
typename Comm>
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++) {
70 getParametersFromXmlFile(
"Tempus_BackwardEuler_CDR.xml");
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;
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);
102 .
set(
"Initial Time Step", dt);
103 integrator = Tempus::createIntegratorBasic<double>(pl, model);
106 bool integratorStatus = integrator->advanceTime();
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();
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);
135 integrator->getSolutionHistory();
136 int nStates = solutionHistory->getNumStates();
137 for (
int i = 0; i < nStates; i++) {
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;
160 writeOrderError(
"Tempus_BackwardEuler_CDR-Error.dat", stepper, StepSize,
161 solutions, xErrorNorm, xSlope, solutionsDot, xDotErrorNorm,
177 getParametersFromXmlFile(
"Tempus_BackwardEuler_CDR.xml");
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;
198 #ifdef TEMPUS_ENABLE_EPETRA_STACK
204 RCP<Epetra_Comm> comm;
205 #ifdef Tempus_ENABLE_MPI
206 comm =
rcp(
new Epetra_MpiComm(MPI_COMM_WORLD));
208 comm =
rcp(
new Epetra_SerialComm);
211 CDR_Test<double, Tempus_Test::CDR_Model<double>>(comm, comm->NumProc(), out,
216 #ifdef TEMPUS_ENABLE_TPETRA_STACK
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;
227 auto comm = Tpetra::getDefaultComm();
229 CDR_Test<SC, Tempus_Test::CDR_Model_Tpetra<SC, LO, GO, Node>>(
230 comm, comm->getSize(), out, success);
T & get(const std::string &name, T def_value)
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)
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
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...
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
Solution state for integrators and steppers.