VertexCFD  0.0-dev
VertexCFD_TempusObserver_IterationOutput_impl.hpp
1 #ifndef VERTEXCFD_TEMPUSOBSERVER_ITERATIONOUTPUT_IMPL_HPP
2 #define VERTEXCFD_TEMPUSOBSERVER_ITERATIONOUTPUT_IMPL_HPP
3 
4 #include <cmath>
5 #include <string>
6 #include <vector>
7 
8 namespace VertexCFD
9 {
10 namespace TempusObserver
11 {
12 //---------------------------------------------------------------------------//
13 template<class Scalar>
14 IterationOutput<Scalar>::IterationOutput(
15  Teuchos::RCP<TempusTimeStepControl::Strategy<Scalar>> dt_strategy)
16  : _ostream(Teuchos::rcp(&std::cout, false))
17  , _dt_strategy(dt_strategy)
18 {
19  _ostream.setShowProcRank(false);
20  _ostream.setOutputToRootOnly(0);
21 }
22 
23 //---------------------------------------------------------------------------//
24 template<class Scalar>
26  const Tempus::Integrator<Scalar>& integrator)
27 {
28  const std::time_t begin = std::time(nullptr);
29  _ostream << "\n==========================================================="
30  "=================\n"
31  << "Time Integration Begin\n"
32  << std::asctime(std::localtime(&begin))
33  << "\n Stepper = " << integrator.getStepper()->description()
34  << "\n Simulation Time Range ["
35  << integrator.getTimeStepControl()->getInitTime() << ", "
36  << integrator.getTimeStepControl()->getFinalTime() << "]"
37  << "\n-----------------------------------------------------------"
38  "-----------------\n";
39 
40  auto tsc = integrator.getTimeStepControl();
41  const auto final_time = tsc->getFinalTime();
42  const auto dt_min = tsc->getMinTimeStep();
43 
44  // Get location of first significant digit of a double
45  auto get_digit = [](double t) -> int {
46  return t > 0.0 ? std::floor(std::log10(t)) : 0;
47  };
48 
49  const int t_digit = get_digit(final_time);
50  const int dt_digit = get_digit(dt_min);
51 
52  // Make sure we always get meaningful time precision
53  _time_precision = t_digit - dt_digit + 2;
54 }
55 
56 //---------------------------------------------------------------------------//
57 template<class Scalar>
59  const Tempus::Integrator<Scalar>& /*integrator*/)
60 {
61 }
62 
63 //---------------------------------------------------------------------------//
64 template<class Scalar>
66  const Tempus::Integrator<Scalar>& /*integrator*/)
67 {
68 }
69 
70 //---------------------------------------------------------------------------//
71 template<class Scalar>
73  const Tempus::Integrator<Scalar>& integrator)
74 {
75  const auto current_time = integrator.getTime();
76  const auto state = integrator.getSolutionHistory()->getWorkingState();
77  _ostream << "\nTime Step = " << state->getIndex();
78  _ostream << "; Order = " << state->getOrder();
79  _ostream << "\nCFL = " << std::setprecision(3) << std::scientific
80  << _dt_strategy->currentCFL();
81  _ostream << "; dt = " << std::setprecision(3) << std::scientific
82  << state->getTimeStep();
83  _ostream << "; Time = " << std::setprecision(_time_precision)
84  << std::scientific << current_time << "\n";
85 }
86 
87 //---------------------------------------------------------------------------//
88 template<class Scalar>
90  const Tempus::Integrator<Scalar>& /*integrator*/)
91 {
92 }
93 
94 //---------------------------------------------------------------------------//
95 template<class Scalar>
97  const Tempus::Integrator<Scalar>& /*integrator*/)
98 {
99 }
100 
101 //---------------------------------------------------------------------------//
102 template<class Scalar>
104  const Tempus::Integrator<Scalar>& integrator)
105 {
106  const auto stepper_time = integrator.getStepperTimer()->totalElapsedTime();
107  integrator.getStepperTimer()->reset();
108  _ostream << "Time step time to completion (s): " << std::setprecision(2)
109  << std::scientific << stepper_time;
110  _ostream << "\n";
111 }
112 
113 //---------------------------------------------------------------------------//
114 template<class Scalar>
116  const Tempus::Integrator<Scalar>& integrator)
117 {
118  std::string exit_status;
119  if (integrator.getSolutionHistory()->getCurrentState()->getSolutionStatus()
120  == Tempus::Status::FAILED
121  or integrator.getStatus() == Tempus::Status::FAILED)
122  {
123  exit_status = "Time integration FAILURE!";
124  }
125  else
126  {
127  exit_status = "Time integration complete.";
128  }
129  const std::time_t end = std::time(nullptr);
130  const auto runtime_sec
131  = integrator.getIntegratorTimer()->totalElapsedTime();
132  const auto runtime_min = runtime_sec / 60;
133  const auto runtime_hr = runtime_min / 60;
134  _ostream << "\n-----------------------------------------------------------"
135  "-----------------\n"
136  << "Total runtime = " << std::setprecision(2) << std::scientific
137  << runtime_sec << " sec\n"
138  << " = " << std::setprecision(2) << std::scientific
139  << runtime_min << " min\n"
140  << " = " << std::setprecision(2) << std::fixed
141  << runtime_hr << " hr\n"
142  << std::asctime(std::localtime(&end)) << exit_status
143  << "\n==========================================================="
144  "=================\n"
145  << "\n";
146 }
147 
148 //---------------------------------------------------------------------------//
149 
150 } // end namespace TempusObserver
151 } // end namespace VertexCFD
152 
153 #endif // end VERTEXCFD_TEMPUSOBSERVER_ITERATIONOUTPUT_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23
VertexCFD::TempusObserver::IterationOutput::observeBeforeTakeStep
void observeBeforeTakeStep(const Tempus::Integrator< Scalar > &integrator) override
Observe before Stepper takes step.
Definition: VertexCFD_TempusObserver_IterationOutput_impl.hpp:72
VertexCFD::TempusObserver::IterationOutput::observeAfterTakeStep
void observeAfterTakeStep(const Tempus::Integrator< Scalar > &integrator) override
Observe after Stepper takes step.
Definition: VertexCFD_TempusObserver_IterationOutput_impl.hpp:89
VertexCFD::TempusObserver::IterationOutput::observeStartIntegrator
void observeStartIntegrator(const Tempus::Integrator< Scalar > &integrator) override
Observe the beginning of the time integrator.
Definition: VertexCFD_TempusObserver_IterationOutput_impl.hpp:25
VertexCFD::TempusObserver::IterationOutput::observeEndTimeStep
void observeEndTimeStep(const Tempus::Integrator< Scalar > &integrator) override
Observe the end of the time step loop.
Definition: VertexCFD_TempusObserver_IterationOutput_impl.hpp:103
VertexCFD::TempusObserver::IterationOutput::observeEndIntegrator
void observeEndIntegrator(const Tempus::Integrator< Scalar > &integrator) override
Observe the end of the time integrator.
Definition: VertexCFD_TempusObserver_IterationOutput_impl.hpp:115
VertexCFD::TempusObserver::IterationOutput::observeStartTimeStep
void observeStartTimeStep(const Tempus::Integrator< Scalar > &integrator) override
Observe the beginning of the time step loop.
Definition: VertexCFD_TempusObserver_IterationOutput_impl.hpp:58
VertexCFD::TempusObserver::IterationOutput::observeNextTimeStep
void observeNextTimeStep(const Tempus::Integrator< Scalar > &integrator) override
Definition: VertexCFD_TempusObserver_IterationOutput_impl.hpp:65
VertexCFD::TempusObserver::IterationOutput::observeAfterCheckTimeStep
void observeAfterCheckTimeStep(const Tempus::Integrator< Scalar > &integrator) override
Definition: VertexCFD_TempusObserver_IterationOutput_impl.hpp:96