libMesh::TimeSolver Class Referenceabstract

#include <time_solver.h>

Inheritance diagram for libMesh::TimeSolver:

Public Types

typedef DifferentiableSystem sys_type
 

Public Member Functions

 TimeSolver (sys_type &s)
 
virtual ~TimeSolver ()
 
virtual void init ()
 
virtual void init_data ()
 
virtual void reinit ()
 
virtual void solve ()
 
virtual void advance_timestep ()
 
virtual void adjoint_advance_timestep ()
 
virtual void retrieve_timestep ()
 
virtual bool element_residual (bool request_jacobian, DiffContext &)=0
 
virtual bool side_residual (bool request_jacobian, DiffContext &)=0
 
virtual bool nonlocal_residual (bool request_jacobian, DiffContext &)=0
 
virtual void before_timestep ()
 
const sys_typesystem () const
 
sys_typesystem ()
 
virtual std::unique_ptr< DiffSolver > & diff_solver ()
 
virtual std::unique_ptr< LinearSolver< Number > > & linear_solver ()
 
virtual Real du (const SystemNorm &norm) const =0
 
virtual bool is_steady () const =0
 
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

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

void increment_constructor_count (const std::string &name)
 
void increment_destructor_count (const std::string &name)
 

Protected Attributes

std::unique_ptr< DiffSolver_diff_solver
 
std::unique_ptr< LinearSolver< Number > > _linear_solver
 
sys_type_system
 
std::unique_ptr< 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
 

Private Attributes

bool _is_adjoint
 

Detailed Description

This is a generic class that defines a solver to handle time integration of DifferentiableSystems.

A user can define a solver by deriving from this class and implementing certain functions.

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
Date
2006

Definition at line 59 of file time_solver.h.

Member Typedef Documentation

◆ Counts

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

◆ ReinitFuncType

typedef void(DiffContext::* libMesh::TimeSolver::ReinitFuncType) (Real)
protected

Definition at line 273 of file time_solver.h.

◆ ResFuncType

typedef bool(DifferentiablePhysics::* libMesh::TimeSolver::ResFuncType) (bool, DiffContext &)
protected

Definitions of argument types for use in refactoring subclasses.

Definition at line 271 of file time_solver.h.

◆ sys_type

The type of system

Definition at line 65 of file time_solver.h.

Constructor & Destructor Documentation

◆ TimeSolver()

libMesh::TimeSolver::TimeSolver ( sys_type s)
explicit

Constructor. Requires a reference to the system to be solved.

Definition at line 29 of file time_solver.C.

30  : quiet (true),
32  _diff_solver (),
33  _linear_solver (),
34  _system (s),
35  solution_history(new NoSolutionHistory()), // Default setting for solution_history
36  _is_adjoint (false)
37 {
38 }
std::unique_ptr< DiffSolver > _diff_solver
Definition: time_solver.h:248
std::unique_ptr< SolutionHistory > solution_history
Definition: time_solver.h:265
sys_type & _system
Definition: time_solver.h:258
unsigned int reduce_deltat_on_diffsolver_failure
Definition: time_solver.h:221
std::unique_ptr< LinearSolver< Number > > _linear_solver
Definition: time_solver.h:253

◆ ~TimeSolver()

libMesh::TimeSolver::~TimeSolver ( )
virtual

Destructor.

Definition at line 42 of file time_solver.C.

43 {
44 }

Member Function Documentation

◆ adjoint_advance_timestep()

void libMesh::TimeSolver::adjoint_advance_timestep ( )
virtual

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 in libMesh::UnsteadySolver, and libMesh::NewmarkSolver.

Definition at line 106 of file time_solver.C.

107 {
108 }

◆ advance_timestep()

void libMesh::TimeSolver::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 in libMesh::EigenTimeSolver, libMesh::UnsteadySolver, libMesh::AdaptiveTimeSolver, and libMesh::NewmarkSolver.

Definition at line 102 of file time_solver.C.

103 {
104 }

◆ before_timestep()

virtual void libMesh::TimeSolver::before_timestep ( )
inlinevirtual

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

Definition at line 167 of file time_solver.h.

167 {}

◆ diff_solver()

virtual std::unique_ptr<DiffSolver>& libMesh::TimeSolver::diff_solver ( )
inlinevirtual

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

Reimplemented in libMesh::AdaptiveTimeSolver.

Definition at line 182 of file time_solver.h.

References _diff_solver.

Referenced by init(), init_data(), reinit(), and solve().

182 { return _diff_solver; }
std::unique_ptr< DiffSolver > _diff_solver
Definition: time_solver.h:248

◆ disable_print_counter_info()

void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited

Definition at line 106 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

Referenced by libMesh::LibMeshInit::LibMeshInit().

107 {
108  _enable_print_counter = false;
109  return;
110 }

◆ du()

virtual Real libMesh::TimeSolver::du ( const SystemNorm norm) const
pure virtual

Computes the size of ||u^{n+1} - u^{n}|| in some norm.

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

Implemented in libMesh::UnsteadySolver, libMesh::EigenTimeSolver, and libMesh::SteadySolver.

◆ element_residual()

virtual bool libMesh::TimeSolver::element_residual ( bool  request_jacobian,
DiffContext  
)
pure virtual

This method uses the DifferentiablePhysics element_time_derivative(), element_constraint(), and mass_residual() to build a full residual on an element. What combination

it uses will depend on the type of solver. See the subclasses for more details.

Implemented in libMesh::EigenTimeSolver, libMesh::NewmarkSolver, libMesh::AdaptiveTimeSolver, libMesh::SteadySolver, libMesh::Euler2Solver, and libMesh::EulerSolver.

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

◆ enable_print_counter_info()

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 }

◆ get_info()

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 (const auto & pr : _counts)
59  {
60  const std::string name(pr.first);
61  const unsigned int creations = pr.second.first;
62  const unsigned int destructions = pr.second.second;
63 
64  oss << "| " << name << " reference count information:\n"
65  << "| Creations: " << creations << '\n'
66  << "| Destructions: " << destructions << '\n';
67  }
68 
69  oss << " ---------------------------------------------------------------------------- \n";
70 
71  return oss.str();
72 
73 #else
74 
75  return "";
76 
77 #endif
78 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:42

◆ increment_constructor_count()

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 181 of file reference_counter.h.

References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

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

182 {
183  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
184  std::pair<unsigned int, unsigned int> & p = _counts[name];
185 
186  p.first++;
187 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:42
spin_mutex spin_mtx
Definition: threads.C:29

◆ increment_destructor_count()

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 194 of file reference_counter.h.

References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

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

195 {
196  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
197  std::pair<unsigned int, unsigned int> & p = _counts[name];
198 
199  p.second++;
200 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:42
spin_mutex spin_mtx
Definition: threads.C:29

◆ init()

void libMesh::TimeSolver::init ( )
virtual

The initialization function. This method is used to initialize internal data structures before a simulation begins.

Reimplemented in libMesh::EigenTimeSolver, libMesh::AdaptiveTimeSolver, libMesh::UnsteadySolver, and libMesh::SecondOrderUnsteadySolver.

Definition at line 64 of file time_solver.C.

References _system, libMesh::LinearSolver< T >::build(), libMesh::DiffSolver::build(), libMesh::ParallelObject::comm(), diff_solver(), and linear_solver().

Referenced by libMesh::UnsteadySolver::init().

65 {
66  // If the user hasn't given us a solver to use,
67  // just build a default solver
68  if (this->diff_solver().get() == nullptr)
70 
71  if (this->linear_solver().get() == nullptr)
73 }
static std::unique_ptr< LinearSolver< T > > build(const libMesh::Parallel::Communicator &comm_in, const SolverPackage solver_package=libMesh::default_solver_package())
Definition: linear_solver.C:57
const Parallel::Communicator & comm() const
static std::unique_ptr< DiffSolver > build(sys_type &s)
Definition: diff_solver.C:53
virtual std::unique_ptr< DiffSolver > & diff_solver()
Definition: time_solver.h:182
sys_type & _system
Definition: time_solver.h:258
virtual std::unique_ptr< LinearSolver< Number > > & linear_solver()
Definition: time_solver.h:187

◆ init_data()

void libMesh::TimeSolver::init_data ( )
virtual

The data initialization function. This method is used to initialize internal data structures after the underlying System has been initialized

Reimplemented in libMesh::UnsteadySolver, and libMesh::SecondOrderUnsteadySolver.

Definition at line 77 of file time_solver.C.

References _system, diff_solver(), linear_solver(), libMesh::System::name(), and libMesh::on_command_line().

Referenced by libMesh::UnsteadySolver::init_data().

78 {
79  this->diff_solver()->init();
80 
81  if (libMesh::on_command_line("--solver-system-names"))
82  this->linear_solver()->init((_system.name()+"_").c_str());
83  else
84  this->linear_solver()->init();
85 }
virtual std::unique_ptr< DiffSolver > & diff_solver()
Definition: time_solver.h:182
sys_type & _system
Definition: time_solver.h:258
virtual std::unique_ptr< LinearSolver< Number > > & linear_solver()
Definition: time_solver.h:187
bool on_command_line(std::string arg)
Definition: libmesh.C:876
const std::string & name() const
Definition: system.h:2017

◆ is_adjoint()

bool libMesh::TimeSolver::is_adjoint ( ) const
inline

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

Definition at line 233 of file time_solver.h.

References _is_adjoint.

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

234  { return _is_adjoint; }

◆ is_steady()

virtual bool libMesh::TimeSolver::is_steady ( ) const
pure virtual

◆ linear_solver()

virtual std::unique_ptr<LinearSolver<Number> >& libMesh::TimeSolver::linear_solver ( )
inlinevirtual

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

Reimplemented in libMesh::AdaptiveTimeSolver.

Definition at line 187 of file time_solver.h.

References _linear_solver.

Referenced by init(), init_data(), and reinit().

187 { return _linear_solver; }
std::unique_ptr< LinearSolver< Number > > _linear_solver
Definition: time_solver.h:253

◆ n_objects()

static unsigned int libMesh::ReferenceCounter::n_objects ( )
inlinestaticinherited

Prints the number of outstanding (created, but not yet destroyed) objects.

Definition at line 83 of file reference_counter.h.

References libMesh::ReferenceCounter::_n_objects.

84  { return _n_objects; }
static Threads::atomic< unsigned int > _n_objects

◆ nonlocal_residual()

virtual bool libMesh::TimeSolver::nonlocal_residual ( bool  request_jacobian,
DiffContext  
)
pure virtual

This method uses the DifferentiablePhysics nonlocal_time_derivative(), nonlocal_constraint(), and nonlocal_mass_residual() to build a full residual of non-local terms.

What combination it uses will depend on the type of solver. See the subclasses for more details.

Implemented in libMesh::NewmarkSolver, libMesh::EigenTimeSolver, libMesh::SteadySolver, libMesh::AdaptiveTimeSolver, libMesh::Euler2Solver, and libMesh::EulerSolver.

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

◆ print_info()

void libMesh::ReferenceCounter::print_info ( std::ostream &  out = libMesh::out)
staticinherited

Prints the reference information, by default to libMesh::out.

Definition at line 87 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter, and libMesh::ReferenceCounter::get_info().

88 {
90  out_stream << ReferenceCounter::get_info();
91 }
static std::string get_info()

◆ reinit()

void libMesh::TimeSolver::reinit ( )
virtual

The reinitialization function. This method is used after changes in the mesh

Reimplemented in libMesh::EigenTimeSolver, libMesh::UnsteadySolver, libMesh::SecondOrderUnsteadySolver, and libMesh::AdaptiveTimeSolver.

Definition at line 48 of file time_solver.C.

References _system, diff_solver(), linear_solver(), libMesh::System::name(), libMesh::on_command_line(), and system().

Referenced by libMesh::UnsteadySolver::reinit().

49 {
50  libmesh_assert(this->diff_solver().get());
51  libmesh_assert_equal_to (&(this->diff_solver()->system()), &(this->system()));
52  this->diff_solver()->reinit();
53 
54  libmesh_assert(this->linear_solver().get());
55  this->linear_solver()->clear();
56  if (libMesh::on_command_line("--solver-system-names"))
57  this->linear_solver()->init((_system.name()+"_").c_str());
58  else
59  this->linear_solver()->init();
60 }
virtual std::unique_ptr< DiffSolver > & diff_solver()
Definition: time_solver.h:182
const sys_type & system() const
Definition: time_solver.h:172
sys_type & _system
Definition: time_solver.h:258
virtual std::unique_ptr< LinearSolver< Number > > & linear_solver()
Definition: time_solver.h:187
bool on_command_line(std::string arg)
Definition: libmesh.C:876
const std::string & name() const
Definition: system.h:2017

◆ retrieve_timestep()

void libMesh::TimeSolver::retrieve_timestep ( )
virtual

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

Reimplemented in libMesh::UnsteadySolver, and libMesh::SecondOrderUnsteadySolver.

Definition at line 110 of file time_solver.C.

111 {
112 }

◆ set_is_adjoint()

void libMesh::TimeSolver::set_is_adjoint ( bool  _is_adjoint_value)
inline

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

Definition at line 240 of file time_solver.h.

References _is_adjoint.

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

241  { _is_adjoint = _is_adjoint_value; }

◆ set_solution_history()

void libMesh::TimeSolver::set_solution_history ( const SolutionHistory _solution_history)

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

98 {
99  solution_history = _solution_history.clone();
100 }
std::unique_ptr< SolutionHistory > solution_history
Definition: time_solver.h:265

◆ side_residual()

virtual bool libMesh::TimeSolver::side_residual ( bool  request_jacobian,
DiffContext  
)
pure virtual

This method uses the DifferentiablePhysics side_time_derivative(), side_constraint(), and side_mass_residual() to build a full residual on an element's side.

What combination it uses will depend on the type of solver. See the subclasses for more details.

Implemented in libMesh::NewmarkSolver, libMesh::EigenTimeSolver, libMesh::AdaptiveTimeSolver, libMesh::SteadySolver, libMesh::Euler2Solver, and libMesh::EulerSolver.

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

◆ solve()

void libMesh::TimeSolver::solve ( )
virtual

This method solves for the solution at the next timestep (or solves for a steady-state solution). Usually we will only need to solve one (non)linear system per timestep, but more complex subclasses may override this.

Reimplemented in libMesh::NewmarkSolver, libMesh::EigenTimeSolver, libMesh::UnsteadySolver, libMesh::AdaptiveTimeSolver, and libMesh::TwostepTimeSolver.

Definition at line 89 of file time_solver.C.

References diff_solver(), and system().

90 {
91  libmesh_assert(this->diff_solver().get());
92  libmesh_assert_equal_to (&(this->diff_solver()->system()), &(this->system()));
93  this->diff_solver()->solve();
94 }
virtual std::unique_ptr< DiffSolver > & diff_solver()
Definition: time_solver.h:182
const sys_type & system() const
Definition: time_solver.h:172

◆ system() [1/2]

const sys_type& libMesh::TimeSolver::system ( ) const
inline
Returns
A constant reference to the system we are solving.

Definition at line 172 of file time_solver.h.

References _system.

Referenced by reinit(), and solve().

172 { return _system; }
sys_type & _system
Definition: time_solver.h:258

◆ system() [2/2]

sys_type& libMesh::TimeSolver::system ( )
inline
Returns
A writable reference to the system we are solving.

Definition at line 177 of file time_solver.h.

References _system.

177 { return _system; }
sys_type & _system
Definition: time_solver.h:258

Member Data Documentation

◆ _counts

ReferenceCounter::Counts libMesh::ReferenceCounter::_counts
staticprotectedinherited

◆ _diff_solver

std::unique_ptr<DiffSolver> libMesh::TimeSolver::_diff_solver
protected

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

Definition at line 248 of file time_solver.h.

Referenced by libMesh::NewmarkSolver::compute_initial_accel(), diff_solver(), and libMesh::UnsteadySolver::solve().

◆ _enable_print_counter

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 141 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::disable_print_counter_info(), libMesh::ReferenceCounter::enable_print_counter_info(), and libMesh::ReferenceCounter::print_info().

◆ _is_adjoint

bool libMesh::TimeSolver::_is_adjoint
private

This boolean tells the TimeSolver whether we are solving a primal or adjoint problem

Definition at line 281 of file time_solver.h.

Referenced by is_adjoint(), and set_is_adjoint().

◆ _linear_solver

std::unique_ptr<LinearSolver<Number> > libMesh::TimeSolver::_linear_solver
protected

An implicit linear solver to use for adjoint problems.

Definition at line 253 of file time_solver.h.

Referenced by linear_solver().

◆ _mutex

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

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 135 of file reference_counter.h.

◆ _n_objects

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 130 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::n_objects(), libMesh::ReferenceCounter::ReferenceCounter(), and libMesh::ReferenceCounter::~ReferenceCounter().

◆ _system

sys_type& libMesh::TimeSolver::_system
protected

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::NewmarkSolver::_general_residual(), libMesh::UnsteadySolver::adjoint_advance_timestep(), libMesh::NewmarkSolver::advance_timestep(), libMesh::AdaptiveTimeSolver::advance_timestep(), libMesh::UnsteadySolver::advance_timestep(), libMesh::NewmarkSolver::compute_initial_accel(), libMesh::FirstOrderUnsteadySolver::compute_second_order_eqns(), libMesh::UnsteadySolver::du(), libMesh::EulerSolver::element_residual(), libMesh::Euler2Solver::element_residual(), libMesh::EigenTimeSolver::element_residual(), libMesh::SecondOrderUnsteadySolver::init(), libMesh::UnsteadySolver::init(), init(), libMesh::EigenTimeSolver::init(), libMesh::SecondOrderUnsteadySolver::init_data(), libMesh::UnsteadySolver::init_data(), init_data(), libMesh::EulerSolver::nonlocal_residual(), libMesh::Euler2Solver::nonlocal_residual(), libMesh::EigenTimeSolver::nonlocal_residual(), libMesh::UnsteadySolver::old_nonlinear_solution(), libMesh::SecondOrderUnsteadySolver::old_solution_accel(), libMesh::SecondOrderUnsteadySolver::old_solution_rate(), libMesh::NewmarkSolver::project_initial_accel(), libMesh::SecondOrderUnsteadySolver::project_initial_rate(), libMesh::SecondOrderUnsteadySolver::reinit(), libMesh::UnsteadySolver::reinit(), reinit(), libMesh::UnsteadySolver::retrieve_timestep(), libMesh::EigenTimeSolver::side_residual(), libMesh::TwostepTimeSolver::solve(), libMesh::UnsteadySolver::solve(), libMesh::EigenTimeSolver::solve(), and system().

◆ quiet

bool libMesh::TimeSolver::quiet

Print extra debugging information if quiet == false.

Definition at line 192 of file time_solver.h.

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

◆ reduce_deltat_on_diffsolver_failure

unsigned int libMesh::TimeSolver::reduce_deltat_on_diffsolver_failure

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
This has no effect for SteadySolvers.
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().

◆ solution_history

std::unique_ptr<SolutionHistory> libMesh::TimeSolver::solution_history
protected

A std::unique_ptr 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 265 of file time_solver.h.

Referenced by libMesh::UnsteadySolver::adjoint_advance_timestep(), libMesh::UnsteadySolver::advance_timestep(), libMesh::UnsteadySolver::retrieve_timestep(), and set_solution_history().


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