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