$extrastylesheet
libMesh::EigenTimeSolver Class Reference

#include <eigen_time_solver.h>

Inheritance diagram for libMesh::EigenTimeSolver:

List of all members.

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 bool nonlocal_residual (bool get_jacobian, DiffContext &)
virtual Real du (const SystemNorm &) const
virtual bool is_steady () const
virtual void init_data ()
virtual void adjoint_advance_timestep ()
virtual void retrieve_timestep ()
virtual void before_timestep ()
const sys_typesystem () const
sys_typesystem ()
virtual UniquePtr< DiffSolver > & diff_solver ()
virtual UniquePtr
< LinearSolver< Number > > & 
linear_solver ()
void set_solution_history (const SolutionHistory &_solution_history)
bool is_adjoint () const
void set_is_adjoint (bool _is_adjoint_value)

Static Public Member Functions

static std::string get_info ()
static void print_info (std::ostream &out=libMesh::out)
static unsigned int n_objects ()
static void enable_print_counter_info ()
static void disable_print_counter_info ()

Public Attributes

UniquePtr< 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 bool(DifferentiablePhysics::* ResFuncType )(bool, DiffContext &)
typedef void(DiffContext::* ReinitFuncType )(Real)
typedef std::map< std::string,
std::pair< unsigned int,
unsigned int > > 
Counts

Protected Member Functions

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

Protected Attributes

UniquePtr< DiffSolver_diff_solver
UniquePtr< LinearSolver< Number > > _linear_solver
sys_type_system
UniquePtr< NumericVector
< Number > > 
old_local_nonlinear_solution
UniquePtr< SolutionHistorysolution_history

Static Protected Attributes

static Counts _counts
static Threads::atomic
< unsigned int > 
_n_objects
static Threads::spin_mutex _mutex
static bool _enable_print_counter = true

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 67 of file eigen_time_solver.h.


Member Typedef Documentation

typedef std::map<std::string, std::pair<unsigned int, unsigned int> > libMesh::ReferenceCounter::Counts [protected, inherited]

Data structure to log the information. The log is identified by the class name.

Definition at line 113 of file reference_counter.h.

The parent class

Definition at line 78 of file eigen_time_solver.h.

typedef void(DiffContext::* libMesh::TimeSolver::ReinitFuncType)(Real) [protected, inherited]

Definition at line 278 of file time_solver.h.

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

Definitions of argument types for use in refactoring subclasses.

Definition at line 276 of file time_solver.h.

The type of system

Reimplemented from libMesh::TimeSolver.

Definition at line 73 of file eigen_time_solver.h.


Member Enumeration Documentation

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 199 of file eigen_time_solver.h.


Constructor & Destructor Documentation

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

Definition at line 32 of file eigen_time_solver.C.

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

  : Parent(s),
    eigen_solver     (EigenSolver<Number>::build(s.comm())),
    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);
}

Destructor.

Definition at line 47 of file eigen_time_solver.C.

{
}

Member Function Documentation

void libMesh::TimeSolver::adjoint_advance_timestep ( ) [virtual, inherited]

This method advances the adjoint solution to the previous timestep, after an adjoint_solve() has been performed. This will be done before every UnsteadySolver::adjoint_solve().

Reimplemented in libMesh::UnsteadySolver.

Definition at line 106 of file time_solver.C.

{
}
virtual void libMesh::EigenTimeSolver::advance_timestep ( ) [inline, virtual]

It doesn't make sense to advance the timestep, so we shouldn't call this.

Reimplemented from libMesh::TimeSolver.

Definition at line 113 of file eigen_time_solver.h.

{ }
virtual void libMesh::TimeSolver::before_timestep ( ) [inline, virtual, inherited]

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

Definition at line 168 of file time_solver.h.

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

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

Reimplemented in libMesh::AdaptiveTimeSolver.

Definition at line 183 of file time_solver.h.

References libMesh::TimeSolver::_diff_solver.

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

{ return _diff_solver; }
virtual Real libMesh::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 libMesh::TimeSolver.

Definition at line 145 of file eigen_time_solver.h.

{ return 0.; }
bool libMesh::EigenTimeSolver::element_residual ( bool  get_jacobian,
DiffContext context 
) [virtual]

Forms either the spatial (Jacobian) or mass matrix part of the operator, depending on which is requested.

Implements libMesh::TimeSolver.

Definition at line 127 of file eigen_time_solver.C.

References libMesh::TimeSolver::_system, libMesh::DiffContext::elem_solution_derivative, libMesh::DiffContext::elem_solution_rate_derivative, libMesh::DifferentiablePhysics::element_constraint(), libMesh::DifferentiablePhysics::element_time_derivative(), libMesh::DiffContext::get_elem_jacobian(), libMesh::libmesh_assert(), libMesh::DifferentiablePhysics::mass_residual(), Matrix_A, Matrix_B, and now_assembling.

{
  // The EigenTimeSolver always computes jacobians!
  libmesh_assert (request_jacobian);

  context.elem_solution_rate_derivative = 1;
  context.elem_solution_derivative = 1;

  // 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 to get positive matrix from new negative
      // mass_residual convention
      context.get_elem_jacobian() *= -1.0;

      return mass_jacobian_computed;
    }

  else
    libmesh_error_msg("Unrecognized value now_assembling = " << now_assembling);
}

Methods to enable/disable the reference counter output from print_info()

Definition at line 100 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

{
  _enable_print_counter = true;
  return;
}
virtual Real libMesh::EigenTimeSolver::error_order ( ) const [inline, virtual]

error convergence order against deltat is not applicable to an eigenvalue problem.

Definition at line 119 of file eigen_time_solver.h.

{ return 0.; }
std::string libMesh::ReferenceCounter::get_info ( ) [static, inherited]

Gets a string containing the reference information.

Definition at line 47 of file reference_counter.C.

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

Referenced by libMesh::ReferenceCounter::print_info().

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

  std::ostringstream oss;

  oss << '\n'
      << " ---------------------------------------------------------------------------- \n"
      << "| Reference count information                                                |\n"
      << " ---------------------------------------------------------------------------- \n";

  for (Counts::iterator it = _counts.begin();
       it != _counts.end(); ++it)
    {
      const std::string name(it->first);
      const unsigned int creations    = it->second.first;
      const unsigned int destructions = it->second.second;

      oss << "| " << name << " reference count information:\n"
          << "|  Creations:    " << creations    << '\n'
          << "|  Destructions: " << destructions << '\n';
    }

  oss << " ---------------------------------------------------------------------------- \n";

  return oss.str();

#else

  return "";

#endif
}
void libMesh::ReferenceCounter::increment_constructor_count ( const std::string &  name) [inline, protected, inherited]

Increments the construction counter. Should be called in the constructor of any derived class that will be reference counted.

Definition at line 163 of file reference_counter.h.

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

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

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

  p.first++;
}
void libMesh::ReferenceCounter::increment_destructor_count ( const std::string &  name) [inline, protected, inherited]

Increments the destruction counter. Should be called in the destructor of any derived class that will be reference counted.

Definition at line 176 of file reference_counter.h.

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

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

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

  p.second++;
}
void libMesh::EigenTimeSolver::init ( ) [virtual]

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

Reimplemented from libMesh::TimeSolver.

Definition at line 56 of file eigen_time_solver.C.

References libMesh::TimeSolver::_system, libMesh::ImplicitSystem::add_matrix(), and libMesh::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");
}
void libMesh::TimeSolver::init_data ( ) [virtual, inherited]

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

Reimplemented in libMesh::UnsteadySolver.

Definition at line 77 of file time_solver.C.

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

{
  this->diff_solver()->init();

  if (libMesh::on_command_line("--solver_system_names"))
    this->linear_solver()->init((_system.name()+"_").c_str());
  else
    this->linear_solver()->init();
}
bool libMesh::TimeSolver::is_adjoint ( ) const [inline, inherited]

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

Definition at line 233 of file time_solver.h.

References libMesh::TimeSolver::_is_adjoint.

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

  { return _is_adjoint; }
virtual bool libMesh::EigenTimeSolver::is_steady ( ) const [inline, virtual]

This is effectively a steady-state solver.

Implements libMesh::TimeSolver.

Definition at line 150 of file eigen_time_solver.h.

{ return true; }
virtual UniquePtr<LinearSolver<Number> >& libMesh::TimeSolver::linear_solver ( ) [inline, virtual, inherited]

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

Reimplemented in libMesh::AdaptiveTimeSolver.

Definition at line 188 of file time_solver.h.

References libMesh::TimeSolver::_linear_solver.

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

{ return _linear_solver; }
static unsigned int libMesh::ReferenceCounter::n_objects ( ) [inline, static, inherited]

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

Definition at line 79 of file reference_counter.h.

References libMesh::ReferenceCounter::_n_objects.

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

  { return _n_objects; }
bool libMesh::EigenTimeSolver::nonlocal_residual ( bool  get_jacobian,
DiffContext context 
) [virtual]

Forms the jacobian of the nonlocal terms.

Implements libMesh::TimeSolver.

Definition at line 221 of file eigen_time_solver.C.

References libMesh::TimeSolver::_system, libMesh::DiffContext::get_elem_jacobian(), libMesh::libmesh_assert(), Matrix_A, Matrix_B, libMesh::DifferentiablePhysics::nonlocal_constraint(), libMesh::DifferentiablePhysics::nonlocal_mass_residual(), libMesh::DifferentiablePhysics::nonlocal_time_derivative(), and now_assembling.

{
  // 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.nonlocal_time_derivative(request_jacobian, context);

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

      bool jacobian_computed2 =
        _system.nonlocal_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.nonlocal_mass_residual(request_jacobian, context);

      // Scale Jacobian by -1 to get positive matrix from new negative
      // mass_residual convention
      context.get_elem_jacobian() *= -1.0;

      return mass_jacobian_computed;
    }

  else
    libmesh_error_msg("Unrecognized value now_assembling = " << now_assembling);
}
void libMesh::ReferenceCounter::print_info ( std::ostream &  out = libMesh::out) [static, inherited]

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

Definition at line 88 of file reference_counter.C.

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

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

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

Reimplemented from libMesh::TimeSolver.

Definition at line 51 of file eigen_time_solver.C.

{
  // empty...
}
void libMesh::TimeSolver::retrieve_timestep ( ) [virtual, inherited]

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

Reimplemented in libMesh::UnsteadySolver.

Definition at line 110 of file time_solver.C.

{
}
void libMesh::TimeSolver::set_is_adjoint ( bool  _is_adjoint_value) [inline, inherited]

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

Definition at line 240 of file time_solver.h.

References libMesh::TimeSolver::_is_adjoint.

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

  { _is_adjoint = _is_adjoint_value; }
void libMesh::TimeSolver::set_solution_history ( const SolutionHistory _solution_history) [inherited]

A setter function users will employ if they need to do something other than save no solution history

Definition at line 97 of file time_solver.C.

References libMesh::SolutionHistory::clone(), and libMesh::TimeSolver::solution_history.

{
  solution_history = _solution_history.clone();
}
bool libMesh::EigenTimeSolver::side_residual ( bool  get_jacobian,
DiffContext context 
) [virtual]

Forms the jacobian of the boundary terms.

Implements libMesh::TimeSolver.

Definition at line 174 of file eigen_time_solver.C.

References libMesh::TimeSolver::_system, libMesh::DiffContext::elem_solution_derivative, libMesh::DiffContext::elem_solution_rate_derivative, libMesh::DiffContext::get_elem_jacobian(), libMesh::libmesh_assert(), Matrix_A, Matrix_B, now_assembling, libMesh::DifferentiablePhysics::side_constraint(), libMesh::DifferentiablePhysics::side_mass_residual(), and libMesh::DifferentiablePhysics::side_time_derivative().

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

  context.elem_solution_rate_derivative = 1;
  context.elem_solution_derivative = 1;

  // 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);

      // Scale Jacobian by -1 to get positive matrix from new negative
      // mass_residual convention
      context.get_elem_jacobian() *= -1.0;

      return mass_jacobian_computed;
    }

  else
    libmesh_error_msg("Unrecognized value now_assembling = " << now_assembling);
}
void libMesh::EigenTimeSolver::solve ( ) [virtual]

Implements the assembly of both matrices A and B, and calls the EigenSolver to compute the eigenvalues.

Reimplemented from libMesh::TimeSolver.

Definition at line 68 of file eigen_time_solver.C.

References libMesh::TimeSolver::_system, libMesh::DifferentiableSystem::assembly(), eigen_solver, libMesh::ImplicitSystem::get_matrix(), libMesh::ImplicitSystem::matrix, Matrix_A, Matrix_B, maxits, n_basis_vectors_to_use, n_converged_eigenpairs, n_eigenpairs_to_compute, n_iterations_reqd, now_assembling, libMesh::out, libMesh::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)
    libMesh::out << "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)
    libMesh::out << "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)
    libMesh::out << "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;
}
const sys_type& libMesh::TimeSolver::system ( ) const [inline, inherited]
Returns:
a constant reference to the system we are solving.

Definition at line 173 of file time_solver.h.

References libMesh::TimeSolver::_system.

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

{ return _system; }
sys_type& libMesh::TimeSolver::system ( ) [inline, inherited]
Returns:
a writeable reference to the system we are solving.

Definition at line 178 of file time_solver.h.

References libMesh::TimeSolver::_system.

{ return _system; }

Member Data Documentation

UniquePtr<DiffSolver> libMesh::TimeSolver::_diff_solver [protected, inherited]

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

Definition at line 248 of file time_solver.h.

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

bool libMesh::ReferenceCounter::_enable_print_counter = true [static, protected, inherited]

Flag to control whether reference count information is printed when print_info is called.

Definition at line 137 of file reference_counter.h.

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

UniquePtr<LinearSolver<Number> > libMesh::TimeSolver::_linear_solver [protected, inherited]

An implicit linear solver to use for adjoint problems.

Definition at line 253 of file time_solver.h.

Referenced by libMesh::TimeSolver::linear_solver().

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 131 of file reference_counter.h.

Threads::atomic< unsigned int > libMesh::ReferenceCounter::_n_objects [static, protected, inherited]

The number of objects. Print the reference count information when the number returns to 0.

Definition at line 126 of file reference_counter.h.

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

The EigenSolver object. This is what actually makes the calls to SLEPc.

Definition at line 156 of file eigen_time_solver.h.

Referenced by EigenTimeSolver(), and solve().

The maximum number of iterations allowed to solve the problem.

Definition at line 167 of file eigen_time_solver.h.

Referenced by solve().

The 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 183 of file eigen_time_solver.h.

Referenced by solve().

After a solve, holds the number of eigenpairs successfully converged.

Definition at line 189 of file eigen_time_solver.h.

Referenced by solve().

The number of eigenvectors/values to be computed.

Definition at line 172 of file eigen_time_solver.h.

Referenced by solve().

After a solve, holds the number of iterations required to converge the requested number of eigenpairs.

Definition at line 195 of file eigen_time_solver.h.

Referenced by solve().

Flag which controls the internals of element_residual() and side_residual().

Definition at line 219 of file eigen_time_solver.h.

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

Serial vector of _system.get_vector("_old_nonlinear_solution")

Reimplemented in libMesh::UnsteadySolver.

Definition at line 263 of file time_solver.h.

bool libMesh::TimeSolver::quiet [inherited]

Print extra debugging information if quiet == false.

Definition at line 193 of file time_solver.h.

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

This value (which defaults to zero) is the number of times the TimeSolver is allowed to halve deltat and let the DiffSolver repeat the latest failed solve with a reduced timestep. Note that this has no effect for SteadySolvers. Note that you must set at least one of the DiffSolver flags "continue_after_max_iterations" or "continue_after_backtrack_failure" to allow the TimeSolver to retry the solve.

Definition at line 221 of file time_solver.h.

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

UniquePtr<SolutionHistory> libMesh::TimeSolver::solution_history [protected, inherited]

An UniquePtr to a SolutionHistory object. Default is NoSolutionHistory, which the user can override by declaring a different kind of SolutionHistory in the application

Definition at line 270 of file time_solver.h.

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

The linear solver tolerance to be used when solving the eigenvalue problem. FIXME: need more info...

Definition at line 162 of file eigen_time_solver.h.

Referenced by solve().


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