1 #ifndef VERTEXCFD_TEMPUSTIMESTEPCONTROL_GLOBALTIMESTEP_IMPL_HPP
2 #define VERTEXCFD_TEMPUSTIMESTEPCONTROL_GLOBALTIMESTEP_IMPL_HPP
4 #include "VertexCFD_TempusTimeStepControl_GlobalTimeStep.hpp"
9 namespace TempusTimeStepControl
12 template<
class Scalar>
13 GlobalTimeStep<Scalar>::GlobalTimeStep(
14 const Teuchos::ParameterList& user_params,
15 Teuchos::RCP<PhysicsManager> physics_manager)
16 : _dt_transition_steps(user_params.isType<int>(
"Time step transition "
18 ? user_params.get<int>(
"Time step transition "
21 , _response_manager(physics_manager)
23 _response_manager.addMinValueResponse(
"global_cfl_time_step",
"local_dt");
28 this->stepType_ =
"Variable";
29 this->strategyType_ =
"Global Time Step Strategy";
30 this->name_ =
"Global Time Step Strategy";
35 template<
class Scalar>
36 void GlobalTimeStep<Scalar>::setNextTimeStep(
37 const Tempus::TimeStepControl<Scalar>& tsc,
38 Teuchos::RCP<Tempus::SolutionHistory<Scalar>> solution_history,
42 auto working_state = solution_history->getWorkingState();
45 _response_manager.evaluateResponses(working_state->getX(),
46 working_state->getXDot());
47 const double dt_cfl1 = _response_manager.value();
50 const auto dt_index = working_state->getIndex() - 1;
52 = _dt_transition_steps > 0
53 ? std::min(
static_cast<double>(dt_index) / _dt_transition_steps,
58 = std::max(tsc.getMinTimeStep(), tsc.getInitTimeStep());
59 const double dt_final = tsc.getMaxTimeStep();
60 const double dt = (1.0 - wt) * dt_init + wt * dt_final;
62 working_state->setTimeStep(dt);
63 working_state->setTime(solution_history->getCurrentState()->getTime() + dt);
66 this->setCurrentCFL(dt / dt_cfl1);
74 #endif // VERTEXCFD_TEMPUSTIMESTEPCONTROL_GLOBALTIMESTEP_IMPL_HPP