Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_BackwardEuler_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 "Tempus_config.hpp"
16 #include "Tempus_IntegratorBasic.hpp"
17 #include "Tempus_IntegratorForwardSensitivity.hpp"
18 #include "Tempus_IntegratorPseudoTransientForwardSensitivity.hpp"
19 
20 #include "Thyra_VectorStdOps.hpp"
21 #include "Thyra_MultiVectorStdOps.hpp"
22 
23 #include "../TestModels/SinCosModel.hpp"
24 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
25 
26 #include "Thyra_DefaultMultiVectorProductVector.hpp"
27 
28 #include <vector>
29 #include <fstream>
30 #include <sstream>
31 #include <limits>
32 
33 namespace Tempus_Test {
34 
35 using Teuchos::getParametersFromXmlFile;
37 using Teuchos::RCP;
38 using Teuchos::sublist;
39 
43 
44 // ************************************************************
45 // ************************************************************
46 void test_sincos_fsa(const bool use_combined_method,
47  const bool use_dfdp_as_tangent, Teuchos::FancyOStream &out,
48  bool &success)
49 {
50  std::vector<double> StepSize;
51  std::vector<double> ErrorNorm;
52  const int nTimeStepSizes = 7;
53  double dt = 0.2;
54  double order = 0.0;
57  for (int n = 0; n < nTimeStepSizes; n++) {
58  // Read params from .xml file
59  RCP<ParameterList> pList =
60  getParametersFromXmlFile("Tempus_BackwardEuler_SinCos.xml");
61 
62  // Setup the SinCosModel
63  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
64  scm_pl->set("Use DfDp as Tangent", use_dfdp_as_tangent);
65  RCP<SinCosModel<double> > model =
66  Teuchos::rcp(new SinCosModel<double>(scm_pl));
67 
68  dt /= 2;
69 
70  // Setup sensitivities
71  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
72  ParameterList &sens_pl = pl->sublist("Sensitivities");
73  if (use_combined_method)
74  sens_pl.set("Sensitivity Method", "Combined");
75  else {
76  sens_pl.set("Sensitivity Method", "Staggered");
77  sens_pl.set("Reuse State Linear Solver", true);
78  }
79  sens_pl.set("Use DfDp as Tangent", use_dfdp_as_tangent);
80 
81  // Setup the Integrator and reset initial time step
82  pl->sublist("Default Integrator")
83  .sublist("Time Step Control")
84  .set("Initial Time Step", dt);
86  Tempus::createIntegratorForwardSensitivity<double>(pl, model);
87  order = integrator->getStepper()->getOrder();
88 
89  // Initial Conditions
90  double t0 = pl->sublist("Default Integrator")
91  .sublist("Time Step Control")
92  .get<double>("Initial Time");
94  model->getExactSolution(t0).get_x();
95  const int num_param = model->get_p_space(0)->dim();
97  Thyra::createMembers(model->get_x_space(), num_param);
98  for (int i = 0; i < num_param; ++i)
99  Thyra::assign(DxDp0->col(i).ptr(),
100  *(model->getExactSensSolution(i, t0).get_x()));
101  integrator->initializeSolutionHistory(t0, x0, Teuchos::null, Teuchos::null,
102  DxDp0, Teuchos::null, Teuchos::null);
103 
104  // Integrate to timeMax
105  bool integratorStatus = integrator->advanceTime();
106  TEST_ASSERT(integratorStatus)
107 
108  // Test if at 'Final Time'
109  double time = integrator->getTime();
110  double timeFinal = pl->sublist("Default Integrator")
111  .sublist("Time Step Control")
112  .get<double>("Final Time");
113  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
114 
115  // Time-integrated solution and the exact solution
116  RCP<const Thyra::VectorBase<double> > x = integrator->getX();
117  RCP<const Thyra::MultiVectorBase<double> > DxDp = integrator->getDxDp();
119  model->getExactSolution(time).get_x();
121  Thyra::createMembers(model->get_x_space(), num_param);
122  for (int i = 0; i < num_param; ++i)
123  Thyra::assign(DxDp_exact->col(i).ptr(),
124  *(model->getExactSensSolution(i, time).get_x()));
125 
126  // Plot sample solution and exact solution
127  if (comm->getRank() == 0 && n == nTimeStepSizes - 1) {
129 
130  std::ofstream ftmp("Tempus_BackwardEuler_SinCos_Sens.dat");
131  RCP<const SolutionHistory<double> > solutionHistory =
132  integrator->getSolutionHistory();
133  RCP<Thyra::MultiVectorBase<double> > DxDp_exact_plot =
134  Thyra::createMembers(model->get_x_space(), num_param);
135  for (int i = 0; i < solutionHistory->getNumStates(); i++) {
136  RCP<const SolutionState<double> > solutionState = (*solutionHistory)[i];
137  double time_i = solutionState->getTime();
138  RCP<const DMVPV> x_prod_plot =
139  Teuchos::rcp_dynamic_cast<const DMVPV>(solutionState->getX());
141  x_prod_plot->getMultiVector()->col(0);
143  x_prod_plot->getMultiVector()->subView(
144  Teuchos::Range1D(1, num_param));
145  RCP<const Thyra::VectorBase<double> > x_exact_plot =
146  model->getExactSolution(time_i).get_x();
147  for (int j = 0; j < num_param; ++j)
148  Thyra::assign(DxDp_exact_plot->col(j).ptr(),
149  *(model->getExactSensSolution(j, time_i).get_x()));
150  ftmp << std::fixed << std::setprecision(7) << time_i << std::setw(11)
151  << get_ele(*(x_plot), 0) << std::setw(11) << get_ele(*(x_plot), 1);
152  for (int j = 0; j < num_param; ++j)
153  ftmp << std::setw(11) << get_ele(*(DxDp_plot->col(j)), 0)
154  << std::setw(11) << get_ele(*(DxDp_plot->col(j)), 1);
155  ftmp << std::setw(11) << get_ele(*(x_exact_plot), 0) << std::setw(11)
156  << get_ele(*(x_exact_plot), 1);
157  for (int j = 0; j < num_param; ++j)
158  ftmp << std::setw(11) << get_ele(*(DxDp_exact_plot->col(j)), 0)
159  << std::setw(11) << get_ele(*(DxDp_exact_plot->col(j)), 1);
160  ftmp << std::endl;
161  }
162  ftmp.close();
163  }
164 
165  // Calculate the error
166  RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
167  RCP<Thyra::MultiVectorBase<double> > DxDpdiff = DxDp->clone_mv();
168  Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
169  Thyra::V_VmV(DxDpdiff.ptr(), *DxDp_exact, *DxDp);
170  StepSize.push_back(dt);
171  double L2norm = Thyra::norm_2(*xdiff);
172  L2norm *= L2norm;
173  Teuchos::Array<double> L2norm_DxDp(num_param);
174  Thyra::norms_2(*DxDpdiff, L2norm_DxDp());
175  for (int i = 0; i < num_param; ++i)
176  L2norm += L2norm_DxDp[i] * L2norm_DxDp[i];
177  L2norm = std::sqrt(L2norm);
178  ErrorNorm.push_back(L2norm);
179 
180  out << " n = " << n << " dt = " << dt << " error = " << L2norm << std::endl;
181  }
182 
183  // Check the order and intercept
184  double slope = computeLinearRegressionLogLog<double>(StepSize, ErrorNorm);
185  out << " Stepper = BackwardEuler" << std::endl;
186  out << " =========================" << std::endl;
187  out << " Expected order: " << order << std::endl;
188  out << " Observed order: " << slope << std::endl;
189  out << " =========================" << std::endl;
190  TEST_FLOATING_EQUALITY(slope, order, 0.015);
191  TEST_FLOATING_EQUALITY(ErrorNorm[0], 0.163653, 1.0e-4);
192 
193  if (comm->getRank() == 0) {
194  std::ofstream ftmp("Tempus_BackwardEuler_SinCos_Sens-Error.dat");
195  double error0 = 0.8 * ErrorNorm[0];
196  for (int n = 0; n < nTimeStepSizes; n++) {
197  ftmp << StepSize[n] << " " << ErrorNorm[n] << " "
198  << error0 * (StepSize[n] / StepSize[0]) << std::endl;
199  }
200  ftmp.close();
201  }
202 }
203 
204 } // namespace Tempus_Test
#define TEST_FLOATING_EQUALITY(v1, v2, tol)
#define TEST_ASSERT(v1)
T * get() const
Sine-Cosine model problem from Rythmos. This is a canonical Sine-Cosine differential equation with a...
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
static Teuchos::RCP< const Comm< OrdinalType > > getComm()
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Ptr< T > ptr() const
void test_sincos_fsa(const bool use_combined_method, const bool use_dfdp_as_tangent, Teuchos::FancyOStream &out, bool &success)
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.