Tempus
Version of the Day
Time Integration
Toggle main menu visibility
Loading...
Searching...
No Matches
src
Tempus_IntegratorObserverBasic_impl.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
#ifndef Tempus_IntegratorObserverBasic_impl_hpp
10
#define Tempus_IntegratorObserverBasic_impl_hpp
11
12
#include "Tempus_Stepper.hpp"
13
14
namespace
Tempus
{
15
16
template
<
class
Scalar>
17
IntegratorObserverBasic<Scalar>::IntegratorObserverBasic
(){}
18
19
template
<
class
Scalar>
20
IntegratorObserverBasic<Scalar>::~IntegratorObserverBasic
(){}
21
22
template
<
class
Scalar>
23
void
IntegratorObserverBasic<Scalar>::
24
observeStartIntegrator
(
const
Integrator<Scalar>
& integrator){
25
26
std::time_t begin = std::time(
nullptr
);
27
const
Teuchos::RCP<Teuchos::FancyOStream> out = integrator.getOStream();
28
out->setOutputToRootOnly(0);
29
Teuchos::OSTab ostab(out,0,
"ScreenOutput"
);
30
*out <<
"\nTempus - IntegratorBasic\n"
31
<< std::asctime(std::localtime(&begin)) <<
"\n"
32
<<
" Stepper = "
<< integrator.
getStepper
()->description() <<
"\n"
33
<<
" Simulation Time Range ["
<< integrator.
getTimeStepControl
()->getInitTime()
34
<<
", "
<< integrator.
getTimeStepControl
()->getFinalTime() <<
"]\n"
35
<<
" Simulation Index Range ["
<< integrator.
getTimeStepControl
()->getInitIndex()
36
<<
", "
<< integrator.
getTimeStepControl
()->getFinalIndex() <<
"]\n"
37
<<
"============================================================================\n"
38
<<
" Step Time dt Abs Error Rel Error Order nFail dCompTime"
39
<< std::endl;
40
}
41
42
template
<
class
Scalar>
43
void
IntegratorObserverBasic<Scalar>::
44
observeStartTimeStep
(
const
Integrator<Scalar>
&
/* integrator */
){}
45
46
template
<
class
Scalar>
47
void
IntegratorObserverBasic<Scalar>::
48
observeNextTimeStep
(
const
Integrator<Scalar>
&
/* integrator */
){}
49
50
template
<
class
Scalar>
51
void
IntegratorObserverBasic<Scalar>::
52
observeBeforeTakeStep
(
const
Integrator<Scalar>
&
/* integrator */
){}
53
54
template
<
class
Scalar>
55
void
IntegratorObserverBasic<Scalar>::
56
observeAfterTakeStep
(
const
Integrator<Scalar>
&
/* integrator */
){}
57
58
template
<
class
Scalar>
59
void
IntegratorObserverBasic<Scalar>::
60
observeAfterCheckTimeStep
(
const
Integrator<Scalar>
&
/* integrator */
){}
61
62
template
<
class
Scalar>
63
void
IntegratorObserverBasic<Scalar>::
64
observeEndTimeStep
(
const
Integrator<Scalar>
& integrator){
65
66
using
Teuchos::RCP;
67
auto
cs = integrator.
getSolutionHistory
()->getCurrentState();
68
69
if
((cs->getOutputScreen() ==
true
) ||
70
(cs->getOutput() ==
true
) ||
71
(cs->getTime() == integrator.
getTimeStepControl
()->getFinalTime())) {
72
73
const
Scalar steppertime = integrator.
getStepperTimer
()->totalElapsedTime();
74
// reset the stepper timer
75
integrator.
getStepperTimer
()->reset();
76
77
const
Teuchos::RCP<Teuchos::FancyOStream> out = integrator.getOStream();
78
out->setOutputToRootOnly(0);
79
Teuchos::OSTab ostab(out, 0,
"ScreenOutput"
);
80
*out<<std::scientific
81
<<std::setw( 6)<<std::setprecision(3)<<cs->getIndex()
82
<<std::setw(11)<<std::setprecision(3)<<cs->getTime()
83
<<std::setw(11)<<std::setprecision(3)<<cs->getTimeStep()
84
<<std::setw(11)<<std::setprecision(3)<<cs->getErrorAbs()
85
<<std::setw(11)<<std::setprecision(3)<<cs->getErrorRel()
86
<<std::fixed <<std::setw( 7)<<std::setprecision(1)<<cs->getOrder()
87
<<std::scientific<<std::setw( 7)<<std::setprecision(3)<<cs->getNFailures()
88
<<std::setw(11)<<std::setprecision(3)<<steppertime
89
<<std::endl;
90
}
91
92
}
93
94
template
<
class
Scalar>
95
void
IntegratorObserverBasic<Scalar>::
96
observeEndIntegrator
(
const
Integrator<Scalar>
& integrator){
97
98
std::string exitStatus;
99
//const Scalar runtime = integrator.getIntegratorTimer()->totalElapsedTime();
100
if
(integrator.
getSolutionHistory
()->getCurrentState()->getSolutionStatus() ==
101
Status::FAILED
|| integrator.
getStatus
() ==
Status::FAILED
) {
102
exitStatus =
"Time integration FAILURE!"
;
103
}
else
{
104
exitStatus =
"Time integration complete."
;
105
}
106
std::time_t end = std::time(
nullptr
);
107
const
Scalar runtime = integrator.
getIntegratorTimer
()->totalElapsedTime();
108
const
Teuchos::RCP<Teuchos::FancyOStream> out = integrator.getOStream();
109
out->setOutputToRootOnly(0);
110
Teuchos::OSTab ostab(out,0,
"ScreenOutput"
);
111
*out <<
"============================================================================\n"
112
<<
" Total runtime = "
<< runtime <<
" sec = "
113
<< runtime/60.0 <<
" min\n"
114
<< std::asctime(std::localtime(&end))
115
<<
"\nNumber of Accepted Steps = "
<< integrator.
getSolutionHistory
()->getCurrentState()->getIndex()
116
<<
"\nNumber of Failures = "
<< integrator.
getSolutionHistory
()->getCurrentState()->getMetaData()->getNRunningFailures()
117
<<
"\n"
118
<< exitStatus <<
"\n"
119
<< std::endl;
120
}
121
122
}
// namespace Tempus
123
#endif
// Tempus_IntegratorObserverBasic_impl_hpp
Tempus::IntegratorObserverBasic::~IntegratorObserverBasic
virtual ~IntegratorObserverBasic()
Destructor.
Definition
Tempus_IntegratorObserverBasic_impl.hpp:20
Tempus::IntegratorObserverBasic::observeEndIntegrator
virtual void observeEndIntegrator(const Integrator< Scalar > &integrator) override
Observe the end of the time integrator.
Definition
Tempus_IntegratorObserverBasic_impl.hpp:96
Tempus::IntegratorObserverBasic::observeEndTimeStep
virtual void observeEndTimeStep(const Integrator< Scalar > &integrator) override
Observe the end of the time step loop.
Definition
Tempus_IntegratorObserverBasic_impl.hpp:64
Tempus::IntegratorObserverBasic::observeStartIntegrator
virtual void observeStartIntegrator(const Integrator< Scalar > &integrator) override
Observe the beginning of the time integrator.
Definition
Tempus_IntegratorObserverBasic_impl.hpp:24
Tempus::IntegratorObserverBasic::observeBeforeTakeStep
virtual void observeBeforeTakeStep(const Integrator< Scalar > &integrator) override
Observe before Stepper takes step.
Definition
Tempus_IntegratorObserverBasic_impl.hpp:52
Tempus::IntegratorObserverBasic::observeNextTimeStep
virtual void observeNextTimeStep(const Integrator< Scalar > &integrator) override
Observe after the next time step size is selected.
Definition
Tempus_IntegratorObserverBasic_impl.hpp:48
Tempus::IntegratorObserverBasic::IntegratorObserverBasic
IntegratorObserverBasic()
Constructor.
Definition
Tempus_IntegratorObserverBasic_impl.hpp:17
Tempus::IntegratorObserverBasic::observeAfterTakeStep
virtual void observeAfterTakeStep(const Integrator< Scalar > &integrator) override
Observe after Stepper takes step.
Definition
Tempus_IntegratorObserverBasic_impl.hpp:56
Tempus::IntegratorObserverBasic::observeStartTimeStep
virtual void observeStartTimeStep(const Integrator< Scalar > &integrator) override
Observe the beginning of the time step loop.
Definition
Tempus_IntegratorObserverBasic_impl.hpp:44
Tempus::IntegratorObserverBasic::observeAfterCheckTimeStep
virtual void observeAfterCheckTimeStep(const Integrator< Scalar > &integrator) override
Observe after checking time step. Observer can still fail the time step here.
Definition
Tempus_IntegratorObserverBasic_impl.hpp:60
Tempus::Integrator
Thyra Base interface for time integrators. Time integrators are designed to advance the solution from...
Definition
Tempus_Integrator.hpp:65
Tempus::Integrator::getIntegratorTimer
virtual Teuchos::RCP< Teuchos::Time > getIntegratorTimer() const =0
Returns the IntegratorTimer_ for this Integrator.
Tempus::Integrator::getSolutionHistory
virtual Teuchos::RCP< const SolutionHistory< Scalar > > getSolutionHistory() const =0
Returns the SolutionHistory for this Integrator.
Tempus::Integrator::getTimeStepControl
virtual Teuchos::RCP< const TimeStepControl< Scalar > > getTimeStepControl() const =0
Returns the TimeStepControl for this Integrator.
Tempus::Integrator::getStepper
virtual Teuchos::RCP< Stepper< Scalar > > getStepper() const =0
Get the stepper.
Tempus::Integrator::getStatus
virtual Tempus::Status getStatus() const =0
Get the Status.
Tempus::Integrator::getStepperTimer
virtual Teuchos::RCP< Teuchos::Time > getStepperTimer() const =0
Tempus
Definition
Tempus_AdjointAuxSensitivityModelEvaluator_decl.hpp:21
Tempus::FAILED
@ FAILED
Definition
Tempus_Types.hpp:21
Generated by
1.17.0