Go to the documentation of this file.
19 #include "libmesh/twostep_time_solver.h"
20 #include "libmesh/diff_system.h"
21 #include "libmesh/euler_solver.h"
22 #include "libmesh/numeric_vector.h"
23 #include "libmesh/auto_ptr.h"
55 bool max_tolerance_met =
false;
58 Real single_norm(0.), double_norm(0.), error_norm(0.),
61 while (!max_tolerance_met)
70 libMesh::out <<
"\n === Computing adaptive timestep === "
80 std::unique_ptr<NumericVector<Number>> double_solution =
82 std::unique_ptr<NumericVector<Number>> old_solution =
88 libMesh::out <<
"Double norm = " << double_norm << std::endl;
112 libMesh::out <<
"Single norm = " << single_norm << std::endl;
130 std::max(double_norm, single_norm);
133 if (!double_norm && !single_norm)
138 libMesh::out <<
"Error norm = " << error_norm << std::endl;
141 std::max(double_norm, single_norm))
145 std::max(double_norm, single_norm))
175 max_tolerance_met =
true;
188 const Real global_shrink_or_growth_factor =
192 const Real local_shrink_or_growth_factor =
194 (error_norm/std::max(double_norm, single_norm)),
200 << global_shrink_or_growth_factor << std::endl;
202 << local_shrink_or_growth_factor << std::endl;
211 Real shrink_or_growth_factor =
213 local_shrink_or_growth_factor;
219 libMesh::out <<
"delta t is constrained by max_growth" << std::endl;
231 libMesh::out <<
"delta t is constrained by maximum-allowable delta t."
242 libMesh::out <<
"delta t is constrained by minimum-allowable delta t."
TwostepTimeSolver(sys_type &s)
Constructor.
Real deltat
For time-dependent problems, this is the amount delta t to advance the solution in time.
virtual Real calculate_norm(System &, NumericVector< Number > &)
A helper function to calculate error norms.
unsigned int reduce_deltat_on_diffsolver_failure
This value (which defaults to zero) is the number of times the TimeSolver is allowed to halve deltat ...
bool quiet
Print extra debugging information if quiet == false.
The libMesh namespace provides an interface to certain functionality in the library.
virtual std::unique_ptr< NumericVector< T > > clone() const =0
Real upper_tolerance
This tolerance is the maximum relative error between an exact time integration and a single time step...
This class provides a specific system class.
~TwostepTimeSolver()
Destructor.
std::unique_ptr< UnsteadySolver > core_time_solver
This object is used to take timesteps.
Real last_deltat
We need to store the value of the last deltat used, so that advance_timestep() will increment the sys...
bool global_tolerance
This flag, which is true by default, grows (shrinks) the timestep based on the expected global accura...
Real max_growth
Do not allow the adaptive time solver to select a new deltat greater than max_growth times the old de...
double pow(double a, int b)
sys_type & _system
A reference to the system we are solving.
Real time
For time-dependent problems, this is the time t at the beginning of the current timestep.
bool first_solve
A bool that will be true the first time solve() is called, and false thereafter.
This class wraps another UnsteadySolver derived class, and compares the results of timestepping with ...
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
virtual void solve() override
This method solves for the solution at the next timestep.
Real max_deltat
Do not allow the adaptive time solver to select deltat > max_deltat.
Real min_deltat
Do not allow the adaptive time solver to select deltat < min_deltat.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const NumericVector< Number > & get_vector(const std::string &vec_name) const
Real target_tolerance
This tolerance is the target relative error between an exact time integration and a single time step ...