Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_ExplicitRK_PseudoTransient_SA.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 
11 #include "Teuchos_TimeMonitor.hpp"
12 
13 #include "Tempus_config.hpp"
14 #include "Tempus_IntegratorPseudoTransientForwardSensitivity.hpp"
15 #include "Tempus_IntegratorPseudoTransientAdjointSensitivity.hpp"
16 
17 #include "Thyra_VectorStdOps.hpp"
18 #include "Thyra_MultiVectorStdOps.hpp"
19 
20 #include "../TestModels/SteadyQuadraticModel.hpp"
21 
22 #include "Thyra_DefaultMultiVectorProductVector.hpp"
23 #include "Thyra_DefaultProductVector.hpp"
24 
25 #include <vector>
26 #include <fstream>
27 
28 namespace Tempus_Test {
29 
30 using Teuchos::getParametersFromXmlFile;
32 using Teuchos::RCP;
33 using Teuchos::sublist;
34 
35 // ************************************************************
36 // ************************************************************
37 void test_pseudotransient_fsa(const bool use_dfdp_as_tangent,
38  Teuchos::FancyOStream &out, bool &success)
39 {
40  // Read params from .xml file
41  RCP<ParameterList> pList =
42  getParametersFromXmlFile("Tempus_ExplicitRK_SteadyQuadratic.xml");
43 
44  // Setup the SteadyQuadraticModel
45  RCP<ParameterList> scm_pl = sublist(pList, "SteadyQuadraticModel", true);
46  scm_pl->set("Use DfDp as Tangent", use_dfdp_as_tangent);
47  RCP<SteadyQuadraticModel<double> > model =
48  Teuchos::rcp(new SteadyQuadraticModel<double>(scm_pl));
49 
50  // Setup sensitivities
51  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
52  ParameterList &sens_pl = pl->sublist("Sensitivities");
53  sens_pl.set("Use DfDp as Tangent", use_dfdp_as_tangent);
54 
55  // Setup the Integrator
56  RCP<Tempus::IntegratorPseudoTransientForwardSensitivity<double> > integrator =
57  Tempus::createIntegratorPseudoTransientForwardSensitivity<double>(pl,
58  model);
59 
60  // Integrate to timeMax
61  bool integratorStatus = integrator->advanceTime();
62  TEST_ASSERT(integratorStatus);
63 
64  // Test if at 'Final Time'
65  double time = integrator->getTime();
66  double timeFinal = pl->sublist("Demo Integrator")
67  .sublist("Time Step Control")
68  .get<double>("Final Time");
69  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
70 
71  // Time-integrated solution and the exact solution
72  RCP<const Thyra::VectorBase<double> > x_vec = integrator->getX();
73  RCP<const Thyra::MultiVectorBase<double> > DxDp_vec = integrator->getDxDp();
74  const double x = Thyra::get_ele(*x_vec, 0);
75  const double dxdb = Thyra::get_ele(*(DxDp_vec->col(0)), 0);
76  const double x_exact = model->getSteadyStateSolution();
77  const double dxdb_exact = model->getSteadyStateSolutionSensitivity();
78 
79  TEST_FLOATING_EQUALITY(x, x_exact, 1.0e-6);
80  TEST_FLOATING_EQUALITY(dxdb, dxdb_exact, 1.0e-6);
81 }
82 
83 TEUCHOS_UNIT_TEST(ExplicitRK, SteadyQuadratic_PseudoTransient_FSA)
84 {
85  test_pseudotransient_fsa(false, out, success);
86 }
87 
88 TEUCHOS_UNIT_TEST(ExplicitRK, SteadyQuadratic_PseudoTransient_FSA_Tangent)
89 {
90  test_pseudotransient_fsa(true, out, success);
91 }
92 
93 // ************************************************************
94 // ************************************************************
95 TEUCHOS_UNIT_TEST(ExplicitRK, SteadyQuadratic_PseudoTransient_ASA)
96 {
97  // Read params from .xml file
98  RCP<ParameterList> pList =
99  getParametersFromXmlFile("Tempus_ExplicitRK_SteadyQuadratic.xml");
100 
101  // Setup the SteadyQuadraticModel
102  RCP<ParameterList> scm_pl = sublist(pList, "SteadyQuadraticModel", true);
105 
106  // Setup sensitivities
107  RCP<ParameterList> pl = sublist(pList, "Tempus", true);
108  // ParameterList& sens_pl = pl->sublist("Sensitivities");
109 
110  // Setup the Integrator
112  Tempus::integratorPseudoTransientAdjointSensitivity<double>(pl, model);
113 
114  // Integrate to timeMax
115  bool integratorStatus = integrator->advanceTime();
116  TEST_ASSERT(integratorStatus);
117 
118  // Test if at 'Final Time'
119  double time = integrator->getTime();
120  double timeFinal = pl->sublist("Demo Integrator")
121  .sublist("Time Step Control")
122  .get<double>("Final Time");
123  TEST_FLOATING_EQUALITY(time, timeFinal, 1.0e-14);
124 
125  // Time-integrated solution and the exact solution using the fact that g = x
126  RCP<const Thyra::VectorBase<double> > x_vec = integrator->getX();
127  RCP<const Thyra::MultiVectorBase<double> > DxDp_vec = integrator->getDgDp();
128  const double x = Thyra::get_ele(*x_vec, 0);
129  const double dxdb = Thyra::get_ele(*(DxDp_vec->col(0)), 0);
130  const double x_exact = model->getSteadyStateSolution();
131  const double dxdb_exact = model->getSteadyStateSolutionSensitivity();
132 
133  TEST_FLOATING_EQUALITY(x, x_exact, 1.0e-6);
134  TEST_FLOATING_EQUALITY(dxdb, dxdb_exact, 1.0e-6);
135 }
136 
137 } // namespace Tempus_Test
void test_pseudotransient_fsa(const bool use_dfdp_as_tangent, Teuchos::FancyOStream &out, bool &success)
#define TEST_FLOATING_EQUALITY(v1, v2, tol)
#define TEST_ASSERT(v1)
T * get() const
TEUCHOS_UNIT_TEST(BackwardEuler, SinCos_ASA)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Simple quadratic equation with a stable steady-state.