Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_IntegratorTest.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_IntegratorBasic.hpp"
14 
15 #include "../TestModels/SinCosModel.hpp"
16 
17 #include <vector>
18 
19 namespace Tempus_Test {
20 
21 using Teuchos::RCP;
23 using Teuchos::sublist;
24 using Teuchos::getParametersFromXmlFile;
25 
29 
30 // Test Integrator construction from ParameterList and ModelEvaluator.
31 TEUCHOS_UNIT_TEST(IntegratorBasic, PL_ME_Construction)
32 {
33  // 1) Setup the ParameterList (here we start with params from .xml file)
34  RCP<ParameterList> pl = getParametersFromXmlFile("Tempus_default.xml");
35 
36  // 2) Setup the ModelEvaluator
38 
39  // 3) Setup the Integrator
40  RCP<ParameterList> tempusPL = sublist(pl, "Tempus", true);
42  Tempus::createIntegratorBasic<double>(tempusPL, model);
43 
44  // Test the ParameterList
45  auto testPL = integrator->getValidParameters();
46  // Write out ParameterList to rebaseline test.
47  //writeParameterListToXmlFile(*testPL, "Tempus_IntegratorBasic_ref-test.xml");
48 
49  // Read params from reference .xml file
50  RCP<ParameterList> referencePL =
51  getParametersFromXmlFile("Tempus_IntegratorBasic_ref.xml");
52 
53  bool pass = haveSameValuesSorted(*testPL, *referencePL, true);
54  if (!pass) {
55  std::cout << std::endl;
56  std::cout << "testPL -------------- \n" << *testPL << std::endl;
57  std::cout << "referencePL -------------- \n" << *referencePL << std::endl;
58  }
59  TEST_ASSERT(pass)
60 }
61 
62 
63 // Test integator construction, and then setParameterList, setStepper, and
64 // initialization.
66 {
67  // 1) Setup the Integrator
69  Tempus::createIntegratorBasic<double>();
70 
71  // 2) Setup the ParameterList
72  // - Start with the default Tempus PL
73  // - Add Stepper PL
74  RCP<ParameterList> tempusPL = Teuchos::rcp_const_cast<ParameterList>(
75  integrator->getValidParameters());
76 
77  tempusPL->sublist("Default Integrator").set("Stepper Name", "Demo Stepper");
78  RCP<ParameterList> stepperPL = Teuchos::parameterList();
79  stepperPL->set("Stepper Type", "Forward Euler");
80  tempusPL->set("Demo Stepper", *stepperPL);
81 
82  // 3) Create integrator
84  integrator = Tempus::createIntegratorBasic<double>(tempusPL, model);
85  integrator->initialize();
86 
87  // Test the ParameterList
88  auto testPL = integrator->getValidParameters();
89  // Write out ParameterList to rebaseline test.
90  //writeParameterListToXmlFile(*testPL,"Tempus_IntegratorBasic_ref2-test.xml");
91 
92  // Read params from reference .xml file
93  RCP<ParameterList> referencePL =
94  getParametersFromXmlFile("Tempus_IntegratorBasic_ref2.xml");
95 
96  bool pass = haveSameValuesSorted(*testPL, *referencePL, true);
97  if (!pass) {
98  std::cout << std::endl;
99  std::cout << "testPL -------------- \n" << *testPL << std::endl;
100  std::cout << "referencePL -------------- \n" << *referencePL << std::endl;
101  }
102  TEST_ASSERT(pass)
103 }
104 
105 } // namespace Tempus_Test
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
#define TEST_ASSERT(v1)
Sine-Cosine model problem from Rythmos. This is a canonical Sine-Cosine differential equation with a...
TEUCHOS_UNIT_TEST(BackwardEuler, SinCos_ASA)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
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. SolutionState contains the metadata for solutions and th...