libMesh::EigenTimeSolver Class Reference
#include <eigen_time_solver.h>

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_type & | system () const |
| sys_type & | system () |
| 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< SolutionHistory > | solution_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.
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.
The type of system
Reimplemented from libMesh::TimeSolver.
Definition at line 73 of file eigen_time_solver.h.
Member Enumeration Documentation
enum libMesh::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 193 of file eigen_time_solver.h.
00193 { 00197 Matrix_A, 00198 00202 Matrix_B, 00203 00207 Invalid_Matrix 00208 };
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.
00033 : Parent(s), 00034 eigen_solver (EigenSolver<Number>::build()), 00035 tol(pow(TOLERANCE, 5./3.)), 00036 maxits(1000), 00037 n_eigenpairs_to_compute(5), 00038 n_basis_vectors_to_use(3*n_eigenpairs_to_compute), 00039 n_converged_eigenpairs(0), 00040 n_iterations_reqd(0) 00041 { 00042 libmesh_experimental(); 00043 eigen_solver->set_eigenproblem_type(GHEP);//or GNHEP 00044 eigen_solver->set_position_of_spectrum(LARGEST_MAGNITUDE); 00045 }
| libMesh::EigenTimeSolver::~EigenTimeSolver | ( | ) | [virtual] |
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 100 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 152 of file time_solver.h.
| virtual AutoPtr<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 167 of file time_solver.h.
References libMesh::TimeSolver::_diff_solver.
00167 { return _diff_solver; }
| void libMesh::ReferenceCounter::disable_print_counter_info | ( | ) | [static, inherited] |
Definition at line 106 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter.
00107 { 00108 _enable_print_counter = false; 00109 return; 00110 }
| 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 139 of file eigen_time_solver.h.
| 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::DifferentiablePhysics::mass_residual(), Matrix_A, Matrix_B, and now_assembling.
00129 { 00130 // The EigenTimeSolver always computes jacobians! 00131 libmesh_assert (request_jacobian); 00132 00133 // Assemble the operator for the spatial part. 00134 if (now_assembling == Matrix_A) 00135 { 00136 bool jacobian_computed = 00137 _system.element_time_derivative(request_jacobian, context); 00138 00139 // The user shouldn't compute a jacobian unless requested 00140 libmesh_assert(request_jacobian || !jacobian_computed); 00141 00142 bool jacobian_computed2 = 00143 _system.element_constraint(jacobian_computed, context); 00144 00145 // The user shouldn't compute a jacobian unless requested 00146 libmesh_assert (jacobian_computed || !jacobian_computed2); 00147 00148 return jacobian_computed && jacobian_computed2; 00149 00150 } 00151 00152 // Assemble the mass matrix operator 00153 else if (now_assembling == Matrix_B) 00154 { 00155 bool mass_jacobian_computed = 00156 _system.mass_residual(request_jacobian, context); 00157 00158 // Scale Jacobian by -1? 00159 //context.elem_jacobian *= -1.0; 00160 00161 return mass_jacobian_computed; 00162 } 00163 00164 else 00165 { 00166 libmesh_error(); 00167 return false; 00168 } 00169 }
| void libMesh::ReferenceCounter::enable_print_counter_info | ( | ) | [static, inherited] |
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.
00101 { 00102 _enable_print_counter = true; 00103 return; 00104 }
| 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.
| 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().
00048 { 00049 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG) 00050 00051 std::ostringstream oss; 00052 00053 oss << '\n' 00054 << " ---------------------------------------------------------------------------- \n" 00055 << "| Reference count information |\n" 00056 << " ---------------------------------------------------------------------------- \n"; 00057 00058 for (Counts::iterator it = _counts.begin(); 00059 it != _counts.end(); ++it) 00060 { 00061 const std::string name(it->first); 00062 const unsigned int creations = it->second.first; 00063 const unsigned int destructions = it->second.second; 00064 00065 oss << "| " << name << " reference count information:\n" 00066 << "| Creations: " << creations << '\n' 00067 << "| Destructions: " << destructions << '\n'; 00068 } 00069 00070 oss << " ---------------------------------------------------------------------------- \n"; 00071 00072 return oss.str(); 00073 00074 #else 00075 00076 return ""; 00077 00078 #endif 00079 }
| 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, and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().
00164 { 00165 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 00166 std::pair<unsigned int, unsigned int>& p = _counts[name]; 00167 00168 p.first++; 00169 }
| 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, and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().
00177 { 00178 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 00179 std::pair<unsigned int, unsigned int>& p = _counts[name]; 00180 00181 p.second++; 00182 }
| 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().
00057 { 00058 // Add matrix "B" to _system if not already there. 00059 // The user may have already added a matrix "B" before 00060 // calling the System initialization. This would be 00061 // necessary if e.g. the System originally started life 00062 // with a different type of TimeSolver and only later 00063 // had its TimeSolver changed to an EigenTimeSolver. 00064 if (!_system.have_matrix("B")) 00065 _system.add_matrix("B"); 00066 }
| 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.
| bool libMesh::TimeSolver::is_adjoint | ( | ) | const [inline, inherited] |
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().
00218 { 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 144 of file eigen_time_solver.h.
| virtual AutoPtr<LinearSolver<Number> >& libMesh::TimeSolver::linear_solver | ( | ) | [inline, virtual, inherited] |
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.
00172 { 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.
00080 { return _n_objects; }
| 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().
00089 { 00090 if( _enable_print_counter ) out_stream << ReferenceCounter::get_info(); 00091 }
| 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.
| 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 104 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 224 of file time_solver.h.
References libMesh::TimeSolver::_is_adjoint.
Referenced by libMesh::DifferentiableSystem::adjoint_solve(), libMesh::FEMSystem::postprocess(), and libMesh::DifferentiableSystem::solve().
00225 { _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.
00092 { 00093 solution_history = _solution_history.clone(); 00094 }
| 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, Matrix_A, Matrix_B, now_assembling, libMesh::DifferentiablePhysics::side_constraint(), libMesh::DifferentiablePhysics::side_mass_residual(), and libMesh::DifferentiablePhysics::side_time_derivative().
00175 { 00176 // The EigenTimeSolver always requests jacobians? 00177 //libmesh_assert (request_jacobian); 00178 00179 // Assemble the operator for the spatial part. 00180 if (now_assembling == Matrix_A) 00181 { 00182 bool jacobian_computed = 00183 _system.side_time_derivative(request_jacobian, context); 00184 00185 // The user shouldn't compute a jacobian unless requested 00186 libmesh_assert (request_jacobian || !jacobian_computed); 00187 00188 bool jacobian_computed2 = 00189 _system.side_constraint(jacobian_computed, context); 00190 00191 // The user shouldn't compute a jacobian unless requested 00192 libmesh_assert (jacobian_computed || !jacobian_computed2); 00193 00194 return jacobian_computed && jacobian_computed2; 00195 00196 } 00197 00198 // There is now a "side" equivalent for the mass matrix 00199 else if (now_assembling == Matrix_B) 00200 { 00201 bool mass_jacobian_computed = 00202 _system.side_mass_residual(request_jacobian, context); 00203 00204 return mass_jacobian_computed; 00205 } 00206 00207 else 00208 { 00209 libmesh_error(); 00210 return false; 00211 } 00212 }
| 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.
00069 { 00070 // The standard implementation is basically to call: 00071 // _diff_solver->solve(); 00072 // which internally assembles (when necessary) matrices and vectors 00073 // and calls linear solver software while also doing Newton steps (see newton_solver.C) 00074 // 00075 // The element_residual and side_residual functions below control 00076 // what happens in the interior of the element assembly loops. 00077 // We have a system reference, so it's possible to call _system.assembly() 00078 // ourselves if we want to... 00079 // 00080 // Interestingly, for the EigenSolver we don't need residuals...just Jacobians. 00081 // The Jacobian should therefore always be requested, and always return 00082 // jacobian_computed as being true. 00083 00084 // The basic plan of attack is: 00085 // .) Construct the Jacobian using _system.assembly(true,true) as we 00086 // would for a steady system. Use a flag in this class to 00087 // control behavior in element_residual and side_residual 00088 // .) Swap _system.matrix to matrix "B" (be sure to add this extra matrix during init) 00089 // .) Call _system.assembly(true,true) again, using the flag in element_residual 00090 // and side_residual to only get the mass matrix terms. 00091 // .) Send A and B to Steffen's EigenSolver interface. 00092 00093 // Assemble the spatial part (matrix A) of the operator 00094 if (!this->quiet) 00095 libMesh::out << "Assembling matrix A." << std::endl; 00096 _system.matrix = &( _system.get_matrix ("System Matrix") ); 00097 this->now_assembling = Matrix_A; 00098 _system.assembly(true, true); 00099 //_system.matrix->print_matlab("matrix_A.m"); 00100 00101 // Point the system's matrix at B, call assembly again. 00102 if (!this->quiet) 00103 libMesh::out << "Assembling matrix B." << std::endl; 00104 _system.matrix = &( _system.get_matrix ("B") ); 00105 this->now_assembling = Matrix_B; 00106 _system.assembly(true, true); 00107 //_system.matrix->print_matlab("matrix_B.m"); 00108 00109 // Send matrices A, B to Steffen's SlepcEigenSolver interface 00110 //libmesh_here(); 00111 if (!this->quiet) 00112 libMesh::out << "Calling the EigenSolver." << std::endl; 00113 std::pair<unsigned int, unsigned int> solve_data = 00114 eigen_solver->solve_generalized (_system.get_matrix ("System Matrix"), 00115 _system.get_matrix ("B"), 00116 n_eigenpairs_to_compute, 00117 n_basis_vectors_to_use, 00118 tol, 00119 maxits); 00120 00121 this->n_converged_eigenpairs = solve_data.first; 00122 this->n_iterations_reqd = solve_data.second; 00123 }
| sys_type& libMesh::TimeSolver::system | ( | ) | [inline, inherited] |
- Returns:
- a writeable reference to the system we are solving.
Definition at line 162 of file time_solver.h.
References libMesh::TimeSolver::_system.
00162 { return _system; }
| const sys_type& libMesh::TimeSolver::system | ( | ) | const [inline, inherited] |
- 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().
00157 { return _system; }
Member Data Documentation
ReferenceCounter::Counts libMesh::ReferenceCounter::_counts [static, protected, inherited] |
Actually holds the data.
Definition at line 118 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::get_info(), libMesh::ReferenceCounter::increment_constructor_count(), and libMesh::ReferenceCounter::increment_destructor_count().
AutoPtr<DiffSolver> libMesh::TimeSolver::_diff_solver [protected, inherited] |
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 [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().
AutoPtr<LinearSolver<Number> > libMesh::TimeSolver::_linear_solver [protected, inherited] |
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 [static, protected, inherited] |
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().
sys_type& libMesh::TimeSolver::_system [protected, inherited] |
A reference to the system we are solving.
Definition at line 242 of file time_solver.h.
Referenced by libMesh::UnsteadySolver::adjoint_advance_timestep(), libMesh::UnsteadySolver::advance_timestep(), libMesh::AdaptiveTimeSolver::advance_timestep(), libMesh::UnsteadySolver::du(), libMesh::SteadySolver::element_residual(), libMesh::EulerSolver::element_residual(), libMesh::Euler2Solver::element_residual(), element_residual(), libMesh::UnsteadySolver::init(), libMesh::TimeSolver::init(), init(), libMesh::UnsteadySolver::init_data(), libMesh::UnsteadySolver::old_nonlinear_solution(), libMesh::UnsteadySolver::reinit(), libMesh::UnsteadySolver::retrieve_timestep(), libMesh::SteadySolver::side_residual(), libMesh::EulerSolver::side_residual(), libMesh::Euler2Solver::side_residual(), side_residual(), libMesh::UnsteadySolver::solve(), libMesh::TwostepTimeSolver::solve(), solve(), and libMesh::TimeSolver::system().
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 [protected, inherited] |
A bool that will be true the first time solve() is called, and false thereafter
Reimplemented in libMesh::UnsteadySolver.
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().
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 [protected, inherited] |
Serial vector of _system.get_vector("_old_nonlinear_solution")
Reimplemented in libMesh::UnsteadySolver.
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::UnsteadySolver::solve(), libMesh::TwostepTimeSolver::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::UnsteadySolver::solve(), and libMesh::TwostepTimeSolver::solve().
AutoPtr<SolutionHistory> libMesh::TimeSolver::solution_history [protected, inherited] |
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().
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 05 2013 19:55:17 UTC
Hosted By: