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: Copyright (2017) Sandia Corporation
4 //
5 // Distributed under BSD 3-clause license (See accompanying file Copyright.txt)
6 // ****************************************************************************
7 // @HEADER
8 
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 "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 void test_sincos_fsa(const bool use_combined_method,
46  const bool use_dfdp_as_tangent, Teuchos::FancyOStream &out,
47  bool &success)
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;
56  for (int n = 0; n < nTimeStepSizes; n++) {
57  // Read params from .xml file
58  RCP<ParameterList> pList =
59  getParametersFromXmlFile("Tempus_BackwardEuler_SinCos.xml");
60 
61  // Setup the SinCosModel
62  RCP<ParameterList> scm_pl = sublist(pList, "SinCosModel", true);
63  scm_pl->set("Use DfDp as Tangent", use_dfdp_as_tangent);
64  RCP<SinCosModel<double> > model =
65  Teuchos::rcp(new SinCosModel<double>(scm_pl));
66 
67  dt /= 2;
68 
69  // Setup sensitivities
70  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
71  ParameterList &sens_pl = pl->sublist("Sensitivities");
72  if (use_combined_method)
73  sens_pl.set("Sensitivity Method", "Combined");
74  else {
75  sens_pl.set("Sensitivity Method", "Staggered");
76  sens_pl.set("Reuse State Linear Solver", true);
77  }
78  sens_pl.set("Use DfDp as Tangent", use_dfdp_as_tangent);
79 
80  // Setup the Integrator and reset initial time step
81  pl->sublist("Default Integrator")
82  .sublist("Time Step Control")
83  .set("Initial Time Step", dt);
85  Tempus::createIntegratorForwardSensitivity<double>(pl, model);
86  order = integrator->getStepper()->getOrder();
87 
88  // Initial Conditions
89  double t0 = pl->sublist("Default Integrator")
90  .sublist("Time Step Control")
91  .get<double>("Initial Time");
93  model->getExactSolution(t0).get_x();
94  const int num_param = model->get_p_space(0)->dim();
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->initializeSolutionHistory(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")
111  .get<double>("Final Time");
112  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
113 
114  // Time-integrated solution and the exact solution
115  RCP<const Thyra::VectorBase<double> > x = integrator->getX();
116  RCP<const Thyra::MultiVectorBase<double> > DxDp = integrator->getDxDp();
118  model->getExactSolution(time).get_x();
120  Thyra::createMembers(model->get_x_space(), num_param);
121  for (int i = 0; i < num_param; ++i)
122  Thyra::assign(DxDp_exact->col(i).ptr(),
123  *(model->getExactSensSolution(i, time).get_x()));
124 
125  // Plot sample solution and exact solution
126  if (comm->getRank() == 0 && n == nTimeStepSizes - 1) {
128 
129  std::ofstream ftmp("Tempus_BackwardEuler_SinCos_Sens.dat");
130  RCP<const SolutionHistory<double> > solutionHistory =
131  integrator->getSolutionHistory();
132  RCP<Thyra::MultiVectorBase<double> > DxDp_exact_plot =
133  Thyra::createMembers(model->get_x_space(), num_param);
134  for (int i = 0; i < solutionHistory->getNumStates(); i++) {
135  RCP<const SolutionState<double> > solutionState = (*solutionHistory)[i];
136  double time_i = solutionState->getTime();
137  RCP<const DMVPV> x_prod_plot =
138  Teuchos::rcp_dynamic_cast<const DMVPV>(solutionState->getX());
140  x_prod_plot->getMultiVector()->col(0);
142  x_prod_plot->getMultiVector()->subView(
143  Teuchos::Range1D(1, num_param));
144  RCP<const Thyra::VectorBase<double> > x_exact_plot =
145  model->getExactSolution(time_i).get_x();
146  for (int j = 0; j < num_param; ++j)
147  Thyra::assign(DxDp_exact_plot->col(j).ptr(),
148  *(model->getExactSensSolution(j, time_i).get_x()));
149  ftmp << std::fixed << std::setprecision(7) << time_i << std::setw(11)
150  << get_ele(*(x_plot), 0) << std::setw(11) << get_ele(*(x_plot), 1);
151  for (int j = 0; j < num_param; ++j)
152  ftmp << std::setw(11) << get_ele(*(DxDp_plot->col(j)), 0)
153  << std::setw(11) << get_ele(*(DxDp_plot->col(j)), 1);
154  ftmp << std::setw(11) << get_ele(*(x_exact_plot), 0) << std::setw(11)
155  << get_ele(*(x_exact_plot), 1);
156  for (int j = 0; j < num_param; ++j)
157  ftmp << std::setw(11) << get_ele(*(DxDp_exact_plot->col(j)), 0)
158  << std::setw(11) << get_ele(*(DxDp_exact_plot->col(j)), 1);
159  ftmp << std::endl;
160  }
161  ftmp.close();
162  }
163 
164  // Calculate the error
165  RCP<Thyra::VectorBase<double> > xdiff = x->clone_v();
166  RCP<Thyra::MultiVectorBase<double> > DxDpdiff = DxDp->clone_mv();
167  Thyra::V_StVpStV(xdiff.ptr(), 1.0, *x_exact, -1.0, *(x));
168  Thyra::V_VmV(DxDpdiff.ptr(), *DxDp_exact, *DxDp);
169  StepSize.push_back(dt);
170  double L2norm = Thyra::norm_2(*xdiff);
171  L2norm *= L2norm;
172  Teuchos::Array<double> L2norm_DxDp(num_param);
173  Thyra::norms_2(*DxDpdiff, L2norm_DxDp());
174  for (int i = 0; i < num_param; ++i)
175  L2norm += L2norm_DxDp[i] * L2norm_DxDp[i];
176  L2norm = std::sqrt(L2norm);
177  ErrorNorm.push_back(L2norm);
178 
179  out << " n = " << n << " dt = " << dt << " error = " << L2norm << std::endl;
180  }
181 
182  // Check the order and intercept
183  double slope = computeLinearRegressionLogLog<double>(StepSize, ErrorNorm);
184  out << " Stepper = BackwardEuler" << std::endl;
185  out << " =========================" << std::endl;
186  out << " Expected order: " << order << std::endl;
187  out << " Observed order: " << slope << std::endl;
188  out << " =========================" << std::endl;
189  TEST_FLOATING_EQUALITY(slope, order, 0.015);
190  TEST_FLOATING_EQUALITY(ErrorNorm[0], 0.163653, 1.0e-4);
191 
192  if (comm->getRank() == 0) {
193  std::ofstream ftmp("Tempus_BackwardEuler_SinCos_Sens-Error.dat");
194  double error0 = 0.8 * ErrorNorm[0];
195  for (int n = 0; n < nTimeStepSizes; n++) {
196  ftmp << StepSize[n] << " " << ErrorNorm[n] << " "
197  << error0 * (StepSize[n] / StepSize[0]) << std::endl;
198  }
199  ftmp.close();
200  }
201 }
202 
203 } // 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_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.