libMesh::EigenTimeSolver Class Reference

#include <eigen_time_solver.h>

Inheritance diagram for libMesh::EigenTimeSolver:

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 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 AutoPtr< DiffSolver > & diff_solver ()
 
virtual AutoPtr< 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

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

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

Protected Attributes

AutoPtr< DiffSolver_diff_solver
 
AutoPtr< LinearSolver< Number > > _linear_solver
 
sys_type_system
 
bool first_solve
 
AutoPtr< NumericVector< Number > > old_local_nonlinear_solution
 
AutoPtr< 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
protectedinherited

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.

The type of system

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

193  {
197  Matrix_A,
198 
202  Matrix_B,
203 
208  };

Constructor & Destructor Documentation

libMesh::EigenTimeSolver::EigenTimeSolver ( sys_type s)
explicit

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.

33  : Parent(s),
35  tol(pow(TOLERANCE, 5./3.)),
36  maxits(1000),
41 {
42  libmesh_experimental();
43  eigen_solver->set_eigenproblem_type(GHEP);//or GNHEP
44  eigen_solver->set_position_of_spectrum(LARGEST_MAGNITUDE);
45 }
libMesh::EigenTimeSolver::~EigenTimeSolver ( )
virtual

Destructor.

Definition at line 47 of file eigen_time_solver.C.

48 {
49 }

Member Function Documentation

void libMesh::TimeSolver::adjoint_advance_timestep ( )
virtualinherited

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 100 of file time_solver.C.

101 {
102 }
virtual void libMesh::EigenTimeSolver::advance_timestep ( )
inlinevirtual

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.

113 { }
virtual void libMesh::TimeSolver::before_timestep ( )
inlinevirtualinherited

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

Definition at line 152 of file time_solver.h.

152 {}
virtual AutoPtr<DiffSolver>& libMesh::TimeSolver::diff_solver ( )
inlinevirtualinherited

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

Reimplemented in libMesh::AdaptiveTimeSolver.

Definition at line 167 of file time_solver.h.

References libMesh::TimeSolver::_diff_solver.

167 { return _diff_solver; }
void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited

Definition at line 106 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

107 {
108  _enable_print_counter = false;
109  return;
110 }
virtual Real libMesh::EigenTimeSolver::du ( const SystemNorm ) const
inlinevirtual

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

139 { 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::DifferentiablePhysics::element_constraint(), libMesh::DifferentiablePhysics::element_time_derivative(), libMesh::libmesh_assert(), libMesh::DifferentiablePhysics::mass_residual(), Matrix_A, Matrix_B, and now_assembling.

129 {
130  // The EigenTimeSolver always computes jacobians!
131  libmesh_assert (request_jacobian);
132 
133  // Assemble the operator for the spatial part.
134  if (now_assembling == Matrix_A)
135  {
136  bool jacobian_computed =
137  _system.element_time_derivative(request_jacobian, context);
138 
139  // The user shouldn't compute a jacobian unless requested
140  libmesh_assert(request_jacobian || !jacobian_computed);
141 
142  bool jacobian_computed2 =
143  _system.element_constraint(jacobian_computed, context);
144 
145  // The user shouldn't compute a jacobian unless requested
146  libmesh_assert (jacobian_computed || !jacobian_computed2);
147 
148  return jacobian_computed && jacobian_computed2;
149 
150  }
151 
152  // Assemble the mass matrix operator
153  else if (now_assembling == Matrix_B)
154  {
155  bool mass_jacobian_computed =
156  _system.mass_residual(request_jacobian, context);
157 
158  // Scale Jacobian by -1?
159  //context.elem_jacobian *= -1.0;
160 
161  return mass_jacobian_computed;
162  }
163 
164  else
165  {
166  libmesh_error();
167  return false;
168  }
169 }
void libMesh::ReferenceCounter::enable_print_counter_info ( )
staticinherited

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.

101 {
102  _enable_print_counter = true;
103  return;
104 }
virtual Real libMesh::EigenTimeSolver::error_order ( ) const
inlinevirtual

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

Definition at line 119 of file eigen_time_solver.h.

119 { return 0.; }
std::string libMesh::ReferenceCounter::get_info ( )
staticinherited

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

48 {
49 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
50 
51  std::ostringstream oss;
52 
53  oss << '\n'
54  << " ---------------------------------------------------------------------------- \n"
55  << "| Reference count information |\n"
56  << " ---------------------------------------------------------------------------- \n";
57 
58  for (Counts::iterator it = _counts.begin();
59  it != _counts.end(); ++it)
60  {
61  const std::string name(it->first);
62  const unsigned int creations = it->second.first;
63  const unsigned int destructions = it->second.second;
64 
65  oss << "| " << name << " reference count information:\n"
66  << "| Creations: " << creations << '\n'
67  << "| Destructions: " << destructions << '\n';
68  }
69 
70  oss << " ---------------------------------------------------------------------------- \n";
71 
72  return oss.str();
73 
74 #else
75 
76  return "";
77 
78 #endif
79 }
void libMesh::ReferenceCounter::increment_constructor_count ( const std::string &  name)
inlineprotectedinherited

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

164 {
165  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
166  std::pair<unsigned int, unsigned int>& p = _counts[name];
167 
168  p.first++;
169 }
void libMesh::ReferenceCounter::increment_destructor_count ( const std::string &  name)
inlineprotectedinherited

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

177 {
178  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
179  std::pair<unsigned int, unsigned int>& p = _counts[name];
180 
181  p.second++;
182 }
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().

57 {
58  // Add matrix "B" to _system if not already there.
59  // The user may have already added a matrix "B" before
60  // calling the System initialization. This would be
61  // necessary if e.g. the System originally started life
62  // with a different type of TimeSolver and only later
63  // had its TimeSolver changed to an EigenTimeSolver.
64  if (!_system.have_matrix("B"))
65  _system.add_matrix("B");
66 }
void libMesh::TimeSolver::init_data ( )
virtualinherited

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.

78 {
79 }
bool libMesh::TimeSolver::is_adjoint ( ) const
inlineinherited

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

Definition at line 217 of file time_solver.h.

References libMesh::TimeSolver::_is_adjoint.

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

218  { return _is_adjoint; }
virtual bool libMesh::EigenTimeSolver::is_steady ( ) const
inlinevirtual

This is effectively a steady-state solver.

Implements libMesh::TimeSolver.

Definition at line 144 of file eigen_time_solver.h.

144 { return true; }
virtual AutoPtr<LinearSolver<Number> >& libMesh::TimeSolver::linear_solver ( )
inlinevirtualinherited

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

Definition at line 172 of file time_solver.h.

References libMesh::TimeSolver::_linear_solver.

172 { return _linear_solver; }
static unsigned int libMesh::ReferenceCounter::n_objects ( )
inlinestaticinherited

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.

80  { return _n_objects; }
void libMesh::ReferenceCounter::print_info ( std::ostream &  out = libMesh::out)
staticinherited

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

89 {
91 }
void libMesh::EigenTimeSolver::reinit ( )
virtual

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.

52 {
53  // empty...
54 }
void libMesh::TimeSolver::retrieve_timestep ( )
virtualinherited

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

Reimplemented in libMesh::UnsteadySolver.

Definition at line 104 of file time_solver.C.

105 {
106 }
void libMesh::TimeSolver::set_is_adjoint ( bool  _is_adjoint_value)
inlineinherited

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

Definition at line 224 of file time_solver.h.

References libMesh::TimeSolver::_is_adjoint.

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

225  { _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 91 of file time_solver.C.

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

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

Forms the jacobian of the boundary terms.

Implements libMesh::TimeSolver.

Definition at line 173 of file eigen_time_solver.C.

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

175 {
176  // The EigenTimeSolver always requests jacobians?
177  //libmesh_assert (request_jacobian);
178 
179  // Assemble the operator for the spatial part.
180  if (now_assembling == Matrix_A)
181  {
182  bool jacobian_computed =
183  _system.side_time_derivative(request_jacobian, context);
184 
185  // The user shouldn't compute a jacobian unless requested
186  libmesh_assert (request_jacobian || !jacobian_computed);
187 
188  bool jacobian_computed2 =
189  _system.side_constraint(jacobian_computed, context);
190 
191  // The user shouldn't compute a jacobian unless requested
192  libmesh_assert (jacobian_computed || !jacobian_computed2);
193 
194  return jacobian_computed && jacobian_computed2;
195 
196  }
197 
198  // There is now a "side" equivalent for the mass matrix
199  else if (now_assembling == Matrix_B)
200  {
201  bool mass_jacobian_computed =
202  _system.side_mass_residual(request_jacobian, context);
203 
204  return mass_jacobian_computed;
205  }
206 
207  else
208  {
209  libmesh_error();
210  return false;
211  }
212 }
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.

69 {
70  // The standard implementation is basically to call:
71  // _diff_solver->solve();
72  // which internally assembles (when necessary) matrices and vectors
73  // and calls linear solver software while also doing Newton steps (see newton_solver.C)
74  //
75  // The element_residual and side_residual functions below control
76  // what happens in the interior of the element assembly loops.
77  // We have a system reference, so it's possible to call _system.assembly()
78  // ourselves if we want to...
79  //
80  // Interestingly, for the EigenSolver we don't need residuals...just Jacobians.
81  // The Jacobian should therefore always be requested, and always return
82  // jacobian_computed as being true.
83 
84  // The basic plan of attack is:
85  // .) Construct the Jacobian using _system.assembly(true,true) as we
86  // would for a steady system. Use a flag in this class to
87  // control behavior in element_residual and side_residual
88  // .) Swap _system.matrix to matrix "B" (be sure to add this extra matrix during init)
89  // .) Call _system.assembly(true,true) again, using the flag in element_residual
90  // and side_residual to only get the mass matrix terms.
91  // .) Send A and B to Steffen's EigenSolver interface.
92 
93  // Assemble the spatial part (matrix A) of the operator
94  if (!this->quiet)
95  libMesh::out << "Assembling matrix A." << std::endl;
96  _system.matrix = &( _system.get_matrix ("System Matrix") );
97  this->now_assembling = Matrix_A;
98  _system.assembly(true, true);
99  //_system.matrix->print_matlab("matrix_A.m");
100 
101  // Point the system's matrix at B, call assembly again.
102  if (!this->quiet)
103  libMesh::out << "Assembling matrix B." << std::endl;
104  _system.matrix = &( _system.get_matrix ("B") );
105  this->now_assembling = Matrix_B;
106  _system.assembly(true, true);
107  //_system.matrix->print_matlab("matrix_B.m");
108 
109  // Send matrices A, B to Steffen's SlepcEigenSolver interface
110  //libmesh_here();
111  if (!this->quiet)
112  libMesh::out << "Calling the EigenSolver." << std::endl;
113  std::pair<unsigned int, unsigned int> solve_data =
114  eigen_solver->solve_generalized (_system.get_matrix ("System Matrix"),
115  _system.get_matrix ("B"),
118  tol,
119  maxits);
120 
121  this->n_converged_eigenpairs = solve_data.first;
122  this->n_iterations_reqd = solve_data.second;
123 }
const sys_type& libMesh::TimeSolver::system ( ) const
inlineinherited
Returns
a constant reference to the system we are solving.

Definition at line 157 of file time_solver.h.

References libMesh::TimeSolver::_system.

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

157 { return _system; }
sys_type& libMesh::TimeSolver::system ( )
inlineinherited
Returns
a writeable reference to the system we are solving.

Definition at line 162 of file time_solver.h.

References libMesh::TimeSolver::_system.

162 { return _system; }

Member Data Documentation

ReferenceCounter::Counts libMesh::ReferenceCounter::_counts
staticprotectedinherited
AutoPtr<DiffSolver> libMesh::TimeSolver::_diff_solver
protectedinherited

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

Definition at line 232 of file time_solver.h.

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

bool libMesh::ReferenceCounter::_enable_print_counter = true
staticprotectedinherited

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

AutoPtr<LinearSolver<Number> > libMesh::TimeSolver::_linear_solver
protectedinherited

An implicit linear solver to use for adjoint problems.

Definition at line 237 of file time_solver.h.

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

Threads::spin_mutex libMesh::ReferenceCounter::_mutex
staticprotectedinherited

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
staticprotectedinherited

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

AutoPtr<EigenSolver<Number> > libMesh::EigenTimeSolver::eigen_solver

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

Definition at line 150 of file eigen_time_solver.h.

Referenced by EigenTimeSolver(), and solve().

bool libMesh::TimeSolver::first_solve
protectedinherited

A bool that will be true the first time solve() is called, and false thereafter

Definition at line 248 of file time_solver.h.

unsigned int libMesh::EigenTimeSolver::maxits

The maximum number of iterations allowed to solve the problem.

Definition at line 161 of file eigen_time_solver.h.

Referenced by solve().

unsigned int libMesh::EigenTimeSolver::n_basis_vectors_to_use

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

Referenced by solve().

unsigned int libMesh::EigenTimeSolver::n_converged_eigenpairs

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

Definition at line 183 of file eigen_time_solver.h.

Referenced by solve().

unsigned int libMesh::EigenTimeSolver::n_eigenpairs_to_compute

The number of eigenvectors/values to be computed.

Definition at line 166 of file eigen_time_solver.h.

Referenced by solve().

unsigned int libMesh::EigenTimeSolver::n_iterations_reqd

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

Definition at line 189 of file eigen_time_solver.h.

Referenced by solve().

NowAssembling libMesh::EigenTimeSolver::now_assembling
private

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

Definition at line 213 of file eigen_time_solver.h.

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

AutoPtr<NumericVector<Number> > libMesh::TimeSolver::old_local_nonlinear_solution
protectedinherited

Serial vector of _system.get_vector("_old_nonlinear_solution")

Definition at line 253 of file time_solver.h.

bool libMesh::TimeSolver::quiet
inherited

Print extra debugging information if quiet == false.

Definition at line 177 of file time_solver.h.

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

unsigned int libMesh::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 205 of file time_solver.h.

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

AutoPtr<SolutionHistory> libMesh::TimeSolver::solution_history
protectedinherited

An AutoPtr 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 260 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().

Real libMesh::EigenTimeSolver::tol

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

Definition at line 156 of file eigen_time_solver.h.

Referenced by solve().


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

Site Created By: libMesh Developers
Last modified: February 07 2014 16:57:59 UTC

Hosted By:
SourceForge.net Logo