Tempus  Version of the Day
Time Integration
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 
9 #include "Teuchos_UnitTestHarness.hpp"
10 #include "Teuchos_XMLParameterListHelpers.hpp"
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::RCP;
33 using Teuchos::ParameterList;
34 using Teuchos::sublist;
35 using Teuchos::getParametersFromXmlFile;
36 
40 
41 
42 // ************************************************************
43 // ************************************************************
44 void test_vdp_fsa(const bool use_combined_method,
45  const bool use_dfdp_as_tangent,
46  Teuchos::FancyOStream &out, 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  std::vector<double> stepperOrders;
55  std::vector<double> stepperErrors;
56  if (use_dfdp_as_tangent) {
57  if (use_combined_method) {
58  stepperOrders.push_back(1.16082);
59  stepperOrders.push_back(1.97231);
60  stepperOrders.push_back(2.5914);
61  stepperOrders.push_back(1.99148);
62 
63  stepperErrors.push_back(0.00820931);
64  stepperErrors.push_back(0.287112);
65  stepperErrors.push_back(0.00646096);
66  stepperErrors.push_back(0.148848);
67  }
68  else {
69  stepperOrders.push_back(1.07932);
70  stepperOrders.push_back(1.97396);
71  stepperOrders.push_back(2.63724);
72  stepperOrders.push_back(1.99133);
73 
74  stepperErrors.push_back(0.055626);
75  stepperErrors.push_back(0.198898);
76  stepperErrors.push_back(0.00614135);
77  stepperErrors.push_back(0.0999881);
78  }
79  }
80  else {
81  if (use_combined_method) {
82  stepperOrders.push_back(1.1198);
83  stepperOrders.push_back(1.98931);
84  stepperOrders.push_back(2.60509);
85  stepperOrders.push_back(1.992);
86 
87  stepperErrors.push_back(0.00619674);
88  stepperErrors.push_back(0.294989);
89  stepperErrors.push_back(0.0062125);
90  stepperErrors.push_back(0.142489);
91  }
92  else {
93  stepperOrders.push_back(1.07932);
94  stepperOrders.push_back(1.97396);
95  stepperOrders.push_back(2.63724);
96  stepperOrders.push_back(1.99133);
97 
98  stepperErrors.push_back(0.055626);
99  stepperErrors.push_back(0.198898);
100  stepperErrors.push_back(0.00614135);
101  stepperErrors.push_back(0.0999881);
102  }
103  }
104 
105  std::vector<double> stepperInitDt;
106  stepperInitDt.push_back(0.0125);
107  stepperInitDt.push_back(0.05);
108  stepperInitDt.push_back(0.05);
109  stepperInitDt.push_back(0.05);
110 
111  Teuchos::RCP<const Teuchos::Comm<int> > comm =
112  Teuchos::DefaultComm<int>::getComm();
113  Teuchos::RCP<Teuchos::FancyOStream> my_out =
114  Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
115  my_out->setProcRankAndSize(comm->getRank(), comm->getSize());
116  my_out->setOutputToRootOnly(0);
117 
118  std::vector<std::string>::size_type m;
119  for(m = 0; m != stepperTypes.size(); m++) {
120 
121  std::string stepperType = stepperTypes[m];
122  std::string stepperName = stepperTypes[m];
123  std::replace(stepperName.begin(), stepperName.end(), ' ', '_');
124  std::replace(stepperName.begin(), stepperName.end(), '/', '.');
125 
126  std::vector<RCP<Thyra::VectorBase<double>>> solutions;
127  std::vector<RCP<Thyra::VectorBase<double>>> sensitivities;
128  std::vector<double> StepSize;
129  std::vector<double> ErrorNorm;
130  const int nTimeStepSizes = 3; // 6 for error plot
131  double dt = stepperInitDt[m];
132  double order = 0.0;
133  for (int n=0; n<nTimeStepSizes; n++) {
134 
135  // Read params from .xml file
136  RCP<ParameterList> pList =
137  getParametersFromXmlFile("Tempus_IMEX_RK_VanDerPol.xml");
138 
139  // Setup the explicit VanDerPol ModelEvaluator
140  RCP<ParameterList> vdpmPL = sublist(pList, "VanDerPolModel", true);
141  vdpmPL->set("Use DfDp as Tangent", use_dfdp_as_tangent);
142  const bool useProductVector = true;
143  RCP<VanDerPol_IMEX_ExplicitModel<double> > explicitModel =
144  Teuchos::rcp(new VanDerPol_IMEX_ExplicitModel<double>(vdpmPL,
145  useProductVector));
146 
147  // Setup the implicit VanDerPol ModelEvaluator (reuse vdpmPL)
148  RCP<VanDerPol_IMEXPart_ImplicitModel<double> > implicitModel =
149  Teuchos::rcp(new VanDerPol_IMEXPart_ImplicitModel<double>(vdpmPL));
150 
151  // Setup the IMEX Pair ModelEvaluator
152  const int numExplicitBlocks = 1;
153  const int parameterIndex = 4;
154  RCP<Tempus::WrapperModelEvaluatorPairPartIMEX_Basic<double> > model =
155  Teuchos::rcp(
157  explicitModel, implicitModel,
158  numExplicitBlocks, parameterIndex));
159 
160  // Setup sensitivities
161  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
162  ParameterList& sens_pl = pl->sublist("Sensitivities");
163  if (use_combined_method)
164  sens_pl.set("Sensitivity Method", "Combined");
165  else {
166  sens_pl.set("Sensitivity Method", "Staggered");
167  sens_pl.set("Reuse State Linear Solver", true);
168  }
169  sens_pl.set("Use DfDp as Tangent", use_dfdp_as_tangent);
170  ParameterList& interp_pl =
171  pl->sublist("Default Integrator").sublist("Solution History").sublist("Interpolator");
172  interp_pl.set("Interpolator Type", "Lagrange");
173  interp_pl.set("Order", 2); // All RK methods here are at most 3rd order
174 
175  // Set the Stepper
176  if (stepperType == "General Partitioned IMEX RK"){
177  // use the appropriate stepper sublist
178  pl->sublist("Default Integrator").set("Stepper Name", "General IMEX RK");
179  } else {
180  pl->sublist("Default Stepper").set("Stepper Type", stepperType);
181  }
182 
183  // Set the step size
184  if (n == nTimeStepSizes-1) dt /= 10.0;
185  else dt /= 2;
186 
187  // Setup the Integrator and reset initial time step
188  pl->sublist("Default Integrator")
189  .sublist("Time Step Control").set("Initial Time Step", dt);
190  RCP<Tempus::IntegratorForwardSensitivity<double> > integrator =
191  Tempus::integratorForwardSensitivity<double>(pl, model);
192  order = integrator->getStepper()->getOrder();
193 
194  // Integrate to timeMax
195  bool integratorStatus = integrator->advanceTime();
196  TEST_ASSERT(integratorStatus)
197 
198  // Test if at 'Final Time'
199  double time = integrator->getTime();
200  double timeFinal =pl->sublist("Default Integrator")
201  .sublist("Time Step Control").get<double>("Final Time");
202  double tol = 100.0 * std::numeric_limits<double>::epsilon();
203  TEST_FLOATING_EQUALITY(time, timeFinal, tol);
204 
205  // Store off the final solution and step size
206  auto solution = Thyra::createMember(model->get_x_space());
207  auto sensitivity = Thyra::createMember(model->get_x_space());
208  Thyra::copy(*(integrator->getX()),solution.ptr());
209  Thyra::copy(*(integrator->getDxDp()->col(0)),sensitivity.ptr());
210  solutions.push_back(solution);
211  sensitivities.push_back(sensitivity);
212  StepSize.push_back(dt);
213 
214  // Output finest temporal solution for plotting
215  if ((n == 0) or (n == nTimeStepSizes-1)) {
216  typedef Thyra::DefaultMultiVectorProductVector<double> DMVPV;
217 
218  std::string fname = "Tempus_"+stepperName+"_VanDerPol_Sens-Ref.dat";
219  if (n == 0) fname = "Tempus_"+stepperName+"_VanDerPol_Sens.dat";
220  std::ofstream ftmp(fname);
221  RCP<const SolutionHistory<double> > solutionHistory =
222  integrator->getSolutionHistory();
223  int nStates = solutionHistory->getNumStates();
224  for (int i=0; i<nStates; i++) {
225  RCP<const SolutionState<double> > solutionState =
226  (*solutionHistory)[i];
227  RCP<const DMVPV> x_prod =
228  Teuchos::rcp_dynamic_cast<const DMVPV>(solutionState->getX());
229  RCP<const Thyra::VectorBase<double> > x =
230  x_prod->getMultiVector()->col(0);
231  RCP<const Thyra::VectorBase<double> > dxdp =
232  x_prod->getMultiVector()->col(1);
233  double ttime = solutionState->getTime();
234  ftmp << std::fixed << std::setprecision(7)
235  << ttime << " "
236  << std::setw(11) << get_ele(*x, 0) << " "
237  << std::setw(11) << get_ele(*x, 1) << " "
238  << std::setw(11) << get_ele(*dxdp, 0) << " "
239  << std::setw(11) << get_ele(*dxdp, 1)
240  << std::endl;
241  }
242  ftmp.close();
243  }
244  }
245 
246  // Calculate the error - use the most temporally refined mesh for
247  // the reference solution.
248  auto ref_solution = solutions[solutions.size()-1];
249  auto ref_sensitivity = sensitivities[solutions.size()-1];
250  std::vector<double> StepSizeCheck;
251  for (std::size_t i=0; i < (solutions.size()-1); ++i) {
252  auto sol = solutions[i];
253  auto sen = sensitivities[i];
254  Thyra::Vp_StV(sol.ptr(), -1.0, *ref_solution);
255  Thyra::Vp_StV(sen.ptr(), -1.0, *ref_sensitivity);
256  const double L2norm_sol = Thyra::norm_2(*sol);
257  const double L2norm_sen = Thyra::norm_2(*sen);
258  const double L2norm =
259  std::sqrt(L2norm_sol*L2norm_sol + L2norm_sen*L2norm_sen);
260  StepSizeCheck.push_back(StepSize[i]);
261  ErrorNorm.push_back(L2norm);
262 
263  *my_out << " n = " << i << " dt = " << StepSize[i]
264  << " error = " << L2norm << std::endl;
265  }
266 
267  // Check the order and intercept
268  double slope = computeLinearRegressionLogLog<double>(StepSizeCheck,ErrorNorm);
269  std::cout << " Stepper = " << stepperType << std::endl;
270  std::cout << " =========================" << std::endl;
271  std::cout << " Expected order: " << order << std::endl;
272  std::cout << " Observed order: " << slope << std::endl;
273  std::cout << " =========================" << std::endl;
274  TEST_FLOATING_EQUALITY( slope, stepperOrders[m], 0.02 );
275  TEST_FLOATING_EQUALITY( ErrorNorm[0], stepperErrors[m], 1.0e-4 );
276 
277  // Write error data
278  {
279  std::ofstream ftmp("Tempus_"+stepperName+"_VanDerPol_Sens-Error.dat");
280  double error0 = 0.8*ErrorNorm[0];
281  for (std::size_t n = 0; n < StepSizeCheck.size(); n++) {
282  ftmp << StepSizeCheck[n] << " " << ErrorNorm[n] << " "
283  << error0*(pow(StepSize[n]/StepSize[0],order)) << std::endl;
284  }
285  ftmp.close();
286  }
287  }
288  Teuchos::TimeMonitor::summarize();
289 }
290 
291 
292 } // namespace Tempus_Test
Tempus_Test
Definition: Tempus_BackwardEuler_ASA.cpp:34
Tempus::WrapperModelEvaluatorPairPartIMEX_Basic
ModelEvaluator pair for implicit and explicit (IMEX) evaulations.
Definition: Tempus_WrapperModelEvaluatorPairPartIMEX_Basic_decl.hpp:37
Tempus::solutionHistory
Teuchos::RCP< SolutionHistory< Scalar > > solutionHistory(Teuchos::RCP< Teuchos::ParameterList > pList=Teuchos::null)
Nonmember constructor.
Definition: Tempus_SolutionHistory_impl.hpp:504
Tempus::IntegratorBasic
Basic time integrator.
Definition: Tempus_IntegratorBasic_decl.hpp:35
Tempus_Test::test_vdp_fsa
void test_vdp_fsa(const bool use_combined_method, const bool use_dfdp_as_tangent, Teuchos::FancyOStream &out, bool &success)
Definition: Tempus_IMEX_RK_FSA.hpp:45
Tempus::SolutionState
Solution state for integrators and steppers. SolutionState contains the metadata for solutions and th...
Definition: Tempus_SolutionState_decl.hpp:56
Tempus::SolutionHistory
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Definition: Tempus_Integrator.hpp:25