$extrastylesheet
#include <adaptive_time_solver.h>

Public Types | |
| typedef UnsteadySolver | Parent |
| typedef DifferentiableSystem | sys_type |
Public Member Functions | |
| AdaptiveTimeSolver (sys_type &s) | |
| virtual | ~AdaptiveTimeSolver () |
| virtual void | init () |
| virtual void | reinit () |
| virtual void | solve ()=0 |
| virtual void | advance_timestep () |
| virtual Real | error_order () const |
| virtual bool | element_residual (bool get_jacobian, DiffContext &) |
| virtual bool | side_residual (bool get_jacobian, DiffContext &) |
| virtual bool | nonlocal_residual (bool get_jacobian, DiffContext &) |
| virtual UniquePtr< DiffSolver > & | diff_solver () |
| virtual UniquePtr < LinearSolver< Number > > & | linear_solver () |
| virtual void | init_data () |
| virtual void | adjoint_advance_timestep () |
| virtual void | retrieve_timestep () |
| Number | old_nonlinear_solution (const dof_id_type global_dof_number) const |
| virtual Real | du (const SystemNorm &norm) const |
| virtual bool | is_steady () const |
| virtual void | before_timestep () |
| const sys_type & | system () const |
| sys_type & | system () |
| void | set_solution_history (const SolutionHistory &_solution_history) |
| bool | is_adjoint () const |
| void | set_is_adjoint (bool _is_adjoint_value) |
Static Public Member Functions | |
| static std::string | get_info () |
| static void | print_info (std::ostream &out=libMesh::out) |
| static unsigned int | n_objects () |
| static void | enable_print_counter_info () |
| static void | disable_print_counter_info () |
Public Attributes | |
| UniquePtr< UnsteadySolver > | core_time_solver |
| SystemNorm | component_norm |
| std::vector< float > | component_scale |
| Real | target_tolerance |
| Real | upper_tolerance |
| Real | max_deltat |
| Real | min_deltat |
| Real | max_growth |
| bool | global_tolerance |
| UniquePtr< NumericVector < Number > > | old_local_nonlinear_solution |
| bool | quiet |
| unsigned int | reduce_deltat_on_diffsolver_failure |
Protected Types | |
| typedef bool(DifferentiablePhysics::* | ResFuncType )(bool, DiffContext &) |
| typedef void(DiffContext::* | ReinitFuncType )(Real) |
| typedef std::map< std::string, std::pair< unsigned int, unsigned int > > | Counts |
Protected Member Functions | |
| virtual Real | calculate_norm (System &, NumericVector< Number > &) |
| void | increment_constructor_count (const std::string &name) |
| void | increment_destructor_count (const std::string &name) |
Protected Attributes | |
| Real | last_deltat |
| bool | first_solve |
| bool | first_adjoint_step |
| UniquePtr< DiffSolver > | _diff_solver |
| UniquePtr< LinearSolver< Number > > | _linear_solver |
| sys_type & | _system |
| UniquePtr< SolutionHistory > | solution_history |
Static Protected Attributes | |
| static Counts | _counts |
| static Threads::atomic < unsigned int > | _n_objects |
| static Threads::spin_mutex | _mutex |
| static bool | _enable_print_counter = true |
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.
Definition at line 51 of file adaptive_time_solver.h.
typedef std::map<std::string, std::pair<unsigned int, unsigned int> > libMesh::ReferenceCounter::Counts [protected, inherited] |
Data structure to log the information. The log is identified by the class name.
Definition at line 113 of file reference_counter.h.
The parent class
Reimplemented in libMesh::TwostepTimeSolver.
Definition at line 57 of file adaptive_time_solver.h.
typedef void(DiffContext::* libMesh::TimeSolver::ReinitFuncType)(Real) [protected, inherited] |
Definition at line 278 of file time_solver.h.
typedef bool(DifferentiablePhysics::* libMesh::TimeSolver::ResFuncType)(bool, DiffContext &) [protected, inherited] |
Definitions of argument types for use in refactoring subclasses.
Definition at line 276 of file time_solver.h.
typedef DifferentiableSystem libMesh::TimeSolver::sys_type [inherited] |
The type of system
Reimplemented in libMesh::EigenTimeSolver, and libMesh::SteadySolver.
Definition at line 66 of file time_solver.h.
| libMesh::AdaptiveTimeSolver::AdaptiveTimeSolver | ( | sys_type & | s | ) | [explicit] |
Constructor. Requires a reference to the system to be solved.
Definition at line 27 of file adaptive_time_solver.C.
References libMesh::UnsteadySolver::old_local_nonlinear_solution.
: UnsteadySolver(s), core_time_solver(), 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(); }
| libMesh::AdaptiveTimeSolver::~AdaptiveTimeSolver | ( | ) | [virtual] |
Destructor.
Definition at line 47 of file adaptive_time_solver.C.
References libMesh::UnsteadySolver::old_local_nonlinear_solution.
{
// 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();
}
| void libMesh::UnsteadySolver::adjoint_advance_timestep | ( | ) | [virtual, inherited] |
This method advances the adjoint solution to the previous timestep, after an adjoint_solve() has been performed. This will be done before every UnsteadySolver::adjoint_solve().
Reimplemented from libMesh::TimeSolver.
Definition at line 178 of file unsteady_solver.C.
References libMesh::TimeSolver::_system, libMesh::DifferentiableSystem::deltat, libMesh::UnsteadySolver::first_adjoint_step, libMesh::System::get_dof_map(), libMesh::DofMap::get_send_list(), libMesh::System::get_vector(), libMesh::NumericVector< T >::localize(), libMesh::UnsteadySolver::old_local_nonlinear_solution, libMesh::TimeSolver::solution_history, and libMesh::System::time.
{
// On the first call of this function, we dont save the adjoint solution or
// decrement the time, we just call the retrieve function below
if(!first_adjoint_step)
{
// Call the store function to store the last adjoint before decrementing the time
solution_history->store();
// Decrement the system time
_system.time -= _system.deltat;
}
else
{
first_adjoint_step = false;
}
// Retrieve the primal solution vectors at this time using the
// solution_history object
solution_history->retrieve();
// Dont forget to localize the old_nonlinear_solution !
_system.get_vector("_old_nonlinear_solution").localize
(*old_local_nonlinear_solution,
_system.get_dof_map().get_send_list());
}
| void libMesh::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 libMesh::UnsteadySolver.
Definition at line 87 of file adaptive_time_solver.C.
References libMesh::TimeSolver::_system, libMesh::UnsteadySolver::first_solve, libMesh::System::get_vector(), last_deltat, libMesh::System::solution, and libMesh::System::time.
{
NumericVector<Number> &old_nonlinear_soln =
_system.get_vector("_old_nonlinear_solution");
NumericVector<Number> &nonlinear_solution =
*(_system.solution);
// _system.get_vector("_nonlinear_solution");
old_nonlinear_soln = nonlinear_solution;
if (!first_solve)
_system.time += last_deltat;
}
| virtual void libMesh::TimeSolver::before_timestep | ( | ) | [inline, virtual, inherited] |
This method is for subclasses or users to override to do arbitrary processing between timesteps
Definition at line 168 of file time_solver.h.
{}
| Real libMesh::AdaptiveTimeSolver::calculate_norm | ( | System & | s, |
| NumericVector< Number > & | v | ||
| ) | [protected, virtual] |
A helper function to calculate error norms
Definition at line 156 of file adaptive_time_solver.C.
References libMesh::System::calculate_norm(), and component_norm.
Referenced by libMesh::TwostepTimeSolver::solve().
{
return s.calculate_norm(v, component_norm);
}
| UniquePtr< DiffSolver > & libMesh::AdaptiveTimeSolver::diff_solver | ( | ) | [virtual] |
An implicit linear or nonlinear solver to use at each timestep.
Reimplemented from libMesh::TimeSolver.
Definition at line 142 of file adaptive_time_solver.C.
References core_time_solver.
{
return core_time_solver->diff_solver();
}
| void libMesh::ReferenceCounter::disable_print_counter_info | ( | ) | [static, inherited] |
Definition at line 106 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter.
Referenced by libMesh::LibMeshInit::LibMeshInit().
{
_enable_print_counter = false;
return;
}
| Real libMesh::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 libMesh::TimeSolver.
Definition at line 227 of file unsteady_solver.C.
References libMesh::TimeSolver::_system, libMesh::System::calculate_norm(), libMesh::System::get_vector(), and libMesh::System::solution.
{
UniquePtr<NumericVector<Number> > solution_copy =
_system.solution->clone();
solution_copy->add(-1., _system.get_vector("_old_nonlinear_solution"));
solution_copy->close();
return _system.calculate_norm(*solution_copy, norm);
}
| bool libMesh::AdaptiveTimeSolver::element_residual | ( | bool | get_jacobian, |
| DiffContext & | context | ||
| ) | [virtual] |
This method is passed on to the core_time_solver
Implements libMesh::TimeSolver.
Definition at line 112 of file adaptive_time_solver.C.
References core_time_solver, and libMesh::libmesh_assert().
{
libmesh_assert(core_time_solver.get());
return core_time_solver->element_residual(request_jacobian, context);
}
| void libMesh::ReferenceCounter::enable_print_counter_info | ( | ) | [static, inherited] |
Methods to enable/disable the reference counter output from print_info()
Definition at line 100 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter.
{
_enable_print_counter = true;
return;
}
| Real libMesh::AdaptiveTimeSolver::error_order | ( | ) | const [virtual] |
This method is passed on to the core_time_solver
Implements libMesh::UnsteadySolver.
Definition at line 103 of file adaptive_time_solver.C.
References core_time_solver, and libMesh::libmesh_assert().
{
libmesh_assert(core_time_solver.get());
return core_time_solver->error_order();
}
| std::string libMesh::ReferenceCounter::get_info | ( | ) | [static, inherited] |
Gets a string containing the reference information.
Definition at line 47 of file reference_counter.C.
References libMesh::ReferenceCounter::_counts, and libMesh::Quality::name().
Referenced by libMesh::ReferenceCounter::print_info().
{
#if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
std::ostringstream oss;
oss << '\n'
<< " ---------------------------------------------------------------------------- \n"
<< "| Reference count information |\n"
<< " ---------------------------------------------------------------------------- \n";
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;
oss << "| " << name << " reference count information:\n"
<< "| Creations: " << creations << '\n'
<< "| Destructions: " << destructions << '\n';
}
oss << " ---------------------------------------------------------------------------- \n";
return oss.str();
#else
return "";
#endif
}
| void libMesh::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 163 of file reference_counter.h.
References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().
{
Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
std::pair<unsigned int, unsigned int>& p = _counts[name];
p.first++;
}
| void libMesh::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 176 of file reference_counter.h.
References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().
{
Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
std::pair<unsigned int, unsigned int>& p = _counts[name];
p.second++;
}
| void libMesh::AdaptiveTimeSolver::init | ( | ) | [virtual] |
The initialization function. This method is used to initialize internal data structures before a simulation begins.
Reimplemented from libMesh::UnsteadySolver.
Definition at line 58 of file adaptive_time_solver.C.
References core_time_solver, libMesh::libmesh_assert(), and libMesh::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 UniquePtrs holding the same pointer - dangerous
// for future use. Replace with shared_ptr?
old_local_nonlinear_solution =
UniquePtr<NumericVector<Number> >(core_time_solver->old_local_nonlinear_solution.get());
}
| void libMesh::UnsteadySolver::init_data | ( | ) | [virtual, inherited] |
The data initialization function. This method is used to initialize internal data structures after the underlying System has been initialized
Reimplemented from libMesh::TimeSolver.
Definition at line 55 of file unsteady_solver.C.
References libMesh::TimeSolver::_system, libMesh::System::get_dof_map(), libMesh::DofMap::get_send_list(), libMesh::GHOSTED, libMesh::System::n_dofs(), libMesh::System::n_local_dofs(), libMesh::UnsteadySolver::old_local_nonlinear_solution, and libMesh::SERIAL.
{
TimeSolver::init_data();
#ifdef LIBMESH_ENABLE_GHOSTED
old_local_nonlinear_solution->init (_system.n_dofs(), _system.n_local_dofs(),
_system.get_dof_map().get_send_list(), false,
GHOSTED);
#else
old_local_nonlinear_solution->init (_system.n_dofs(), false, SERIAL);
#endif
}
| bool libMesh::TimeSolver::is_adjoint | ( | ) | const [inline, inherited] |
Accessor for querying whether we need to do a primal or adjoint solve
Definition at line 233 of file time_solver.h.
References libMesh::TimeSolver::_is_adjoint.
Referenced by libMesh::FEMSystem::build_context().
{ return _is_adjoint; }
| virtual bool libMesh::UnsteadySolver::is_steady | ( | ) | const [inline, virtual, inherited] |
This is not a steady-state solver.
Implements libMesh::TimeSolver.
Definition at line 146 of file unsteady_solver.h.
{ return false; }
| UniquePtr< LinearSolver< Number > > & libMesh::AdaptiveTimeSolver::linear_solver | ( | ) | [virtual] |
An implicit linear solver to use for adjoint and sensitivity problems.
Reimplemented from libMesh::TimeSolver.
Definition at line 149 of file adaptive_time_solver.C.
References core_time_solver.
{
return core_time_solver->linear_solver();
}
| static unsigned int libMesh::ReferenceCounter::n_objects | ( | ) | [inline, static, inherited] |
Prints the number of outstanding (created, but not yet destroyed) objects.
Definition at line 79 of file reference_counter.h.
References libMesh::ReferenceCounter::_n_objects.
Referenced by libMesh::LibMeshInit::~LibMeshInit().
{ return _n_objects; }
| bool libMesh::AdaptiveTimeSolver::nonlocal_residual | ( | bool | get_jacobian, |
| DiffContext & | context | ||
| ) | [virtual] |
This method is passed on to the core_time_solver
Implements libMesh::TimeSolver.
Definition at line 132 of file adaptive_time_solver.C.
References core_time_solver, and libMesh::libmesh_assert().
{
libmesh_assert(core_time_solver.get());
return core_time_solver->nonlocal_residual(request_jacobian, context);
}
| Number libMesh::UnsteadySolver::old_nonlinear_solution | ( | const dof_id_type | global_dof_number | ) | const [inherited] |
Definition at line 216 of file unsteady_solver.C.
References libMesh::TimeSolver::_system, libMesh::System::get_dof_map(), libMesh::DofMap::n_dofs(), and libMesh::UnsteadySolver::old_local_nonlinear_solution.
Referenced by libMesh::EulerSolver::_general_residual(), and libMesh::Euler2Solver::_general_residual().
{
libmesh_assert_less (global_dof_number, _system.get_dof_map().n_dofs());
libmesh_assert_less (global_dof_number, old_local_nonlinear_solution->size());
return (*old_local_nonlinear_solution)(global_dof_number);
}
| void libMesh::ReferenceCounter::print_info | ( | std::ostream & | out = libMesh::out | ) | [static, inherited] |
Prints the reference information, by default to libMesh::out.
Definition at line 88 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter, and libMesh::ReferenceCounter::get_info().
Referenced by libMesh::LibMeshInit::~LibMeshInit().
{
if( _enable_print_counter ) out_stream << ReferenceCounter::get_info();
}
| void libMesh::AdaptiveTimeSolver::reinit | ( | ) | [virtual] |
The reinitialization function. This method is used to resize internal data vectors after a mesh change.
Reimplemented from libMesh::UnsteadySolver.
Definition at line 77 of file adaptive_time_solver.C.
References core_time_solver, and libMesh::libmesh_assert().
{
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();
}
| void libMesh::UnsteadySolver::retrieve_timestep | ( | ) | [virtual, inherited] |
This method retrieves all the stored solutions at the current system.time
Reimplemented from libMesh::TimeSolver.
Definition at line 204 of file unsteady_solver.C.
References libMesh::TimeSolver::_system, libMesh::System::get_dof_map(), libMesh::DofMap::get_send_list(), libMesh::System::get_vector(), libMesh::NumericVector< T >::localize(), libMesh::UnsteadySolver::old_local_nonlinear_solution, and libMesh::TimeSolver::solution_history.
{
// Retrieve all the stored vectors at the current time
solution_history->retrieve();
// Dont forget to localize the old_nonlinear_solution !
_system.get_vector("_old_nonlinear_solution").localize
(*old_local_nonlinear_solution,
_system.get_dof_map().get_send_list());
}
| void libMesh::TimeSolver::set_is_adjoint | ( | bool | _is_adjoint_value | ) | [inline, inherited] |
Accessor for setting whether we need to do a primal or adjoint solve
Definition at line 240 of file time_solver.h.
References libMesh::TimeSolver::_is_adjoint.
Referenced by libMesh::DifferentiableSystem::adjoint_solve(), libMesh::FEMSystem::postprocess(), and libMesh::DifferentiableSystem::solve().
{ _is_adjoint = _is_adjoint_value; }
| void libMesh::TimeSolver::set_solution_history | ( | const SolutionHistory & | _solution_history | ) | [inherited] |
A setter function users will employ if they need to do something other than save no solution history
Definition at line 97 of file time_solver.C.
References libMesh::SolutionHistory::clone(), and libMesh::TimeSolver::solution_history.
{
solution_history = _solution_history.clone();
}
| bool libMesh::AdaptiveTimeSolver::side_residual | ( | bool | get_jacobian, |
| DiffContext & | context | ||
| ) | [virtual] |
This method is passed on to the core_time_solver
Implements libMesh::TimeSolver.
Definition at line 122 of file adaptive_time_solver.C.
References core_time_solver, and libMesh::libmesh_assert().
{
libmesh_assert(core_time_solver.get());
return core_time_solver->side_residual(request_jacobian, context);
}
| virtual void libMesh::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.
Reimplemented from libMesh::UnsteadySolver.
Implemented in libMesh::TwostepTimeSolver.
| const sys_type& libMesh::TimeSolver::system | ( | ) | const [inline, inherited] |
Definition at line 173 of file time_solver.h.
References libMesh::TimeSolver::_system.
Referenced by libMesh::TimeSolver::reinit(), and libMesh::TimeSolver::solve().
{ return _system; }
| sys_type& libMesh::TimeSolver::system | ( | ) | [inline, inherited] |
Definition at line 178 of file time_solver.h.
References libMesh::TimeSolver::_system.
{ return _system; }
ReferenceCounter::Counts libMesh::ReferenceCounter::_counts [static, protected, inherited] |
Actually holds the data.
Definition at line 118 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::get_info(), libMesh::ReferenceCounter::increment_constructor_count(), and libMesh::ReferenceCounter::increment_destructor_count().
UniquePtr<DiffSolver> libMesh::TimeSolver::_diff_solver [protected, inherited] |
An implicit linear or nonlinear solver to use at each timestep.
Definition at line 248 of file time_solver.h.
Referenced by libMesh::TimeSolver::diff_solver(), and libMesh::UnsteadySolver::solve().
bool libMesh::ReferenceCounter::_enable_print_counter = true [static, protected, inherited] |
Flag to control whether reference count information is printed when print_info is called.
Definition at line 137 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::disable_print_counter_info(), libMesh::ReferenceCounter::enable_print_counter_info(), and libMesh::ReferenceCounter::print_info().
UniquePtr<LinearSolver<Number> > libMesh::TimeSolver::_linear_solver [protected, inherited] |
An implicit linear solver to use for adjoint problems.
Definition at line 253 of file time_solver.h.
Referenced by libMesh::TimeSolver::linear_solver().
Threads::spin_mutex libMesh::ReferenceCounter::_mutex [static, protected, inherited] |
Mutual exclusion object to enable thread-safe reference counting.
Definition at line 131 of file reference_counter.h.
Threads::atomic< unsigned int > libMesh::ReferenceCounter::_n_objects [static, protected, inherited] |
The number of objects. Print the reference count information when the number returns to 0.
Definition at line 126 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::n_objects(), libMesh::ReferenceCounter::ReferenceCounter(), and libMesh::ReferenceCounter::~ReferenceCounter().
sys_type& libMesh::TimeSolver::_system [protected, inherited] |
A reference to the system we are solving.
Definition at line 258 of file time_solver.h.
Referenced by libMesh::EulerSolver::_general_residual(), libMesh::Euler2Solver::_general_residual(), libMesh::SteadySolver::_general_residual(), libMesh::UnsteadySolver::adjoint_advance_timestep(), advance_timestep(), libMesh::UnsteadySolver::advance_timestep(), libMesh::UnsteadySolver::du(), libMesh::EigenTimeSolver::element_residual(), libMesh::UnsteadySolver::init(), libMesh::TimeSolver::init(), libMesh::EigenTimeSolver::init(), libMesh::UnsteadySolver::init_data(), libMesh::TimeSolver::init_data(), libMesh::EigenTimeSolver::nonlocal_residual(), libMesh::UnsteadySolver::old_nonlinear_solution(), libMesh::UnsteadySolver::reinit(), libMesh::TimeSolver::reinit(), libMesh::UnsteadySolver::retrieve_timestep(), libMesh::EigenTimeSolver::side_residual(), libMesh::TwostepTimeSolver::solve(), libMesh::UnsteadySolver::solve(), libMesh::EigenTimeSolver::solve(), and libMesh::TimeSolver::system().
Error calculations are done in this norm, DISCRETE_L2 by default.
Definition at line 121 of file adaptive_time_solver.h.
Referenced by calculate_norm().
| std::vector<float> libMesh::AdaptiveTimeSolver::component_scale |
If 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 129 of file adaptive_time_solver.h.
This object is used to take timesteps
Definition at line 116 of file adaptive_time_solver.h.
Referenced by diff_solver(), element_residual(), error_order(), init(), linear_solver(), nonlocal_residual(), reinit(), side_residual(), libMesh::TwostepTimeSolver::solve(), and libMesh::TwostepTimeSolver::TwostepTimeSolver().
bool libMesh::UnsteadySolver::first_adjoint_step [protected, inherited] |
A bool that will be true the first time adjoint_advance_timestep() is called, (when the primal solution is to be used to set adjoint boundary conditions) and false thereafter
Definition at line 160 of file unsteady_solver.h.
Referenced by libMesh::UnsteadySolver::adjoint_advance_timestep().
bool libMesh::UnsteadySolver::first_solve [protected, inherited] |
A bool that will be true the first time solve() is called, and false thereafter
Definition at line 154 of file unsteady_solver.h.
Referenced by advance_timestep(), libMesh::UnsteadySolver::advance_timestep(), libMesh::TwostepTimeSolver::solve(), and libMesh::UnsteadySolver::solve().
This 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 199 of file adaptive_time_solver.h.
Referenced by libMesh::TwostepTimeSolver::solve().
Real libMesh::AdaptiveTimeSolver::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 208 of file adaptive_time_solver.h.
Referenced by advance_timestep(), and libMesh::TwostepTimeSolver::solve().
Do 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 169 of file adaptive_time_solver.h.
Referenced by libMesh::TwostepTimeSolver::solve().
Do 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 183 of file adaptive_time_solver.h.
Referenced by libMesh::TwostepTimeSolver::solve().
Do not allow the adaptive time solver to select deltat < min_deltat. The default value is 0.0.
Definition at line 175 of file adaptive_time_solver.h.
Referenced by libMesh::TwostepTimeSolver::solve().
UniquePtr<NumericVector<Number> > libMesh::UnsteadySolver::old_local_nonlinear_solution [inherited] |
Serial vector of _system.get_vector("_old_nonlinear_solution")
Reimplemented from libMesh::TimeSolver.
Definition at line 130 of file unsteady_solver.h.
Referenced by AdaptiveTimeSolver(), libMesh::UnsteadySolver::adjoint_advance_timestep(), libMesh::UnsteadySolver::advance_timestep(), init(), libMesh::UnsteadySolver::init_data(), libMesh::UnsteadySolver::old_nonlinear_solution(), libMesh::UnsteadySolver::reinit(), libMesh::UnsteadySolver::retrieve_timestep(), and ~AdaptiveTimeSolver().
bool libMesh::TimeSolver::quiet [inherited] |
Print extra debugging information if quiet == false.
Definition at line 193 of file time_solver.h.
Referenced by libMesh::TwostepTimeSolver::solve(), libMesh::UnsteadySolver::solve(), and libMesh::EigenTimeSolver::solve().
unsigned int libMesh::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 221 of file time_solver.h.
Referenced by libMesh::TwostepTimeSolver::solve(), and libMesh::UnsteadySolver::solve().
UniquePtr<SolutionHistory> libMesh::TimeSolver::solution_history [protected, inherited] |
An UniquePtr to a SolutionHistory object. Default is NoSolutionHistory, which the user can override by declaring a different kind of SolutionHistory in the application
Definition at line 270 of file time_solver.h.
Referenced by libMesh::UnsteadySolver::adjoint_advance_timestep(), libMesh::UnsteadySolver::advance_timestep(), libMesh::UnsteadySolver::retrieve_timestep(), and libMesh::TimeSolver::set_solution_history().
This tolerance is the target relative error between an exact time integration and a single time step output, scaled by deltat. integrator, scaled by deltat. If the estimated error exceeds or undershoots the target error tolerance, 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.
If a *negative* target_tolerance is specified, then its absolute value is used to scale the estimated error from the first simulation time step, and this becomes the target tolerance of all future time steps.
Definition at line 146 of file adaptive_time_solver.h.
Referenced by libMesh::TwostepTimeSolver::solve().
This tolerance is the maximum relative error between an exact time integration and a single time step output, scaled by deltat. If this error tolerance is exceeded by the estimated error of the current time step, that time step will be repeated with a smaller deltat.
If you use the default upper_tolerance=0.0, then the current time step will not be repeated regardless of estimated error.
If a *negative* upper_tolerance is specified, then its absolute value is used to scale the estimated error from the first simulation time step, and this becomes the upper tolerance of all future time steps.
Definition at line 163 of file adaptive_time_solver.h.
Referenced by libMesh::TwostepTimeSolver::solve().