VertexCFD  0.0-dev
VertexCFD_TempusTimeStepControl_GlobalTimeStep_impl.hpp
1 #ifndef VERTEXCFD_TEMPUSTIMESTEPCONTROL_GLOBALTIMESTEP_IMPL_HPP
2 #define VERTEXCFD_TEMPUSTIMESTEPCONTROL_GLOBALTIMESTEP_IMPL_HPP
3 
4 #include "VertexCFD_TempusTimeStepControl_GlobalTimeStep.hpp"
5 #include <algorithm>
6 
7 namespace VertexCFD
8 {
9 namespace TempusTimeStepControl
10 {
11 //---------------------------------------------------------------------------//
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 "
17  "steps")
18  ? user_params.get<int>("Time step transition "
19  "steps")
20  : 0)
21  , _response_manager(physics_manager)
22 {
23  _response_manager.addMinValueResponse("global_cfl_time_step", "local_dt");
24 
25  // For the new Tempus::TimeStepControlStrategy interface, we need to
26  // set a few base class member variables. In particular, incorrect
27  // behavior may result if "stepType_" is not set to "Variable".
28  this->stepType_ = "Variable";
29  this->strategyType_ = "Global Time Step Strategy";
30  this->name_ = "Global Time Step Strategy";
31 }
32 
33 //---------------------------------------------------------------------------//
34 // Determine the time step size.
35 template<class Scalar>
36 void GlobalTimeStep<Scalar>::setNextTimeStep(
37  const Tempus::TimeStepControl<Scalar>& tsc,
38  Teuchos::RCP<Tempus::SolutionHistory<Scalar>> solution_history,
39  Tempus::Status&)
40 {
41  // Get the working state.
42  auto working_state = solution_history->getWorkingState();
43 
44  // Get minimum time step that ensures CFL <= 1
45  _response_manager.evaluateResponses(working_state->getX(),
46  working_state->getXDot());
47  const double dt_cfl1 = _response_manager.value();
48 
49  // Get time step index (1-based) and compute linear weight
50  const auto dt_index = working_state->getIndex() - 1;
51  const auto wt
52  = _dt_transition_steps > 0
53  ? std::min(static_cast<double>(dt_index) / _dt_transition_steps,
54  1.0)
55  : 1.0;
56 
57  const double dt_init
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;
61 
62  working_state->setTimeStep(dt);
63  working_state->setTime(solution_history->getCurrentState()->getTime() + dt);
64 
65  // Save current CFL so it may be accessed elsewhere
66  this->setCurrentCFL(dt / dt_cfl1);
67 }
68 
69 //---------------------------------------------------------------------------//
70 
71 } // namespace TempusTimeStepControl
72 } // namespace VertexCFD
73 
74 #endif // VERTEXCFD_TEMPUSTIMESTEPCONTROL_GLOBALTIMESTEP_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23