Go to the documentation of this file.
20 #ifndef LIBMESH_RB_SCM_CONSTRUCTION_H
21 #define LIBMESH_RB_SCM_CONSTRUCTION_H
24 #include "libmesh/libmesh_config.h"
28 #if defined(LIBMESH_HAVE_SLEPC) && (LIBMESH_HAVE_GLPK)
31 #include "libmesh/rb_construction_base.h"
34 #include "libmesh/condensed_eigen_system.h"
42 class RBSCMEvaluation;
62 const std::string & name_in,
63 const unsigned int number_in);
84 virtual void clear ()
override;
189 virtual void enrich_C_J(
unsigned int new_C_J_index);
243 #endif // LIBMESH_HAVE_SLEPC && LIBMESH_HAVE_GLPK
245 #endif // LIBMESH_RB_SCM_CONSTRUCTION_H
virtual void perform_SCM_greedy()
Perform the SCM greedy algorithm to develop a lower bound over the training set.
virtual std::pair< unsigned int, Real > compute_SCM_bounds_on_training_set()
Compute upper and lower bounds for each SCM training point.
virtual void resize_SCM_vectors()
Clear and resize the SCM data vectors.
virtual void load_matrix_B()
Copy over the matrix to store in matrix_B, usually this is the mass or inner-product matrix,...
Number B_inner_product(const NumericVector< Number > &v, const NumericVector< Number > &w) const
Compute the inner product between two vectors using the system's matrix_B.
void set_rb_scm_evaluation(RBSCMEvaluation &rb_scm_eval_in)
Set the RBSCMEvaluation object.
Real SCM_training_tolerance
Tolerance which controls when to terminate the SCM Greedy.
RBSCMConstruction sys_type
The type of system.
The libMesh namespace provides an interface to certain functionality in the library.
virtual void clear() override
Clear all the data structures associated with the system.
virtual void set_eigensolver_properties(int)
This function is called before truth eigensolves in compute_SCM_bounding_box and evaluate_stability_c...
This class is part of the rbOOmit framework.
void set_RB_system_name(const std::string &new_name)
Set the name of the associated RB system — we need this to load the (symmetrized) affine operators.
RBSCMConstruction(EquationSystems &es, const std::string &name_in, const unsigned int number_in)
Constructor.
RBThetaExpansion & get_rb_theta_expansion()
Get a reference to the RBThetaExpansion object.
void set_SCM_training_tolerance(Real SCM_training_tolerance_in)
Number Aq_inner_product(unsigned int q, const NumericVector< Number > &v, const NumericVector< Number > &w)
Compute the inner product between two vectors using matrix Aq.
virtual void process_parameters_file(const std::string ¶meters_filename)
Read in the parameters from file specified by parameters_filename and set the this system's member va...
RBSCMEvaluation & get_rb_scm_evaluation()
Get a reference to the RBSCMEvaluation object.
virtual void compute_SCM_bounding_box()
Compute the SCM bounding box.
Real get_SCM_training_tolerance() const
Get/set SCM_training_tolerance: tolerance for SCM greedy.
This class is part of the rbOOmit framework.
virtual void enrich_C_J(unsigned int new_C_J_index)
Enrich C_J by adding the element of SCM_training_samples that has the largest gap between alpha_LB an...
virtual ~RBSCMConstruction()
Destructor.
This is the EquationSystems class.
virtual void evaluate_stability_constant()
Compute the stability constant for current_parameters by solving a generalized eigenvalue problem ove...
virtual Real SCM_greedy_error_indicator(Real LB, Real UB)
Helper function which provides an error indicator to be used in the SCM greedy.
RBSCMEvaluation * rb_scm_eval
The current RBSCMEvaluation object we are using to perform the Evaluation stage of the SCM.
std::string RB_system_name
The name of the associated RB system.
virtual void print_info()
Print out info that describes the current setup of this RBSCMConstruction.
This class is part of the rbOOmit framework.
virtual void attach_deflation_space()
Attach the deflation space defined by the specified vector, can be useful in solving constrained eige...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
This class stores the set of RBTheta functor objects that define the "parameter-dependent expansion" ...
virtual void add_scaled_symm_Aq(unsigned int q_a, Number scalar)
Add the scaled symmetrized affine matrix from the associated RBSystem to matrix_A.
RBConstructionBase< CondensedEigenSystem > Parent
The type of the parent.