Go to the documentation of this file.
20 #ifndef LIBMESH_RB_CONSTRUCTION_H
21 #define LIBMESH_RB_CONSTRUCTION_H
24 #include "libmesh/rb_construction_base.h"
27 #include "libmesh/linear_implicit_system.h"
28 #include "libmesh/dense_vector.h"
29 #include "libmesh/dense_matrix.h"
30 #include "libmesh/dg_fem_context.h"
31 #include "libmesh/dirichlet_boundaries.h"
38 class RBThetaExpansion;
39 class RBAssemblyExpansion;
62 const std::string &
name,
63 const unsigned int number);
128 virtual void clear ()
override;
274 bool skip_vector_assembly=
false);
323 bool skip_vector_assembly);
354 const bool write_binary_residual_representors);
363 const bool write_binary_residual_representors);
386 bool deterministic_training_in,
387 unsigned int training_parameters_random_seed_in,
389 unsigned int Nmax_in,
390 Real rel_training_tolerance_in,
391 Real abs_training_tolerance_in,
392 bool normalize_rb_error_bound_in_greedy_in,
395 std::map<std::string, std::vector<Real>> discrete_parameter_values_in,
396 std::map<std::string,bool> log_scaling);
441 #ifdef LIBMESH_ENABLE_DIRICHLET
609 Real initial_greedy_error,
629 bool symmetrize=
false,
630 bool apply_dof_constraints=
true);
826 std::vector<std::unique_ptr<SparseMatrix<Number>>>
Aq_vector;
832 std::vector<std::unique_ptr<NumericVector<Number>>>
Fq_vector;
867 #endif // LIBMESH_RB_CONSTRUCTION_H
std::unique_ptr< SparseMatrix< Number > > inner_product_matrix
The inner product matrix.
virtual void enrich_RB_space()
Add a new basis function to the RB space.
virtual void post_process_truth_solution()
Similarly, provide an opportunity to post-process the truth solution after the solve is complete.
virtual Real get_RB_error_bound()
NumericVector< Number > * get_non_dirichlet_Fq(unsigned int q)
Get a pointer to non-Dirichlet Fq.
SparseMatrix< Number > * get_Aq(unsigned int q)
Get a pointer to Aq.
unsigned int Nmax
Maximum number of reduced basis functions we are willing to use.
virtual void process_parameters_file(const std::string ¶meters_filename)
Read in from the file specified by parameters_filename and set the this system's member variables acc...
virtual void get_all_vectors(std::map< std::string, NumericVector< Number > * > &all_vectors)
Get a map that stores pointers to all of the vectors.
NumericVector< Number > * get_non_dirichlet_Fq_if_avail(unsigned int q)
Get a pointer to non_dirichlet_Fq if it's available, otherwise get Fq.
virtual void initialize_rb_construction(bool skip_matrix_assembly=false, bool skip_vector_assembly=false)
Allocate all the data structures necessary for the construction stage of the RB method.
void set_rb_evaluation(RBEvaluation &rb_eval_in)
Set the RBEvaluation object.
This class extends FEMContext in order to provide extra data required to perform local element residu...
Real abs_training_tolerance
SparseMatrix< Number > * get_non_dirichlet_inner_product_matrix()
Get the non-Dirichlet (or more generally no-constraints) version of the inner-product matrix.
Real rel_training_tolerance
Relative and absolute tolerances for training reduced basis using the Greedy scheme.
The libMesh namespace provides an interface to certain functionality in the library.
std::vector< Number > Fq_representor_innerprods
Vectors storing the residual representor inner products to be used in computing the residuals online.
virtual void compute_output_dual_innerprods()
Compute and store the dual norm of each output functional.
virtual ~RBConstruction()
Destructor.
virtual void post_process_elem_matrix_and_vector(DGFEMContext &)
This function is called from add_scaled_matrix_and_vector() before each element matrix and vector are...
std::vector< std::unique_ptr< NumericVector< Number > > > Fq_representor
Vector storing the residual representors associated with the right-hand side.
bool assert_convergence
A boolean flag to indicate whether to check for proper convergence after each solve.
This class is part of the rbOOmit framework.
virtual SparseMatrix< Number > & get_matrix_for_output_dual_solves()
Return the matrix for the output residual dual norm solves.
static std::unique_ptr< DirichletBoundary > build_zero_dirichlet_boundary_object()
It's helpful to be able to generate a DirichletBoundary that stores a ZeroFunction in order to impose...
unsigned int number() const
std::vector< std::unique_ptr< SparseMatrix< Number > > > Aq_vector
Vector storing the Q_a matrices from the affine expansion.
std::vector< std::unique_ptr< NumericVector< Number > > > non_dirichlet_Fq_vector
virtual void init_context(FEMContext &)
Initialize the FEMContext prior to performing an element loop.
std::vector< std::vector< std::unique_ptr< NumericVector< Number > > > > outputs_vector
The libMesh vectors that define the output functionals.
RBEvaluation * rb_eval
The current RBEvaluation object we are using to perform the Evaluation stage of the reduced basis met...
void set_inner_product_assembly(ElemAssembly &inner_product_assembly_in)
Set the rb_assembly_expansion object.
bool use_empty_rb_solve_in_greedy
A boolean flag to indicate whether or not we initialize the Greedy algorithm by performing rb_solves ...
virtual void recompute_all_residual_terms(const bool compute_inner_products=true)
This function computes all of the residual representors, can be useful when restarting a basis traini...
void zero_constrained_dofs_on_vector(NumericVector< Number > &vector)
It is sometimes useful to be able to zero vector entries that correspond to constrained dofs.
virtual bool greedy_termination_test(Real abs_greedy_error, Real initial_greedy_error, int count)
Function that indicates when to terminate the Greedy basis training.
virtual void assemble_misc_matrices()
Assemble and store all the inner-product matrix, the constraint matrix (for constrained problems) and...
void update_greedy_param_list()
Update the list of Greedily chosen parameters with current_parameters.
void print_basis_function_orthogonality()
Print out a matrix that shows the orthogonality of the RB basis functions.
This class stores the set of ElemAssembly functor objects that define the "parameter-independent expa...
bool get_convergence_assertion_flag() const
Getter for the flag determining if convergence should be checked after each solve.
SparseMatrix< Number > * get_inner_product_matrix()
Get a pointer to inner_product_matrix.
ElemAssembly * inner_product_assembly
Pointer to inner product assembly.
bool skip_degenerate_sides
In some cases meshes are intentionally created with degenerate sides as a way to represent,...
virtual void update_residual_terms(bool compute_inner_products=true)
Compute the terms that are combined ‘online’ to determine the dual norm of the residual.
std::unique_ptr< LinearSolver< Number > > inner_product_solver
We store an extra linear solver object which we can optionally use for solving all systems in which t...
virtual void clear() override
Clear all the data structures associated with the system.
virtual void allocate_data_structures()
Helper function that actually allocates all the data structures required by this class.
void assemble_Fq_vector(unsigned int q, NumericVector< Number > *input_vector, bool apply_dof_constraints=true)
Assemble the q^th affine vector and store it in input_matrix.
virtual void update_system()
Update the system after enriching the RB space; this calls a series of functions to update the system...
void set_convergence_assertion_flag(bool flag)
Setter for the flag determining if convergence should be checked after each solve.
RBConstruction(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
void set_rel_training_tolerance(Real new_training_tolerance)
Get/set the relative tolerance for the basis training.
SparseMatrix< Number > * get_non_dirichlet_Aq_if_avail(unsigned int q)
Get a pointer to non_dirichlet_Aq if it's available, otherwise get Aq.
bool is_rb_eval_initialized() const
virtual void write_riesz_representors_to_files(const std::string &riesz_representors_dir, const bool write_binary_residual_representors)
Write out all the Riesz representor data to files.
virtual void get_output_vectors(std::map< std::string, NumericVector< Number > * > &all_vectors)
Get a map that stores pointers to all of the vectors.
std::vector< Number > energy_inner_product_coeffs
We may optionally want to use the "energy inner-product" rather than the inner-product assembly speci...
std::vector< std::vector< std::unique_ptr< NumericVector< Number > > > > non_dirichlet_outputs_vector
NumericVector< Number > * get_output_vector(unsigned int n, unsigned int q_l)
Get a pointer to the n^th output.
virtual Real compute_max_error_bound()
(i) Compute the a posteriori error bound for each set of parameters in the training set,...
LinearSolver< Number > * extra_linear_solver
Also, we store a pointer to an extra linear solver.
void assemble_Aq_matrix(unsigned int q, SparseMatrix< Number > *input_matrix, bool apply_dof_constraints=true)
Assemble the q^th affine matrix and store it in input_matrix.
This class is part of the rbOOmit framework.
virtual void solve_for_matrix_and_rhs(LinearSolver< Number > &input_solver, SparseMatrix< Number > &input_matrix, NumericVector< Number > &input_rhs)
Assembles & solves the linear system A*x=b for the specified matrix input_matrix and right-hand side ...
RBEvaluation & get_rb_evaluation()
Get a reference to the RBEvaluation object.
RBThetaExpansion & get_rb_theta_expansion()
Get a reference to the RBThetaExpansion object that that belongs to rb_eval.
void enrich_basis_from_rhs_terms(const bool resize_rb_eval_data=true)
This function computes one basis function for each rhs term.
This class is part of the rbOOmit framework.
RBConstructionBase< LinearImplicitSystem > Parent
The type of the parent.
unsigned int get_Nmax() const
Get/set Nmax, the maximum number of RB functions we are willing to compute.
void set_rb_assembly_expansion(RBAssemblyExpansion &rb_assembly_expansion_in)
Set the rb_assembly_expansion object.
void add_scaled_Aq(Number scalar, unsigned int q_a, SparseMatrix< Number > *input_matrix, bool symmetrize)
Add the scaled q^th affine matrix to input_matrix.
const RBParameters & get_greedy_parameter(unsigned int i)
Return the parameters chosen during the i^th step of the Greedy algorithm.
RBConstruction sys_type
The type of system.
void check_convergence(LinearSolver< Number > &input_solver)
Check if the linear solver reports convergence.
This is the EquationSystems class.
virtual void print_info()
Print out info that describes the current setup of this RBConstruction.
SparseMatrix< Number > * get_non_dirichlet_inner_product_matrix_if_avail()
Get the non-Dirichlet inner-product matrix if it's available, otherwise get the inner-product matrix ...
virtual void assemble_affine_expansion(bool skip_matrix_assembly, bool skip_vector_assembly)
Assemble the matrices and vectors for this system.
bool get_normalize_rb_bound_in_greedy()
void set_normalize_rb_bound_in_greedy(bool normalize_rb_bound_in_greedy_in)
Get/set the boolean to indicate if we normalize the RB error in the greedy.
virtual void update_RB_system_matrices()
Compute the reduced basis matrices for the current basis.
std::vector< std::unique_ptr< NumericVector< Number > > > Fq_vector
Vector storing the Q_f vectors in the affine decomposition of the right-hand side.
virtual void read_riesz_representors_from_files(const std::string &riesz_representors_dir, const bool write_binary_residual_representors)
Read in all the Riesz representor data from files.
std::vector< Number > truth_outputs
Vector storing the truth output values from the most recent truth solve.
Real compute_residual_dual_norm_slow(const unsigned int N)
The slow (but simple, non-error prone) way to compute the residual dual norm.
bool skip_residual_in_train_reduced_basis
Boolean flag to indicate if we skip residual calculations in train_reduced_basis.
bool compute_RB_inner_product
Boolean flag to indicate whether we compute the RB_inner_product_matrix.
SparseMatrix< Number > * get_non_dirichlet_Aq(unsigned int q)
Get a pointer to non_dirichlet_Aq.
unsigned int delta_N
The number of basis functions that we add at each greedy step.
unsigned int get_delta_N() const
Get delta_N, the number of basis functions we add to the RB space per iteration of the greedy algorit...
ElemAssembly provides a per-element (interior and boundary) assembly functionality.
std::vector< std::unique_ptr< SparseMatrix< Number > > > non_dirichlet_Aq_vector
We may also need a second set of matrices/vectors that do not have the Dirichlet boundary conditions ...
bool normalize_rb_bound_in_greedy
This boolean indicates if we normalize the RB error in the greedy using RBEvaluation::get_error_bound...
void add_scaled_matrix_and_vector(Number scalar, ElemAssembly *elem_assembly, SparseMatrix< Number > *input_matrix, NumericVector< Number > *input_vector, bool symmetrize=false, bool apply_dof_constraints=true)
This function loops over the mesh and applies the specified interior and/or boundary assembly routine...
RBAssemblyExpansion * rb_assembly_expansion
This member holds the (parameter independent) assembly functors that define the "affine expansion" of...
virtual void load_basis_function(unsigned int i)
Load the i^th RB function into the RBConstruction solution vector.
std::vector< std::vector< Number > > output_dual_innerprods
The vector storing the dual norm inner product terms for each output.
const std::string & name() const
bool impose_internal_fluxes
Boolean flag to indicate whether we impose "fluxes" (i.e.
This class is part of the rbOOmit framework.
Real get_abs_training_tolerance()
virtual void load_rb_solution()
Load the RB solution from the most recent solve with rb_eval into this system's solution vector.
Real get_rel_training_tolerance()
bool use_energy_inner_product
Boolean to indicate whether we're using the energy inner-product.
bool exit_on_repeated_greedy_parameters
Boolean flag to indicate whether we exit the greedy if we select the same parameters twice in a row.
std::vector< Real > training_error_bounds
Vector storing the values of the error bound for each parameter in the training set — the parameter g...
virtual void truth_assembly()
Assemble the truth matrix and right-hand side for current_parameters.
void set_abs_training_tolerance(Real new_training_tolerance)
Get/set the absolute tolerance for the basis training.
virtual void get_all_matrices(std::map< std::string, SparseMatrix< Number > * > &all_matrices)
Get a map that stores pointers to all of the matrices.
NumericVector< Number > * get_Fq(unsigned int q)
Get a pointer to Fq.
bool output_dual_innerprods_computed
A boolean flag to indicate whether or not the output dual norms have already been computed — used to ...
std::unique_ptr< SparseMatrix< Number > > non_dirichlet_inner_product_matrix
virtual void assemble_all_affine_vectors()
Assemble and store the affine RHS vectors.
bool store_non_dirichlet_operators
Boolean flag to indicate whether we store a second copy of each affine operator and vector which does...
virtual Real truth_solve(int plot_solution)
Perform a "truth" solve, i.e.
NumericVector< Number > * get_non_dirichlet_output_vector(unsigned int n, unsigned int q_l)
Get a pointer to non-Dirichlet output vector.
virtual std::string system_type() const override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void assemble_inner_product_matrix(SparseMatrix< Number > *input_matrix, bool apply_dof_constraints=true)
Assemble the inner product matrix and store it in input_matrix.
virtual void compute_Fq_representor_innerprods(bool compute_inner_products=true)
Compute the terms that are combined ‘online’ to determine the dual norm of the residual.
virtual void set_context_solution_vec(NumericVector< Number > &vec)
Set current_local_solution = vec so that we can access vec from FEMContext during assembly.
virtual void assemble_all_output_vectors()
Assemble and store the output vectors.
virtual void assemble_all_affine_operators()
Assemble and store all Q_a affine operators as well as the inner-product matrix.
This class stores the set of RBTheta functor objects that define the "parameter-dependent expansion" ...
bool Fq_representor_innerprods_computed
A boolean flag to indicate whether or not the Fq representor norms have already been computed — used ...
virtual void set_Nmax(unsigned int Nmax)
ElemAssembly & get_inner_product_assembly()
void set_energy_inner_product(const std::vector< Number > &energy_inner_product_coeffs_in)
Specify the coefficients of the A_q operators to be used in the energy inner-product.
virtual Real train_reduced_basis(const bool resize_rb_eval_data=true)
Train the reduced basis.
RBAssemblyExpansion & get_rb_assembly_expansion()
void set_rb_construction_parameters(unsigned int n_training_samples_in, bool deterministic_training_in, unsigned int training_parameters_random_seed_in, bool quiet_mode_in, unsigned int Nmax_in, Real rel_training_tolerance_in, Real abs_training_tolerance_in, bool normalize_rb_error_bound_in_greedy_in, RBParameters mu_min_in, RBParameters mu_max_in, std::map< std::string, std::vector< Real >> discrete_parameter_values_in, std::map< std::string, bool > log_scaling)
Set the state of this RBConstruction object based on the arguments to this function.
This class provides all data required for a physics package (e.g.
virtual std::unique_ptr< DGFEMContext > build_context()
Builds a DGFEMContext object with enough information to do evaluations on each element.