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