Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_BackwardEuler_ASA.cpp
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_IntegratorAdjointSensitivity.hpp"
18 
19 #include "Thyra_VectorStdOps.hpp"
20 #include "Thyra_MultiVectorStdOps.hpp"
21 
22 #include "../TestModels/SinCosModel.hpp"
23 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
24 
25 #include "Thyra_DefaultMultiVectorProductVector.hpp"
26 
27 #include <vector>
28 #include <fstream>
29 #include <sstream>
30 #include <limits>
31 
32 namespace Tempus_Test {
33 
34 using Teuchos::getParametersFromXmlFile;
36 using Teuchos::RCP;
37 using Teuchos::sublist;
38 
42 
43 // ************************************************************
44 // ************************************************************
45 TEUCHOS_UNIT_TEST(BackwardEuler, SinCos_ASA)
46 {
47  std::vector<double> StepSize;
48  std::vector<double> ErrorNorm;
49  const int nTimeStepSizes = 7;
50  double dt = 0.2;
51  double order = 0.0;
54  for (int n = 0; n < nTimeStepSizes; n++) {
55  // Read params from .xml file
56  RCP<ParameterList> pList =
57  getParametersFromXmlFile("Tempus_BackwardEuler_SinCos_ASA.xml");
58 
59  // Setup the SinCosModel
60  // Here we test using an explicit adjoint model for adjoint sensitivities
61  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
62  RCP<SinCosModel<double> > model =
63  Teuchos::rcp(new SinCosModel<double>(scm_pl));
64  RCP<SinCosModelAdjoint<double> > adjoint_model =
66 
67  dt /= 2;
68 
69  // Setup sensitivities
70  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
71  ParameterList& sens_pl = pl->sublist("Sensitivities");
72  sens_pl.set("Mass Matrix Is Identity", false); // Just for testing
73  ParameterList& interp_pl = pl->sublist("Default Integrator")
74  .sublist("Solution History")
75  .sublist("Interpolator");
76  interp_pl.set("Interpolator Type", "Lagrange");
77  interp_pl.set("Order", 0);
78 
79  // Set FSAL to false, because it is not currently setup for ASA.
80  pl->sublist("Default Stepper").set("Use FSAL", false);
81 
82  // Set IC consistency check to false, because it is not currently
83  // setup for ASA.
84  pl->sublist("Default Stepper")
85  .set("Initial Condition Consistency Check", false);
86 
87  // Setup the Integrator and reset initial time step
88  pl->sublist("Default Integrator")
89  .sublist("Time Step Control")
90  .set("Initial Time Step", dt);
92  Tempus::createIntegratorAdjointSensitivity<double>(pl, model,
93  adjoint_model);
94  order = integrator->getStepper()->getOrder();
95 
96  // Initial Conditions
97  double t0 = pl->sublist("Default Integrator")
98  .sublist("Time Step Control")
99  .get<double>("Initial Time");
101  model->getExactSolution(t0).get_x();
102  const int num_param = model->get_p_space(0)->dim();
104  Thyra::createMembers(model->get_x_space(), num_param);
105  for (int i = 0; i < num_param; ++i)
106  Thyra::assign(DxDp0->col(i).ptr(),
107  *(model->getExactSensSolution(i, t0).get_x()));
108  integrator->initializeSolutionHistory(t0, x0, Teuchos::null, Teuchos::null,
109  DxDp0, Teuchos::null, Teuchos::null);
110 
111  // Integrate to timeMax
112  bool integratorStatus = integrator->advanceTime();
113  TEST_ASSERT(integratorStatus)
114 
115  // Test if at 'Final Time'
116  double time = integrator->getTime();
117  double timeFinal = pl->sublist("Default Integrator")
118  .sublist("Time Step Control")
119  .get<double>("Final Time");
120  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
121 
122  // Time-integrated solution and the exact solution along with
123  // sensitivities (relying on response g(x) = x). Note we must transpose
124  // dg/dp since the integrator returns it in gradient form.
125  RCP<const Thyra::VectorBase<double> > x = integrator->getX();
126  RCP<const Thyra::MultiVectorBase<double> > DgDp = integrator->getDgDp();
128  Thyra::createMembers(model->get_x_space(), num_param);
129  {
132  const int num_g = DgDp->domain()->dim();
133  for (int i = 0; i < num_g; ++i)
134  for (int j = 0; j < num_param; ++j) dxdp_view(i, j) = dgdp_view(j, i);
135  }
137  model->getExactSolution(time).get_x();
139  Thyra::createMembers(model->get_x_space(), num_param);
140  for (int i = 0; i < num_param; ++i)
141  Thyra::assign(DxDp_exact->col(i).ptr(),
142  *(model->getExactSensSolution(i, time).get_x()));
143 
144  // Plot sample solution, exact solution, and adjoint solution
145  if (comm->getRank() == 0 && n == nTimeStepSizes - 1) {
148 
149  std::ofstream ftmp("Tempus_BackwardEuler_SinCos_AdjSens.dat");
150  RCP<const SolutionHistory<double> > solutionHistory =
151  integrator->getSolutionHistory();
152  for (int i = 0; i < solutionHistory->getNumStates(); i++) {
153  RCP<const SolutionState<double> > solutionState = (*solutionHistory)[i];
154  const double time_i = solutionState->getTime();
155  RCP<const DPV> x_prod_plot =
156  Teuchos::rcp_dynamic_cast<const DPV>(solutionState->getX());
158  x_prod_plot->getVectorBlock(0);
159  RCP<const DMVPV> adjoint_prod_plot =
160  Teuchos::rcp_dynamic_cast<const DMVPV>(
161  x_prod_plot->getVectorBlock(1));
163  adjoint_prod_plot->getMultiVector();
164  RCP<const Thyra::VectorBase<double> > x_exact_plot =
165  model->getExactSolution(time_i).get_x();
166  ftmp << std::fixed << std::setprecision(7) << time_i << std::setw(11)
167  << get_ele(*(x_plot), 0) << std::setw(11) << get_ele(*(x_plot), 1)
168  << std::setw(11) << get_ele(*(adjoint_plot->col(0)), 0)
169  << std::setw(11) << get_ele(*(adjoint_plot->col(0)), 1)
170  << std::setw(11) << get_ele(*(adjoint_plot->col(1)), 0)
171  << std::setw(11) << get_ele(*(adjoint_plot->col(1)), 1)
172  << std::setw(11) << get_ele(*(x_exact_plot), 0) << std::setw(11)
173  << get_ele(*(x_exact_plot), 1) << std::endl;
174  }
175  ftmp.close();
176  }
177 
178  // Calculate the error
179  RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
180  RCP<Thyra::MultiVectorBase<double> > DxDpdiff = DxDp->clone_mv();
181  Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
182  Thyra::V_VmV(DxDpdiff.ptr(), *DxDp_exact, *DxDp);
183  StepSize.push_back(dt);
184  double L2norm = Thyra::norm_2(*xdiff);
185  L2norm *= L2norm;
186  Teuchos::Array<double> L2norm_DxDp(num_param);
187  Thyra::norms_2(*DxDpdiff, L2norm_DxDp());
188  for (int i = 0; i < num_param; ++i)
189  L2norm += L2norm_DxDp[i] * L2norm_DxDp[i];
190  L2norm = std::sqrt(L2norm);
191  ErrorNorm.push_back(L2norm);
192 
193  // out << " n = " << n << " dt = " << dt << " error = " << L2norm
194  // << std::endl;
195  }
196 
197  // Check the order and intercept
198  double slope = computeLinearRegressionLogLog<double>(StepSize, ErrorNorm);
199  out << " Stepper = BackwardEuler" << std::endl;
200  out << " =========================" << std::endl;
201  out << " Expected order: " << order << std::endl;
202  out << " Observed order: " << slope << std::endl;
203  out << " =========================" << std::endl;
204  TEST_FLOATING_EQUALITY(slope, order, 0.015);
205  TEST_FLOATING_EQUALITY(ErrorNorm[0], 0.142525, 1.0e-4);
206 
207  if (comm->getRank() == 0) {
208  std::ofstream ftmp("Tempus_BackwardEuler_SinCos_AdjSens-Error.dat");
209  double error0 = 0.8 * ErrorNorm[0];
210  for (int n = 0; n < nTimeStepSizes; n++) {
211  ftmp << StepSize[n] << " " << ErrorNorm[n] << " "
212  << error0 * (StepSize[n] / StepSize[0]) << std::endl;
213  }
214  ftmp.close();
215  }
216 }
217 
218 } // namespace Tempus_Test
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
#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...
static Teuchos::RCP< const Comm< OrdinalType > > getComm()
TEUCHOS_UNIT_TEST(BackwardEuler, SinCos_ASA)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Ptr< T > ptr() const
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.