Tempus  Version of the Day
Time Integration
Tempus_BackwardEuler_ASA.cpp
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 "Tempus_config.hpp"
15 #include "Tempus_IntegratorBasic.hpp"
16 #include "Tempus_IntegratorAdjointSensitivity.hpp"
17 
18 #include "Thyra_VectorStdOps.hpp"
19 #include "Thyra_MultiVectorStdOps.hpp"
20 
21 #include "../TestModels/SinCosModel.hpp"
22 #include "../TestUtils/Tempus_ConvergenceTestUtils.hpp"
23 
24 #include "Stratimikos_DefaultLinearSolverBuilder.hpp"
25 #include "Thyra_LinearOpWithSolveFactoryHelpers.hpp"
26 #include "Thyra_DefaultMultiVectorProductVector.hpp"
27 #include "Thyra_DefaultProductVector.hpp"
28 
29 #include <vector>
30 #include <fstream>
31 #include <sstream>
32 #include <limits>
33 
34 namespace Tempus_Test {
35 
36 using Teuchos::RCP;
37 using Teuchos::ParameterList;
38 using Teuchos::sublist;
39 using Teuchos::getParametersFromXmlFile;
40 
44 
45 // ************************************************************
46 // ************************************************************
47 TEUCHOS_UNIT_TEST(BackwardEuler, SinCos_ASA)
48 {
49  std::vector<double> StepSize;
50  std::vector<double> ErrorNorm;
51  const int nTimeStepSizes = 7;
52  double dt = 0.2;
53  double order = 0.0;
54  Teuchos::RCP<const Teuchos::Comm<int> > comm =
55  Teuchos::DefaultComm<int>::getComm();
56  Teuchos::RCP<Teuchos::FancyOStream> my_out =
57  Teuchos::fancyOStream(Teuchos::rcpFromRef(std::cout));
58  my_out->setProcRankAndSize(comm->getRank(), comm->getSize());
59  my_out->setOutputToRootOnly(0);
60  for (int n=0; n<nTimeStepSizes; n++) {
61 
62  // Read params from .xml file
63  RCP<ParameterList> pList =
64  getParametersFromXmlFile("Tempus_BackwardEuler_SinCos.xml");
65 
66  // Setup the SinCosModel
67  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
68  RCP<SinCosModel<double> > model =
69  Teuchos::rcp(new SinCosModel<double>(scm_pl));
70 
71  dt /= 2;
72 
73  // Setup sensitivities
74  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
75  ParameterList& sens_pl = pl->sublist("Sensitivities");
76  sens_pl.set("Mass Matrix Is Identity", false); // Just for testing
77  ParameterList& interp_pl =
78  pl->sublist("Default Integrator").sublist("Solution History").sublist("Interpolator");
79  interp_pl.set("Interpolator Type", "Lagrange");
80  interp_pl.set("Order", 0);
81 
82  // Setup the Integrator and reset initial time step
83  pl->sublist("Default Integrator")
84  .sublist("Time Step Control").set("Initial Time Step", dt);
85  RCP<Tempus::IntegratorAdjointSensitivity<double> > integrator =
86  Tempus::integratorAdjointSensitivity<double>(pl, model);
87  order = integrator->getStepper()->getOrder();
88 
89  // Initial Conditions
90  double t0 = pl->sublist("Default Integrator")
91  .sublist("Time Step Control").get<double>("Initial Time");
92  RCP<const Thyra::VectorBase<double> > x0 =
93  model->getExactSolution(t0).get_x();
94  const int num_param = model->get_p_space(0)->dim();
95  RCP<Thyra::MultiVectorBase<double> > DxDp0 =
96  Thyra::createMembers(model->get_x_space(), num_param);
97  for (int i=0; i<num_param; ++i)
98  Thyra::assign(DxDp0->col(i).ptr(),
99  *(model->getExactSensSolution(i, t0).get_x()));
100  integrator->setInitialState(t0, x0, Teuchos::null, Teuchos::null,
101  DxDp0, Teuchos::null, Teuchos::null);
102 
103  // Integrate to timeMax
104  bool integratorStatus = integrator->advanceTime();
105  TEST_ASSERT(integratorStatus)
106 
107  // Test if at 'Final Time'
108  double time = integrator->getTime();
109  double timeFinal =pl->sublist("Default Integrator")
110  .sublist("Time Step Control").get<double>("Final Time");
111  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
112 
113  // Time-integrated solution and the exact solution along with
114  // sensitivities (relying on response g(x) = x). Note we must transpose
115  // dg/dp since the integrator returns it in gradient form.
116  RCP<const Thyra::VectorBase<double> > x = integrator->getX();
117  RCP<const Thyra::MultiVectorBase<double> > DgDp = integrator->getDgDp();
118  RCP<Thyra::MultiVectorBase<double> > DxDp =
119  Thyra::createMembers(model->get_x_space(), num_param);
120  {
121  Thyra::ConstDetachedMultiVectorView<double> dgdp_view(*DgDp);
122  Thyra::DetachedMultiVectorView<double> dxdp_view(*DxDp);
123  const int num_g = DgDp->domain()->dim();
124  for (int i=0; i<num_g; ++i)
125  for (int j=0; j<num_param; ++j)
126  dxdp_view(i,j) = dgdp_view(j,i);
127  }
128  RCP<const Thyra::VectorBase<double> > x_exact =
129  model->getExactSolution(time).get_x();
130  RCP<Thyra::MultiVectorBase<double> > DxDp_exact =
131  Thyra::createMembers(model->get_x_space(), num_param);
132  for (int i=0; i<num_param; ++i)
133  Thyra::assign(DxDp_exact->col(i).ptr(),
134  *(model->getExactSensSolution(i, time).get_x()));
135 
136  // Plot sample solution, exact solution, and adjoint solution
137  if (comm->getRank() == 0 && n == nTimeStepSizes-1) {
138  typedef Thyra::DefaultProductVector<double> DPV;
139  typedef Thyra::DefaultMultiVectorProductVector<double> DMVPV;
140 
141  std::ofstream ftmp("Tempus_BackwardEuler_SinCos_AdjSens.dat");
142  RCP<const SolutionHistory<double> > solutionHistory =
143  integrator->getSolutionHistory();
144  for (int i=0; i<solutionHistory->getNumStates(); i++) {
145  RCP<const SolutionState<double> > solutionState = (*solutionHistory)[i];
146  const double time = solutionState->getTime();
147  RCP<const DPV> x_prod_plot =
148  Teuchos::rcp_dynamic_cast<const DPV>(solutionState->getX());
149  RCP<const Thyra::VectorBase<double> > x_plot =
150  x_prod_plot->getVectorBlock(0);
151  RCP<const DMVPV > adjoint_prod_plot =
152  Teuchos::rcp_dynamic_cast<const DMVPV>(x_prod_plot->getVectorBlock(1));
153  RCP<const Thyra::MultiVectorBase<double> > adjoint_plot =
154  adjoint_prod_plot->getMultiVector();
155  RCP<const Thyra::VectorBase<double> > x_exact_plot =
156  model->getExactSolution(time).get_x();
157  ftmp << std::fixed << std::setprecision(7)
158  << time
159  << std::setw(11) << get_ele(*(x_plot), 0)
160  << std::setw(11) << get_ele(*(x_plot), 1)
161  << std::setw(11) << get_ele(*(adjoint_plot->col(0)), 0)
162  << std::setw(11) << get_ele(*(adjoint_plot->col(0)), 1)
163  << std::setw(11) << get_ele(*(adjoint_plot->col(1)), 0)
164  << std::setw(11) << get_ele(*(adjoint_plot->col(1)), 1)
165  << std::setw(11) << get_ele(*(x_exact_plot), 0)
166  << std::setw(11) << get_ele(*(x_exact_plot), 1)
167  << std::endl;
168  }
169  ftmp.close();
170  }
171 
172  // Calculate the error
173  RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
174  RCP<Thyra::MultiVectorBase<double> > DxDpdiff = DxDp->clone_mv();
175  Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
176  Thyra::V_VmV(DxDpdiff.ptr(), *DxDp_exact, *DxDp);
177  StepSize.push_back(dt);
178  double L2norm = Thyra::norm_2(*xdiff);
179  L2norm *= L2norm;
180  Teuchos::Array<double> L2norm_DxDp(num_param);
181  Thyra::norms_2(*DxDpdiff, L2norm_DxDp());
182  for (int i=0; i<num_param; ++i)
183  L2norm += L2norm_DxDp[i]*L2norm_DxDp[i];
184  L2norm = std::sqrt(L2norm);
185  ErrorNorm.push_back(L2norm);
186 
187  *my_out << " n = " << n << " dt = " << dt << " error = " << L2norm
188  << std::endl;
189 
190  }
191 
192  // Check the order and intercept
193  double slope = computeLinearRegressionLogLog<double>(StepSize, ErrorNorm);
194  *my_out << " Stepper = BackwardEuler" << std::endl;
195  *my_out << " =========================" << std::endl;
196  *my_out << " Expected order: " << order << std::endl;
197  *my_out << " Observed order: " << slope << std::endl;
198  *my_out << " =========================" << std::endl;
199  TEST_FLOATING_EQUALITY( slope, order, 0.015 );
200  TEST_FLOATING_EQUALITY( ErrorNorm[0], 0.151746, 1.0e-4 );
201 
202  if (comm->getRank() == 0) {
203  std::ofstream ftmp("Tempus_BackwardEuler_SinCos_AdjSens-Error.dat");
204  double error0 = 0.8*ErrorNorm[0];
205  for (int n=0; n<nTimeStepSizes; n++) {
206  ftmp << StepSize[n] << " " << ErrorNorm[n] << " "
207  << error0*(StepSize[n]/StepSize[0]) << std::endl;
208  }
209  ftmp.close();
210  }
211 
212 }
213 
214 } // namespace Tempus_Test
Tempus_Test
Definition: Tempus_BackwardEuler_ASA.cpp:34
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::SinCosModel
Sine-Cosine model problem from Rythmos. This is a canonical Sine-Cosine differential equation.
Definition: SinCosModel_decl.hpp:91
Tempus::SolutionState
Solution state for integrators and steppers. SolutionState contains the metadata for solutions and th...
Definition: Tempus_SolutionState_decl.hpp:56
Tempus_Test::TEUCHOS_UNIT_TEST
TEUCHOS_UNIT_TEST(BackwardEuler, SinCos_ASA)
Definition: Tempus_BackwardEuler_ASA.cpp:47
Tempus::SolutionHistory
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Definition: Tempus_Integrator.hpp:25