libMesh::AdaptiveTimeSolver Class Referenceabstract

#include <adaptive_time_solver.h>

Inheritance diagram for libMesh::AdaptiveTimeSolver:

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
protectedinherited

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

Definition at line 57 of file adaptive_time_solver.h.

The type of system

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.

28  : UnsteadySolver(s),
29  core_time_solver(NULL),
30  target_tolerance(1.e-3), upper_tolerance(0.0),
31  max_deltat(0.),
32  min_deltat(0.),
33  max_growth(0.),
34  global_tolerance(true)
35 {
36  // the child class must populate core_time_solver
37  // with whatever actual time solver is to be used
38 
39  // As an UnsteadySolver, we have an old_local_nonlinear_solution, but we're
40  // going to drop it and use our core time solver's instead.
42 }
libMesh::AdaptiveTimeSolver::~AdaptiveTimeSolver ( )
virtual

Destructor.

Definition at line 46 of file adaptive_time_solver.C.

References libMesh::UnsteadySolver::old_local_nonlinear_solution.

47 {
48  // As an UnsteadySolver, we have an old_local_nonlinear_solution, but it
49  // is being managed by our core_time_solver. Make sure we don't delete
50  // it out from under them, in case the user wants to keep using the core
51  // solver after they're done with us.
53 }

Member Function Documentation

void libMesh::UnsteadySolver::adjoint_advance_timestep ( )
virtualinherited

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

177 {
178  // On the first call of this function, we dont save the adjoint solution or
179  // decrement the time, we just call the retrieve function below
180  if(!first_adjoint_step)
181  {
182  // Call the store function to store the last adjoint before decrementing the time
183  solution_history->store();
184  // Decrement the system time
186  }
187  else
188  {
189  first_adjoint_step = false;
190  }
191 
192  // Retrieve the primal solution vectors at this time using the
193  // solution_history object
194  solution_history->retrieve();
195 
196  // Dont forget to localize the old_nonlinear_solution !
197  _system.get_vector("_old_nonlinear_solution").localize
200 }
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.

87 {
88  NumericVector<Number> &old_nonlinear_soln =
89  _system.get_vector("_old_nonlinear_solution");
90  NumericVector<Number> &nonlinear_solution =
91  *(_system.solution);
92 // _system.get_vector("_nonlinear_solution");
93 
94  old_nonlinear_soln = nonlinear_solution;
95 
96  if (!first_solve)
98 }
virtual void libMesh::TimeSolver::before_timestep ( )
inlinevirtualinherited

This method is for subclasses or users to override to do arbitrary processing between timesteps

Definition at line 152 of file time_solver.h.

152 {}
Real libMesh::AdaptiveTimeSolver::calculate_norm ( System s,
NumericVector< Number > &  v 
)
protectedvirtual

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

140 {
141  return s.calculate_norm(v, component_norm);
142 }
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.

132 {
133  return core_time_solver->diff_solver();
134 }
void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited

Definition at line 106 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

107 {
108  _enable_print_counter = false;
109  return;
110 }
Real libMesh::UnsteadySolver::du ( const SystemNorm norm) const
virtualinherited

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

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

226 {
227 
228  AutoPtr<NumericVector<Number> > solution_copy =
229  _system.solution->clone();
230 
231  solution_copy->add(-1., _system.get_vector("_old_nonlinear_solution"));
232 
233  solution_copy->close();
234 
235  return _system.calculate_norm(*solution_copy, norm);
236 }
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::libmesh_assert().

113 {
115 
116  return core_time_solver->element_residual(request_jacobian, context);
117 }
void libMesh::ReferenceCounter::enable_print_counter_info ( )
staticinherited

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.

101 {
102  _enable_print_counter = true;
103  return;
104 }
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::libmesh_assert().

103 {
105 
106  return core_time_solver->error_order();
107 }
std::string libMesh::ReferenceCounter::get_info ( )
staticinherited

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

48 {
49 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
50 
51  std::ostringstream oss;
52 
53  oss << '\n'
54  << " ---------------------------------------------------------------------------- \n"
55  << "| Reference count information |\n"
56  << " ---------------------------------------------------------------------------- \n";
57 
58  for (Counts::iterator it = _counts.begin();
59  it != _counts.end(); ++it)
60  {
61  const std::string name(it->first);
62  const unsigned int creations = it->second.first;
63  const unsigned int destructions = it->second.second;
64 
65  oss << "| " << name << " reference count information:\n"
66  << "| Creations: " << creations << '\n'
67  << "| Destructions: " << destructions << '\n';
68  }
69 
70  oss << " ---------------------------------------------------------------------------- \n";
71 
72  return oss.str();
73 
74 #else
75 
76  return "";
77 
78 #endif
79 }
void libMesh::ReferenceCounter::increment_constructor_count ( const std::string &  name)
inlineprotectedinherited

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

164 {
165  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
166  std::pair<unsigned int, unsigned int>& p = _counts[name];
167 
168  p.first++;
169 }
void libMesh::ReferenceCounter::increment_destructor_count ( const std::string &  name)
inlineprotectedinherited

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

177 {
178  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
179  std::pair<unsigned int, unsigned int>& p = _counts[name];
180 
181  p.second++;
182 }
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::libmesh_assert(), and libMesh::UnsteadySolver::old_local_nonlinear_solution.

58 {
60 
61  // We override this because our core_time_solver is the one that
62  // needs to handle new vectors, diff_solver->init(), etc
63  core_time_solver->init();
64 
65  // As an UnsteadySolver, we have an old_local_nonlinear_solution, but it
66  // isn't pointing to the right place - fix it
67  //
68  // This leaves us with two AutoPtrs holding the same pointer - dangerous
69  // for future use. Replace with shared_ptr?
71  AutoPtr<NumericVector<Number> >(core_time_solver->old_local_nonlinear_solution.get());
72 }
void libMesh::UnsteadySolver::init_data ( )
virtualinherited

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.

56 {
57 #ifdef LIBMESH_ENABLE_GHOSTED
60  GHOSTED);
61 #else
63 #endif
64 }
bool libMesh::TimeSolver::is_adjoint ( ) const
inlineinherited

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

218  { return _is_adjoint; }
virtual bool libMesh::UnsteadySolver::is_steady ( ) const
inlinevirtualinherited

This is not a steady-state solver.

Implements libMesh::TimeSolver.

Definition at line 149 of file unsteady_solver.h.

149 { return false; }
virtual AutoPtr<LinearSolver<Number> >& libMesh::TimeSolver::linear_solver ( )
inlinevirtualinherited

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.

172 { return _linear_solver; }
static unsigned int libMesh::ReferenceCounter::n_objects ( )
inlinestaticinherited

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.

80  { 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 214 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().

216 {
217  libmesh_assert_less (global_dof_number, _system.get_dof_map().n_dofs());
218  libmesh_assert_less (global_dof_number, old_local_nonlinear_solution->size());
219 
220  return (*old_local_nonlinear_solution)(global_dof_number);
221 }
void libMesh::ReferenceCounter::print_info ( std::ostream &  out = libMesh::out)
staticinherited

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

89 {
91 }
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::libmesh_assert().

77 {
79 
80  // We override this because our core_time_solver is the one that
81  // needs to handle new vectors, diff_solver->reinit(), etc
82  core_time_solver->reinit();
83 }
void libMesh::UnsteadySolver::retrieve_timestep ( )
virtualinherited

This method retrieves all the stored solutions at the current system.time

Reimplemented from libMesh::TimeSolver.

Definition at line 202 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.

203  {
204  // Retrieve all the stored vectors at the current time
205  solution_history->retrieve();
206 
207  // Dont forget to localize the old_nonlinear_solution !
208  _system.get_vector("_old_nonlinear_solution").localize
211  }
void libMesh::TimeSolver::set_is_adjoint ( bool  _is_adjoint_value)
inlineinherited

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

225  { _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.

92  {
93  solution_history = _solution_history.clone();
94  }
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::libmesh_assert().

123 {
125 
126  return core_time_solver->side_residual(request_jacobian, context);
127 }
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
inlineinherited
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().

157 { return _system; }
sys_type& libMesh::TimeSolver::system ( )
inlineinherited
Returns
a writeable reference to the system we are solving.

Definition at line 162 of file time_solver.h.

References libMesh::TimeSolver::_system.

162 { return _system; }

Member Data Documentation

ReferenceCounter::Counts libMesh::ReferenceCounter::_counts
staticprotectedinherited
AutoPtr<DiffSolver> libMesh::TimeSolver::_diff_solver
protectedinherited

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
staticprotectedinherited

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

AutoPtr<LinearSolver<Number> > libMesh::TimeSolver::_linear_solver
protectedinherited

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

Threads::spin_mutex libMesh::ReferenceCounter::_mutex
staticprotectedinherited

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
staticprotectedinherited

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

SystemNorm libMesh::AdaptiveTimeSolver::component_norm

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

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

AutoPtr<UnsteadySolver> libMesh::AdaptiveTimeSolver::core_time_solver
bool libMesh::UnsteadySolver::first_adjoint_step
protectedinherited

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
protectedinherited

A bool that will be true the first time solve() is called, and false thereafter

Definition at line 157 of file unsteady_solver.h.

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

bool libMesh::AdaptiveTimeSolver::global_tolerance

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

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

Real libMesh::AdaptiveTimeSolver::max_deltat

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

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

Real libMesh::AdaptiveTimeSolver::max_growth

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

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

Real libMesh::AdaptiveTimeSolver::min_deltat

Do not allow the adaptive time solver to select deltat < min_deltat. The default value is 0.0.

Definition at line 163 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::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 205 of file time_solver.h.

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

AutoPtr<SolutionHistory> libMesh::TimeSolver::solution_history
protectedinherited

An AutoPtr 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 260 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().

Real libMesh::AdaptiveTimeSolver::target_tolerance

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

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

Real libMesh::AdaptiveTimeSolver::upper_tolerance

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 151 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 07 2014 16:57:58 UTC

Hosted By:
SourceForge.net Logo