20 #include "libmesh/libmesh_config.h" 21 #ifdef LIBMESH_HAVE_SLEPC 23 #include "libmesh/diff_system.h" 24 #include "libmesh/eigen_time_solver.h" 25 #include "libmesh/eigen_solver.h" 26 #include "libmesh/sparse_matrix.h" 27 #include "libmesh/enum_eigen_solver_type.h" 38 n_eigenpairs_to_compute(5),
39 n_basis_vectors_to_use(3*n_eigenpairs_to_compute),
40 n_converged_eigenpairs(0),
43 libmesh_experimental();
111 libMesh::out <<
"Calling the EigenSolver." << std::endl;
135 bool jacobian_computed =
141 bool jacobian_computed2 =
147 return jacobian_computed && jacobian_computed2;
154 bool mass_jacobian_computed =
161 return mass_jacobian_computed;
165 libmesh_error_msg(
"Unrecognized value now_assembling = " <<
now_assembling);
182 bool jacobian_computed =
188 bool jacobian_computed2 =
194 return jacobian_computed && jacobian_computed2;
201 bool mass_jacobian_computed =
208 return mass_jacobian_computed;
212 libmesh_error_msg(
"Unrecognized value now_assembling = " <<
now_assembling);
226 bool jacobian_computed =
232 bool jacobian_computed2 =
238 return jacobian_computed && jacobian_computed2;
245 bool mass_jacobian_computed =
252 return mass_jacobian_computed;
256 libmesh_error_msg(
"Unrecognized value now_assembling = " <<
now_assembling);
261 #endif // LIBMESH_HAVE_SLEPC bool quiet
Print extra debugging information if quiet == false.
virtual bool side_constraint(bool request_jacobian, DiffContext &)
Adds the constraint contribution on side of elem to elem_residual.
const DenseMatrix< Number > & get_elem_jacobian() const
Const accessor for element Jacobian.
virtual void init() override
The initialization function.
This class provides all data required for a physics package (e.g.
virtual bool element_time_derivative(bool request_jacobian, DiffContext &)
Adds the time derivative contribution on elem to elem_residual.
virtual void assembly(bool get_residual, bool get_jacobian, bool apply_heterogeneous_constraints=false, bool apply_no_constraints=false) override=0
Assembles a residual in rhs and/or a jacobian in matrix, as requested.
This is a generic class that defines a solver to handle time integration of DifferentiableSystems.
static constexpr Real TOLERANCE
virtual void solve() override
Implements the assembly of both matrices A and B, and calls the EigenSolver to compute the eigenvalue...
virtual bool nonlocal_residual(bool get_jacobian, DiffContext &) override
Forms the jacobian of the nonlocal terms.
virtual bool nonlocal_time_derivative(bool request_jacobian, DiffContext &)
Adds any nonlocal time derivative contributions (e.g.
The libMesh namespace provides an interface to certain functionality in the library.
virtual bool nonlocal_constraint(bool request_jacobian, DiffContext &)
Adds any nonlocal constraint contributions (e.g.
virtual bool mass_residual(bool request_jacobian, DiffContext &)
Subtracts a mass vector contribution on elem from elem_residual.
bool have_matrix(std::string_view mat_name) const
This class provides a specific system class.
sys_type & _system
A reference to the system we are solving.
Real elem_solution_rate_derivative
The derivative of elem_solution_rate with respect to the current nonlinear solution, for use by systems with non default mass_residual terms.
Real elem_solution_derivative
The derivative of elem_solution with respect to the current nonlinear solution.
unsigned int maxits
The maximum number of iterations allowed to solve the problem.
unsigned int n_iterations_reqd
After a solve, holds the number of iterations required to converge the requested number of eigenpairs...
const DifferentiablePhysics * get_physics() const
virtual ~EigenTimeSolver()
Destructor.
NowAssembling now_assembling
Flag which controls the internals of element_residual() and side_residual().
double tol
The linear solver tolerance to be used when solving the eigenvalue problem.
SparseMatrix< Number > * matrix
The system matrix.
unsigned int n_basis_vectors_to_use
The number of basis vectors to use in the computation.
EigenTimeSolver(sys_type &s)
Constructor.
virtual bool side_residual(bool get_jacobian, DiffContext &) override
Forms the jacobian of the boundary terms.
virtual void reinit() override
The reinitialization function.
virtual bool side_mass_residual(bool request_jacobian, DiffContext &)
Subtracts a mass vector contribution on side of elem from elem_residual.
virtual bool side_time_derivative(bool request_jacobian, DiffContext &)
Adds the time derivative contribution on side of elem to elem_residual.
virtual bool element_residual(bool get_jacobian, DiffContext &) override
Forms either the spatial (Jacobian) or mass matrix part of the operator, depending on which is reques...
unsigned int n_converged_eigenpairs
After a solve, holds the number of eigenpairs successfully converged.
virtual bool element_constraint(bool request_jacobian, DiffContext &)
Adds the constraint contribution on elem to elem_residual.
SparseMatrix< Number > & add_matrix(std::string_view mat_name, ParallelType type=PARALLEL, MatrixBuildType mat_build_type=MatrixBuildType::AUTOMATIC)
Adds the additional matrix mat_name to this system.
const SparseMatrix< Number > & get_matrix(std::string_view mat_name) const
unsigned int n_eigenpairs_to_compute
The number of eigenvectors/values to be computed.
The matrix associated with the spatial part of the operator.
const SparseMatrix< Number > & get_system_matrix() const
The matrix associated with the time derivative (mass matrix).
virtual bool nonlocal_mass_residual(bool request_jacobian, DiffContext &c)
Subtracts any nonlocal mass vector contributions (e.g.
std::unique_ptr< EigenSolver< Number > > eigen_solver
The EigenSolver object.
This class provides an interface to solvers for eigenvalue problems.