Poster of Linux kernelThe best gift for a Linux geek
TwostepTimeSolver

TwostepTimeSolver

Section: C Library Functions (3) Updated: Thu Apr 7 2011
Local index Up
 

NAME

TwostepTimeSolver -  

SYNOPSIS


#include <twostep_time_solver.h>

Inherits AdaptiveTimeSolver.  

Public Types


typedef AdaptiveTimeSolver Parent

typedef DifferentiableSystem sys_type
 

Public Member Functions


TwostepTimeSolver (sys_type &s)

~TwostepTimeSolver ()

void solve ()

virtual void init ()

virtual void reinit ()

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

Number old_nonlinear_solution (const unsigned int global_dof_number) const

virtual Real du (const SystemNorm &norm) const

virtual void adjoint_solve ()

virtual void adjoint_recede_timestep ()

virtual void before_timestep ()

const sys_type & system () const

virtual AutoPtr< DiffSolver > & adjoint_solver ()
 

Static Public Member Functions


static std::string get_info ()

static void print_info ()

static unsigned int n_objects ()
 

Public Attributes


AutoPtr< UnsteadySolver > core_time_solver

SystemNorm component_norm

std::vector< float > component_scale

Real target_tolerance

Real upper_tolerance

Real max_deltat

Real min_deltat

Real max_growth

bool global_tolerance

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 > &)

sys_type & system ()

void increment_constructor_count (const std::string &name)

void increment_destructor_count (const std::string &name)
 

Protected Attributes


Real last_deltat

bool first_solve

AutoPtr< DiffSolver > _diff_solver

AutoPtr< DiffSolver > _adjoint_solver

sys_type & _system
 

Static Protected Attributes


static Counts _counts

static Threads::atomic< unsigned int > _n_objects

static Threads::spin_mutex _mutex
 

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 twostep_time_solver.h.  

Member Typedef Documentation

 

typedef std::map<std::string, std::pair<unsigned int, unsigned int> > ReferenceCounter::Counts [protected, inherited]Data structure to log the information. The log is identified by the class name.

Definition at line 105 of file reference_counter.h.  

typedef AdaptiveTimeSolver TwostepTimeSolver::ParentThe parent class

Reimplemented from AdaptiveTimeSolver.

Definition at line 57 of file twostep_time_solver.h.  

typedef DifferentiableSystem TimeSolver::sys_type [inherited]The type of system

Reimplemented in EigenTimeSolver, and SteadySolver.

Definition at line 62 of file time_solver.h.  

Constructor & Destructor Documentation

 

TwostepTimeSolver::TwostepTimeSolver (sys_type &s)Constructor. Requires a reference to the system to be solved.

Definition at line 28 of file twostep_time_solver.C.

References AdaptiveTimeSolver::core_time_solver, and AutoPtr< Tp >::reset().

  : AdaptiveTimeSolver(s)

{
  // We start with a reasonable time solver: implicit Euler
  core_time_solver.reset(new EulerSolver(s));
}
 

TwostepTimeSolver::~TwostepTimeSolver ()Destructor.

Definition at line 38 of file twostep_time_solver.C.

{
}
 

Member Function Documentation

 

void TimeSolver::adjoint_recede_timestep () [virtual, inherited]This method advances the adjoint solution to the previous timestep, after an adjoint_solve() has been performed. This will probably be done after every UnsteadySolver::adjoint_solve().

Definition at line 90 of file time_solver.C.

{
}
 

void TimeSolver::adjoint_solve () [virtual, inherited]This method solves for the adjoint solution at the previous timestep (or solves a steady adjoint problem). Usually we will only need to solve one linear system per timestep, but more complex subclasses may override this.

Definition at line 77 of file time_solver.C.

References TimeSolver::_adjoint_solver.

{
  _adjoint_solver->adjoint_solve();
}
 

virtual AutoPtr<DiffSolver>& TimeSolver::adjoint_solver () [inline, virtual, inherited]An implicit linear solver to use for adjoint problems.

Definition at line 157 of file time_solver.h.

References TimeSolver::_adjoint_solver.

{ return _adjoint_solver; }
 

void AdaptiveTimeSolver::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 UnsteadySolver.

Definition at line 86 of file adaptive_time_solver.C.

References TimeSolver::_system, UnsteadySolver::first_solve, System::get_vector(), AdaptiveTimeSolver::last_deltat, UnsteadySolver::old_nonlinear_solution(), System::solution, and DifferentiableSystem::time.

{
  NumericVector<Number> &old_nonlinear_solution =
  _system.get_vector('_old_nonlinear_solution');
  NumericVector<Number> &nonlinear_solution =
    *(_system.solution);
//    _system.get_vector('_nonlinear_solution');

  old_nonlinear_solution = nonlinear_solution;

  if (!first_solve)
    _system.time += last_deltat;
}
 

virtual void TimeSolver::before_timestep () [inline, virtual, inherited]This method is for subclasses or users to override to do arbitrary processing between timesteps

Definition at line 142 of file time_solver.h.

{}
 

Real AdaptiveTimeSolver::calculate_norm (System &s, NumericVector< Number > &v) [protected, virtual, inherited]A helper function to calculate error norms

Definition at line 138 of file adaptive_time_solver.C.

References System::calculate_norm(), and AdaptiveTimeSolver::component_norm.

Referenced by solve().

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

AutoPtr< DiffSolver > & AdaptiveTimeSolver::diff_solver () [virtual, inherited]An implicit linear or nonlinear solver to use at each timestep.

Reimplemented from TimeSolver.

Definition at line 131 of file adaptive_time_solver.C.

References AdaptiveTimeSolver::core_time_solver.

{
  return core_time_solver->diff_solver();
}
 

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

Definition at line 150 of file unsteady_solver.C.

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

{

  AutoPtr<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 AdaptiveTimeSolver::element_residual (boolget_jacobian, DiffContext &context) [virtual, inherited]This method is passed on to the core_time_solver

Implements TimeSolver.

Definition at line 111 of file adaptive_time_solver.C.

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

{
  libmesh_assert(core_time_solver.get());

  return core_time_solver->element_residual(request_jacobian, context);
}
 

Real AdaptiveTimeSolver::error_order () const [virtual, inherited]This method is passed on to the core_time_solver

Implements UnsteadySolver.

Definition at line 102 of file adaptive_time_solver.C.

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

{
  libmesh_assert(core_time_solver.get());

  return core_time_solver->error_order();
}
 

std::string ReferenceCounter::get_info () [static, inherited]Gets a string containing the reference information.

Definition at line 45 of file reference_counter.C.

References ReferenceCounter::_counts, and Quality::name().

Referenced by ReferenceCounter::print_info().

{
#if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)

  std::ostringstream out;
  
  out << '
      << ' ---------------------------------------------------------------------------- 
      << '| Reference count information                                                |
      << ' ---------------------------------------------------------------------------- ;
  
  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;

      out << '| ' << name << ' reference count information:
          << '|  Creations:    ' << creations    << '
          << '|  Destructions: ' << destructions << ';
    }
  
  out << ' ---------------------------------------------------------------------------- ;

  return out.str();

#else

  return '';
  
#endif
}
 

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

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

Referenced by ReferenceCountedObject< Value >::ReferenceCountedObject().

{
  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
  std::pair<unsigned int, unsigned int>& p = _counts[name];

  p.first++;
}
 

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

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

Referenced by ReferenceCountedObject< Value >::~ReferenceCountedObject().

{
  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
  std::pair<unsigned int, unsigned int>& p = _counts[name];

  p.second++;
}
 

void AdaptiveTimeSolver::init () [virtual, inherited]The initialization function. This method is used to initialize internal data structures before a simulation begins.

Reimplemented from UnsteadySolver.

Definition at line 57 of file adaptive_time_solver.C.

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

{
  libmesh_assert(core_time_solver.get());

  // We override this because our core_time_solver is the one that
  // needs to handle new vectors, diff_solver->init(), etc
  core_time_solver->init();

  // As an UnsteadySolver, we have an old_local_nonlinear_solution, but it
  // isn't pointing to the right place - fix it
  //
  // This leaves us with two AutoPtrs holding the same pointer - dangerous
  // for future use.  Replace with shared_ptr?
  old_local_nonlinear_solution = 
    AutoPtr<NumericVector<Number> >(core_time_solver->old_local_nonlinear_solution.get());
}
 

static unsigned int ReferenceCounter::n_objects () [inline, static, inherited]Prints the number of outstanding (created, but not yet destroyed) objects.

Definition at line 76 of file reference_counter.h.

References ReferenceCounter::_n_objects.

Referenced by System::read_serialized_blocked_dof_objects(), and System::write_serialized_blocked_dof_objects().

  { return _n_objects; }
 

Number UnsteadySolver::old_nonlinear_solution (const unsigned intglobal_dof_number) const [inherited]Returns:

the old nonlinear solution for the specified global DOF.

Definition at line 139 of file unsteady_solver.C.

References TimeSolver::_system, System::get_dof_map(), DofMap::n_dofs(), and UnsteadySolver::old_local_nonlinear_solution.

Referenced by UnsteadySolver::advance_timestep(), AdaptiveTimeSolver::advance_timestep(), EulerSolver::element_residual(), Euler2Solver::element_residual(), FEMSystem::eulerian_residual(), EulerSolver::side_residual(), and Euler2Solver::side_residual().

{
  libmesh_assert (global_dof_number < _system.get_dof_map().n_dofs());
  libmesh_assert (global_dof_number < old_local_nonlinear_solution->size());

  return (*old_local_nonlinear_solution)(global_dof_number);
}
 

void ReferenceCounter::print_info () [static, inherited]Prints the reference information to std::cout.

Definition at line 83 of file reference_counter.C.

References ReferenceCounter::get_info().

{
#if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
  
  std::cout << ReferenceCounter::get_info();
  
#endif
}
 

void AdaptiveTimeSolver::reinit () [virtual, inherited]The reinitialization function. This method is used after changes in the mesh

Reimplemented from TimeSolver.

Definition at line 76 of file adaptive_time_solver.C.

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

{
  libmesh_assert(core_time_solver.get());

  // We override this because our core_time_solver is the one that
  // needs to handle new vectors, diff_solver->reinit(), etc
  core_time_solver->reinit();
}
 

bool AdaptiveTimeSolver::side_residual (boolget_jacobian, DiffContext &context) [virtual, inherited]This method is passed on to the core_time_solver

Implements TimeSolver.

Definition at line 121 of file adaptive_time_solver.C.

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

{
  libmesh_assert(core_time_solver.get());

  return core_time_solver->side_residual(request_jacobian, context);
}
 

void TwostepTimeSolver::solve () [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.

Implements AdaptiveTimeSolver.

Definition at line 44 of file twostep_time_solver.C.

References TimeSolver::_system, AdaptiveTimeSolver::calculate_norm(), NumericVector< T >::clone(), AdaptiveTimeSolver::core_time_solver, DifferentiableSystem::deltat, UnsteadySolver::first_solve, AutoPtr< Tp >::get(), System::get_vector(), AdaptiveTimeSolver::global_tolerance, AdaptiveTimeSolver::last_deltat, std::max(), AdaptiveTimeSolver::max_deltat, AdaptiveTimeSolver::max_growth, AdaptiveTimeSolver::min_deltat, Utility::pow(), TimeSolver::quiet, TimeSolver::reduce_deltat_on_diffsolver_failure, System::solution, AdaptiveTimeSolver::target_tolerance, DifferentiableSystem::time, and AdaptiveTimeSolver::upper_tolerance.

{
  libmesh_assert(core_time_solver.get());

  // The core_time_solver will handle any first_solve actions
  first_solve = false;

  // We may have to repeat timesteps entirely if our error is bad
  // enough
  bool max_tolerance_met = false;

  // Calculating error values each time
  Real single_norm(0.), double_norm(0.), error_norm(0.),
       relative_error(0.);
 
  while (!max_tolerance_met)
    {
      // If we've been asked to reduce deltat if necessary, make sure
      // the core timesolver does so
      core_time_solver->reduce_deltat_on_diffsolver_failure =
        this->reduce_deltat_on_diffsolver_failure;

      if (!quiet)
        {
          std::cout << '=== Computing adaptive timestep === ' 
                    << std::endl;
        }

      // Use the double-length timestep first (so the
      // old_nonlinear_solution won't have to change)
      core_time_solver->solve();

      // Save a copy of the double-length nonlinear solution
      // and the old nonlinear solution
      AutoPtr<NumericVector<Number> > double_solution =
        _system.solution->clone();
      AutoPtr<NumericVector<Number> > old_solution =
        _system.get_vector('_old_nonlinear_solution').clone();

      double_norm = calculate_norm(_system, *double_solution);
      if (!quiet)
        {
          std::cout << 'Double norm = ' << double_norm << std::endl;
        }

      // Then reset the initial guess for our single-length calcs
      *(_system.solution) = _system.get_vector('_old_nonlinear_solution');

      // Call two single-length timesteps
      // Be sure that the core_time_solver does not change the
      // timestep here.  (This is unlikely because it just succeeded
      // with a timestep twice as large!)
      // FIXME: even if diffsolver failure is unlikely, we ought to
      // do *something* if it happens
      core_time_solver->reduce_deltat_on_diffsolver_failure = 0;
  
      Real old_time = _system.time;
      Real old_deltat = _system.deltat;
      _system.deltat *= 0.5;
      core_time_solver->solve();
      core_time_solver->advance_timestep();
      core_time_solver->solve();

      single_norm = calculate_norm(_system, *_system.solution);
      if (!quiet)
        {
          std::cout << 'Single norm = ' << single_norm << std::endl;
        }

      // Reset the core_time_solver's reduce_deltat... value.
      core_time_solver->reduce_deltat_on_diffsolver_failure =
        this->reduce_deltat_on_diffsolver_failure;
  
      // But then back off just in case our advance_timestep() isn't
      // called.
      // FIXME: this probably doesn't work with multistep methods
      _system.get_vector('_old_nonlinear_solution') = *old_solution;
      _system.time = old_time;
      _system.deltat = old_deltat;

      // Find the relative error
      *double_solution -= *(_system.solution);
      error_norm  = calculate_norm(_system, *double_solution);
      relative_error = error_norm / _system.deltat /
        std::max(double_norm, single_norm);

      // If the relative error makes no sense, we're done
      if (!double_norm && !single_norm)
        return;

      if (!quiet)
        {
          std::cout << 'Error norm = ' << error_norm << std::endl;
          std::cout << 'Local relative error = '
                    << (error_norm /
                        std::max(double_norm, single_norm))
                    << std::endl;
          std::cout << 'Global relative error = '
                    << (error_norm / _system.deltat / 
                        std::max(double_norm, single_norm)) 
                    << std::endl;
          std::cout << 'old delta t = ' << _system.deltat << std::endl;
        }

      // If we haven't met our upper error tolerance, we'll have to
      // repeat this timestep entirely
      if (this->upper_tolerance && relative_error > this->upper_tolerance)
        {
          // Reset the initial guess for our next try
          *(_system.solution) =
            _system.get_vector('_old_nonlinear_solution');

          // Chop delta t in half
          _system.deltat /= 2.;

          if (!quiet)
            {
              std::cout << 'Failed to meet upper error tolerance' 
                        << std::endl;
              std::cout << 'Retrying with delta t = '
                        << _system.deltat << std::endl;
            }
        }
      else
        max_tolerance_met = true;
    }

  
  // Otherwise, compare the relative error to the tolerance
  // and adjust deltat
  last_deltat = _system.deltat;

  const Real global_shrink_or_growth_factor =
    std::pow(this->target_tolerance / relative_error,
             1. / core_time_solver->error_order());

  const Real local_shrink_or_growth_factor =
    std::pow(this->target_tolerance /
             (error_norm/std::max(double_norm, single_norm)),
             1. / (core_time_solver->error_order()+1.));

  if (!quiet)
    {
      std::cout << 'The global growth/shrink factor is: '
                << global_shrink_or_growth_factor << std::endl;
      std::cout << 'The local growth/shrink factor is: '
                << local_shrink_or_growth_factor << std::endl;
    }

  // The local s.o.g. factor is based on the expected **local**
  // truncation error for the timestepping method, the global
  // s.o.g. factor is based on the method's **global** truncation
  // error.  You can shrink/grow the timestep to attempt to satisfy
  // either a global or local time-discretization error tolerance.
 
  Real shrink_or_growth_factor =
    this->global_tolerance ? global_shrink_or_growth_factor :
                             local_shrink_or_growth_factor;

  if (this->max_growth && this->max_growth < shrink_or_growth_factor)
    {
      if (!quiet && this->global_tolerance)
        {
          std::cout << 'delta t is constrained by max_growth' << std::endl;
        }
        shrink_or_growth_factor = this->max_growth;
    }

  _system.deltat *= shrink_or_growth_factor;
  
  // Restrict deltat to max-allowable value if necessary
  if ((this->max_deltat != 0.0) && (_system.deltat > this->max_deltat))
    {
      if (!quiet)
        {
          std::cout << 'delta t is constrained by maximum-allowable delta t.'
                    << std::endl;
        }
      _system.deltat = this->max_deltat;
    }

  // Restrict deltat to min-allowable value if necessary
  if ((this->min_deltat != 0.0) && (_system.deltat < this->min_deltat))
    {
      if (!quiet)
        {
          std::cout << 'delta t is constrained by minimum-allowable delta t.'
                    << std::endl;
        }
      _system.deltat = this->min_deltat;
    }
  
  if (!quiet)
    {
      std::cout << 'new delta t = ' << _system.deltat << std::endl;
    }
}
 

const sys_type& TimeSolver::system () const [inline, inherited]Returns:

a constant reference to the system we are solving.

Definition at line 147 of file time_solver.h.

References TimeSolver::_system.

{ return _system; }
 

sys_type& TimeSolver::system () [inline, protected, inherited]Returns:

a writeable reference to the system we are solving.

Definition at line 202 of file time_solver.h.

References TimeSolver::_system.

{ return _system; }
 

Member Data Documentation

 

AutoPtr<DiffSolver> TimeSolver::_adjoint_solver [protected, inherited]An implicit linear solver to use for adjoint problems.

Definition at line 197 of file time_solver.h.

Referenced by TimeSolver::adjoint_solve(), TimeSolver::adjoint_solver(), TimeSolver::init(), and TimeSolver::reinit().  

ReferenceCounter::Counts ReferenceCounter::_counts [static, protected, inherited]Actually holds the data.

Definition at line 110 of file reference_counter.h.

Referenced by ReferenceCounter::get_info(), ReferenceCounter::increment_constructor_count(), and ReferenceCounter::increment_destructor_count().  

AutoPtr<DiffSolver> TimeSolver::_diff_solver [protected, inherited]An implicit linear or nonlinear solver to use at each timestep.

Definition at line 192 of file time_solver.h.

Referenced by TimeSolver::diff_solver(), TimeSolver::init(), TimeSolver::reinit(), UnsteadySolver::solve(), and TimeSolver::solve().  

Threads::spin_mutex ReferenceCounter::_mutex [static, protected, inherited]Mutual exclusion object to enable thread-safe reference counting.

Definition at line 123 of file reference_counter.h.  

Threads::atomic< unsigned int > ReferenceCounter::_n_objects [static, protected, inherited]The number of objects. Print the reference count information when the number returns to 0.

Definition at line 118 of file reference_counter.h.

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

sys_type& TimeSolver::_system [protected, inherited]A reference to the system we are solving.

Definition at line 207 of file time_solver.h.

Referenced by UnsteadySolver::advance_timestep(), AdaptiveTimeSolver::advance_timestep(), UnsteadySolver::du(), SteadySolver::element_residual(), EulerSolver::element_residual(), Euler2Solver::element_residual(), EigenTimeSolver::element_residual(), UnsteadySolver::init(), TimeSolver::init(), EigenTimeSolver::init(), UnsteadySolver::old_nonlinear_solution(), SteadySolver::side_residual(), EulerSolver::side_residual(), Euler2Solver::side_residual(), EigenTimeSolver::side_residual(), UnsteadySolver::solve(), solve(), EigenTimeSolver::solve(), and TimeSolver::system().  

SystemNorm AdaptiveTimeSolver::component_norm [inherited]Error calculations are done in this norm, DISCRETE_L2 by default.

Definition at line 106 of file adaptive_time_solver.h.

Referenced by AdaptiveTimeSolver::calculate_norm().  

std::vector<float> AdaptiveTimeSolver::component_scale [inherited]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 114 of file adaptive_time_solver.h.  

AutoPtr<UnsteadySolver> AdaptiveTimeSolver::core_time_solver [inherited]This object is used to take timesteps

Definition at line 101 of file adaptive_time_solver.h.

Referenced by AdaptiveTimeSolver::diff_solver(), AdaptiveTimeSolver::element_residual(), AdaptiveTimeSolver::error_order(), AdaptiveTimeSolver::init(), AdaptiveTimeSolver::reinit(), AdaptiveTimeSolver::side_residual(), solve(), and TwostepTimeSolver().  

bool UnsteadySolver::first_solve [protected, inherited]A bool that will be true the first time solve() is called, and false thereafter

Reimplemented from TimeSolver.

Definition at line 124 of file unsteady_solver.h.

Referenced by UnsteadySolver::advance_timestep(), AdaptiveTimeSolver::advance_timestep(), UnsteadySolver::solve(), and solve().  

bool AdaptiveTimeSolver::global_tolerance [inherited]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 170 of file adaptive_time_solver.h.

Referenced by solve().  

Real AdaptiveTimeSolver::last_deltat [protected, inherited]We need to store the value of the last deltat used, so that advance_timestep() will increment the system time correctly.

Definition at line 179 of file adaptive_time_solver.h.

Referenced by AdaptiveTimeSolver::advance_timestep(), and solve().  

Real AdaptiveTimeSolver::max_deltat [inherited]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 140 of file adaptive_time_solver.h.

Referenced by solve().  

Real AdaptiveTimeSolver::max_growth [inherited]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 154 of file adaptive_time_solver.h.

Referenced by solve().  

Real AdaptiveTimeSolver::min_deltat [inherited]Do not allow the adaptive time solver to select deltat < min_deltat. The default value is 0.0.

Definition at line 146 of file adaptive_time_solver.h.

Referenced by solve().  

AutoPtr<NumericVector<Number> > UnsteadySolver::old_local_nonlinear_solution [inherited]Serial vector of _system.get_vector('_old_nonlinear_solution')

Reimplemented from TimeSolver.

Definition at line 105 of file unsteady_solver.h.

Referenced by AdaptiveTimeSolver::AdaptiveTimeSolver(), AdaptiveTimeSolver::init(), UnsteadySolver::old_nonlinear_solution(), UnsteadySolver::solve(), and AdaptiveTimeSolver::~AdaptiveTimeSolver().  

bool TimeSolver::quiet [inherited]Print extra debugging information if quiet == false.

Definition at line 162 of file time_solver.h.

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

unsigned int 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 185 of file time_solver.h.

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

Real AdaptiveTimeSolver::target_tolerance [inherited]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 125 of file adaptive_time_solver.h.

Referenced by solve().  

Real AdaptiveTimeSolver::upper_tolerance [inherited]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 134 of file adaptive_time_solver.h.

Referenced by solve().

 

Author

Generated automatically by Doxygen for libMesh from the source code.


 

Index

NAME
SYNOPSIS
Public Types
Public Member Functions
Static Public Member Functions
Public Attributes
Protected Types
Protected Member Functions
Protected Attributes
Static Protected Attributes
Detailed Description
Member Typedef Documentation
typedef std::map<std::string, std::pair<unsigned int, unsigned int> > ReferenceCounter::Counts [protected, inherited]Data structure to log the information. The log is identified by the class name.
typedef AdaptiveTimeSolver TwostepTimeSolver::ParentThe parent class
typedef DifferentiableSystem TimeSolver::sys_type [inherited]The type of system
Constructor & Destructor Documentation
TwostepTimeSolver::TwostepTimeSolver (sys_type &s)Constructor. Requires a reference to the system to be solved.
TwostepTimeSolver::~TwostepTimeSolver ()Destructor.
Member Function Documentation
void TimeSolver::adjoint_recede_timestep () [virtual, inherited]This method advances the adjoint solution to the previous timestep, after an adjoint_solve() has been performed. This will probably be done after every UnsteadySolver::adjoint_solve().
void TimeSolver::adjoint_solve () [virtual, inherited]This method solves for the adjoint solution at the previous timestep (or solves a steady adjoint problem). Usually we will only need to solve one linear system per timestep, but more complex subclasses may override this.
virtual AutoPtr<DiffSolver>& TimeSolver::adjoint_solver () [inline, virtual, inherited]An implicit linear solver to use for adjoint problems.
void AdaptiveTimeSolver::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.
virtual void TimeSolver::before_timestep () [inline, virtual, inherited]This method is for subclasses or users to override to do arbitrary processing between timesteps
Real AdaptiveTimeSolver::calculate_norm (System &s, NumericVector< Number > &v) [protected, virtual, inherited]A helper function to calculate error norms
AutoPtr< DiffSolver > & AdaptiveTimeSolver::diff_solver () [virtual, inherited]An implicit linear or nonlinear solver to use at each timestep.
Real UnsteadySolver::du (const SystemNorm &norm) const [virtual, inherited]Computes the size of ||u^{n+1} - u^{n}|| in some norm.
bool AdaptiveTimeSolver::element_residual (boolget_jacobian, DiffContext &context) [virtual, inherited]This method is passed on to the core_time_solver
Real AdaptiveTimeSolver::error_order () const [virtual, inherited]This method is passed on to the core_time_solver
std::string ReferenceCounter::get_info () [static, inherited]Gets a string containing the reference information.
void 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.
void 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.
void AdaptiveTimeSolver::init () [virtual, inherited]The initialization function. This method is used to initialize internal data structures before a simulation begins.
static unsigned int ReferenceCounter::n_objects () [inline, static, inherited]Prints the number of outstanding (created, but not yet destroyed) objects.
Number UnsteadySolver::old_nonlinear_solution (const unsigned intglobal_dof_number) const [inherited]Returns:
void ReferenceCounter::print_info () [static, inherited]Prints the reference information to std::cout.
void AdaptiveTimeSolver::reinit () [virtual, inherited]The reinitialization function. This method is used after changes in the mesh
bool AdaptiveTimeSolver::side_residual (boolget_jacobian, DiffContext &context) [virtual, inherited]This method is passed on to the core_time_solver
void TwostepTimeSolver::solve () [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.
const sys_type& TimeSolver::system () const [inline, inherited]Returns:
sys_type& TimeSolver::system () [inline, protected, inherited]Returns:
Member Data Documentation
AutoPtr<DiffSolver> TimeSolver::_adjoint_solver [protected, inherited]An implicit linear solver to use for adjoint problems.
ReferenceCounter::Counts ReferenceCounter::_counts [static, protected, inherited]Actually holds the data.
AutoPtr<DiffSolver> TimeSolver::_diff_solver [protected, inherited]An implicit linear or nonlinear solver to use at each timestep.
Threads::spin_mutex ReferenceCounter::_mutex [static, protected, inherited]Mutual exclusion object to enable thread-safe reference counting.
Threads::atomic< unsigned int > ReferenceCounter::_n_objects [static, protected, inherited]The number of objects. Print the reference count information when the number returns to 0.
sys_type& TimeSolver::_system [protected, inherited]A reference to the system we are solving.
SystemNorm AdaptiveTimeSolver::component_norm [inherited]Error calculations are done in this norm, DISCRETE_L2 by default.
std::vector<float> AdaptiveTimeSolver::component_scale [inherited]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.
AutoPtr<UnsteadySolver> AdaptiveTimeSolver::core_time_solver [inherited]This object is used to take timesteps
bool UnsteadySolver::first_solve [protected, inherited]A bool that will be true the first time solve() is called, and false thereafter
bool AdaptiveTimeSolver::global_tolerance [inherited]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.
Real AdaptiveTimeSolver::last_deltat [protected, inherited]We need to store the value of the last deltat used, so that advance_timestep() will increment the system time correctly.
Real AdaptiveTimeSolver::max_deltat [inherited]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.
Real AdaptiveTimeSolver::max_growth [inherited]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.
Real AdaptiveTimeSolver::min_deltat [inherited]Do not allow the adaptive time solver to select deltat < min_deltat. The default value is 0.0.
AutoPtr<NumericVector<Number> > UnsteadySolver::old_local_nonlinear_solution [inherited]Serial vector of _system.get_vector('_old_nonlinear_solution')
bool TimeSolver::quiet [inherited]Print extra debugging information if quiet == false.
unsigned int 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.
Real AdaptiveTimeSolver::target_tolerance [inherited]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.
Real AdaptiveTimeSolver::upper_tolerance [inherited]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.
Author

This document was created by man2html, using the manual pages.
Time: 21:57:52 GMT, April 16, 2011