Poster of Linux kernelThe best gift for a Linux geek
EigenTimeSolver

EigenTimeSolver

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

NAME

EigenTimeSolver -  

SYNOPSIS


#include <eigen_time_solver.h>

Inherits TimeSolver.  

Public Types


typedef DifferentiableSystem sys_type

typedef TimeSolver Parent
 

Public Member Functions


EigenTimeSolver (sys_type &s)

virtual ~EigenTimeSolver ()

virtual void init ()

virtual void reinit ()

virtual void solve ()

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 Real du (const SystemNorm &) const

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


AutoPtr< EigenSolver< Number > > eigen_solver

Real tol

unsigned int maxits

unsigned int n_eigenpairs_to_compute

unsigned int n_basis_vectors_to_use

unsigned int n_converged_eigenpairs

unsigned int n_iterations_reqd

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


AutoPtr< DiffSolver > _diff_solver

AutoPtr< DiffSolver > _adjoint_solver

sys_type & _system

bool first_solve

AutoPtr< NumericVector< Number > > old_local_nonlinear_solution
 

Static Protected Attributes


static Counts _counts

static Threads::atomic< unsigned int > _n_objects

static Threads::spin_mutex _mutex
 

Private Types


enum NowAssembling { Matrix_A, Matrix_B, Invalid_Matrix }
 

Private Attributes


NowAssembling now_assembling
 

Detailed Description

The name of this class is confusing...it's meant to refer to the base class (TimeSolver) while still telling one that it's for solving (generalized) EigenValue problems that arise from finite element discretizations. For a time-dependent problem du/dt=F(u), with a steady solution 0=F(u_0), we look at the time evolution of a small perturbation, p=u-u_0, for which the (linearized) governing equation is

dp/dt = F'(u_0)p

where F'(u_0) is the Jacobian. The generalized eigenvalue problem arises by considering perturbations of the general form p = exp(lambda*t)x, which leads to

Ax = lambda*Bx

where A is the (discretized by FEM) Jacobian matrix and B is the (discretized by FEM) mass matrix.

The EigenSystem class (by Steffen Petersen) is related but does not fall under the FEMSystem paradigm invented by Roy Stogner. The EigenSolver class (also by Steffen) is meant to provide a generic 'linear solver' interface for EigenValue software. The only current concrete implementation is a SLEPc-based eigensolver class, which we make use of here as well.

Author:

John W. Peterson 2007

Definition at line 66 of file eigen_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 TimeSolver EigenTimeSolver::ParentThe parent class

Definition at line 77 of file eigen_time_solver.h.  

typedef DifferentiableSystem EigenTimeSolver::sys_typeThe type of system

Reimplemented from TimeSolver.

Definition at line 72 of file eigen_time_solver.h.  

Member Enumeration Documentation

 

enum EigenTimeSolver::NowAssembling [private]

Enumerator:

Matrix_A
The matrix associated with the spatial part of the operator.
Matrix_B
The matrix associated with the time derivative (mass matrix).
Invalid_Matrix
The enum is in an invalid state.

Definition at line 186 of file eigen_time_solver.h.

                     {
    Matrix_A,

    Matrix_B,

    Invalid_Matrix
  };
 

Constructor & Destructor Documentation

 

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

Definition at line 32 of file eigen_time_solver.C.

References eigen_solver, libMeshEnums::GHEP, and libMeshEnums::LARGEST_MAGNITUDE.

  : Parent(s),
    eigen_solver     (EigenSolver<Number>::build()),
    tol(pow(TOLERANCE, 5./3.)),
    maxits(1000),
    n_eigenpairs_to_compute(5),
    n_basis_vectors_to_use(3*n_eigenpairs_to_compute),
    n_converged_eigenpairs(0),    
    n_iterations_reqd(0)
{
  libmesh_experimental();
  eigen_solver->set_eigenproblem_type(GHEP);//or GNHEP
  eigen_solver->set_position_of_spectrum(LARGEST_MAGNITUDE);
}
 

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

Definition at line 47 of file eigen_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; }
 

virtual void EigenTimeSolver::advance_timestep () [inline, virtual]It doesn't make sense to advance the timestep, so we shouldn't call this.

Reimplemented from TimeSolver.

Definition at line 111 of file eigen_time_solver.h.

{ }
 

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

virtual Real EigenTimeSolver::du (const SystemNorm &) const [inline, virtual]Nominally computes the size of the difference between successive solution iterates ||u^{n+1} - u^{n}|| in some norm, but for this class just returns 0.

Implements TimeSolver.

Definition at line 137 of file eigen_time_solver.h.

{ return 0.; }
 

bool EigenTimeSolver::element_residual (boolget_jacobian, DiffContext &context) [virtual]Forms either the spatial (Jacobian) or mass matrix part of the operator, depending on which is requested.

Implements TimeSolver.

Definition at line 127 of file eigen_time_solver.C.

References TimeSolver::_system, DifferentiableSystem::element_constraint(), DifferentiableSystem::element_time_derivative(), DifferentiableSystem::mass_residual(), Matrix_A, Matrix_B, and now_assembling.

{
  // The EigenTimeSolver always computes jacobians!
  libmesh_assert (request_jacobian);
  
  // Assemble the operator for the spatial part.
  if (now_assembling == Matrix_A)
    {
      bool jacobian_computed =
        _system.element_time_derivative(request_jacobian, context);

      // The user shouldn't compute a jacobian unless requested
      libmesh_assert(request_jacobian || !jacobian_computed);

      bool jacobian_computed2 =
        _system.element_constraint(jacobian_computed, context);

      // The user shouldn't compute a jacobian unless requested
      libmesh_assert (jacobian_computed || !jacobian_computed2);

      return jacobian_computed && jacobian_computed2;
  
    }

  // Assemble the mass matrix operator
  else if (now_assembling == Matrix_B)
    {
      bool mass_jacobian_computed =
        _system.mass_residual(request_jacobian, context);

      // Scale Jacobian by -1?
      //context.elem_jacobian *= -1.0;

      return mass_jacobian_computed;
    }

  else
    {
      libmesh_error();
      return false;
    }
}
 

virtual Real EigenTimeSolver::error_order () const [inline, virtual]error convergence order against deltat is not applicable to an eigenvalue problem.

Definition at line 117 of file eigen_time_solver.h.

{ return 0.; }
 

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

Reimplemented from TimeSolver.

Definition at line 56 of file eigen_time_solver.C.

References TimeSolver::_system, ImplicitSystem::add_matrix(), and ImplicitSystem::have_matrix().

{
  // Add matrix 'B' to _system if not already there.
  // The user may have already added a matrix 'B' before
  // calling the System initialization.  This would be
  // necessary if e.g. the System originally started life
  // with a different type of TimeSolver and only later
  // had its TimeSolver changed to an EigenTimeSolver.
  if (!_system.have_matrix('B'))
    _system.add_matrix('B');
}
 

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

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 EigenTimeSolver::reinit () [virtual]The reinitialization function. This method is used after changes in the mesh

Reimplemented from TimeSolver.

Definition at line 51 of file eigen_time_solver.C.

{
  // empty...
}
 

bool EigenTimeSolver::side_residual (boolget_jacobian, DiffContext &context) [virtual]Forms the jacobian of the boundary terms.

Implements TimeSolver.

Definition at line 173 of file eigen_time_solver.C.

References TimeSolver::_system, Matrix_A, Matrix_B, now_assembling, DifferentiableSystem::side_constraint(), DifferentiableSystem::side_mass_residual(), and DifferentiableSystem::side_time_derivative().

{
  // The EigenTimeSolver always requests jacobians?
  //libmesh_assert (request_jacobian);

  // Assemble the operator for the spatial part.
  if (now_assembling == Matrix_A)
    {
      bool jacobian_computed =
        _system.side_time_derivative(request_jacobian, context);

      // The user shouldn't compute a jacobian unless requested
      libmesh_assert (request_jacobian || !jacobian_computed);

      bool jacobian_computed2 =
        _system.side_constraint(jacobian_computed, context);

      // The user shouldn't compute a jacobian unless requested
      libmesh_assert (jacobian_computed || !jacobian_computed2);
  
      return jacobian_computed && jacobian_computed2;
  
    }

  // There is now a 'side' equivalent for the mass matrix
  else if (now_assembling == Matrix_B)
    {
      bool mass_jacobian_computed =
        _system.side_mass_residual(request_jacobian, context);

      return mass_jacobian_computed;
    }

  else
    {
      libmesh_error();
      return false;
    }
}
 

void EigenTimeSolver::solve () [virtual]Implements the assembly of both matrices A and B, and calls the EigenSolver to compute the eigenvalues.

Reimplemented from TimeSolver.

Definition at line 68 of file eigen_time_solver.C.

References TimeSolver::_system, DifferentiableSystem::assembly(), eigen_solver, ImplicitSystem::get_matrix(), ImplicitSystem::matrix, Matrix_A, Matrix_B, maxits, n_basis_vectors_to_use, n_converged_eigenpairs, n_eigenpairs_to_compute, n_iterations_reqd, now_assembling, TimeSolver::quiet, and tol.

{
  // The standard implementation is basically to call:
  // _diff_solver->solve();
  // which internally assembles (when necessary) matrices and vectors
  // and calls linear solver software while also doing Newton steps (see newton_solver.C)
  //
  // The element_residual and side_residual functions below control
  // what happens in the interior of the element assembly loops.
  // We have a system reference, so it's possible to call _system.assembly()
  // ourselves if we want to...
  //
  // Interestingly, for the EigenSolver we don't need residuals...just Jacobians.
  // The Jacobian should therefore always be requested, and always return
  // jacobian_computed as being true.

  // The basic plan of attack is:
  // .) Construct the Jacobian using _system.assembly(true,true) as we
  //    would for a steady system.  Use a flag in this class to
  //    control behavior in element_residual and side_residual
  // .) Swap _system.matrix to matrix 'B' (be sure to add this extra matrix during init)
  // .) Call _system.assembly(true,true) again, using the flag in element_residual
  //    and side_residual to only get the mass matrix terms.
  // .) Send A and B to Steffen's EigenSolver interface.

  // Assemble the spatial part (matrix A) of the operator
  if (!this->quiet)
    std::cout << 'Assembling matrix A.' << std::endl;
  _system.matrix =   &( _system.get_matrix ('System Matrix') );
  this->now_assembling = Matrix_A;
  _system.assembly(true, true);
  //_system.matrix->print_matlab('matrix_A.m');
  
  // Point the system's matrix at B, call assembly again.
  if (!this->quiet)
    std::cout << 'Assembling matrix B.' << std::endl;
  _system.matrix =   &( _system.get_matrix ('B') );
  this->now_assembling = Matrix_B;
  _system.assembly(true, true);
  //_system.matrix->print_matlab('matrix_B.m');
  
  // Send matrices A, B to Steffen's SlepcEigenSolver interface
  //libmesh_here();
  if (!this->quiet)
    std::cout << 'Calling the EigenSolver.' << std::endl;
  std::pair<unsigned int, unsigned int> solve_data =
    eigen_solver->solve_generalized (_system.get_matrix ('System Matrix'),
                                     _system.get_matrix ('B'),
                                     n_eigenpairs_to_compute,
                                     n_basis_vectors_to_use,
                                     tol,
                                     maxits);
  
  this->n_converged_eigenpairs = solve_data.first;
  this->n_iterations_reqd      = solve_data.second;
}
 

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

AutoPtr<EigenSolver<Number> > EigenTimeSolver::eigen_solverThe EigenSolver object. This is what actually makes the calls to SLEPc.

Definition at line 143 of file eigen_time_solver.h.

Referenced by EigenTimeSolver(), and solve().  

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

Reimplemented in UnsteadySolver.

Definition at line 213 of file time_solver.h.  

unsigned int EigenTimeSolver::maxitsThe maximum number of iterations allowed to solve the problem.

Definition at line 154 of file eigen_time_solver.h.

Referenced by solve().  

unsigned int EigenTimeSolver::n_basis_vectors_to_useThe number of basis vectors to use in the computation. According to ex16, the number of basis vectors must be >= the number of eigenpairs requested, and ncv >= 2*nev is recommended. Increasing this number, even by a little bit, can _greatly_ reduce the number of (EigenSolver) iterations required to compute the desired number of eigenpairs, but the _cost per iteration_ goes up drastically as well.

Definition at line 170 of file eigen_time_solver.h.

Referenced by solve().  

unsigned int EigenTimeSolver::n_converged_eigenpairsAfter a solve, holds the number of eigenpairs successfully converged.

Definition at line 176 of file eigen_time_solver.h.

Referenced by solve().  

unsigned int EigenTimeSolver::n_eigenpairs_to_computeThe number of eigenvectors/values to be computed.

Definition at line 159 of file eigen_time_solver.h.

Referenced by solve().  

unsigned int EigenTimeSolver::n_iterations_reqdAfter a solve, holds the number of iterations required to converge the requested number of eigenpairs.

Definition at line 182 of file eigen_time_solver.h.

Referenced by solve().  

NowAssembling EigenTimeSolver::now_assembling [private]Flag which controls the internals of element_residual() and side_residual().

Definition at line 206 of file eigen_time_solver.h.

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

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

Reimplemented in UnsteadySolver.

Definition at line 218 of file time_solver.h.  

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 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 EigenTimeSolver::tolThe linear solver tolerance to be used when solving the eigenvalue problem. FIXME: need more info...

Definition at line 149 of file eigen_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
Private Types
Private 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 TimeSolver EigenTimeSolver::ParentThe parent class
typedef DifferentiableSystem EigenTimeSolver::sys_typeThe type of system
Member Enumeration Documentation
enum EigenTimeSolver::NowAssembling [private]
Constructor & Destructor Documentation
EigenTimeSolver::EigenTimeSolver (sys_type &s)Constructor. Requires a reference to the system to be solved.
EigenTimeSolver::~EigenTimeSolver () [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.
virtual void EigenTimeSolver::advance_timestep () [inline, virtual]It doesn't make sense to advance the timestep, so we shouldn't call this.
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.
virtual Real EigenTimeSolver::du (const SystemNorm &) const [inline, virtual]Nominally computes the size of the difference between successive solution iterates ||u^{n+1} - u^{n}|| in some norm, but for this class just returns 0.
bool EigenTimeSolver::element_residual (boolget_jacobian, DiffContext &context) [virtual]Forms either the spatial (Jacobian) or mass matrix part of the operator, depending on which is requested.
virtual Real EigenTimeSolver::error_order () const [inline, virtual]error convergence order against deltat is not applicable to an eigenvalue problem.
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 EigenTimeSolver::init () [virtual]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.
void ReferenceCounter::print_info () [static, inherited]Prints the reference information to std::cout.
void EigenTimeSolver::reinit () [virtual]The reinitialization function. This method is used after changes in the mesh
bool EigenTimeSolver::side_residual (boolget_jacobian, DiffContext &context) [virtual]Forms the jacobian of the boundary terms.
void EigenTimeSolver::solve () [virtual]Implements the assembly of both matrices A and B, and calls the EigenSolver to compute the eigenvalues.
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.
AutoPtr<EigenSolver<Number> > EigenTimeSolver::eigen_solverThe EigenSolver object. This is what actually makes the calls to SLEPc.
bool TimeSolver::first_solve [protected, inherited]A bool that will be true the first time solve() is called, and false thereafter
unsigned int EigenTimeSolver::maxitsThe maximum number of iterations allowed to solve the problem.
unsigned int EigenTimeSolver::n_basis_vectors_to_useThe number of basis vectors to use in the computation. According to ex16, the number of basis vectors must be >= the number of eigenpairs requested, and ncv >= 2*nev is recommended. Increasing this number, even by a little bit, can _greatly_ reduce the number of (EigenSolver) iterations required to compute the desired number of eigenpairs, but the _cost per iteration_ goes up drastically as well.
unsigned int EigenTimeSolver::n_converged_eigenpairsAfter a solve, holds the number of eigenpairs successfully converged.
unsigned int EigenTimeSolver::n_eigenpairs_to_computeThe number of eigenvectors/values to be computed.
unsigned int EigenTimeSolver::n_iterations_reqdAfter a solve, holds the number of iterations required to converge the requested number of eigenpairs.
NowAssembling EigenTimeSolver::now_assembling [private]Flag which controls the internals of element_residual() and side_residual().
AutoPtr<NumericVector<Number> > TimeSolver::old_local_nonlinear_solution [protected, 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 EigenTimeSolver::tolThe linear solver tolerance to be used when solving the eigenvalue problem. FIXME: need more info...
Author

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