20 #ifndef LIBMESH_TRANSIENT_RB_CONSTRUCTION_H 21 #define LIBMESH_TRANSIENT_RB_CONSTRUCTION_H 24 #include "libmesh/rb_construction.h" 25 #include "libmesh/transient_rb_evaluation.h" 26 #include "libmesh/rb_temporal_discretization.h" 29 #include "libmesh/transient_system.h" 57 const std::string &
name,
58 const unsigned int number);
85 virtual void clear ()
override;
96 bool skip_vector_assembly=
false)
override;
128 Real initial_greedy_error,
146 bool apply_dirichlet_bc=
true);
184 bool apply_dirichlet_bc=
true);
255 const bool write_binary_residual_representors)
override;
262 const bool write_binary_residual_representors)
override;
281 std::vector<std::unique_ptr<SparseMatrix<Number>>>
M_q_vector;
328 bool skip_vector_assembly)
override;
419 #endif // LIBMESH_TRANSIENT_RB_CONSTRUCTION_H std::unique_ptr< SparseMatrix< Number > > L2_matrix
The L2 matrix.
virtual void allocate_data_structures() override
Helper function that actually allocates all the data structures required by this class.
virtual void load_rb_solution() override
Load the RB solution from the current time-level into the libMesh solution vector.
This is the EquationSystems class.
int max_truth_solves
Maximum number of truth solves in the POD-Greedy.
Real get_POD_tol() const
Get/set POD_tol.
void set_max_truth_solves(int max_truth_solves_in)
void assemble_mass_matrix(SparseMatrix< Number > *input_matrix)
Assemble the mass matrix at the current parameter and store it in input_matrix.
virtual void assemble_all_affine_operators() override
Assemble and store all the affine operators.
virtual SparseMatrix< Number > & get_matrix_for_output_dual_solves() override
Override to return the L2 product matrix for output dual norm solves for transient state problems...
SparseMatrix< Number > * get_M_q(unsigned int q)
Get a pointer to M_q.
virtual Real train_reduced_basis(const bool resize_rb_eval_data=true) override
Train the reduced basis.
virtual void truth_assembly() override
Assemble the truth system in the transient linear case.
std::vector< std::unique_ptr< SparseMatrix< Number > > > non_dirichlet_M_q_vector
We sometimes also need a second set of M_q matrices that do not have the Dirichlet boundary condition...
virtual void initialize_truth()
This function imposes a truth initial condition, defaults to zero initial condition if the flag nonze...
virtual void get_all_matrices(std::map< std::string, SparseMatrix< Number > *> &all_matrices) override
Get a map that stores pointers to all of the matrices.
TransientSystem< RBConstruction > Parent
The type of the parent.
virtual Real truth_solve(int write_interval) override
Perform a truth solve at the current parameter.
Manages storage and variables for transient systems.
Real POD_tol
If positive, this tolerance determines the number of POD modes we add to the space on a call to enric...
std::vector< std::vector< Number > > truth_outputs_all_k
The truth outputs for all time-levels from the most recent truth_solve.
void update_RB_initial_condition_all_N()
Compute the L2 projection of the initial condition onto the RB space for 1 <= N <= RB_size and store ...
The libMesh namespace provides an interface to certain functionality in the library.
std::string init_filename
The filename of the file containing the initial condition projected onto the truth mesh...
virtual void update_RB_system_matrices() override
Compute the reduced basis matrices for the current basis.
virtual void read_riesz_representors_from_files(const std::string &riesz_representors_dir, const bool write_binary_residual_representors) override
Write out all the Riesz representor data to files.
void assemble_L2_matrix(SparseMatrix< Number > *input_matrix, bool apply_dirichlet_bc=true)
Assemble the L2 matrix.
unsigned int number() const
virtual void write_riesz_representors_to_files(const std::string &riesz_representors_dir, const bool write_binary_residual_representors) override
Write out all the Riesz representor data to files.
std::unique_ptr< SparseMatrix< Number > > non_dirichlet_L2_matrix
The L2 matrix without Dirichlet conditions enforced.
void set_L2_assembly(ElemAssembly &L2_assembly_in)
Set the L2 object.
virtual void enrich_RB_space() override
Add a new basis functions to the RB space.
int get_max_truth_solves() const
Get/set max_truth_solves, the maximum number of RB truth solves we are willing to compute in the tran...
Number set_error_temporal_data()
Set column k (i.e.
virtual void clear() override
Clear all the data structures associated with the system.
bool nonzero_initialization
Boolean flag to indicate whether we are using a non-zero initialization.
This class is part of the rbOOmit framework.
ElemAssembly * L2_assembly
Function pointer for assembling the L2 matrix.
ElemAssembly & get_L2_assembly()
unsigned int delta_N
The number of basis functions that we add at each greedy step.
virtual void update_residual_terms(bool compute_inner_products) override
Compute the terms that are combined ‘online’ to determine the dual norm of the residual.
virtual bool greedy_termination_test(Real abs_greedy_error, Real initial_greedy_error, int count) override
Function that indicates when to terminate the Greedy basis training.
TransientRBConstruction sys_type
The type of system.
TransientRBConstruction(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
void set_POD_tol(const Real POD_tol_in)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void assemble_affine_expansion(bool skip_matrix_assembly, bool skip_vector_assembly) override
Override assemble_affine_expansion to also initialize RB_ic_proj_rhs_all_N, if necessary.
std::vector< std::unique_ptr< NumericVector< Number > > > temporal_data
Dense matrix to store the data that we use for the temporal POD.
virtual void process_parameters_file(const std::string ¶meters_filename) override
Read in the parameters from file and set up the system accordingly.
ElemAssembly provides a per-element (interior and boundary) assembly functionality.
void assemble_Mq_matrix(unsigned int q, SparseMatrix< Number > *input_matrix, bool apply_dirichlet_bc=true)
Assemble the q^th affine term of the mass matrix and store it in input_matrix.
void add_scaled_mass_matrix(Number scalar, SparseMatrix< Number > *input_matrix)
Add the scaled mass matrix (assembled for the current parameter) to input_matrix. ...
virtual ~TransientRBConstruction()
const NumericVector< Number > & get_error_temporal_data()
Get the column of temporal_data corresponding to the current time level.
Define a class that encapsulates the details of a "generalized Euler" temporal discretization to be u...
bool compute_truth_projection_error
Boolean flag that indicates whether we will compute the projection error for the truth solution into ...
void mass_matrix_scaled_matvec(Number scalar, NumericVector< Number > &dest, NumericVector< Number > &arg)
Perform a matrix-vector multiplication with the current mass matrix and store the result in dest...
const std::string & name() const
void set_delta_N(const unsigned int new_delta_N)
Set delta_N, the number of basis functions we add to the RB space from each POD.
TransientRBConstruction & operator=(const TransientRBConstruction &)=delete
void add_IC_to_RB_space()
Initialize RB space by adding the truth initial condition as the first RB basis function.
SparseMatrix< Number > * get_non_dirichlet_M_q(unsigned int q)
Get a pointer to non_dirichlet_M_q.
std::vector< std::unique_ptr< SparseMatrix< Number > > > M_q_vector
Vector storing the Q_m matrices from the mass operator.
virtual void assemble_misc_matrices() override
Override to assemble the L2 matrix as well.
virtual void update_system() override
Update the system after enriching the RB space.
virtual void print_info() const override
Print out info that describes the current setup of this RBConstruction.
DenseVector< Number > RB_ic_proj_rhs_all_N
The vector that stores the right-hand side for the initial condition projections. ...
virtual void initialize_rb_construction(bool skip_matrix_assembly=false, bool skip_vector_assembly=false) override
Allocate all the data structures necessary for the construction stage of the RB method.