28 virtual std::string
solverTypeString(
unsigned int solver_sys_num = 0)
override;
30 #ifdef LIBMESH_HAVE_SLEPC 31 virtual void solve(
const unsigned int nl_sys_num)
override;
33 virtual void init()
override;
129 const std::set<TagID> & tags);
142 const unsigned int nl_sys_num)
override;
296 #ifdef LIBMESH_HAVE_SLEPC 302 mooseError(
"eigen problems do not currently support multiple nonlinear eigen systems");
Nonlinear eigenvalue system to be solved.
bool isGeneralizedEigenvalueProblem() const
bool constantMatrices() const
Whether or not matrices are constant.
bool wereMatricesFormed() const
Whether or not constant matrices were already formed.
virtual void computeResidualTag(const NumericVector< Number > &soln, NumericVector< Number > &residual, TagID tag) override
Form a vector for all kernels and BCs with a given tag.
static InputParameters validParams()
void onLinearSolver(bool ols)
Set a flag to indicate whether or not we are in a linear solver iteration.
bool _matrices_formed
Whether or not matrices had been formed.
virtual void initPetscOutputAndSomeSolverSettings() override
Hook up monitors for SNES and KSP.
PostprocessorName _normalization
Postprocessor used to compute a factor from eigenvector.
unsigned int activeEigenvalueIndex() const
Which eigenvalue is active.
virtual void computeJacobianBlocks(std::vector< JacobianBlock *> &blocks, const unsigned int nl_sys_num) override
Computes several Jacobian blocks simultaneously, summing their contributions into smaller preconditio...
bool outputInverseEigenvalue() const
Whether or not to output eigenvalue inverse.
void setBxNorm(const PostprocessorName &bx_norm)
Set the Bx norm postprocessor programatically.
bool _on_linear_solver
Whether or not we are in linear solver.
void initEigenvector(const Real initial_value)
For nonlinear eigen solver, a good initial value can help convergence.
bool _has_normalization
Whether or not we normalize eigenvector.
EigenProblem(const InputParameters ¶meters)
virtual void execute(const ExecFlagType &exec_type) override
Convenience function for performing execution of MOOSE systems.
unsigned int _n_eigen_pairs_required
bool negativeSignEigenKernel() const
A flag indicates if a negative sign is used in eigen kernels.
void adjustEigenVector(const Real value, bool scaling)
Adjust eigen vector by either scaling the existing values or setting new values The operations are ap...
void doFreeNonlinearPowerIterations(unsigned int free_power_iterations)
Do some free/extra power iterations.
bool _output_inverse_eigenvalue
Whether or not output eigenvalue as its inverse. By default, we output regular eigenvalue.
bool _negative_sign_eigen_kernel
Whether or not use negative sign for Bx.
bool _constant_matrices
Whether or not require constant matrices.
void setNEigenPairsRequired(unsigned int n_eigen_pairs)
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
void preScaleEigenVector(const std::pair< Real, Real > &eig)
Eigenvector need to be scaled back if it was scaled in an earlier stage Scaling eigen vector does not...
bool isNonlinearEigenvalueSolver(unsigned int eigen_sys_num) const
auto max(const L &left, const R &right)
std::shared_ptr< NonlinearEigenSystem > _nl_eigen
bool _generalized_eigenvalue_problem
bool _do_free_power_iteration
Whether or not we are doing free power iteration.
Real formNorm()
Form the Bx norm.
bool & _first_solve
A flag to indicate if it is the first time calling the solve.
unsigned int getNEigenPairsRequired() const
std::vector< std::shared_ptr< NonlinearSystemBase > > _nl
The nonlinear systems.
void constantMatrices(bool cm)
Set a flag to indicate whether or not we use constant matrices.
virtual Real computeResidualL2Norm() override
Compute the residual of Ax - Bx.
void setNormalization(const PostprocessorName &pp, const Real value=std::numeric_limits< Real >::max())
Set postprocessor and normalization factor 'Postprocessor' is often used to compute an integral of ph...
virtual void init() override
virtual void computeJacobian(const NumericVector< libMesh::Number > &soln, libMesh::SparseMatrix< libMesh::Number > &jacobian, const unsigned int nl_sys_num)
Form a Jacobian matrix with the default tag (system).
void computeResidualAB(const NumericVector< Number > &soln, NumericVector< Number > &residualA, NumericVector< Number > &residualB, TagID tagA, TagID tagB)
Form two vetors, where each is associated with one tag, through one element-loop. ...
void computeMatricesTags(const NumericVector< Number > &soln, const std::vector< SparseMatrix< Number > *> &jacobians, const std::set< TagID > &tags)
Form several matrices simultaneously.
Real _initial_eigenvalue
A value used for initial normalization.
unsigned int _active_eigen_index
Which eigenvalue is used to compute residual.
NonlinearEigenSystem & getCurrentNonlinearEigenSystem()
void outputInverseEigenvalue(bool inverse)
Set a flag to indicate whether or not to output eigenvalue inverse.
NonlinearEigenSystem & getNonlinearEigenSystem(const unsigned int nl_sys_num)
void scaleEigenvector(const Real scaling_factor)
Scale eigenvector.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
bool onLinearSolver() const
Whether or not we are in a linear solver iteration.
Class for containing MooseEnum item information.
void postScaleEigenVector()
Normalize eigen vector.
virtual std::string solverTypeString(unsigned int solver_sys_num=0) override
Return solver type as a human readable string.
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.
const InputParameters & parameters() const
Get the parameters of the object.
Real _normal_factor
Postprocessor target value.
void doFreePowerIteration(bool do_power)
Set a flag to indicate whether or not we are doing free power iterations.
void setInitialEigenvalue(const Real initial_eigenvalue)
Set an initial eigenvalue for initial normalization.
std::optional< PostprocessorName > _bx_norm_name
The name of the Postprocessor providing the Bx norm.
Problem for solving eigenvalue problems.
bool bxNormProvided() const
Whether a Bx norm postprocessor has been provided.
void wereMatricesFormed(bool mf)
Set a flag to indicate whether or not constant matrices were already formed.
virtual void checkProblemIntegrity() override
Method called to perform a series of sanity checks before a simulation is run.
EigenProblemType
Type of the eigen problem.
bool doFreePowerIteration() const
Whether or not we are doing free power iteration.
void setEigenproblemType(Moose::EigenProblemType eigen_problem_type)
Set eigen problem type.
virtual void solve(const unsigned int nl_sys_num) override
void computeJacobianAB(const NumericVector< Number > &soln, SparseMatrix< Number > &jacobianA, SparseMatrix< Number > &jacobianB, TagID tagA, TagID tagB)
Form two Jacobian matrices, where each is associated with one tag, through one element-loop.
virtual void computeJacobianTag(const NumericVector< Number > &soln, SparseMatrix< Number > &jacobian, TagID tag) override
Form a Jacobian matrix for all kernels and BCs with a given tag.
virtual bool solverSystemConverged(const unsigned int solver_sys_num) override