This class wraps another UnsteadySolver derived class, and compares the results of timestepping with deltat and timestepping with 2*deltat to adjust future timestep lengths.
Currently this class only works on fully coupled Systems
This class is part of the new DifferentiableSystem framework, which is still experimental. Users of this framework should beware of bugs and future API changes.
Author:
Roy H. Stogner 2007
Definition at line 49 of file adaptive_time_solver.h.
Member Typedef Documentation
typedef std::map<std::string, std::pair<unsigned int, unsigned int> > ReferenceCounter::Counts [protected, inherited]Data structure to log the information. The log is identified by the class name.
Definition at line 105 of file reference_counter.h.
typedef UnsteadySolverAdaptiveTimeSolver::ParentThe parent class
Reimplemented in TwostepTimeSolver.
Definition at line 55 of file adaptive_time_solver.h.
typedef DifferentiableSystemTimeSolver::sys_type [inherited]The type of system
Reimplemented in EigenTimeSolver, and SteadySolver.
Definition at line 62 of file time_solver.h.
Constructor & Destructor Documentation
AdaptiveTimeSolver::AdaptiveTimeSolver (sys_type &s)Constructor. Requires a reference to the system to be solved.
Definition at line 27 of file adaptive_time_solver.C.
References UnsteadySolver::old_local_nonlinear_solution, and AutoPtr< Tp >::reset().
: UnsteadySolver(s),
core_time_solver(NULL),
target_tolerance(1.e-3), upper_tolerance(0.0),
max_deltat(0.),
min_deltat(0.),
max_growth(0.),
global_tolerance(true)
{
// the child class must populate core_time_solver
// with whatever actual time solver is to be used
// As an UnsteadySolver, we have an old_local_nonlinear_solution, but we're
// going to drop it and use our core time solver's instead.
old_local_nonlinear_solution.reset();
}
Definition at line 46 of file adaptive_time_solver.C.
References UnsteadySolver::old_local_nonlinear_solution, and AutoPtr< Tp >::release().
{
// As an UnsteadySolver, we have an old_local_nonlinear_solution, but it
// is being managed by our core_time_solver. Make sure we don't delete
// it out from under them, in case the user wants to keep using the core
// solver after they're done with us.
old_local_nonlinear_solution.release();
}
Member Function Documentation
void TimeSolver::adjoint_recede_timestep () [virtual, inherited]This method advances the adjoint solution to the previous timestep, after an adjoint_solve() has been performed. This will probably be done after every UnsteadySolver::adjoint_solve().
Definition at line 90 of file time_solver.C.
{
}
void TimeSolver::adjoint_solve () [virtual, inherited]This method solves for the adjoint solution at the previous timestep (or solves a steady adjoint problem). Usually we will only need to solve one linear system per timestep, but more complex subclasses may override this.
Definition at line 77 of file time_solver.C.
References TimeSolver::_adjoint_solver.
{
_adjoint_solver->adjoint_solve();
}
virtual AutoPtr<DiffSolver>& TimeSolver::adjoint_solver () [inline, virtual, inherited]An implicit linear solver to use for adjoint problems.
Definition at line 157 of file time_solver.h.
References TimeSolver::_adjoint_solver.
{ return _adjoint_solver; }
void AdaptiveTimeSolver::advance_timestep () [virtual]This method advances the solution to the next timestep, after a solve() has been performed. Often this will be done after every UnsteadySolver::solve(), but adaptive mesh refinement and/or adaptive time step selection may require some solve() steps to be repeated.
Reimplemented from UnsteadySolver.
Definition at line 86 of file adaptive_time_solver.C.
References TimeSolver::_system, UnsteadySolver::first_solve, System::get_vector(), last_deltat, UnsteadySolver::old_nonlinear_solution(), System::solution, and DifferentiableSystem::time.
virtual void TimeSolver::before_timestep () [inline, virtual, inherited]This method is for subclasses or users to override to do arbitrary processing between timesteps
Definition at line 142 of file time_solver.h.
{}
Real AdaptiveTimeSolver::calculate_norm (System &s, NumericVector< Number > &v) [protected, virtual]A helper function to calculate error norms
Definition at line 138 of file adaptive_time_solver.C.
References System::calculate_norm(), and component_norm.
Referenced by TwostepTimeSolver::solve().
{
return s.calculate_norm(v, component_norm);
}
AutoPtr< DiffSolver > & AdaptiveTimeSolver::diff_solver () [virtual]An implicit linear or nonlinear solver to use at each timestep.
Reimplemented from TimeSolver.
Definition at line 131 of file adaptive_time_solver.C.
References core_time_solver.
{
return core_time_solver->diff_solver();
}
Real UnsteadySolver::du (const SystemNorm &norm) const [virtual, inherited]Computes the size of ||u^{n+1} - u^{n}|| in some norm.
Note that, while you can always call this function, its result may or may not be very meaningful. For example, if you call this function right after calling advance_timestep() then you'll get a result of zero since old_nonlinear_solution is set equal to nonlinear_solution in this function.
Implements TimeSolver.
Definition at line 150 of file unsteady_solver.C.
References TimeSolver::_system, System::calculate_norm(), System::get_vector(), and System::solution.
std::string ReferenceCounter::get_info () [static, inherited]Gets a string containing the reference information.
Definition at line 45 of file reference_counter.C.
References ReferenceCounter::_counts, and Quality::name().
Referenced by ReferenceCounter::print_info().
{
#if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
std::ostringstream out;
out << '
<< ' ----------------------------------------------------------------------------
<< '| Reference count information |
<< ' ---------------------------------------------------------------------------- ;
for (Counts::iterator it = _counts.begin();
it != _counts.end(); ++it)
{
const std::string name(it->first);
const unsigned int creations = it->second.first;
const unsigned int destructions = it->second.second;
out << '| ' << name << ' reference count information:
<< '| Creations: ' << creations << '
<< '| Destructions: ' << destructions << ';
}
out << ' ---------------------------------------------------------------------------- ;
return out.str();
#else
return '';
#endif
}
void ReferenceCounter::increment_constructor_count (const std::string &name) [inline, protected, inherited]Increments the construction counter. Should be called in the constructor of any derived class that will be reference counted.
Definition at line 149 of file reference_counter.h.
References ReferenceCounter::_counts, Quality::name(), and Threads::spin_mtx.
Referenced by ReferenceCountedObject< Value >::ReferenceCountedObject().
void ReferenceCounter::increment_destructor_count (const std::string &name) [inline, protected, inherited]Increments the destruction counter. Should be called in the destructor of any derived class that will be reference counted.
Definition at line 167 of file reference_counter.h.
References ReferenceCounter::_counts, Quality::name(), and Threads::spin_mtx.
Referenced by ReferenceCountedObject< Value >::~ReferenceCountedObject().
void AdaptiveTimeSolver::init () [virtual]The initialization function. This method is used to initialize internal data structures before a simulation begins.
Reimplemented from UnsteadySolver.
Definition at line 57 of file adaptive_time_solver.C.
References core_time_solver, AutoPtr< Tp >::get(), and UnsteadySolver::old_local_nonlinear_solution.
{
libmesh_assert(core_time_solver.get());
// We override this because our core_time_solver is the one that
// needs to handle new vectors, diff_solver->init(), etc
core_time_solver->init();
// As an UnsteadySolver, we have an old_local_nonlinear_solution, but it
// isn't pointing to the right place - fix it
//
// This leaves us with two AutoPtrs holding the same pointer - dangerous
// for future use. Replace with shared_ptr?
old_local_nonlinear_solution =
AutoPtr<NumericVector<Number> >(core_time_solver->old_local_nonlinear_solution.get());
}
static unsigned int ReferenceCounter::n_objects () [inline, static, inherited]Prints the number of outstanding (created, but not yet destroyed) objects.
Definition at line 76 of file reference_counter.h.
References ReferenceCounter::_n_objects.
Referenced by System::read_serialized_blocked_dof_objects(), and System::write_serialized_blocked_dof_objects().
{ return _n_objects; }
Number UnsteadySolver::old_nonlinear_solution (const unsigned intglobal_dof_number) const [inherited]Returns:
the old nonlinear solution for the specified global DOF.
Definition at line 139 of file unsteady_solver.C.
References TimeSolver::_system, System::get_dof_map(), DofMap::n_dofs(), and UnsteadySolver::old_local_nonlinear_solution.
Referenced by UnsteadySolver::advance_timestep(), advance_timestep(), EulerSolver::element_residual(), Euler2Solver::element_residual(), FEMSystem::eulerian_residual(), EulerSolver::side_residual(), and Euler2Solver::side_residual().
void AdaptiveTimeSolver::reinit () [virtual]The reinitialization function. This method is used after changes in the mesh
Reimplemented from TimeSolver.
Definition at line 76 of file adaptive_time_solver.C.
References core_time_solver, and AutoPtr< Tp >::get().
{
libmesh_assert(core_time_solver.get());
// We override this because our core_time_solver is the one that
// needs to handle new vectors, diff_solver->reinit(), etc
core_time_solver->reinit();
}
bool AdaptiveTimeSolver::side_residual (boolget_jacobian, DiffContext &context) [virtual]This method is passed on to the core_time_solver
Implements TimeSolver.
Definition at line 121 of file adaptive_time_solver.C.
References core_time_solver, and AutoPtr< Tp >::get().
virtual void AdaptiveTimeSolver::solve () [pure virtual]This method solves for the solution at the next timestep. Usually we will only need to solve one (non)linear system per timestep, but more complex subclasses may override this.
Definition at line 123 of file reference_counter.h.
Threads::atomic< unsigned int > ReferenceCounter::_n_objects [static, protected, inherited]The number of objects. Print the reference count information when the number returns to 0.
Definition at line 118 of file reference_counter.h.
Referenced by ReferenceCounter::n_objects(), ReferenceCounter::ReferenceCounter(), and ReferenceCounter::~ReferenceCounter().
sys_type& TimeSolver::_system [protected, inherited]A reference to the system we are solving.
SystemNormAdaptiveTimeSolver::component_normError calculations are done in this norm, DISCRETE_L2 by default.
Definition at line 106 of file adaptive_time_solver.h.
Referenced by calculate_norm().
std::vector<float> AdaptiveTimeSolver::component_scaleIf component_norms is non-empty, each variable's contribution to the error of a system will also be scaled by component_scale[var], unless component_scale is empty in which case all variables will be weighted equally.
Definition at line 114 of file adaptive_time_solver.h.
AutoPtr<UnsteadySolver> AdaptiveTimeSolver::core_time_solverThis object is used to take timesteps
Definition at line 101 of file adaptive_time_solver.h.
Referenced by diff_solver(), element_residual(), error_order(), init(), reinit(), side_residual(), TwostepTimeSolver::solve(), and TwostepTimeSolver::TwostepTimeSolver().
bool UnsteadySolver::first_solve [protected, inherited]A bool that will be true the first time solve() is called, and false thereafter
Reimplemented from TimeSolver.
Definition at line 124 of file unsteady_solver.h.
Referenced by UnsteadySolver::advance_timestep(), advance_timestep(), UnsteadySolver::solve(), and TwostepTimeSolver::solve().
bool AdaptiveTimeSolver::global_toleranceThis flag, which is true by default, grows (shrinks) the timestep based on the expected global accuracy of the timestepping scheme. Global in this sense means the cumulative final-time accuracy of the scheme. For example, the backward Euler scheme's truncation error is locally of order 2, so that after N timesteps of size deltat, the result is first-order accurate. If you set this to false, you can grow (shrink) your timestep based on the local accuracy rather than the global accuracy of the core TimeSolver. Note that by setting this value to false you may fail to achieve the predicted convergence in time of the underlying method, however it may be possible to get more fine-grained control over step sizes as well.
Definition at line 170 of file adaptive_time_solver.h.
Referenced by TwostepTimeSolver::solve().
RealAdaptiveTimeSolver::last_deltat [protected]We need to store the value of the last deltat used, so that advance_timestep() will increment the system time correctly.
Definition at line 179 of file adaptive_time_solver.h.
Referenced by advance_timestep(), and TwostepTimeSolver::solve().
RealAdaptiveTimeSolver::max_deltatDo not allow the adaptive time solver to select deltat > max_deltat. If you use the default max_deltat=0.0, then deltat is unlimited.
Definition at line 140 of file adaptive_time_solver.h.
Referenced by TwostepTimeSolver::solve().
RealAdaptiveTimeSolver::max_growthDo not allow the adaptive time solver to select a new deltat greater than max_growth times the old deltat. If you use the default max_growth=0.0, then the deltat growth is unlimited.
Definition at line 154 of file adaptive_time_solver.h.
Referenced by TwostepTimeSolver::solve().
RealAdaptiveTimeSolver::min_deltatDo not allow the adaptive time solver to select deltat < min_deltat. The default value is 0.0.
Definition at line 146 of file adaptive_time_solver.h.
Referenced by TwostepTimeSolver::solve().
AutoPtr<NumericVector<Number> > UnsteadySolver::old_local_nonlinear_solution [inherited]Serial vector of _system.get_vector('_old_nonlinear_solution')
Reimplemented from TimeSolver.
Definition at line 105 of file unsteady_solver.h.
Referenced by AdaptiveTimeSolver(), init(), UnsteadySolver::old_nonlinear_solution(), UnsteadySolver::solve(), and ~AdaptiveTimeSolver().
bool TimeSolver::quiet [inherited]Print extra debugging information if quiet == false.
Definition at line 162 of file time_solver.h.
Referenced by UnsteadySolver::solve(), TwostepTimeSolver::solve(), and EigenTimeSolver::solve().
unsigned int TimeSolver::reduce_deltat_on_diffsolver_failure [inherited]This value (which defaults to zero) is the number of times the TimeSolver is allowed to halve deltat and let the DiffSolver repeat the latest failed solve with a reduced timestep. Note that this has no effect for SteadySolvers. Note that you must set at least one of the DiffSolver flags 'continue_after_max_iterations' or 'continue_after_backtrack_failure' to allow the TimeSolver to retry the solve.
Definition at line 185 of file time_solver.h.
Referenced by UnsteadySolver::solve(), and TwostepTimeSolver::solve().
RealAdaptiveTimeSolver::target_toleranceThis tolerance is the target relative error between double-deltat and single-deltat timesteps, scaled by deltat. If this error tolerance is exceeded or not met, future timesteps will be run with deltat shrunk or grown to compensate.
The default value is 1.0e-2; obviously users should select their own tolerance.
Definition at line 125 of file adaptive_time_solver.h.
Referenced by TwostepTimeSolver::solve().
RealAdaptiveTimeSolver::upper_toleranceThis tolerance is the maximum relative error between double-deltat and single-deltat timesteps, scaled by deltat. If this error tolerance is exceeded, the current timestep will be repeated with a smaller deltat.
If you use the default upper_tolerance=0.0,
Definition at line 134 of file adaptive_time_solver.h.
Referenced by TwostepTimeSolver::solve().
Author
Generated automatically by Doxygen for libMesh from the source code.