$extrastylesheet
libMesh::EulerSolver Class Reference

#include <euler_solver.h>

Inheritance diagram for libMesh::EulerSolver:

List of all members.

Public Types

typedef UnsteadySolver Parent
typedef DifferentiableSystem sys_type

Public Member Functions

 EulerSolver (sys_type &s)
virtual ~EulerSolver ()
virtual Real error_order () const
virtual bool element_residual (bool request_jacobian, DiffContext &)
virtual bool side_residual (bool request_jacobian, DiffContext &)
virtual bool nonlocal_residual (bool request_jacobian, DiffContext &)
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 ()
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 UniquePtr< DiffSolver > & diff_solver ()
virtual UniquePtr
< 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

Real theta
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 bool _general_residual (bool request_jacobian, DiffContext &, ResFuncType mass, ResFuncType time_deriv, ResFuncType constraint, ReinitFuncType reinit)
void increment_constructor_count (const std::string &name)
void increment_destructor_count (const std::string &name)

Protected Attributes

bool first_solve
bool first_adjoint_step
UniquePtr< DiffSolver_diff_solver
UniquePtr< LinearSolver< Number > > _linear_solver
sys_type_system
UniquePtr< 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 defines a theta-method Euler (defaulting to Backward Euler with theta = 1.0) solver to handle time integration of DifferentiableSystems.

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 2006

Definition at line 45 of file euler_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

Definition at line 51 of file euler_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.

The type of system

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

Definition at line 66 of file time_solver.h.


Constructor & Destructor Documentation

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

Definition at line 27 of file euler_solver.C.

  : UnsteadySolver(s), theta(1.)
{
}

Destructor.

Definition at line 34 of file euler_solver.C.

{
}

Member Function Documentation

bool libMesh::EulerSolver::_general_residual ( bool  request_jacobian,
DiffContext context,
ResFuncType  mass,
ResFuncType  time_deriv,
ResFuncType  constraint,
ReinitFuncType  reinit 
) [protected, virtual]

This method is the underlying implementation of the public residual methods.

Definition at line 89 of file euler_solver.C.

References libMesh::TimeSolver::_system, libMesh::DenseVector< T >::add(), libMesh::DifferentiableSystem::deltat, libMesh::DiffContext::elem_solution_derivative, libMesh::DiffContext::elem_solution_rate_derivative, libMesh::DiffContext::fixed_solution_derivative, libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_fixed_solution(), libMesh::DiffContext::get_elem_jacobian(), libMesh::DiffContext::get_elem_solution(), libMesh::DiffContext::get_elem_solution_rate(), libMesh::UnsteadySolver::old_nonlinear_solution(), libMesh::DenseVector< T >::size(), libMesh::DenseVector< T >::swap(), libMesh::DenseMatrix< T >::swap(), theta, and libMesh::System::use_fixed_solution.

Referenced by element_residual(), nonlocal_residual(), and side_residual().

{
  unsigned int n_dofs = context.get_elem_solution().size();

  // We might need to save the old jacobian in case one of our physics
  // terms later is unable to update it analytically.
  DenseMatrix<Number> old_elem_jacobian(n_dofs, n_dofs);
  if (request_jacobian)
    old_elem_jacobian.swap(context.get_elem_jacobian());

  // Local nonlinear solution at old timestep
  DenseVector<Number> old_elem_solution(n_dofs);
  for (unsigned int i=0; i != n_dofs; ++i)
    old_elem_solution(i) =
      old_nonlinear_solution(context.get_dof_indices()[i]);

  // Local time derivative of solution
  context.get_elem_solution_rate() = context.get_elem_solution();
  context.get_elem_solution_rate() -= old_elem_solution;
  context.elem_solution_rate_derivative = 1 / _system.deltat;
  context.get_elem_solution_rate() *=
    context.elem_solution_rate_derivative;

  // Local nonlinear solution at time t_theta
  DenseVector<Number> theta_solution(context.get_elem_solution());
  theta_solution *= theta;
  theta_solution.add(1. - theta, old_elem_solution);

  context.elem_solution_derivative = theta;
  context.fixed_solution_derivative = theta;

  // If a fixed solution is requested, we'll use theta_solution
  if (_system.use_fixed_solution)
    context.get_elem_fixed_solution() = theta_solution;

  // Move theta_->elem_, elem_->theta_
  context.get_elem_solution().swap(theta_solution);

  // Move the mesh into place first if necessary
  (context.*reinit_func)(theta);

  // Get the time derivative at t_theta
  bool jacobian_computed =
    (_system.*time_deriv)(request_jacobian, context);

  jacobian_computed = (_system.*mass)(jacobian_computed, context) &&
    jacobian_computed;

  // Restore the elem position if necessary
  (context.*reinit_func)(1);

  // Move elem_->elem_, theta_->theta_
  context.get_elem_solution().swap(theta_solution);
  context.elem_solution_derivative = 1;

  // Add the constraint term
  jacobian_computed = (_system.*constraint)(jacobian_computed, context) &&
    jacobian_computed;

  // Add back (or restore) the old jacobian
  if (request_jacobian)
    {
      if (jacobian_computed)
        context.get_elem_jacobian() += old_elem_jacobian;
      else
        context.get_elem_jacobian().swap(old_elem_jacobian);
    }

  return jacobian_computed;
}

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::UnsteadySolver::advance_timestep ( ) [virtual, inherited]

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::TimeSolver.

Reimplemented in libMesh::AdaptiveTimeSolver.

Definition at line 152 of file unsteady_solver.C.

References libMesh::TimeSolver::_system, libMesh::DifferentiableSystem::deltat, libMesh::UnsteadySolver::first_solve, 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::System::solution, libMesh::TimeSolver::solution_history, and libMesh::System::time.

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

{
  if (!first_solve)
    {
      // Store the solution, does nothing by default
      // User has to attach appropriate solution_history object for this to
      // actually store anything anywhere
      solution_history->store();

      _system.time += _system.deltat;
    }

  NumericVector<Number> &old_nonlinear_soln =
    _system.get_vector("_old_nonlinear_solution");
  NumericVector<Number> &nonlinear_solution =
    *(_system.solution);

  old_nonlinear_soln = nonlinear_solution;

  old_nonlinear_soln.localize
    (*old_local_nonlinear_solution,
     _system.get_dof_map().get_send_list());
}
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.

{}
virtual UniquePtr<DiffSolver>& libMesh::TimeSolver::diff_solver ( ) [inline, virtual, inherited]

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

Reimplemented in libMesh::AdaptiveTimeSolver.

Definition at line 183 of file time_solver.h.

References libMesh::TimeSolver::_diff_solver.

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

{ return _diff_solver; }
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::EulerSolver::element_residual ( bool  request_jacobian,
DiffContext context 
) [virtual]

This method uses the DifferentiablePhysics' element_time_derivative() and element_constraint() to build a full residual on an element. What combination it uses will depend on theta.

Implements libMesh::TimeSolver.

Definition at line 50 of file euler_solver.C.

References libMesh::DifferentiablePhysics::_eulerian_time_deriv(), _general_residual(), libMesh::DiffContext::elem_reinit(), libMesh::DifferentiablePhysics::element_constraint(), and libMesh::DifferentiablePhysics::mass_residual().

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;
}

Error convergence order: 2 for Crank-Nicolson, 1 otherwise

Implements libMesh::UnsteadySolver.

Definition at line 40 of file euler_solver.C.

References theta.

{
  if (theta == 0.5)
    return 2.;
  return 1.;
}
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::UnsteadySolver::init ( ) [virtual, inherited]

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

Reimplemented from libMesh::TimeSolver.

Reimplemented in libMesh::AdaptiveTimeSolver.

Definition at line 46 of file unsteady_solver.C.

References libMesh::TimeSolver::_system, and libMesh::System::add_vector().

{
  TimeSolver::init();

  _system.add_vector("_old_nonlinear_solution");
}
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.

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; }
virtual UniquePtr<LinearSolver<Number> >& libMesh::TimeSolver::linear_solver ( ) [inline, virtual, inherited]

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

Reimplemented in libMesh::AdaptiveTimeSolver.

Definition at line 188 of file time_solver.h.

References libMesh::TimeSolver::_linear_solver.

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

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

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

  { return _n_objects; }
bool libMesh::EulerSolver::nonlocal_residual ( bool  request_jacobian,
DiffContext context 
) [virtual]
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 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 _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().

void libMesh::UnsteadySolver::reinit ( ) [virtual, inherited]

The reinitialization function. This method is used to resize internal data vectors after a mesh change.

Reimplemented from libMesh::TimeSolver.

Reimplemented in libMesh::AdaptiveTimeSolver.

Definition at line 70 of file unsteady_solver.C.

References libMesh::TimeSolver::_system, libMesh::System::get_dof_map(), libMesh::DofMap::get_send_list(), libMesh::System::get_vector(), libMesh::GHOSTED, libMesh::NumericVector< T >::localize(), libMesh::System::n_dofs(), libMesh::System::n_local_dofs(), libMesh::UnsteadySolver::old_local_nonlinear_solution, and libMesh::SERIAL.

{
  TimeSolver::reinit();

#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

  // localize the old solution
  NumericVector<Number> &old_nonlinear_soln =
    _system.get_vector("_old_nonlinear_solution");

  old_nonlinear_soln.localize
    (*old_local_nonlinear_solution,
     _system.get_dof_map().get_send_list());
}
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::EulerSolver::side_residual ( bool  request_jacobian,
DiffContext context 
) [virtual]
void libMesh::UnsteadySolver::solve ( ) [virtual, inherited]

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::TimeSolver.

Reimplemented in libMesh::AdaptiveTimeSolver, and libMesh::TwostepTimeSolver.

Definition at line 93 of file unsteady_solver.C.

References libMesh::TimeSolver::_diff_solver, libMesh::TimeSolver::_system, libMesh::UnsteadySolver::advance_timestep(), libMesh::DifferentiableSystem::deltat, libMesh::DiffSolver::DIVERGED_BACKTRACKING_FAILURE, libMesh::DiffSolver::DIVERGED_MAX_NONLINEAR_ITERATIONS, libMesh::UnsteadySolver::first_solve, libMesh::out, libMesh::TimeSolver::quiet, and libMesh::TimeSolver::reduce_deltat_on_diffsolver_failure.

{
  if (first_solve)
    {
      advance_timestep();
      first_solve = false;
    }

  unsigned int solve_result = _diff_solver->solve();

  // If we requested the UnsteadySolver to attempt reducing dt after a
  // failed DiffSolver solve, check the results of the solve now.
  if (reduce_deltat_on_diffsolver_failure)
    {
      bool backtracking_failed =
        solve_result & DiffSolver::DIVERGED_BACKTRACKING_FAILURE;

      bool max_iterations =
        solve_result & DiffSolver::DIVERGED_MAX_NONLINEAR_ITERATIONS;

      if (backtracking_failed || max_iterations)
        {
          // Cut timestep in half
          for (unsigned int nr=0; nr<reduce_deltat_on_diffsolver_failure; ++nr)
            {
              _system.deltat *= 0.5;
              libMesh::out << "Newton backtracking failed.  Trying with smaller timestep, dt="
                           << _system.deltat << std::endl;

              solve_result = _diff_solver->solve();

              // Check solve results with reduced timestep
              bool backtracking_still_failed =
                solve_result & DiffSolver::DIVERGED_BACKTRACKING_FAILURE;

              bool backtracking_max_iterations =
                solve_result & DiffSolver::DIVERGED_MAX_NONLINEAR_ITERATIONS;

              if (!backtracking_still_failed && !backtracking_max_iterations)
                {
                  if (!quiet)
                    libMesh::out << "Reduced dt solve succeeded." << std::endl;
                  return;
                }
            }

          // If we made it here, we still couldn't converge the solve after
          // reducing deltat
          libMesh::out << "DiffSolver::solve() did not succeed after "
                       << reduce_deltat_on_diffsolver_failure
                       << " attempts." << std::endl;
          libmesh_convergence_failure();

        } // end if (backtracking_failed || max_iterations)
    } // end if (reduce_deltat_on_diffsolver_failure)
}
const sys_type& libMesh::TimeSolver::system ( ) const [inline, inherited]
Returns:
a constant reference to the system we are solving.

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]
Returns:
a writeable reference to the system we are solving.

Definition at line 178 of file time_solver.h.

References libMesh::TimeSolver::_system.

{ return _system; }

Member Data Documentation

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

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

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 libMesh::AdaptiveTimeSolver::advance_timestep(), libMesh::UnsteadySolver::advance_timestep(), libMesh::TwostepTimeSolver::solve(), and libMesh::UnsteadySolver::solve().

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

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

The value for the theta method to employ: 1.0 corresponds to backwards Euler, 0.0 corresponds to forwards Euler, 0.5 corresponds to Crank-Nicolson.

Definition at line 102 of file euler_solver.h.

Referenced by _general_residual(), and error_order().


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