1 #ifndef VERTEXCFD_TEMPUSOBSERVER_ITERATIONOUTPUT_IMPL_HPP
2 #define VERTEXCFD_TEMPUSOBSERVER_ITERATIONOUTPUT_IMPL_HPP
10 namespace TempusObserver
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)
19 _ostream.setShowProcRank(
false);
20 _ostream.setOutputToRootOnly(0);
24 template<
class Scalar>
26 const Tempus::Integrator<Scalar>& integrator)
28 const std::time_t begin = std::time(
nullptr);
29 _ostream <<
"\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";
40 auto tsc = integrator.getTimeStepControl();
41 const auto final_time = tsc->getFinalTime();
42 const auto dt_min = tsc->getMinTimeStep();
45 auto get_digit = [](
double t) ->
int {
46 return t > 0.0 ? std::floor(std::log10(t)) : 0;
49 const int t_digit = get_digit(final_time);
50 const int dt_digit = get_digit(dt_min);
53 _time_precision = t_digit - dt_digit + 2;
57 template<
class Scalar>
59 const Tempus::Integrator<Scalar>& )
64 template<
class Scalar>
66 const Tempus::Integrator<Scalar>& )
71 template<
class Scalar>
73 const Tempus::Integrator<Scalar>& integrator)
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";
88 template<
class Scalar>
90 const Tempus::Integrator<Scalar>& )
95 template<
class Scalar>
97 const Tempus::Integrator<Scalar>& )
102 template<
class Scalar>
104 const Tempus::Integrator<Scalar>& integrator)
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;
114 template<
class Scalar>
116 const Tempus::Integrator<Scalar>& integrator)
118 std::string exit_status;
119 if (integrator.getSolutionHistory()->getCurrentState()->getSolutionStatus()
120 == Tempus::Status::FAILED
121 or integrator.getStatus() == Tempus::Status::FAILED)
123 exit_status =
"Time integration FAILURE!";
127 exit_status =
"Time integration complete.";
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"
153 #endif // end VERTEXCFD_TEMPUSOBSERVER_ITERATIONOUTPUT_IMPL_HPP