Poster of Linux kernelThe best gift for a Linux geek
EulerSolver

EulerSolver

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

NAME

EulerSolver -  

SYNOPSIS


#include <euler_solver.h>

Inherits UnsteadySolver.  

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

virtual void solve ()

virtual void advance_timestep ()

Number old_nonlinear_solution (const unsigned int global_dof_number) const

virtual Real du (const SystemNorm &norm) const

virtual void reinit ()

virtual void adjoint_solve ()

virtual void adjoint_recede_timestep ()

virtual void before_timestep ()

const sys_type & system () const

virtual AutoPtr< DiffSolver > & diff_solver ()

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


Real theta

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


sys_type & system ()

void increment_constructor_count (const std::string &name)

void increment_destructor_count (const std::string &name)
 

Protected Attributes


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 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 44 of file euler_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 UnsteadySolver EulerSolver::ParentThe parent class

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

 

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

Definition at line 27 of file euler_solver.C.

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

EulerSolver::~EulerSolver () [virtual]Destructor.

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

Reimplemented in AdaptiveTimeSolver.

Definition at line 124 of file unsteady_solver.C.

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

Referenced by UnsteadySolver::solve().

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

  old_nonlinear_solution = nonlinear_solution;

  if (!first_solve)
    _system.time += _system.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.

{}
 

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

Reimplemented in AdaptiveTimeSolver.

Definition at line 152 of file time_solver.h.

References TimeSolver::_diff_solver.

{ return _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 EulerSolver::element_residual (boolrequest_jacobian, DiffContext &context) [virtual]This method uses the DifferentiableSystem's element_time_derivative() and element_constraint() to build a full residual on an element. What combination it uses will depend on theta.

Implements TimeSolver.

Definition at line 50 of file euler_solver.C.

References TimeSolver::_system, DenseVector< T >::add(), DifferentiableSystem::deltat, DiffContext::dof_indices, DiffContext::elem_fixed_solution, DiffContext::elem_jacobian, DiffContext::elem_reinit(), DiffContext::elem_residual, DiffContext::elem_solution, DiffContext::elem_solution_derivative, DifferentiableSystem::element_constraint(), DifferentiableSystem::element_time_derivative(), DifferentiableSystem::eulerian_residual(), DiffContext::fixed_solution_derivative, DifferentiableSystem::mass_residual(), UnsteadySolver::old_nonlinear_solution(), DenseVector< T >::size(), DenseMatrix< T >::swap(), DenseVector< T >::swap(), theta, DifferentiableSystem::use_fixed_solution, DenseVector< T >::zero(), and DenseMatrix< T >::zero().

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

  // 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.dof_indices[i]);

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

  // Technically the elem_solution_derivative is either theta
  // or -1.0 in this implementation, but we scale the former part
  // ourselves
  context.elem_solution_derivative = 1.0;

// Technically the fixed_solution_derivative is always theta,
// but we're scaling a whole jacobian by theta after these first
// evaluations
  context.fixed_solution_derivative = 1.0;

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

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

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

  // We're going to compute just the change in elem_residual
  // (and possibly elem_jacobian), then add back the old values
  DenseVector<Number> old_elem_residual(context.elem_residual);
  DenseMatrix<Number> old_elem_jacobian;
  if (request_jacobian)
    {
      old_elem_jacobian = context.elem_jacobian;
      context.elem_jacobian.zero();
    }
  context.elem_residual.zero();

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

  // For a moving mesh problem we may need the pseudoconvection term too
  jacobian_computed =
    _system.eulerian_residual(jacobian_computed, context) && jacobian_computed;

  // Scale the time-dependent residual and jacobian correctly
  context.elem_residual *= _system.deltat;
  if (jacobian_computed)
    context.elem_jacobian *= (theta * _system.deltat);

  // The fixed_solution_derivative is always theta,
  // and now we're done scaling jacobians
  context.fixed_solution_derivative = theta;

  // We evaluate mass_residual with the change in solution
  // to get the mass matrix, reusing old_elem_solution to hold that
  // delta_solution.  We're solving dt*F(u) - du = 0, so our
  // delta_solution is old_solution - new_solution.
  // We're still keeping elem_solution in theta_solution for now
  old_elem_solution -= theta_solution;

  // Move old_->elem_, theta_->old_
  context.elem_solution.swap(old_elem_solution);

  // We do a trick here to avoid using a non-1
  // elem_solution_derivative:
  context.elem_jacobian *= -1.0;
  jacobian_computed = _system.mass_residual(jacobian_computed, context) &&
    jacobian_computed;
  context.elem_jacobian *= -1.0;

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

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

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

  // Add back the old residual and jacobian
  context.elem_residual += old_elem_residual;
  if (request_jacobian)
    {
      if (jacobian_computed)
        context.elem_jacobian += old_elem_jacobian;
      else
        context.elem_jacobian.swap(old_elem_jacobian);
    }

  return jacobian_computed;
}
 

Real EulerSolver::error_order () const [virtual]Error convergence order: 2 for Crank-Nicolson, 1 otherwise

Implements UnsteadySolver.

Definition at line 40 of file euler_solver.C.

References theta.

{
  if (theta == 0.5)
    return 2.;
  return 1.;
}
 

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 UnsteadySolver::init () [virtual, inherited]The initialization function. This method is used to initialize internal data structures before a simulation begins.

Reimplemented from TimeSolver.

Reimplemented in AdaptiveTimeSolver.

Definition at line 44 of file unsteady_solver.C.

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

{
  TimeSolver::init();

  _system.add_vector('_old_nonlinear_solution');
}
 

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(), element_residual(), Euler2Solver::element_residual(), FEMSystem::eulerian_residual(), 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 TimeSolver::reinit () [virtual, inherited]The reinitialization function. This method is used after changes in the mesh

Reimplemented in AdaptiveTimeSolver, and EigenTimeSolver.

Definition at line 46 of file time_solver.C.

References TimeSolver::_adjoint_solver, and TimeSolver::_diff_solver.

{
  _diff_solver->reinit();
  _adjoint_solver->reinit();
}
 

bool EulerSolver::side_residual (boolrequest_jacobian, DiffContext &context) [virtual]This method uses the DifferentiableSystem's side_time_derivative() and side_constraint() to build a full residual on an element's side. What combination it uses will depend on theta.

Implements TimeSolver.

Definition at line 156 of file euler_solver.C.

References TimeSolver::_system, DenseVector< T >::add(), DifferentiableSystem::deltat, DiffContext::dof_indices, DiffContext::elem_fixed_solution, DiffContext::elem_jacobian, DiffContext::elem_residual, DiffContext::elem_side_reinit(), DiffContext::elem_solution, DiffContext::elem_solution_derivative, DiffContext::fixed_solution_derivative, UnsteadySolver::old_nonlinear_solution(), DifferentiableSystem::side_constraint(), DifferentiableSystem::side_mass_residual(), DifferentiableSystem::side_time_derivative(), DenseVector< T >::size(), DenseMatrix< T >::swap(), DenseVector< T >::swap(), theta, DifferentiableSystem::use_fixed_solution, DenseVector< T >::zero(), and DenseMatrix< T >::zero().

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

  // 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.dof_indices[i]);

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

  // Technically the elem_solution_derivative is either theta
  // or 1.0 in this implementation, but we scale the former part
  // ourselves
  context.elem_solution_derivative = 1.0;

// Technically the fixed_solution_derivative is always theta,
// but we're scaling a whole jacobian by theta after these first
// evaluations
  context.fixed_solution_derivative = 1.0;

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

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

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

  // We're going to compute just the change in elem_residual
  // (and possibly elem_jacobian), then add back the old values
  DenseVector<Number> old_elem_residual(context.elem_residual);
  DenseMatrix<Number> old_elem_jacobian;
  if (request_jacobian)
    {
      old_elem_jacobian = context.elem_jacobian;
      context.elem_jacobian.zero();
    }
  context.elem_residual.zero();

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

  // Scale the time-dependent residual and jacobian correctly
  context.elem_residual *= _system.deltat;
  if (jacobian_computed)
    context.elem_jacobian *= (theta * _system.deltat);

  // The fixed_solution_derivative is always theta,
  // and now we're done scaling jacobians
  context.fixed_solution_derivative = theta;

  // We evaluate side_mass_residual with the change in solution
  // to get the mass matrix, reusing old_elem_solution to hold that
  // delta_solution.  We're solving dt*F(u) - du = 0, so our
  // delta_solution is old_solution - new_solution.
  // We're still keeping elem_solution in theta_solution for now
  old_elem_solution -= theta_solution;

  // Move old_->elem_, theta_->old_
  context.elem_solution.swap(old_elem_solution);

  // We do a trick here to avoid using a non-1
  // elem_solution_derivative:
  context.elem_jacobian *= -1.0;
  jacobian_computed = _system.side_mass_residual(jacobian_computed, context) &&
    jacobian_computed;
  context.elem_jacobian *= -1.0;

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

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

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

  // Add back the old residual and jacobian
  context.elem_residual += old_elem_residual;
  if (request_jacobian)
    {
      if (jacobian_computed)
        context.elem_jacobian += old_elem_jacobian;
      else
        context.elem_jacobian.swap(old_elem_jacobian);
    }

  return jacobian_computed;
}
 

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

Reimplemented in AdaptiveTimeSolver, and TwostepTimeSolver.

Definition at line 53 of file unsteady_solver.C.

References TimeSolver::_diff_solver, TimeSolver::_system, UnsteadySolver::advance_timestep(), DifferentiableSystem::deltat, DiffSolver::DIVERGED_BACKTRACKING_FAILURE, DiffSolver::DIVERGED_MAX_NONLINEAR_ITERATIONS, UnsteadySolver::first_solve, System::get_dof_map(), DofMap::get_send_list(), System::get_vector(), libMeshEnums::GHOSTED, NumericVector< T >::localize(), System::n_dofs(), System::n_local_dofs(), UnsteadySolver::old_local_nonlinear_solution, TimeSolver::quiet, TimeSolver::reduce_deltat_on_diffsolver_failure, and libMeshEnums::SERIAL.

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

#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

  _system.get_vector('_old_nonlinear_solution').localize
    (*old_local_nonlinear_solution,
     _system.get_dof_map().get_send_list());

  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;
              std::cout << '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_failed =
                solve_result & DiffSolver::DIVERGED_BACKTRACKING_FAILURE;
              
              bool max_iterations =
                solve_result & DiffSolver::DIVERGED_MAX_NONLINEAR_ITERATIONS;

              if (!backtracking_failed && !max_iterations)
                {
                  if (!quiet)
                    std::cout << 'Reduced dt solve succeeded.' << std::endl;
                  return;
                }
            }

          // If we made it here, we still couldn't converge the solve after
          // reducing deltat
          std::cout << '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)
}
 

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

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

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(), element_residual(), Euler2Solver::element_residual(), EigenTimeSolver::element_residual(), UnsteadySolver::init(), TimeSolver::init(), EigenTimeSolver::init(), UnsteadySolver::old_nonlinear_solution(), SteadySolver::side_residual(), side_residual(), Euler2Solver::side_residual(), EigenTimeSolver::side_residual(), UnsteadySolver::solve(), TwostepTimeSolver::solve(), EigenTimeSolver::solve(), and TimeSolver::system().  

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

Real EulerSolver::thetaThe 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 91 of file euler_solver.h.

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

 

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 UnsteadySolver EulerSolver::ParentThe parent class
typedef DifferentiableSystem TimeSolver::sys_type [inherited]The type of system
Constructor & Destructor Documentation
EulerSolver::EulerSolver (sys_type &s)Constructor. Requires a reference to the system to be solved.
EulerSolver::~EulerSolver () [virtual]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 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.
virtual void TimeSolver::before_timestep () [inline, virtual, inherited]This method is for subclasses or users to override to do arbitrary processing between timesteps
virtual AutoPtr<DiffSolver>& TimeSolver::diff_solver () [inline, 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 EulerSolver::element_residual (boolrequest_jacobian, DiffContext &context) [virtual]This method uses the DifferentiableSystem's element_time_derivative() and element_constraint() to build a full residual on an element. What combination it uses will depend on theta.
Real EulerSolver::error_order () const [virtual]Error convergence order: 2 for Crank-Nicolson, 1 otherwise
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 UnsteadySolver::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 TimeSolver::reinit () [virtual, inherited]The reinitialization function. This method is used after changes in the mesh
bool EulerSolver::side_residual (boolrequest_jacobian, DiffContext &context) [virtual]This method uses the DifferentiableSystem's side_time_derivative() and side_constraint() to build a full residual on an element's side. What combination it uses will depend on theta.
void 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.
sys_type& TimeSolver::system () [inline, protected, inherited]Returns:
const sys_type& TimeSolver::system () const [inline, 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.
bool UnsteadySolver::first_solve [protected, inherited]A bool that will be true the first time solve() is called, and false thereafter
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 EulerSolver::thetaThe 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.
Author

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