Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_OperatorSplitTest.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 
12 #include "Teuchos_TimeMonitor.hpp"
13 
14 #include "Thyra_VectorStdOps.hpp"
15 
16 #include "Tempus_IntegratorBasic.hpp"
17 #include "Tempus_StepperFactory.hpp"
18 #include "Tempus_StepperOperatorSplit.hpp"
19 #include "Tempus_StepperForwardEuler.hpp"
20 #include "Tempus_StepperBackwardEuler.hpp"
21 
22 #include "../TestModels/VanDerPol_IMEX_ExplicitModel.hpp"
23 #include "../TestModels/VanDerPol_IMEX_ImplicitModel.hpp"
24 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
25 
26 #include <fstream>
27 #include <vector>
28 
29 namespace Tempus_Test {
30 
31 using Teuchos::getParametersFromXmlFile;
33 using Teuchos::RCP;
34 using Teuchos::rcp;
35 using Teuchos::rcp_const_cast;
36 using Teuchos::sublist;
37 
41 
42 // ************************************************************
43 // ************************************************************
44 TEUCHOS_UNIT_TEST(OperatorSplit, ConstructingFromDefaults)
45 {
46  double dt = 0.05;
47 
48  // Read params from .xml file
49  RCP<ParameterList> pList =
50  getParametersFromXmlFile("Tempus_OperatorSplit_VanDerPol.xml");
51  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
52 
53  // Setup the explicit VanDerPol ModelEvaluator
54  RCP<ParameterList> vdpmPL = sublist(pList, "VanDerPolModel", true);
57 
58  // Setup the implicit VanDerPol ModelEvaluator (reuse vdpmPL)
61 
62  // Setup Stepper for field solve ----------------------------
63  auto stepper = rcp(new Tempus::StepperOperatorSplit<double>());
64 
65  auto subStepper1 =
66  Tempus::createStepperForwardEuler(explicitModel, Teuchos::null);
67  auto subStepper2 =
68  Tempus::createStepperBackwardEuler(implicitModel, Teuchos::null);
69 
70  stepper->addStepper(subStepper1);
71  stepper->addStepper(subStepper2);
72  stepper->initialize();
73 
74  // Setup TimeStepControl ------------------------------------
75  auto timeStepControl = rcp(new Tempus::TimeStepControl<double>());
76  ParameterList tscPL =
77  pl->sublist("Demo Integrator").sublist("Time Step Control");
78  timeStepControl->setInitIndex(tscPL.get<int>("Initial Time Index"));
79  timeStepControl->setInitTime(tscPL.get<double>("Initial Time"));
80  timeStepControl->setFinalTime(tscPL.get<double>("Final Time"));
81  timeStepControl->setInitTimeStep(dt);
82  timeStepControl->initialize();
83 
84  // Setup initial condition SolutionState --------------------
85  auto inArgsIC = stepper->getModel()->getNominalValues();
86  auto icX = rcp_const_cast<Thyra::VectorBase<double>>(inArgsIC.get_x());
87  auto icXDot = rcp_const_cast<Thyra::VectorBase<double>>(inArgsIC.get_x_dot());
88  auto icState = Tempus::createSolutionStateX(icX, icXDot);
89  icState->setTime(timeStepControl->getInitTime());
90  icState->setIndex(timeStepControl->getInitIndex());
91  icState->setTimeStep(0.0);
92  icState->setOrder(stepper->getOrder());
93  icState->setSolutionStatus(Tempus::Status::PASSED); // ICs are passing.
94 
95  // Setup SolutionHistory ------------------------------------
96  auto solutionHistory = rcp(new Tempus::SolutionHistory<double>());
97  solutionHistory->setName("Forward States");
98  solutionHistory->setStorageType(Tempus::STORAGE_TYPE_STATIC);
99  solutionHistory->setStorageLimit(2);
100  solutionHistory->addState(icState);
101 
102  // Setup Integrator -----------------------------------------
104  Tempus::createIntegratorBasic<double>();
105  integrator->setStepper(stepper);
106  integrator->setTimeStepControl(timeStepControl);
107  integrator->setSolutionHistory(solutionHistory);
108  // integrator->setObserver(...);
109  integrator->initialize();
110 
111  // Integrate to timeMax
112  bool integratorStatus = integrator->advanceTime();
113  TEST_ASSERT(integratorStatus)
114 
115  // Test if at 'Final Time'
116  double time = integrator->getTime();
117  double timeFinal = pl->sublist("Demo Integrator")
118  .sublist("Time Step Control")
119  .get<double>("Final Time");
120  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
121 
122  // Time-integrated solution and the exact solution
123  RCP<Thyra::VectorBase<double>> x = integrator->getX();
124 
125  // Check the order and intercept
126  out << " Stepper = " << stepper->description() << std::endl;
127  out << " =========================" << std::endl;
128  out << " Computed solution: " << get_ele(*(x), 0) << " "
129  << get_ele(*(x), 1) << std::endl;
130  out << " =========================" << std::endl;
131  TEST_FLOATING_EQUALITY(get_ele(*(x), 0), -2.223910, 1.0e-4);
132  TEST_FLOATING_EQUALITY(get_ele(*(x), 1), 0.565441, 1.0e-4);
133 }
134 
135 // ************************************************************
136 // ************************************************************
137 TEUCHOS_UNIT_TEST(OperatorSplit, VanDerPol)
138 {
140  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
141  std::vector<RCP<Thyra::VectorBase<double>>> solutionsDot;
142  std::vector<double> StepSize;
143  std::vector<double> xErrorNorm;
144  std::vector<double> xDotErrorNorm;
145  const int nTimeStepSizes = 4; // 8 for Error plot
146  double dt = 0.1;
147  double time = 0.0;
148  for (int n = 0; n < nTimeStepSizes; n++) {
149  // Read params from .xml file
150  RCP<ParameterList> pList =
151  getParametersFromXmlFile("Tempus_OperatorSplit_VanDerPol.xml");
152 
153  // Setup the explicit VanDerPol ModelEvaluator
154  RCP<ParameterList> vdpmPL = sublist(pList, "VanDerPolModel", true);
155  auto explicitModel = rcp(new VanDerPol_IMEX_ExplicitModel<double>(vdpmPL));
156 
157  // Setup the implicit VanDerPol ModelEvaluator (reuse vdpmPL)
158  auto implicitModel = rcp(new VanDerPol_IMEX_ImplicitModel<double>(vdpmPL));
159 
160  // Setup vector of models
161  std::vector<RCP<const Thyra::ModelEvaluator<double>>> models;
162  models.push_back(explicitModel);
163  models.push_back(implicitModel);
164 
165  // Set the step size
166  dt /= 2;
167  if (n == nTimeStepSizes - 1) dt /= 10.0;
168 
169  // Setup the Integrator and reset initial time step
170  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
171  pl->sublist("Demo Integrator")
172  .sublist("Time Step Control")
173  .set("Initial Time Step", dt);
174  integrator = Tempus::createIntegratorBasic<double>(pl, models);
175 
176  // Integrate to timeMax
177  bool integratorStatus = integrator->advanceTime();
178  TEST_ASSERT(integratorStatus)
179 
180  // Test if at 'Final Time'
181  time = integrator->getTime();
182  double timeFinal = pl->sublist("Demo Integrator")
183  .sublist("Time Step Control")
184  .get<double>("Final Time");
185  double tol = 100.0 * std::numeric_limits<double>::epsilon();
186  TEST_FLOATING_EQUALITY(time, timeFinal, tol);
187 
188  // Store off the final solution and step size
189  StepSize.push_back(dt);
190  auto solution = Thyra::createMember(implicitModel->get_x_space());
191  Thyra::copy(*(integrator->getX()), solution.ptr());
192  solutions.push_back(solution);
193  auto solutionDot = Thyra::createMember(implicitModel->get_x_space());
194  Thyra::copy(*(integrator->getXDot()), solutionDot.ptr());
195  solutionsDot.push_back(solutionDot);
196 
197  // Output finest temporal solution for plotting
198  // This only works for ONE MPI process
199  if ((n == 0) || (n == nTimeStepSizes - 1)) {
200  std::string fname = "Tempus_OperatorSplit_VanDerPol-Ref.dat";
201  if (n == 0) fname = "Tempus_OperatorSplit_VanDerPol.dat";
202  RCP<const SolutionHistory<double>> solutionHistory =
203  integrator->getSolutionHistory();
204  writeSolution(fname, solutionHistory);
205  // solutionHistory->printHistory("medium");
206  }
207  }
208 
209  // Check the order and intercept
210  double xSlope = 0.0;
211  double xDotSlope = 0.0;
212  RCP<Tempus::Stepper<double>> stepper = integrator->getStepper();
213  double order = stepper->getOrder();
214  writeOrderError("Tempus_OperatorSplit_VanDerPol-Error.dat", stepper, StepSize,
215  solutions, xErrorNorm, xSlope, solutionsDot, xDotErrorNorm,
216  xDotSlope, out);
217 
218  TEST_FLOATING_EQUALITY(xSlope, order, 0.05);
219  TEST_FLOATING_EQUALITY(xDotSlope, order, 0.05); //=order at small dt
220  TEST_FLOATING_EQUALITY(xErrorNorm[0], 1.27294, 1.0e-4);
221  TEST_FLOATING_EQUALITY(xDotErrorNorm[0], 12.7102, 1.0e-4);
222 
224 }
225 
226 } // namespace Tempus_Test
Teuchos::RCP< SolutionState< Scalar > > createSolutionStateX(const Teuchos::RCP< Thyra::VectorBase< Scalar > > &x, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &xdot=Teuchos::null, const Teuchos::RCP< Thyra::VectorBase< Scalar > > &xdotdot=Teuchos::null)
Nonmember constructor from non-const solution vectors, x.
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)
#define TEST_ASSERT(v1)
T * get() const
void writeSolution(const std::string filename, Teuchos::RCP< const Tempus::SolutionHistory< Scalar >> solutionHistory)
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
OperatorSplit stepper loops through the Stepper list.
Teuchos::RCP< StepperForwardEuler< Scalar > > createStepperForwardEuler(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model, Teuchos::RCP< Teuchos::ParameterList > pl)
Nonmember constructor - ModelEvaluator and ParameterList.
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)
TimeStepControl manages the time step size. There several mechanisms that effect the time step size a...
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Keep a fix number of states.
Teuchos::RCP< StepperBackwardEuler< Scalar > > createStepperBackwardEuler(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model, Teuchos::RCP< Teuchos::ParameterList > pl)
Nonmember constructor - ModelEvaluator and ParameterList.
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
Solution state for integrators and steppers.