20 #ifndef __mast__nonlinear_system_h__ 21 #define __mast__nonlinear_system_h__ 30 #include "libmesh/nonlinear_implicit_system.h" 31 #include "libmesh/enum_eigen_solver_type.h" 32 #include "libmesh/eigen_system.h" 40 class SlepcEigenSolver;
42 class PhysicsDisciplineBase;
43 class AssemblyElemOperations;
44 class OutputAssemblyElemOperations;
46 class EigenproblemAssembly;
57 public libMesh::NonlinearImplicitSystem {
65 const std::string& name,
66 const unsigned int number);
111 virtual void clear () libmesh_override;
118 virtual void reinit () libmesh_override;
125 virtual std::pair<unsigned int, Real>
143 bool if_localize_sol,
147 bool if_assemble_jacobian =
true);
155 virtual void adjoint_solve(
const libMesh::NumericVector<Real>& X,
156 bool if_localize_sol,
160 bool if_assemble_jacobian =
true);
181 std::vector<Real>& sens,
182 const std::vector<unsigned int>* indices=
nullptr);
212 libMesh::NumericVector<Real>& vec_re,
213 libMesh::NumericVector<Real>* vec_im =
nullptr);
254 libMesh::EigenProblemType
303 const std::string & directory_name,
304 const std::string & data_name,
305 const bool write_binary_vectors);
312 const std::string & directory_name,
313 const std::string & data_name,
314 const bool read_binary_vectors);
318 libMesh::FunctionBase<Real>& f)
const;
327 virtual void init_data () libmesh_override;
405 #endif // __mast__nonlinear_system_h__
unsigned int _n_iterations
The number of iterations of the eigen solver algorithm.
void initialize_condensed_dofs(MAST::PhysicsDisciplineBase &physics)
Loop over the dofs on each processor to initialize the list of non-condensed dofs.
void read_in_vector(libMesh::NumericVector< Real > &vec, const std::string &directory_name, const std::string &data_name, const bool read_binary_vectors)
reads the specified vector with the specified name in a directory.
This class implements a system for solution of nonlinear systems.
bool _is_generalized_eigenproblem
A boolean flag to indicate whether we are dealing with a generalized eigenvalue problem.
virtual void clear() libmesh_override
Clear all the data structures associated with the system.
void set_n_requested_eigenvalues(unsigned int n)
sets the number of eigenvalues requested
void write_out_vector(libMesh::NumericVector< Real > &vec, const std::string &directory_name, const std::string &data_name, const bool write_binary_vectors)
writes the specified vector with the specified name in a directory.
unsigned int _n_requested_eigenpairs
The number of requested eigenpairs.
This provides the base class for definitin of element level contribution of output quantity in an ana...
std::unique_ptr< MAST::SlepcEigenSolver > eigen_solver
The EigenSolver, definig which interface, i.e solver package to use.
bool _initialize_B_matrix
initialize the B matrix in addition to A, which might be needed for solution of complex system of equ...
void set_eigenproblem_type(libMesh::EigenProblemType ept)
Sets the type of the current eigen problem.
virtual void reinit() libmesh_override
Reinitializes the member data fields associated with the system, so that, e.g., assemble() may be use...
MAST::NonlinearSystem::Operation _operation
current operation of the system
virtual void eigenproblem_sensitivity_solve(MAST::AssemblyElemOperations &elem_ops, MAST::EigenproblemAssembly &assembly, const MAST::FunctionBase &f, std::vector< Real > &sens, const std::vector< unsigned int > *indices=nullptr)
Solves the sensitivity system, for the provided parameters.
unsigned int n_global_non_condensed_dofs() const
Assembles the system of equations for an eigenproblem of type .
unsigned int get_n_requested_eigenvalues() const
libMesh::SparseMatrix< Real > * matrix_A
The system matrix for standard eigenvalue problems.
unsigned int get_n_converged_eigenvalues() const
void project_vector_without_dirichlet(libMesh::NumericVector< Real > &new_vector, libMesh::FunctionBase< Real > &f) const
bool _condensed_dofs_initialized
A private flag to indicate whether the condensed dofs have been initialized.
void set_n_iterations(unsigned int its)
Set the _n_iterations member, useful for subclasses of EigenSystem.
virtual void eigenproblem_solve(MAST::AssemblyElemOperations &elem_ops, MAST::EigenproblemAssembly &assembly)
Assembles & solves the eigen system.
libMesh::EigenProblemType get_eigenproblem_type() const
virtual std::pair< unsigned int, Real > get_linear_solve_parameters()
calls NonlinearImplicitSystem::set_solver_parameters() before accessing the values.
virtual void init_data() libmesh_override
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used...
libMesh::SparseMatrix< Real > * matrix_B
A second system matrix for generalized eigenvalue problems.
libMesh::EigenProblemType _eigen_problem_type
The type of the eigenvalue problem.
void set_init_B_matrix()
flag to also initialize the B matrix.
void set_operation(MAST::NonlinearSystem::Operation op)
sets the current operation of the system
virtual void get_eigenvalue(unsigned int i, Real &re, Real &im)
gets the real and imaginary parts of the ith eigenvalue for the eigenproblem , and the associated eig...
virtual ~NonlinearSystem()
std::unique_ptr< libMesh::LinearSolver< Real > > linear_solver
The LinearSolver for solution of the linear equations.
unsigned int _n_converged_eigenpairs
The number of converged eigenpairs.
std::vector< libMesh::dof_id_type > _local_non_condensed_dofs_vector
Vector storing the local dof indices that will not be condensed.
virtual void get_eigenpair(unsigned int i, Real &re, Real &im, libMesh::NumericVector< Real > &vec_re, libMesh::NumericVector< Real > *vec_im=nullptr)
gets the real and imaginary parts of the ith eigenvalue for the eigenproblem , and the associated eig...
MAST::NonlinearSystem::Operation operation()
bool _exchange_A_and_B
flag to exchange the A and B matrices in the eigenproblem solution
unsigned int get_n_iterations() const
void set_n_converged(unsigned int nconv)
Set the _n_converged_eigenpairs member, useful for subclasses of EigenSystem.
virtual void sensitivity_solve(const libMesh::NumericVector< Real > &X, bool if_localize_sol, MAST::AssemblyElemOperations &elem_ops, MAST::AssemblyBase &assembly, const MAST::FunctionBase &p, bool if_assemble_jacobian=true)
Solves the sensitivity problem for the provided parameter.
void set_exchange_A_and_B(bool flag)
sets the flag to exchange the A and B matrices for a generalized eigenvalue problem.
NonlinearSystem(libMesh::EquationSystems &es, const std::string &name, const unsigned int number)
Default constructor.
virtual void adjoint_solve(const libMesh::NumericVector< Real > &X, bool if_localize_sol, MAST::AssemblyElemOperations &elem_ops, MAST::OutputAssemblyElemOperations &output, MAST::AssemblyBase &assembly, bool if_assemble_jacobian=true)
solves the adjoint problem for the provided output function.
virtual void solve(MAST::AssemblyElemOperations &elem_ops, MAST::AssemblyBase &assembly)
solves the nonlinear problem with the specified assembly operation object