transient_rb_construction.h
Go to the documentation of this file.
1 // rbOOmit: An implementation of the Certified Reduced Basis method.
2 // Copyright (C) 2009, 2010 David J. Knezevic
3 
4 // This file is part of rbOOmit.
5 
6 // rbOOmit is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License, or (at your option) any later version.
10 
11 // rbOOmit is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
15 
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19 
20 #ifndef LIBMESH_TRANSIENT_RB_CONSTRUCTION_H
21 #define LIBMESH_TRANSIENT_RB_CONSTRUCTION_H
22 
23 // rbOOmit includes
27 
28 // libMesh includes
30 
31 // C++ includes
32 
33 namespace libMesh
34 {
35 
49 {
50 public:
51 
57  const std::string & name,
58  const unsigned int number);
59 
63  virtual ~TransientRBConstruction ();
64 
69 
74 
79  virtual void clear () override;
80 
89  virtual void initialize_rb_construction(bool skip_matrix_assembly=false,
90  bool skip_vector_assembly=false) override;
91 
95  virtual Real truth_solve(int write_interval) override;
96 
104  virtual Real train_reduced_basis(const bool resize_rb_eval_data=true) override;
105 
110  virtual void process_parameters_file (const std::string & parameters_filename) override;
111 
115  virtual void print_info() override;
116 
121  virtual bool greedy_termination_test(Real abs_greedy_error,
122  Real initial_greedy_error,
123  int count) override;
124 
129  virtual void assemble_all_affine_operators() override;
130 
134  virtual void assemble_misc_matrices() override;
135 
139  void assemble_L2_matrix(SparseMatrix<Number> * input_matrix,
140  bool apply_dirichlet_bc=true);
141 
146  void assemble_mass_matrix(SparseMatrix<Number> * input_matrix);
147 
152  void add_scaled_mass_matrix(Number scalar,
153  SparseMatrix<Number> * input_matrix);
154 
159  void mass_matrix_scaled_matvec(Number scalar,
160  NumericVector<Number> & dest,
161  NumericVector<Number> & arg);
162 
166  void set_L2_assembly(ElemAssembly & L2_assembly_in);
167 
172 
176  void assemble_Mq_matrix(unsigned int q,
177  SparseMatrix<Number> * input_matrix,
178  bool apply_dirichlet_bc=true);
179 
183  SparseMatrix<Number> * get_M_q(unsigned int q);
184 
189 
193  virtual void get_all_matrices(std::map<std::string, SparseMatrix<Number> *> & all_matrices) override;
194 
198  virtual void truth_assembly() override;
199 
208  int get_max_truth_solves() const { return max_truth_solves; }
209  void set_max_truth_solves(int max_truth_solves_in) { this->max_truth_solves = max_truth_solves_in; }
210 
214  Real get_POD_tol() const { return POD_tol; }
215  void set_POD_tol(const Real POD_tol_in) { this->POD_tol = POD_tol_in; }
216 
221  void set_delta_N(const unsigned int new_delta_N) { this->delta_N = new_delta_N; }
222 
227  virtual void load_rb_solution() override;
228 
236 
243 
248  virtual void write_riesz_representors_to_files(const std::string & riesz_representors_dir,
249  const bool write_binary_residual_representors) override;
250 
255  virtual void read_riesz_representors_from_files(const std::string & riesz_representors_dir,
256  const bool write_binary_residual_representors) override;
257 
258 
259  //----------- PUBLIC DATA MEMBERS -----------//
260 
264  std::unique_ptr<SparseMatrix<Number>> L2_matrix;
265 
270  std::unique_ptr<SparseMatrix<Number>> non_dirichlet_L2_matrix;
271 
275  std::vector<std::unique_ptr<SparseMatrix<Number>>> M_q_vector;
276 
282  std::vector<std::unique_ptr<SparseMatrix<Number>>> non_dirichlet_M_q_vector;
283 
288  std::vector<std::vector<Number>> truth_outputs_all_k;
289 
295 
302 
307  std::string init_filename;
308 
309 protected:
310 
315  virtual void allocate_data_structures() override;
316 
321  virtual void assemble_affine_expansion(bool skip_matrix_assembly,
322  bool skip_vector_assembly) override;
323 
329  virtual void initialize_truth();
330 
336 
341  void add_IC_to_RB_space();
342 
348  virtual void enrich_RB_space() override;
349 
353  virtual void update_system() override;
354 
358  virtual void update_RB_system_matrices() override;
359 
364  virtual void update_residual_terms(bool compute_inner_products) override;
365 
372 
373  //----------- PROTECTED DATA MEMBERS -----------//
374 
381 
389 
394 
400 
401 private:
402 
403  //----------- PRIVATE DATA MEMBERS -----------//
404 
408  std::vector<std::unique_ptr<NumericVector<Number>>> temporal_data;
409 };
410 
411 } // namespace libMesh
412 
413 #endif // LIBMESH_TRANSIENT_RB_CONSTRUCTION_H
virtual void print_info() override
SparseMatrix< Number > * get_non_dirichlet_M_q(unsigned int q)
std::unique_ptr< SparseMatrix< Number > > L2_matrix
virtual void read_riesz_representors_from_files(const std::string &riesz_representors_dir, const bool write_binary_residual_representors) override
virtual void allocate_data_structures() override
virtual void initialize_rb_construction(bool skip_matrix_assembly=false, bool skip_vector_assembly=false) override
Manages multiples systems of equations.
virtual void truth_assembly() override
virtual void assemble_all_affine_operators() override
void set_max_truth_solves(int max_truth_solves_in)
void mass_matrix_scaled_matvec(Number scalar, NumericVector< Number > &dest, NumericVector< Number > &arg)
virtual void get_all_matrices(std::map< std::string, SparseMatrix< Number > *> &all_matrices) override
std::vector< std::unique_ptr< SparseMatrix< Number > > > non_dirichlet_M_q_vector
TransientSystem< RBConstruction > Parent
void assemble_Mq_matrix(unsigned int q, SparseMatrix< Number > *input_matrix, bool apply_dirichlet_bc=true)
Manages storage and variables for transient systems.
std::vector< std::vector< Number > > truth_outputs_all_k
virtual Real truth_solve(int write_interval) override
TransientRBConstruction(EquationSystems &es, const std::string &name, const unsigned int number)
void assemble_mass_matrix(SparseMatrix< Number > *input_matrix)
virtual bool greedy_termination_test(Real abs_greedy_error, Real initial_greedy_error, int count) override
virtual void load_rb_solution() override
virtual void assemble_misc_matrices() override
unsigned int number() const
Definition: system.h:2025
std::unique_ptr< SparseMatrix< Number > > non_dirichlet_L2_matrix
void assemble_L2_matrix(SparseMatrix< Number > *input_matrix, bool apply_dirichlet_bc=true)
SparseMatrix< Number > * get_M_q(unsigned int q)
virtual void update_residual_terms(bool compute_inner_products) override
const NumericVector< Number > & get_error_temporal_data()
virtual Real train_reduced_basis(const bool resize_rb_eval_data=true) override
void set_POD_tol(const Real POD_tol_in)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void update_RB_system_matrices() override
std::vector< std::unique_ptr< NumericVector< Number > > > temporal_data
virtual void enrich_RB_space() override
virtual void clear() override
virtual void write_riesz_representors_to_files(const std::string &riesz_representors_dir, const bool write_binary_residual_representors) override
virtual void assemble_affine_expansion(bool skip_matrix_assembly, bool skip_vector_assembly) override
Encapsulates the details of the generalized Euler discretization.
void add_scaled_mass_matrix(Number scalar, SparseMatrix< Number > *input_matrix)
void set_L2_assembly(ElemAssembly &L2_assembly_in)
const std::string & name() const
Definition: system.h:2017
void set_delta_N(const unsigned int new_delta_N)
virtual SparseMatrix< Number > & get_matrix_for_output_dual_solves() override
virtual void update_system() override
std::vector< std::unique_ptr< SparseMatrix< Number > > > M_q_vector
virtual void process_parameters_file(const std::string &parameters_filename) override