libMesh
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_RB_CONSTRUCTION_H
21 #define LIBMESH_RB_CONSTRUCTION_H
22 
23 // rbOOmit includes
24 #include "libmesh/rb_construction_base.h"
25 
26 // libMesh includes
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"
32 
33 // C++ includes
34 
35 namespace libMesh
36 {
37 
38 class RBThetaExpansion;
39 class RBAssemblyExpansion;
40 class RBEvaluation;
41 class ElemAssembly;
42 
53 class RBConstruction : public RBConstructionBase<LinearImplicitSystem>
54 {
55 public:
56 
62  const std::string & name,
63  const unsigned int number);
64 
68  virtual ~RBConstruction ();
69 
74 
79  virtual void solve_for_matrix_and_rhs (LinearSolver<Number> & input_solver,
80  SparseMatrix<Number> & input_matrix,
81  NumericVector<Number> & input_rhs);
82 
86  void set_rb_evaluation(RBEvaluation & rb_eval_in);
87 
92 
96  bool is_rb_eval_initialized() const;
97 
103 
107  void set_rb_assembly_expansion(RBAssemblyExpansion & rb_assembly_expansion_in);
108 
113 
117  sys_type & system () { return *this; }
118 
123 
128  virtual void clear () override;
129 
133  virtual std::string system_type () const override;
134 
141  virtual Real truth_solve(int plot_solution);
142 
159  virtual Real train_reduced_basis(const bool resize_rb_eval_data=true);
160 
166  void enrich_basis_from_rhs_terms(const bool resize_rb_eval_data=true);
167 
173  virtual Real compute_max_error_bound();
174 
179  const RBParameters & get_greedy_parameter(unsigned int i);
180 
184  void set_rel_training_tolerance(Real new_training_tolerance)
185  {this->rel_training_tolerance = new_training_tolerance; }
187 
191  void set_abs_training_tolerance(Real new_training_tolerance)
192  {this->abs_training_tolerance = new_training_tolerance; }
194 
198  void set_normalize_rb_bound_in_greedy(bool normalize_rb_bound_in_greedy_in)
199  {this->normalize_rb_bound_in_greedy = normalize_rb_bound_in_greedy_in; }
201 
206  unsigned int get_Nmax() const { return Nmax; }
207  virtual void set_Nmax(unsigned int Nmax);
208 
213  virtual void load_basis_function(unsigned int i);
214 
219  virtual void load_rb_solution();
220 
225  Real compute_residual_dual_norm_slow(const unsigned int N);
226 
234 
241 
247 
251  SparseMatrix<Number> * get_Aq(unsigned int q);
252 
256  SparseMatrix<Number> * get_non_dirichlet_Aq(unsigned int q);
257 
263 
273  virtual void initialize_rb_construction(bool skip_matrix_assembly=false,
274  bool skip_vector_assembly=false);
275 
279  NumericVector<Number> * get_Fq(unsigned int q);
280 
285 
291 
295  NumericVector<Number> * get_output_vector(unsigned int n, unsigned int q_l);
296 
300  NumericVector<Number> * get_non_dirichlet_output_vector(unsigned int n, unsigned int q_l);
301 
305  virtual void get_all_matrices(std::map<std::string, SparseMatrix<Number> *> & all_matrices);
306 
310  virtual void get_all_vectors(std::map<std::string, NumericVector<Number> *> & all_vectors);
311 
315  virtual void get_output_vectors(std::map<std::string, NumericVector<Number> *> & all_vectors);
316 
322  virtual void assemble_affine_expansion(bool skip_matrix_assembly,
323  bool skip_vector_assembly);
324 
328  void assemble_inner_product_matrix(SparseMatrix<Number> * input_matrix, bool apply_dof_constraints=true);
329 
333  void assemble_Aq_matrix(unsigned int q, SparseMatrix<Number> * input_matrix, bool apply_dof_constraints=true);
334 
338  void assemble_Fq_vector(unsigned int q, NumericVector<Number> * input_vector, bool apply_dof_constraints=true);
339 
344  void add_scaled_Aq(Number scalar, unsigned int q_a,
345  SparseMatrix<Number> * input_matrix,
346  bool symmetrize);
347 
353  virtual void write_riesz_representors_to_files(const std::string & riesz_representors_dir,
354  const bool write_binary_residual_representors);
355 
362  virtual void read_riesz_representors_from_files(const std::string & riesz_representors_dir,
363  const bool write_binary_residual_representors);
364 
365 
373  virtual void recompute_all_residual_terms(const bool compute_inner_products=true);
374 
379  virtual void process_parameters_file(const std::string & parameters_filename);
380 
385  void set_rb_construction_parameters(unsigned int n_training_samples_in,
386  bool deterministic_training_in,
387  unsigned int training_parameters_random_seed_in,
388  bool quiet_mode_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,
393  RBParameters mu_min_in,
394  RBParameters mu_max_in,
395  std::map<std::string, std::vector<Real>> discrete_parameter_values_in,
396  std::map<std::string,bool> log_scaling);
397 
401  virtual void print_info();
402 
409 
417  unsigned int get_delta_N() const { return delta_N; }
418 
422  void set_inner_product_assembly(ElemAssembly & inner_product_assembly_in);
423 
428 
433  void set_energy_inner_product(const std::vector<Number> & energy_inner_product_coeffs_in);
434 
440 
441 #ifdef LIBMESH_ENABLE_DIRICHLET
442 
446  static std::unique_ptr<DirichletBoundary> build_zero_dirichlet_boundary_object();
447 #endif
448 
453  void set_convergence_assertion_flag(bool flag);
454 
455  //----------- PUBLIC DATA MEMBERS -----------//
456 
464  std::vector<Real> training_error_bounds;
465 
471  std::unique_ptr<LinearSolver<Number>> inner_product_solver;
472 
480 
484  std::unique_ptr<SparseMatrix<Number>> inner_product_matrix;
485 
490  std::vector<Number > truth_outputs;
491 
496  std::vector<std::vector<Number >> output_dual_innerprods;
497 
504  std::vector<std::unique_ptr<NumericVector<Number>>> Fq_representor;
505 
513  std::vector<Number> Fq_representor_innerprods;
514 
523 
531 
538 
546 
554 
561 
568 
575 
576 protected:
577 
582  virtual void allocate_data_structures();
583 
588  virtual void truth_assembly();
589 
595  virtual std::unique_ptr<DGFEMContext> build_context();
596 
603 
608  virtual bool greedy_termination_test(Real abs_greedy_error,
609  Real initial_greedy_error,
610  int count);
611 
617 
626  ElemAssembly * elem_assembly,
627  SparseMatrix<Number> * input_matrix,
628  NumericVector<Number> * input_vector,
629  bool symmetrize=false,
630  bool apply_dof_constraints=true);
631 
641  (DGFEMContext & /*context*/) {}
642 
651  virtual void post_process_truth_solution() {}
652 
659 
666  virtual void assemble_misc_matrices();
667 
672  virtual void assemble_all_affine_operators();
673 
677  virtual void assemble_all_affine_vectors();
678 
682  virtual void assemble_all_output_vectors();
683 
687  virtual void compute_output_dual_innerprods();
688 
701  virtual void compute_Fq_representor_innerprods(bool compute_inner_products=true);
702 
707  virtual void enrich_RB_space();
708 
713  virtual void update_system();
714 
720  virtual Real get_RB_error_bound();
721 
725  virtual void update_RB_system_matrices();
726 
734  virtual void update_residual_terms(bool compute_inner_products=true);
735 
742  virtual void init_context(FEMContext &) {}
743 
748  bool get_convergence_assertion_flag() const;
749 
754  void check_convergence(LinearSolver<Number> & input_solver);
755 
756  //----------- PROTECTED DATA MEMBERS -----------//
757 
762  unsigned int Nmax;
763 
768  unsigned int delta_N;
769 
776 
782 
783 private:
784 
785  //----------- PRIVATE DATA MEMBERS -----------//
786 
793 
799 
804 
810 
821  std::vector<Number> energy_inner_product_coeffs;
822 
826  std::vector<std::unique_ptr<SparseMatrix<Number>>> Aq_vector;
827 
832  std::vector<std::unique_ptr<NumericVector<Number>>> Fq_vector;
833 
838  std::vector<std::vector<std::unique_ptr<NumericVector<Number>>>> outputs_vector;
839 
845  std::vector<std::unique_ptr<SparseMatrix<Number>>> non_dirichlet_Aq_vector;
846  std::vector<std::unique_ptr<NumericVector<Number>>> non_dirichlet_Fq_vector;
847  std::vector<std::vector<std::unique_ptr<NumericVector<Number>>>> non_dirichlet_outputs_vector;
848  std::unique_ptr<SparseMatrix<Number>> non_dirichlet_inner_product_matrix;
849 
856 
862 
863 };
864 
865 } // namespace libMesh
866 
867 #endif // LIBMESH_RB_CONSTRUCTION_H
libMesh::RBConstruction::inner_product_matrix
std::unique_ptr< SparseMatrix< Number > > inner_product_matrix
The inner product matrix.
Definition: rb_construction.h:484
libMesh::RBConstruction::enrich_RB_space
virtual void enrich_RB_space()
Add a new basis function to the RB space.
Definition: rb_construction.C:1316
libMesh::RBConstruction::post_process_truth_solution
virtual void post_process_truth_solution()
Similarly, provide an opportunity to post-process the truth solution after the solve is complete.
Definition: rb_construction.h:651
libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::RBConstruction::get_RB_error_bound
virtual Real get_RB_error_bound()
Definition: rb_construction.C:1356
libMesh::RBConstruction::get_non_dirichlet_Fq
NumericVector< Number > * get_non_dirichlet_Fq(unsigned int q)
Get a pointer to non-Dirichlet Fq.
Definition: rb_construction.C:1913
libMesh::RBConstruction::get_Aq
SparseMatrix< Number > * get_Aq(unsigned int q)
Get a pointer to Aq.
Definition: rb_construction.C:1876
libMesh::RBConstruction::Nmax
unsigned int Nmax
Maximum number of reduced basis functions we are willing to use.
Definition: rb_construction.h:762
libMesh::RBConstruction::process_parameters_file
virtual void process_parameters_file(const std::string &parameters_filename)
Read in from the file specified by parameters_filename and set the this system's member variables acc...
Definition: rb_construction.C:193
libMesh::RBConstruction::get_all_vectors
virtual void get_all_vectors(std::map< std::string, NumericVector< Number > * > &all_vectors)
Get a map that stores pointers to all of the vectors.
Definition: rb_construction.C:1977
libMesh::RBConstruction::get_non_dirichlet_Fq_if_avail
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.
Definition: rb_construction.C:1924
libMesh::RBConstruction::initialize_rb_construction
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.
Definition: rb_construction.C:419
libMesh::RBConstruction::set_rb_evaluation
void set_rb_evaluation(RBEvaluation &rb_eval_in)
Set the RBEvaluation object.
Definition: rb_construction.C:170
libMesh::DGFEMContext
This class extends FEMContext in order to provide extra data required to perform local element residu...
Definition: dg_fem_context.h:39
libMesh::RBConstruction::abs_training_tolerance
Real abs_training_tolerance
Definition: rb_construction.h:855
libMesh::RBConstruction::get_non_dirichlet_inner_product_matrix
SparseMatrix< Number > * get_non_dirichlet_inner_product_matrix()
Get the non-Dirichlet (or more generally no-constraints) version of the inner-product matrix.
Definition: rb_construction.C:1858
libMesh::RBConstruction::rel_training_tolerance
Real rel_training_tolerance
Relative and absolute tolerances for training reduced basis using the Greedy scheme.
Definition: rb_construction.h:854
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::RBConstruction::Fq_representor_innerprods
std::vector< Number > Fq_representor_innerprods
Vectors storing the residual representor inner products to be used in computing the residuals online.
Definition: rb_construction.h:513
libMesh::RBConstruction::compute_output_dual_innerprods
virtual void compute_output_dual_innerprods()
Compute and store the dual norm of each output functional.
Definition: rb_construction.C:1635
libMesh::RBConstruction::~RBConstruction
virtual ~RBConstruction()
Destructor.
Definition: rb_construction.C:94
libMesh::RBConstruction::post_process_elem_matrix_and_vector
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...
Definition: rb_construction.h:641
libMesh::RBConstruction::Fq_representor
std::vector< std::unique_ptr< NumericVector< Number > > > Fq_representor
Vector storing the residual representors associated with the right-hand side.
Definition: rb_construction.h:504
libMesh::RBConstruction::assert_convergence
bool assert_convergence
A boolean flag to indicate whether to check for proper convergence after each solve.
Definition: rb_construction.h:781
libMesh::RBParameters
This class is part of the rbOOmit framework.
Definition: rb_parameters.h:42
libMesh::RBConstruction::get_matrix_for_output_dual_solves
virtual SparseMatrix< Number > & get_matrix_for_output_dual_solves()
Return the matrix for the output residual dual norm solves.
Definition: rb_construction.C:1630
libMesh::RBConstruction::build_zero_dirichlet_boundary_object
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...
Definition: rb_construction.C:2018
libMesh::System::number
unsigned int number() const
Definition: system.h:2075
libMesh::RBConstruction::Aq_vector
std::vector< std::unique_ptr< SparseMatrix< Number > > > Aq_vector
Vector storing the Q_a matrices from the affine expansion.
Definition: rb_construction.h:826
libMesh::RBConstruction::non_dirichlet_Fq_vector
std::vector< std::unique_ptr< NumericVector< Number > > > non_dirichlet_Fq_vector
Definition: rb_construction.h:846
libMesh::RBConstruction::init_context
virtual void init_context(FEMContext &)
Initialize the FEMContext prior to performing an element loop.
Definition: rb_construction.h:742
libMesh::RBConstruction::outputs_vector
std::vector< std::vector< std::unique_ptr< NumericVector< Number > > > > outputs_vector
The libMesh vectors that define the output functionals.
Definition: rb_construction.h:838
libMesh::RBConstruction::rb_eval
RBEvaluation * rb_eval
The current RBEvaluation object we are using to perform the Evaluation stage of the reduced basis met...
Definition: rb_construction.h:792
libMesh::RBConstruction::set_inner_product_assembly
void set_inner_product_assembly(ElemAssembly &inner_product_assembly_in)
Set the rb_assembly_expansion object.
Definition: rb_construction.C:379
libMesh::RBConstruction::use_empty_rb_solve_in_greedy
bool use_empty_rb_solve_in_greedy
A boolean flag to indicate whether or not we initialize the Greedy algorithm by performing rb_solves ...
Definition: rb_construction.h:567
libMesh::RBConstruction::recompute_all_residual_terms
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...
Definition: rb_construction.C:1380
libMesh::RBConstruction::zero_constrained_dofs_on_vector
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.
Definition: rb_construction.C:402
libMesh::SparseMatrix< Number >
libMesh::RBConstruction::greedy_termination_test
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.
Definition: rb_construction.C:1192
libMesh::RBConstruction::assemble_misc_matrices
virtual void assemble_misc_matrices()
Assemble and store all the inner-product matrix, the constraint matrix (for constrained problems) and...
Definition: rb_construction.C:918
libMesh::RBConstruction::update_greedy_param_list
void update_greedy_param_list()
Update the list of Greedily chosen parameters with current_parameters.
Definition: rb_construction.C:1229
libMesh::NumericVector< Number >
libMesh::RBConstruction::print_basis_function_orthogonality
void print_basis_function_orthogonality()
Print out a matrix that shows the orthogonality of the RB basis functions.
Definition: rb_construction.C:348
libMesh::RBAssemblyExpansion
This class stores the set of ElemAssembly functor objects that define the "parameter-independent expa...
Definition: rb_assembly_expansion.h:44
libMesh::RBConstruction::get_convergence_assertion_flag
bool get_convergence_assertion_flag() const
Getter for the flag determining if convergence should be checked after each solve.
Definition: rb_construction.C:2245
libMesh::RBConstruction::get_inner_product_matrix
SparseMatrix< Number > * get_inner_product_matrix()
Get a pointer to inner_product_matrix.
Definition: rb_construction.C:1853
libMesh::RBConstruction::inner_product_assembly
ElemAssembly * inner_product_assembly
Pointer to inner product assembly.
Definition: rb_construction.h:803
libMesh::RBConstruction::skip_degenerate_sides
bool skip_degenerate_sides
In some cases meshes are intentionally created with degenerate sides as a way to represent,...
Definition: rb_construction.h:545
libMesh::RBConstruction::update_residual_terms
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.
Definition: rb_construction.C:1535
libMesh::RBConstruction::inner_product_solver
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...
Definition: rb_construction.h:471
libMesh::RBConstruction::clear
virtual void clear() override
Clear all the data structures associated with the system.
Definition: rb_construction.C:100
libMesh::RBConstruction::allocate_data_structures
virtual void allocate_data_structures()
Helper function that actually allocates all the data structures required by this class.
Definition: rb_construction.C:467
libMesh::RBConstruction::assemble_Fq_vector
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.
Definition: rb_construction.C:971
libMesh::RBConstruction::update_system
virtual void update_system()
Update the system after enriching the RB space; this calls a series of functions to update the system...
Definition: rb_construction.C:1350
libMesh::RBConstruction::set_convergence_assertion_flag
void set_convergence_assertion_flag(bool flag)
Setter for the flag determining if convergence should be checked after each solve.
Definition: rb_construction.C:2250
libMesh::RBConstruction::RBConstruction
RBConstruction(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
Definition: rb_construction.C:62
libMesh::RBConstruction::set_rel_training_tolerance
void set_rel_training_tolerance(Real new_training_tolerance)
Get/set the relative tolerance for the basis training.
Definition: rb_construction.h:184
libMesh::RBConstruction::get_non_dirichlet_Aq_if_avail
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.
Definition: rb_construction.C:1895
libMesh::RBConstruction::is_rb_eval_initialized
bool is_rb_eval_initialized() const
Definition: rb_construction.C:183
libMesh::RBConstruction::write_riesz_representors_to_files
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.
Definition: rb_construction.C:2031
libMesh::RBConstruction::get_output_vectors
virtual void get_output_vectors(std::map< std::string, NumericVector< Number > * > &all_vectors)
Get a map that stores pointers to all of the vectors.
Definition: rb_construction.C:1997
libMesh::RBConstruction::energy_inner_product_coeffs
std::vector< Number > energy_inner_product_coeffs
We may optionally want to use the "energy inner-product" rather than the inner-product assembly speci...
Definition: rb_construction.h:821
libMesh::RBConstruction::non_dirichlet_outputs_vector
std::vector< std::vector< std::unique_ptr< NumericVector< Number > > > > non_dirichlet_outputs_vector
Definition: rb_construction.h:847
libMesh::RBConstruction::get_output_vector
NumericVector< Number > * get_output_vector(unsigned int n, unsigned int q_l)
Get a pointer to the n^th output.
Definition: rb_construction.C:1934
libMesh::RBConstruction::compute_max_error_bound
virtual Real compute_max_error_bound()
(i) Compute the a posteriori error bound for each set of parameters in the training set,...
Definition: rb_construction.C:1395
libMesh::RBConstruction::extra_linear_solver
LinearSolver< Number > * extra_linear_solver
Also, we store a pointer to an extra linear solver.
Definition: rb_construction.h:479
libMesh::RBConstruction::assemble_Aq_matrix
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.
Definition: rb_construction.C:876
libMesh::RBConstructionBase
This class is part of the rbOOmit framework.
Definition: rb_construction_base.h:54
libMesh::RBConstruction::solve_for_matrix_and_rhs
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 ...
Definition: rb_construction.C:130
libMesh::RBConstruction::get_rb_evaluation
RBEvaluation & get_rb_evaluation()
Get a reference to the RBEvaluation object.
Definition: rb_construction.C:175
libMesh::RBConstruction::get_rb_theta_expansion
RBThetaExpansion & get_rb_theta_expansion()
Get a reference to the RBThetaExpansion object that that belongs to rb_eval.
Definition: rb_construction.C:188
libMesh::RBConstruction::enrich_basis_from_rhs_terms
void enrich_basis_from_rhs_terms(const bool resize_rb_eval_data=true)
This function computes one basis function for each rhs term.
Definition: rb_construction.C:1134
libMesh::RBConstruction
This class is part of the rbOOmit framework.
Definition: rb_construction.h:53
libMesh::RBConstruction::system
sys_type & system()
Definition: rb_construction.h:117
libMesh::RBConstruction::Parent
RBConstructionBase< LinearImplicitSystem > Parent
The type of the parent.
Definition: rb_construction.h:122
libMesh::RBConstruction::get_Nmax
unsigned int get_Nmax() const
Get/set Nmax, the maximum number of RB functions we are willing to compute.
Definition: rb_construction.h:206
libMesh::RBConstruction::set_rb_assembly_expansion
void set_rb_assembly_expansion(RBAssemblyExpansion &rb_assembly_expansion_in)
Set the rb_assembly_expansion object.
Definition: rb_construction.C:366
libMesh::RBConstruction::add_scaled_Aq
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.
Definition: rb_construction.C:893
libMesh::RBConstruction::get_greedy_parameter
const RBParameters & get_greedy_parameter(unsigned int i)
Return the parameters chosen during the i^th step of the Greedy algorithm.
Definition: rb_construction.C:1234
libMesh::RBConstruction::sys_type
RBConstruction sys_type
The type of system.
Definition: rb_construction.h:73
libMesh::RBConstruction::check_convergence
void check_convergence(LinearSolver< Number > &input_solver)
Check if the linear solver reports convergence.
Definition: rb_construction.C:2231
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::RBConstruction::print_info
virtual void print_info()
Print out info that describes the current setup of this RBConstruction.
Definition: rb_construction.C:312
libMesh::RBConstruction::get_non_dirichlet_inner_product_matrix_if_avail
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 ...
Definition: rb_construction.C:1866
libMesh::RBConstruction::assemble_affine_expansion
virtual void assemble_affine_expansion(bool skip_matrix_assembly, bool skip_vector_assembly)
Assemble the matrices and vectors for this system.
Definition: rb_construction.C:449
libMesh::RBConstruction::get_normalize_rb_bound_in_greedy
bool get_normalize_rb_bound_in_greedy()
Definition: rb_construction.h:200
libMesh::RBConstruction::set_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.
Definition: rb_construction.h:198
libMesh::RBConstruction::update_RB_system_matrices
virtual void update_RB_system_matrices()
Compute the reduced basis matrices for the current basis.
Definition: rb_construction.C:1466
libMesh::RBConstruction::Fq_vector
std::vector< std::unique_ptr< NumericVector< Number > > > Fq_vector
Vector storing the Q_f vectors in the affine decomposition of the right-hand side.
Definition: rb_construction.h:832
libMesh::RBConstruction::read_riesz_representors_from_files
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.
Definition: rb_construction.C:2140
libMesh::RBConstruction::truth_outputs
std::vector< Number > truth_outputs
Vector storing the truth output values from the most recent truth solve.
Definition: rb_construction.h:490
libMesh::RBConstruction::compute_residual_dual_norm_slow
Real compute_residual_dual_norm_slow(const unsigned int N)
The slow (but simple, non-error prone) way to compute the residual dual norm.
Definition: rb_construction.C:1819
libMesh::RBConstruction::skip_residual_in_train_reduced_basis
bool skip_residual_in_train_reduced_basis
Boolean flag to indicate if we skip residual calculations in train_reduced_basis.
Definition: rb_construction.h:522
libMesh::RBConstruction::compute_RB_inner_product
bool compute_RB_inner_product
Boolean flag to indicate whether we compute the RB_inner_product_matrix.
Definition: rb_construction.h:553
libMesh::RBConstruction::get_non_dirichlet_Aq
SparseMatrix< Number > * get_non_dirichlet_Aq(unsigned int q)
Get a pointer to non_dirichlet_Aq.
Definition: rb_construction.C:1884
libMesh::RBConstruction::delta_N
unsigned int delta_N
The number of basis functions that we add at each greedy step.
Definition: rb_construction.h:768
libMesh::RBConstruction::get_delta_N
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...
Definition: rb_construction.h:417
libMesh::ElemAssembly
ElemAssembly provides a per-element (interior and boundary) assembly functionality.
Definition: elem_assembly.h:38
libMesh::RBConstruction::non_dirichlet_Aq_vector
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 ...
Definition: rb_construction.h:845
libMesh::RBConstruction::normalize_rb_bound_in_greedy
bool normalize_rb_bound_in_greedy
This boolean indicates if we normalize the RB error in the greedy using RBEvaluation::get_error_bound...
Definition: rb_construction.h:861
libMesh::RBConstruction::add_scaled_matrix_and_vector
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...
Definition: rb_construction.C:588
libMesh::RBConstruction::rb_assembly_expansion
RBAssemblyExpansion * rb_assembly_expansion
This member holds the (parameter independent) assembly functors that define the "affine expansion" of...
Definition: rb_construction.h:798
libMesh::RBConstruction::load_basis_function
virtual void load_basis_function(unsigned int i)
Load the i^th RB function into the RBConstruction solution vector.
Definition: rb_construction.C:1305
libMesh::RBConstruction::output_dual_innerprods
std::vector< std::vector< Number > > output_dual_innerprods
The vector storing the dual norm inner product terms for each output.
Definition: rb_construction.h:496
libMesh::System::name
const std::string & name() const
Definition: system.h:2067
libMesh::RBConstruction::impose_internal_fluxes
bool impose_internal_fluxes
Boolean flag to indicate whether we impose "fluxes" (i.e.
Definition: rb_construction.h:537
libMesh::RBEvaluation
This class is part of the rbOOmit framework.
Definition: rb_evaluation.h:50
libMesh::RBConstruction::get_abs_training_tolerance
Real get_abs_training_tolerance()
Definition: rb_construction.h:193
libMesh::LinearSolver< Number >
libMesh::RBConstruction::load_rb_solution
virtual void load_rb_solution()
Load the RB solution from the most recent solve with rb_eval into this system's solution vector.
Definition: rb_construction.C:1800
libMesh::RBConstruction::get_rel_training_tolerance
Real get_rel_training_tolerance()
Definition: rb_construction.h:186
libMesh::RBConstruction::use_energy_inner_product
bool use_energy_inner_product
Boolean to indicate whether we're using the energy inner-product.
Definition: rb_construction.h:809
libMesh::RBConstruction::exit_on_repeated_greedy_parameters
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.
Definition: rb_construction.h:530
libMesh::RBConstruction::training_error_bounds
std::vector< Real > training_error_bounds
Vector storing the values of the error bound for each parameter in the training set — the parameter g...
Definition: rb_construction.h:464
libMesh::RBConstruction::truth_assembly
virtual void truth_assembly()
Assemble the truth matrix and right-hand side for current_parameters.
Definition: rb_construction.C:805
libMesh::RBConstruction::set_abs_training_tolerance
void set_abs_training_tolerance(Real new_training_tolerance)
Get/set the absolute tolerance for the basis training.
Definition: rb_construction.h:191
libMesh::RBConstruction::get_all_matrices
virtual void get_all_matrices(std::map< std::string, SparseMatrix< Number > * > &all_matrices)
Get a map that stores pointers to all of the matrices.
Definition: rb_construction.C:1952
libMesh::RBConstruction::get_Fq
NumericVector< Number > * get_Fq(unsigned int q)
Get a pointer to Fq.
Definition: rb_construction.C:1905
libMesh::RBConstruction::output_dual_innerprods_computed
bool output_dual_innerprods_computed
A boolean flag to indicate whether or not the output dual norms have already been computed — used to ...
Definition: rb_construction.h:775
libMesh::RBConstruction::non_dirichlet_inner_product_matrix
std::unique_ptr< SparseMatrix< Number > > non_dirichlet_inner_product_matrix
Definition: rb_construction.h:848
libMesh::RBConstruction::assemble_all_affine_vectors
virtual void assemble_all_affine_vectors()
Assemble and store the affine RHS vectors.
Definition: rb_construction.C:950
libMesh::RBConstruction::store_non_dirichlet_operators
bool store_non_dirichlet_operators
Boolean flag to indicate whether we store a second copy of each affine operator and vector which does...
Definition: rb_construction.h:560
libMesh::RBConstruction::truth_solve
virtual Real truth_solve(int plot_solution)
Perform a "truth" solve, i.e.
Definition: rb_construction.C:1242
libMesh::RBConstruction::get_non_dirichlet_output_vector
NumericVector< Number > * get_non_dirichlet_output_vector(unsigned int n, unsigned int q_l)
Get a pointer to non-Dirichlet output vector.
Definition: rb_construction.C:1943
libMesh::RBConstruction::system_type
virtual std::string system_type() const override
Definition: rb_construction.C:125
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::RBConstruction::assemble_inner_product_matrix
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.
Definition: rb_construction.C:841
libMesh::RBConstruction::compute_Fq_representor_innerprods
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.
Definition: rb_construction.C:1730
libMesh::RBConstruction::set_context_solution_vec
virtual void set_context_solution_vec(NumericVector< Number > &vec)
Set current_local_solution = vec so that we can access vec from FEMContext during assembly.
Definition: rb_construction.C:797
libMesh::RBConstruction::assemble_all_output_vectors
virtual void assemble_all_output_vectors()
Assemble and store the output vectors.
Definition: rb_construction.C:988
libMesh::RBConstruction::assemble_all_affine_operators
virtual void assemble_all_affine_operators()
Assemble and store all Q_a affine operators as well as the inner-product matrix.
Definition: rb_construction.C:930
libMesh::RBThetaExpansion
This class stores the set of RBTheta functor objects that define the "parameter-dependent expansion" ...
Definition: rb_theta_expansion.h:44
libMesh::RBConstruction::Fq_representor_innerprods_computed
bool Fq_representor_innerprods_computed
A boolean flag to indicate whether or not the Fq representor norms have already been computed — used ...
Definition: rb_construction.h:574
libMesh::RBConstruction::set_Nmax
virtual void set_Nmax(unsigned int Nmax)
Definition: rb_construction.C:1300
libMesh::RBConstruction::get_inner_product_assembly
ElemAssembly & get_inner_product_assembly()
Definition: rb_construction.C:385
libMesh::RBConstruction::set_energy_inner_product
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.
Definition: rb_construction.C:396
libMesh::RBConstruction::train_reduced_basis
virtual Real train_reduced_basis(const bool resize_rb_eval_data=true)
Train the reduced basis.
Definition: rb_construction.C:1024
libMesh::RBConstruction::get_rb_assembly_expansion
RBAssemblyExpansion & get_rb_assembly_expansion()
Definition: rb_construction.C:371
libMesh::RBConstruction::set_rb_construction_parameters
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.
Definition: rb_construction.C:272
libMesh::FEMContext
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:62
libMesh::RBConstruction::build_context
virtual std::unique_ptr< DGFEMContext > build_context()
Builds a DGFEMContext object with enough information to do evaluations on each element.
Definition: rb_construction.C:583