Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_IMEX_RK_FSA.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ****************************************************************************
3 // Tempus: Copyright (2017) Sandia Corporation
4 //
5 // Distributed under BSD 3-clause license (See accompanying file Copyright.txt)
6 // ****************************************************************************
7 // @HEADER
8 
11 #include "Teuchos_TimeMonitor.hpp"
12 #include "Teuchos_DefaultComm.hpp"
13 
14 #include "Thyra_VectorStdOps.hpp"
15 #include "Thyra_MultiVectorStdOps.hpp"
16 #include "Thyra_DefaultMultiVectorProductVector.hpp"
17 #include "Thyra_DefaultProductVector.hpp"
18 
19 #include "Tempus_IntegratorBasic.hpp"
20 #include "Tempus_IntegratorForwardSensitivity.hpp"
21 #include "Tempus_WrapperModelEvaluatorPairIMEX_Basic.hpp"
22 
23 #include "../TestModels/VanDerPolModel.hpp"
24 #include "../TestModels/VanDerPol_IMEX_ExplicitModel.hpp"
25 #include "../TestModels/VanDerPol_IMEX_ImplicitModel.hpp"
26 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
27 
28 #include <fstream>
29 #include <vector>
30 
31 namespace Tempus_Test {
32 
33 using Teuchos::getParametersFromXmlFile;
35 using Teuchos::RCP;
36 using Teuchos::sublist;
37 
41 
42 // ************************************************************
43 // ************************************************************
44 void test_vdp_fsa(const bool use_combined_method,
45  const bool use_dfdp_as_tangent, Teuchos::FancyOStream& out,
46  bool& success)
47 {
48  std::vector<std::string> stepperTypes;
49  stepperTypes.push_back("IMEX RK 1st order");
50  stepperTypes.push_back("IMEX RK SSP2");
51  stepperTypes.push_back("IMEX RK ARS 233");
52  stepperTypes.push_back("General IMEX RK");
53 
54  std::vector<double> stepperOrders;
55  std::vector<double> stepperErrors;
56  if (use_combined_method) {
57  stepperOrders.push_back(1.1198);
58  stepperOrders.push_back(1.98931);
59  stepperOrders.push_back(2.60509);
60  stepperOrders.push_back(1.992);
61 
62  stepperErrors.push_back(0.00619674);
63  stepperErrors.push_back(0.294989);
64  stepperErrors.push_back(0.0062125);
65  stepperErrors.push_back(0.142489);
66  }
67  else {
68  stepperOrders.push_back(1.1198);
69  stepperOrders.push_back(2.05232);
70  stepperOrders.push_back(2.43013);
71  stepperOrders.push_back(2.07975);
72 
73  stepperErrors.push_back(0.00619674);
74  stepperErrors.push_back(0.0534172);
75  stepperErrors.push_back(0.00224533);
76  stepperErrors.push_back(0.032632);
77  }
78  std::vector<double> stepperInitDt;
79  stepperInitDt.push_back(0.0125);
80  stepperInitDt.push_back(0.05);
81  stepperInitDt.push_back(0.05);
82  stepperInitDt.push_back(0.05);
83 
86 
87  std::vector<std::string>::size_type m;
88  for (m = 0; m != stepperTypes.size(); m++) {
89  std::string stepperType = stepperTypes[m];
90  std::string stepperName = stepperTypes[m];
91  std::replace(stepperName.begin(), stepperName.end(), ' ', '_');
92  std::replace(stepperName.begin(), stepperName.end(), '/', '.');
93 
94  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
95  std::vector<RCP<Thyra::VectorBase<double>>> sensitivities;
96  std::vector<double> StepSize;
97  std::vector<double> ErrorNorm;
98  const int nTimeStepSizes = 3; // 6 for error plot
99  double dt = stepperInitDt[m];
100  double order = 0.0;
101  for (int n = 0; n < nTimeStepSizes; n++) {
102  // Read params from .xml file
103  RCP<ParameterList> pList =
104  getParametersFromXmlFile("Tempus_IMEX_RK_VanDerPol.xml");
105 
106  // Setup the explicit VanDerPol ModelEvaluator
107  RCP<ParameterList> vdpmPL = sublist(pList, "VanDerPolModel", true);
108  vdpmPL->set("Use DfDp as Tangent", use_dfdp_as_tangent);
111 
112  // Setup the implicit VanDerPol ModelEvaluator (reuse vdpmPL)
115 
116  // Setup the IMEX Pair ModelEvaluator
119  explicitModel, implicitModel));
120 
121  // Setup sensitivities
122  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
123  ParameterList& sens_pl = pl->sublist("Sensitivities");
124  if (use_combined_method)
125  sens_pl.set("Sensitivity Method", "Combined");
126  else {
127  sens_pl.set("Sensitivity Method", "Staggered");
128  sens_pl.set("Reuse State Linear Solver", true);
129  }
130  sens_pl.set("Use DfDp as Tangent", use_dfdp_as_tangent);
131  ParameterList& interp_pl = pl->sublist("Default Integrator")
132  .sublist("Solution History")
133  .sublist("Interpolator");
134  interp_pl.set("Interpolator Type", "Lagrange");
135  interp_pl.set("Order", 2); // All RK methods here are at most 3rd order
136 
137  // Set the Stepper
138  if (stepperType == "General IMEX RK") {
139  // use the appropriate stepper sublist
140  pl->sublist("Default Integrator")
141  .set("Stepper Name", "General IMEX RK");
142  }
143  else {
144  pl->sublist("Default Stepper").set("Stepper Type", stepperType);
145  }
146 
147  // Set the step size
148  if (n == nTimeStepSizes - 1)
149  dt /= 10.0;
150  else
151  dt /= 2;
152 
153  // Setup the Integrator and reset initial time step
154  pl->sublist("Default Integrator")
155  .sublist("Time Step Control")
156  .set("Initial Time Step", dt);
157  pl->sublist("Default Integrator")
158  .sublist("Time Step Control")
159  .remove("Time Step Control Strategy");
161  Tempus::createIntegratorForwardSensitivity<double>(pl, model);
162  order = integrator->getStepper()->getOrder();
163 
164  // Integrate to timeMax
165  bool integratorStatus = integrator->advanceTime();
166  TEST_ASSERT(integratorStatus)
167 
168  // Test if at 'Final Time'
169  double time = integrator->getTime();
170  double timeFinal = pl->sublist("Default Integrator")
171  .sublist("Time Step Control")
172  .get<double>("Final Time");
173  double tol = 100.0 * std::numeric_limits<double>::epsilon();
174  TEST_FLOATING_EQUALITY(time, timeFinal, tol);
175 
176  // Store off the final solution and step size
177  auto solution = Thyra::createMember(model->get_x_space());
178  auto sensitivity = Thyra::createMember(model->get_x_space());
179  Thyra::copy(*(integrator->getX()), solution.ptr());
180  Thyra::copy(*(integrator->getDxDp()->col(0)), sensitivity.ptr());
181  solutions.push_back(solution);
182  sensitivities.push_back(sensitivity);
183  StepSize.push_back(dt);
184 
185  // Output finest temporal solution for plotting
186  if (comm->getRank() == 0 && ((n == 0) || (n == nTimeStepSizes - 1))) {
188 
189  std::string fname = "Tempus_" + stepperName + "_VanDerPol_Sens-Ref.dat";
190  if (n == 0) fname = "Tempus_" + stepperName + "_VanDerPol_Sens.dat";
191  std::ofstream ftmp(fname);
192  RCP<const SolutionHistory<double>> solutionHistory =
193  integrator->getSolutionHistory();
194  int nStates = solutionHistory->getNumStates();
195  for (int i = 0; i < nStates; i++) {
196  RCP<const SolutionState<double>> solutionState =
197  (*solutionHistory)[i];
198  RCP<const DMVPV> x_prod =
199  Teuchos::rcp_dynamic_cast<const DMVPV>(solutionState->getX());
201  x_prod->getMultiVector()->col(0);
203  x_prod->getMultiVector()->col(1);
204  double ttime = solutionState->getTime();
205  ftmp << std::fixed << std::setprecision(7) << ttime << " "
206  << std::setw(11) << get_ele(*x, 0) << " " << std::setw(11)
207  << get_ele(*x, 1) << " " << std::setw(11) << get_ele(*dxdp, 0)
208  << " " << std::setw(11) << get_ele(*dxdp, 1) << std::endl;
209  }
210  ftmp.close();
211  }
212  }
213 
214  // Calculate the error - use the most temporally refined mesh for
215  // the reference solution.
216  auto ref_solution = solutions[solutions.size() - 1];
217  auto ref_sensitivity = sensitivities[solutions.size() - 1];
218  std::vector<double> StepSizeCheck;
219  for (std::size_t i = 0; i < (solutions.size() - 1); ++i) {
220  auto sol = solutions[i];
221  auto sen = sensitivities[i];
222  Thyra::Vp_StV(sol.ptr(), -1.0, *ref_solution);
223  Thyra::Vp_StV(sen.ptr(), -1.0, *ref_sensitivity);
224  const double L2norm_sol = Thyra::norm_2(*sol);
225  const double L2norm_sen = Thyra::norm_2(*sen);
226  const double L2norm =
227  std::sqrt(L2norm_sol * L2norm_sol + L2norm_sen * L2norm_sen);
228  StepSizeCheck.push_back(StepSize[i]);
229  ErrorNorm.push_back(L2norm);
230 
231  out << " n = " << i << " dt = " << StepSize[i] << " error = " << L2norm
232  << std::endl;
233  }
234 
235  // Check the order and intercept
236  double slope =
237  computeLinearRegressionLogLog<double>(StepSizeCheck, ErrorNorm);
238  out << " Stepper = " << stepperType << std::endl;
239  out << " =========================" << std::endl;
240  out << " Expected order: " << order << std::endl;
241  out << " Observed order: " << slope << std::endl;
242  out << " =========================" << std::endl;
243  TEST_FLOATING_EQUALITY(slope, stepperOrders[m], 0.02);
244  TEST_FLOATING_EQUALITY(ErrorNorm[0], stepperErrors[m], 1.0e-4);
245 
246  // Write error data
247  {
248  std::ofstream ftmp("Tempus_" + stepperName + "_VanDerPol_Sens-Error.dat");
249  double error0 = 0.8 * ErrorNorm[0];
250  for (std::size_t n = 0; n < StepSizeCheck.size(); n++) {
251  ftmp << StepSizeCheck[n] << " " << ErrorNorm[n] << " "
252  << error0 * (pow(StepSize[n] / StepSize[0], order)) << std::endl;
253  }
254  ftmp.close();
255  }
256  }
258 }
259 
260 } // namespace Tempus_Test
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
ModelEvaluator pair for implicit and explicit (IMEX) evaulations.
#define TEST_FLOATING_EQUALITY(v1, v2, tol)
#define TEST_ASSERT(v1)
T * get() const
void test_vdp_fsa(const bool use_combined_method, const bool use_dfdp_as_tangent, Teuchos::FancyOStream &out, bool &success)
static Teuchos::RCP< const Comm< OrdinalType > > getComm()
bool remove(std::string const &name, bool throwIfNotExists=true)
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)
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.