1 #ifndef VERTEXCFD_TEMPUSOBSERVER_RESPONSEOUTPUT_IMPL_HPP
2 #define VERTEXCFD_TEMPUSOBSERVER_RESPONSEOUTPUT_IMPL_HPP
10 namespace TempusObserver
13 template<
class Scalar>
14 ResponseOutput<Scalar>::ResponseOutput(
15 Teuchos::RCP<Response::ResponseManager> response_manager,
16 std::vector<int> output_freq)
17 : _ostream(Teuchos::rcp(&std::cout, false))
18 , _response_manager(response_manager)
19 , _output_freq(std::move(output_freq))
21 _ostream.setShowProcRank(
false);
22 _ostream.setOutputToRootOnly(0);
26 template<
class Scalar>
28 const Tempus::Integrator<Scalar>& integrator)
32 outputResponses(integrator, integrator.getIndex());
36 template<
class Scalar>
38 const Tempus::Integrator<Scalar>& )
43 template<
class Scalar>
45 const Tempus::Integrator<Scalar>& )
50 template<
class Scalar>
52 const Tempus::Integrator<Scalar>& )
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)
76 outputResponses(integrator, integrator.getIndex());
80 template<
class Scalar>
82 const Tempus::Integrator<Scalar>& integrator)
85 outputResponses(integrator);
89 template<
class Scalar>
91 const Tempus::Integrator<Scalar>& integrator,
const int current_index)
93 const int num_resp = _response_manager->numResponses();
99 _response_manager->deactivateAll();
103 for (
int i = 0; i < num_resp; ++i)
105 if (0 == current_index % _output_freq[i])
107 _response_manager->activateResponse(i);
113 if (num_outputs == 0)
117 const auto state = integrator.getSolutionHistory()->getCurrentState();
118 _response_manager->evaluateResponses(state->getX(), state->getXDot());
121 _ostream <<
"Scalar Responses:\n";
122 for (
int i = 0; i < num_resp; ++i)
124 if (0 == current_index % _output_freq[i])
126 const auto& name = _response_manager->name(i);
127 const auto value = _response_manager->value(i);
129 constexpr
int prec = std::numeric_limits<double>::digits10 + 1;
131 _ostream <<
" " << name <<
" = " << std::setprecision(prec)
141 #endif // VERTEXCFD_TEMPUSOBSERVER_RESPONSEOUTPUT_IMPL_HPP