Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_IMEX_RK_Partitioned_FSA.hpp
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 #include "Teuchos_DefaultComm.hpp"
14 
15 #include "Thyra_VectorStdOps.hpp"
16 #include "Thyra_MultiVectorStdOps.hpp"
17 #include "Thyra_DefaultMultiVectorProductVector.hpp"
18 #include "Thyra_DefaultProductVector.hpp"
19 
20 #include "Tempus_IntegratorBasic.hpp"
21 #include "Tempus_IntegratorForwardSensitivity.hpp"
22 #include "Tempus_WrapperModelEvaluatorPairPartIMEX_Basic.hpp"
23 
24 #include "../TestModels/VanDerPol_IMEX_ExplicitModel.hpp"
25 #include "../TestModels/VanDerPol_IMEXPart_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 std::string& method_name,
45  const bool use_combined_method,
46  const bool use_dfdp_as_tangent, Teuchos::FancyOStream& out,
47  bool& success)
48 {
49  std::vector<std::string> stepperTypes;
50  stepperTypes.push_back("Partitioned IMEX RK 1st order");
51  stepperTypes.push_back("Partitioned IMEX RK SSP2");
52  stepperTypes.push_back("Partitioned IMEX RK ARS 233");
53  stepperTypes.push_back("General Partitioned IMEX RK");
54 
55  // Check that method_name is valid
56  if (method_name != "") {
57  auto it = std::find(stepperTypes.begin(), stepperTypes.end(), method_name);
58  TEUCHOS_TEST_FOR_EXCEPTION(it == stepperTypes.end(), std::logic_error,
59  "Invalid stepper type " << method_name);
60  }
61 
62  std::vector<double> stepperOrders;
63  std::vector<double> stepperErrors;
64  if (use_dfdp_as_tangent) {
65  if (use_combined_method) {
66  stepperOrders.push_back(1.16082);
67  stepperOrders.push_back(1.97231);
68  stepperOrders.push_back(2.5914);
69  stepperOrders.push_back(1.99148);
70 
71  stepperErrors.push_back(0.00820931);
72  stepperErrors.push_back(0.287112);
73  stepperErrors.push_back(0.00646096);
74  stepperErrors.push_back(0.148848);
75  }
76  else {
77  stepperOrders.push_back(1.07932);
78  stepperOrders.push_back(1.97396);
79  stepperOrders.push_back(2.63724);
80  stepperOrders.push_back(1.99133);
81 
82  stepperErrors.push_back(0.055626);
83  stepperErrors.push_back(0.198898);
84  stepperErrors.push_back(0.00614135);
85  stepperErrors.push_back(0.0999881);
86  }
87  }
88  else {
89  if (use_combined_method) {
90  stepperOrders.push_back(1.1198);
91  stepperOrders.push_back(1.98931);
92  stepperOrders.push_back(2.60509);
93  stepperOrders.push_back(1.992);
94 
95  stepperErrors.push_back(0.00619674);
96  stepperErrors.push_back(0.294989);
97  stepperErrors.push_back(0.0062125);
98  stepperErrors.push_back(0.142489);
99  }
100  else {
101  stepperOrders.push_back(1.07932);
102  stepperOrders.push_back(1.97396);
103  stepperOrders.push_back(2.63724);
104  stepperOrders.push_back(1.99133);
105 
106  stepperErrors.push_back(0.055626);
107  stepperErrors.push_back(0.198898);
108  stepperErrors.push_back(0.00614135);
109  stepperErrors.push_back(0.0999881);
110  }
111  }
112 
113  std::vector<double> stepperInitDt;
114  stepperInitDt.push_back(0.0125);
115  stepperInitDt.push_back(0.05);
116  stepperInitDt.push_back(0.05);
117  stepperInitDt.push_back(0.05);
118 
121 
122  std::vector<std::string>::size_type m;
123  for (m = 0; m != stepperTypes.size(); m++) {
124  // If we were given a method to run, skip this method if it doesn't match
125  if (method_name != "" && stepperTypes[m] != method_name) continue;
126 
127  std::string stepperType = stepperTypes[m];
128  std::string stepperName = stepperTypes[m];
129  std::replace(stepperName.begin(), stepperName.end(), ' ', '_');
130  std::replace(stepperName.begin(), stepperName.end(), '/', '.');
131 
132  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
133  std::vector<RCP<Thyra::VectorBase<double>>> sensitivities;
134  std::vector<double> StepSize;
135  std::vector<double> ErrorNorm;
136  const int nTimeStepSizes = 3; // 6 for error plot
137  double dt = stepperInitDt[m];
138  double order = 0.0;
139  for (int n = 0; n < nTimeStepSizes; n++) {
140  // Read params from .xml file
141  RCP<ParameterList> pList =
142  getParametersFromXmlFile("Tempus_IMEX_RK_VanDerPol.xml");
143 
144  // Setup the explicit VanDerPol ModelEvaluator
145  RCP<ParameterList> vdpmPL = sublist(pList, "VanDerPolModel", true);
146  vdpmPL->set("Use DfDp as Tangent", use_dfdp_as_tangent);
147  const bool useProductVector = true;
149  new VanDerPol_IMEX_ExplicitModel<double>(vdpmPL, useProductVector));
150 
151  // Setup the implicit VanDerPol ModelEvaluator (reuse vdpmPL)
154 
155  // Setup the IMEX Pair ModelEvaluator
156  const int numExplicitBlocks = 1;
157  const int parameterIndex = 4;
159  Teuchos::rcp(
161  explicitModel, implicitModel, numExplicitBlocks,
162  parameterIndex));
163 
164  // Setup sensitivities
165  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
166  ParameterList& sens_pl = pl->sublist("Sensitivities");
167  if (use_combined_method)
168  sens_pl.set("Sensitivity Method", "Combined");
169  else {
170  sens_pl.set("Sensitivity Method", "Staggered");
171  sens_pl.set("Reuse State Linear Solver", true);
172  }
173  sens_pl.set("Use DfDp as Tangent", use_dfdp_as_tangent);
174  ParameterList& interp_pl = pl->sublist("Default Integrator")
175  .sublist("Solution History")
176  .sublist("Interpolator");
177  interp_pl.set("Interpolator Type", "Lagrange");
178  interp_pl.set("Order", 2); // All RK methods here are at most 3rd order
179 
180  // Set the Stepper
181  if (stepperType == "General Partitioned IMEX RK") {
182  // use the appropriate stepper sublist
183  pl->sublist("Default Integrator")
184  .set("Stepper Name", "General IMEX RK");
185  }
186  else {
187  pl->sublist("Default Stepper").set("Stepper Type", stepperType);
188  }
189 
190  // Set the step size
191  if (n == nTimeStepSizes - 1)
192  dt /= 10.0;
193  else
194  dt /= 2;
195 
196  // Setup the Integrator and reset initial time step
197  pl->sublist("Default Integrator")
198  .sublist("Time Step Control")
199  .set("Initial Time Step", dt);
200  pl->sublist("Default Integrator")
201  .sublist("Time Step Control")
202  .remove("Time Step Control Strategy");
204  Tempus::createIntegratorForwardSensitivity<double>(pl, model);
205  order = integrator->getStepper()->getOrder();
206 
207  // Integrate to timeMax
208  bool integratorStatus = integrator->advanceTime();
209  TEST_ASSERT(integratorStatus)
210 
211  // Test if at 'Final Time'
212  double time = integrator->getTime();
213  double timeFinal = pl->sublist("Default Integrator")
214  .sublist("Time Step Control")
215  .get<double>("Final Time");
216  double tol = 100.0 * std::numeric_limits<double>::epsilon();
217  TEST_FLOATING_EQUALITY(time, timeFinal, tol);
218 
219  // Store off the final solution and step size
220  auto solution = Thyra::createMember(model->get_x_space());
221  auto sensitivity = Thyra::createMember(model->get_x_space());
222  Thyra::copy(*(integrator->getX()), solution.ptr());
223  Thyra::copy(*(integrator->getDxDp()->col(0)), sensitivity.ptr());
224  solutions.push_back(solution);
225  sensitivities.push_back(sensitivity);
226  StepSize.push_back(dt);
227 
228  // Output finest temporal solution for plotting
229  if ((n == 0) || (n == nTimeStepSizes - 1)) {
231 
232  std::string fname = "Tempus_" + stepperName + "_VanDerPol_Sens-Ref.dat";
233  if (n == 0) fname = "Tempus_" + stepperName + "_VanDerPol_Sens.dat";
234  std::ofstream ftmp(fname);
235  RCP<const SolutionHistory<double>> solutionHistory =
236  integrator->getSolutionHistory();
237  int nStates = solutionHistory->getNumStates();
238  for (int i = 0; i < nStates; i++) {
239  RCP<const SolutionState<double>> solutionState =
240  (*solutionHistory)[i];
241  RCP<const DMVPV> x_prod =
242  Teuchos::rcp_dynamic_cast<const DMVPV>(solutionState->getX());
244  x_prod->getMultiVector()->col(0);
246  x_prod->getMultiVector()->col(1);
247  double ttime = solutionState->getTime();
248  ftmp << std::fixed << std::setprecision(7) << ttime << " "
249  << std::setw(11) << get_ele(*x, 0) << " " << std::setw(11)
250  << get_ele(*x, 1) << " " << std::setw(11) << get_ele(*dxdp, 0)
251  << " " << std::setw(11) << get_ele(*dxdp, 1) << std::endl;
252  }
253  ftmp.close();
254  }
255  }
256 
257  // Calculate the error - use the most temporally refined mesh for
258  // the reference solution.
259  auto ref_solution = solutions[solutions.size() - 1];
260  auto ref_sensitivity = sensitivities[solutions.size() - 1];
261  std::vector<double> StepSizeCheck;
262  for (std::size_t i = 0; i < (solutions.size() - 1); ++i) {
263  auto sol = solutions[i];
264  auto sen = sensitivities[i];
265  Thyra::Vp_StV(sol.ptr(), -1.0, *ref_solution);
266  Thyra::Vp_StV(sen.ptr(), -1.0, *ref_sensitivity);
267  const double L2norm_sol = Thyra::norm_2(*sol);
268  const double L2norm_sen = Thyra::norm_2(*sen);
269  const double L2norm =
270  std::sqrt(L2norm_sol * L2norm_sol + L2norm_sen * L2norm_sen);
271  StepSizeCheck.push_back(StepSize[i]);
272  ErrorNorm.push_back(L2norm);
273 
274  // out << " n = " << i << " dt = " << StepSize[i]
275  // << " error = " << L2norm << std::endl;
276  }
277 
278  // Check the order and intercept
279  double slope =
280  computeLinearRegressionLogLog<double>(StepSizeCheck, ErrorNorm);
281  out << " Stepper = " << stepperType << std::endl;
282  out << " =========================" << std::endl;
283  out << " Expected order: " << order << std::endl;
284  out << " Observed order: " << slope << std::endl;
285  out << " =========================" << std::endl;
286  TEST_FLOATING_EQUALITY(slope, stepperOrders[m], 0.02);
287  TEST_FLOATING_EQUALITY(ErrorNorm[0], stepperErrors[m], 1.0e-4);
288 
289  // Write error data
290  {
291  std::ofstream ftmp("Tempus_" + stepperName + "_VanDerPol_Sens-Error.dat");
292  double error0 = 0.8 * ErrorNorm[0];
293  for (std::size_t n = 0; n < StepSizeCheck.size(); n++) {
294  ftmp << StepSizeCheck[n] << " " << ErrorNorm[n] << " "
295  << error0 * (pow(StepSize[n] / StepSize[0], order)) << std::endl;
296  }
297  ftmp.close();
298  }
299  }
301 }
302 
303 } // namespace Tempus_Test
van der Pol model formulated for the partitioned IMEX-RK.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
#define TEST_FLOATING_EQUALITY(v1, v2, tol)
#define TEST_ASSERT(v1)
T * get() const
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
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)
ModelEvaluator pair for implicit and explicit (IMEX) evaulations.
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...
std::string method_name
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
Solution state for integrators and steppers.