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