libMesh::AdaptiveTimeSolver Class Reference

#include <adaptive_time_solver.h>

Inheritance diagram for libMesh::AdaptiveTimeSolver:

List of all members.

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 AutoPtr< DiffSolver > & diff_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_typesystem () const
sys_typesystem ()
virtual AutoPtr< LinearSolver
< Number > > & 
linear_solver ()
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

AutoPtr< UnsteadySolvercore_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
AutoPtr< NumericVector< Number > > old_local_nonlinear_solution
bool quiet
unsigned int reduce_deltat_on_diffsolver_failure

Protected Types

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
AutoPtr< DiffSolver_diff_solver
AutoPtr< LinearSolver< Number > > _linear_solver
sys_type_system
AutoPtr< SolutionHistorysolution_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

Detailed Description

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 51 of file adaptive_time_solver.h.


Member Typedef Documentation

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.

The type of system

Reimplemented in libMesh::EigenTimeSolver, and libMesh::SteadySolver.

Definition at line 66 of file time_solver.h.


Constructor & Destructor Documentation

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, and libMesh::AutoPtr< Tp >::reset().

00028  : UnsteadySolver(s),
00029    core_time_solver(NULL),
00030    target_tolerance(1.e-3), upper_tolerance(0.0),
00031    max_deltat(0.),
00032    min_deltat(0.),
00033    max_growth(0.),
00034    global_tolerance(true)
00035 {
00036   // the child class must populate core_time_solver
00037   // with whatever actual time solver is to be used
00038 
00039   // As an UnsteadySolver, we have an old_local_nonlinear_solution, but we're
00040   // going to drop it and use our core time solver's instead.
00041   old_local_nonlinear_solution.reset();
00042 }

libMesh::AdaptiveTimeSolver::~AdaptiveTimeSolver (  )  [virtual]

Destructor.

Definition at line 46 of file adaptive_time_solver.C.

References libMesh::UnsteadySolver::old_local_nonlinear_solution, and libMesh::AutoPtr< Tp >::release().

00047 {
00048   // As an UnsteadySolver, we have an old_local_nonlinear_solution, but it
00049   // is being managed by our core_time_solver.  Make sure we don't delete
00050   // it out from under them, in case the user wants to keep using the core
00051   // solver after they're done with us.
00052   old_local_nonlinear_solution.release();
00053 }


Member Function Documentation

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 169 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.

00170 {
00171   // On the first call of this function, we dont save the adjoint solution or
00172   // decrement the time, we just call the retrieve function below
00173   if(!first_adjoint_step)
00174     {
00175       // Call the store function to store the last adjoint before decrementing the time
00176       solution_history->store();
00177       // Decrement the system time
00178       _system.time -= _system.deltat;
00179     }
00180   else
00181     {
00182       first_adjoint_step = false;
00183     }
00184 
00185   // Retrieve the primal solution vectors at this time using the
00186   // solution_history object
00187   solution_history->retrieve();
00188 
00189   // Dont forget to localize the old_nonlinear_solution !
00190   _system.get_vector("_old_nonlinear_solution").localize
00191     (*old_local_nonlinear_solution,
00192      _system.get_dof_map().get_send_list());
00193 }

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 86 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.

00087 {
00088   NumericVector<Number> &old_nonlinear_soln =
00089   _system.get_vector("_old_nonlinear_solution");
00090   NumericVector<Number> &nonlinear_solution =
00091     *(_system.solution);
00092 //    _system.get_vector("_nonlinear_solution");
00093 
00094   old_nonlinear_soln = nonlinear_solution;
00095 
00096   if (!first_solve)
00097     _system.time += last_deltat;
00098 }

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 152 of file time_solver.h.

00152 {}

Real libMesh::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 libMesh::System::calculate_norm(), and component_norm.

Referenced by libMesh::TwostepTimeSolver::solve().

00140 {
00141   return s.calculate_norm(v, component_norm);
00142 }

AutoPtr< DiffSolver > & libMesh::AdaptiveTimeSolver::diff_solver (  )  [virtual]

An implicit linear or nonlinear solver to use at each timestep.

Reimplemented from libMesh::TimeSolver.

Definition at line 131 of file adaptive_time_solver.C.

References core_time_solver.

00132 {
00133   return core_time_solver->diff_solver();
00134 }

void libMesh::ReferenceCounter::disable_print_counter_info (  )  [static, inherited]

Definition at line 106 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

00107 {
00108   _enable_print_counter = false;
00109   return;
00110 }

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 218 of file unsteady_solver.C.

References libMesh::TimeSolver::_system, libMesh::System::calculate_norm(), libMesh::System::get_vector(), and libMesh::System::solution.

00219 {
00220 
00221   AutoPtr<NumericVector<Number> > solution_copy =
00222     _system.solution->clone();
00223 
00224   solution_copy->add(-1., _system.get_vector("_old_nonlinear_solution"));
00225 
00226   solution_copy->close();
00227 
00228   return _system.calculate_norm(*solution_copy, norm);
00229 }

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 111 of file adaptive_time_solver.C.

References core_time_solver, and libMesh::AutoPtr< Tp >::get().

00113 {
00114   libmesh_assert(core_time_solver.get());
00115 
00116   return core_time_solver->element_residual(request_jacobian, context);
00117 }

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.

00101 {
00102   _enable_print_counter = true;
00103   return;
00104 }

Real libMesh::AdaptiveTimeSolver::error_order (  )  const [virtual]

This method is passed on to the core_time_solver

Implements libMesh::UnsteadySolver.

Definition at line 102 of file adaptive_time_solver.C.

References core_time_solver, and libMesh::AutoPtr< Tp >::get().

00103 {
00104   libmesh_assert(core_time_solver.get());
00105 
00106   return core_time_solver->error_order();
00107 }

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().

00048 {
00049 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
00050 
00051   std::ostringstream oss;
00052 
00053   oss << '\n'
00054       << " ---------------------------------------------------------------------------- \n"
00055       << "| Reference count information                                                |\n"
00056       << " ---------------------------------------------------------------------------- \n";
00057 
00058   for (Counts::iterator it = _counts.begin();
00059        it != _counts.end(); ++it)
00060     {
00061       const std::string name(it->first);
00062       const unsigned int creations    = it->second.first;
00063       const unsigned int destructions = it->second.second;
00064 
00065       oss << "| " << name << " reference count information:\n"
00066           << "|  Creations:    " << creations    << '\n'
00067           << "|  Destructions: " << destructions << '\n';
00068     }
00069 
00070   oss << " ---------------------------------------------------------------------------- \n";
00071 
00072   return oss.str();
00073 
00074 #else
00075 
00076   return "";
00077 
00078 #endif
00079 }

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, and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().

00164 {
00165   Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
00166   std::pair<unsigned int, unsigned int>& p = _counts[name];
00167 
00168   p.first++;
00169 }

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, and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().

00177 {
00178   Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
00179   std::pair<unsigned int, unsigned int>& p = _counts[name];
00180 
00181   p.second++;
00182 }

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 57 of file adaptive_time_solver.C.

References core_time_solver, libMesh::AutoPtr< Tp >::get(), and libMesh::UnsteadySolver::old_local_nonlinear_solution.

00058 {
00059   libmesh_assert(core_time_solver.get());
00060 
00061   // We override this because our core_time_solver is the one that
00062   // needs to handle new vectors, diff_solver->init(), etc
00063   core_time_solver->init();
00064 
00065   // As an UnsteadySolver, we have an old_local_nonlinear_solution, but it
00066   // isn't pointing to the right place - fix it
00067   //
00068   // This leaves us with two AutoPtrs holding the same pointer - dangerous
00069   // for future use.  Replace with shared_ptr?
00070   old_local_nonlinear_solution =
00071     AutoPtr<NumericVector<Number> >(core_time_solver->old_local_nonlinear_solution.get());
00072 }

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(), libMeshEnums::GHOSTED, libMesh::System::n_dofs(), libMesh::System::n_local_dofs(), libMesh::UnsteadySolver::old_local_nonlinear_solution, and libMeshEnums::SERIAL.

00056 {
00057 #ifdef LIBMESH_ENABLE_GHOSTED
00058   old_local_nonlinear_solution->init (_system.n_dofs(), _system.n_local_dofs(),
00059                                       _system.get_dof_map().get_send_list(), false,
00060                                       GHOSTED);
00061 #else
00062   old_local_nonlinear_solution->init (_system.n_dofs(), false, SERIAL);
00063 #endif
00064 }

bool libMesh::TimeSolver::is_adjoint (  )  const [inline, inherited]

Accessor for querying whether we need to do a primal or adjoint solve

Definition at line 217 of file time_solver.h.

References libMesh::TimeSolver::_is_adjoint.

Referenced by libMesh::FEMSystem::build_context().

00218   { 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 149 of file unsteady_solver.h.

00149 { return false; }

virtual AutoPtr<LinearSolver<Number> >& libMesh::TimeSolver::linear_solver (  )  [inline, virtual, inherited]

An implicit linear solver to use for adjoint and sensitivity problems.

Definition at line 172 of file time_solver.h.

References libMesh::TimeSolver::_linear_solver.

00172 { return _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.

00080   { return _n_objects; }

Number libMesh::UnsteadySolver::old_nonlinear_solution ( const dof_id_type  global_dof_number  )  const [inherited]
Returns:
the old nonlinear solution for the specified global DOF.

Definition at line 207 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::element_residual(), libMesh::Euler2Solver::element_residual(), libMesh::EulerSolver::side_residual(), and libMesh::Euler2Solver::side_residual().

00209 {
00210   libmesh_assert_less (global_dof_number, _system.get_dof_map().n_dofs());
00211   libmesh_assert_less (global_dof_number, old_local_nonlinear_solution->size());
00212 
00213   return (*old_local_nonlinear_solution)(global_dof_number);
00214 }

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().

00089 {
00090   if( _enable_print_counter ) out_stream << ReferenceCounter::get_info();
00091 }

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 76 of file adaptive_time_solver.C.

References core_time_solver, and libMesh::AutoPtr< Tp >::get().

00077 {
00078   libmesh_assert(core_time_solver.get());
00079 
00080   // We override this because our core_time_solver is the one that
00081   // needs to handle new vectors, diff_solver->reinit(), etc
00082   core_time_solver->reinit();
00083 }

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 195 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.

00196   {
00197     // Retrieve all the stored vectors at the current time
00198     solution_history->retrieve();
00199 
00200     // Dont forget to localize the old_nonlinear_solution !
00201     _system.get_vector("_old_nonlinear_solution").localize
00202     (*old_local_nonlinear_solution,
00203      _system.get_dof_map().get_send_list());
00204   }

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 224 of file time_solver.h.

References libMesh::TimeSolver::_is_adjoint.

Referenced by libMesh::DifferentiableSystem::adjoint_solve(), libMesh::FEMSystem::postprocess(), and libMesh::DifferentiableSystem::solve().

00225   { _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 91 of file time_solver.C.

References libMesh::SolutionHistory::clone(), and libMesh::TimeSolver::solution_history.

00092  {
00093    solution_history = _solution_history.clone();
00094  }

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 121 of file adaptive_time_solver.C.

References core_time_solver, and libMesh::AutoPtr< Tp >::get().

00123 {
00124   libmesh_assert(core_time_solver.get());
00125 
00126   return core_time_solver->side_residual(request_jacobian, context);
00127 }

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.

sys_type& libMesh::TimeSolver::system (  )  [inline, inherited]
Returns:
a writeable reference to the system we are solving.

Definition at line 162 of file time_solver.h.

References libMesh::TimeSolver::_system.

00162 { return _system; }

const sys_type& libMesh::TimeSolver::system (  )  const [inline, inherited]
Returns:
a constant reference to the system we are solving.

Definition at line 157 of file time_solver.h.

References libMesh::TimeSolver::_system.

Referenced by libMesh::TimeSolver::reinit(), and libMesh::TimeSolver::solve().

00157 { return _system; }


Member Data Documentation

An implicit linear or nonlinear solver to use at each timestep.

Definition at line 232 of file time_solver.h.

Referenced by libMesh::TimeSolver::diff_solver(), libMesh::TimeSolver::init(), libMesh::TimeSolver::reinit(), libMesh::UnsteadySolver::solve(), and libMesh::TimeSolver::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().

An implicit linear solver to use for adjoint problems.

Definition at line 237 of file time_solver.h.

Referenced by libMesh::TimeSolver::init(), libMesh::TimeSolver::linear_solver(), and libMesh::TimeSolver::reinit().

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().

Error calculations are done in this norm, DISCRETE_L2 by default.

Definition at line 109 of file adaptive_time_solver.h.

Referenced by calculate_norm().

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 117 of file adaptive_time_solver.h.

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 163 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

Reimplemented from libMesh::TimeSolver.

Definition at line 157 of file unsteady_solver.h.

Referenced by libMesh::UnsteadySolver::advance_timestep(), advance_timestep(), libMesh::UnsteadySolver::solve(), and libMesh::TwostepTimeSolver::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 173 of file adaptive_time_solver.h.

Referenced by libMesh::TwostepTimeSolver::solve().

We need to store the value of the last deltat used, so that advance_timestep() will increment the system time correctly.

Definition at line 182 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 143 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 157 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 149 of file adaptive_time_solver.h.

Referenced by libMesh::TwostepTimeSolver::solve().

bool libMesh::TimeSolver::quiet [inherited]

Print extra debugging information if quiet == false.

Definition at line 177 of file time_solver.h.

Referenced by libMesh::UnsteadySolver::solve(), libMesh::TwostepTimeSolver::solve(), and libMesh::EigenTimeSolver::solve().

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 205 of file time_solver.h.

Referenced by libMesh::UnsteadySolver::solve(), and libMesh::TwostepTimeSolver::solve().

This 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 128 of file adaptive_time_solver.h.

Referenced by libMesh::TwostepTimeSolver::solve().

This 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 137 of file adaptive_time_solver.h.

Referenced by libMesh::TwostepTimeSolver::solve().


The documentation for this class was generated from the following files:

Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:06 UTC

Hosted By:
SourceForge.net Logo