libMesh
rb_scm_construction.C
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 // Configuration data
21 #include "libmesh/libmesh_config.h"
22 
23 // Currently, the RBSCMConstruction should only be available
24 // if SLEPc support is enabled.
25 #if defined(LIBMESH_HAVE_SLEPC) && (LIBMESH_HAVE_GLPK)
26 
27 #include "libmesh/rb_scm_construction.h"
28 #include "libmesh/rb_construction.h"
29 #include "libmesh/rb_scm_evaluation.h"
30 
31 #include "libmesh/libmesh_logging.h"
32 #include "libmesh/numeric_vector.h"
33 #include "libmesh/sparse_matrix.h"
34 #include "libmesh/equation_systems.h"
35 #include "libmesh/getpot.h"
36 #include "libmesh/dof_map.h"
37 #include "libmesh/enum_eigen_solver_type.h"
38 
39 // For creating a directory
40 #include <sys/types.h>
41 #include <sys/stat.h>
42 #include <errno.h>
43 
44 namespace libMesh
45 {
46 
48  const std::string & name_in,
49  const unsigned int number_in)
50  : Parent(es, name_in, number_in),
51  SCM_training_tolerance(0.5),
52  RB_system_name(""),
53  rb_scm_eval(nullptr)
54 {
55 
56  // set assemble_before_solve flag to false
57  // so that we control matrix assembly.
58  assemble_before_solve = false;
59 
60  // We symmetrize all operators hence use symmetric solvers
62 
63 }
64 
66 {
67  this->clear();
68 }
69 
70 
72 {
73  Parent::clear();
74 }
75 
77 {
78  rb_scm_eval = &rb_scm_eval_in;
79 }
80 
82 {
83  if (!rb_scm_eval)
84  libmesh_error_msg("Error: RBSCMEvaluation object hasn't been initialized yet");
85 
86  return *rb_scm_eval;
87 }
88 
90 {
92 }
93 
94 void RBSCMConstruction::process_parameters_file(const std::string & parameters_filename)
95 {
96  // First read in data from parameters_filename
97  GetPot infile(parameters_filename);
98  const unsigned int n_training_samples = infile("n_training_samples",1);
99  const bool deterministic_training = infile("deterministic_training",false);
100 
101  // Read in training_parameters_random_seed value. This is used to
102  // seed the RNG when picking the training parameters. By default the
103  // value is -1, which means use std::time to seed the RNG.
104  unsigned int training_parameters_random_seed_in = static_cast<int>(-1);
105  training_parameters_random_seed_in = infile("training_parameters_random_seed",
106  training_parameters_random_seed_in);
107  set_training_random_seed(training_parameters_random_seed_in);
108 
109  // SCM Greedy termination tolerance
110  const Real SCM_training_tolerance_in = infile("SCM_training_tolerance", SCM_training_tolerance);
111  set_SCM_training_tolerance(SCM_training_tolerance_in);
112 
113  // Initialize the parameter ranges and the parameters themselves
114  unsigned int n_continuous_parameters = infile.vector_variable_size("parameter_names");
115  RBParameters mu_min_in;
116  RBParameters mu_max_in;
117  for (unsigned int i=0; i<n_continuous_parameters; i++)
118  {
119  // Read in the parameter names
120  std::string param_name = infile("parameter_names", "NONE", i);
121 
122  {
123  Real min_val = infile(param_name, 0., 0);
124  mu_min_in.set_value(param_name, min_val);
125  }
126 
127  {
128  Real max_val = infile(param_name, 0., 1);
129  mu_max_in.set_value(param_name, max_val);
130  }
131  }
132 
133  std::map<std::string, std::vector<Real>> discrete_parameter_values_in;
134 
135  unsigned int n_discrete_parameters = infile.vector_variable_size("discrete_parameter_names");
136  for (unsigned int i=0; i<n_discrete_parameters; i++)
137  {
138  std::string param_name = infile("discrete_parameter_names", "NONE", i);
139 
140  unsigned int n_vals_for_param = infile.vector_variable_size(param_name);
141  std::vector<Real> vals_for_param(n_vals_for_param);
142  for (unsigned int j=0; j != n_vals_for_param; j++)
143  vals_for_param[j] = infile(param_name, 0., j);
144 
145  discrete_parameter_values_in[param_name] = vals_for_param;
146  }
147 
148  initialize_parameters(mu_min_in, mu_max_in, discrete_parameter_values_in);
149 
150  std::map<std::string,bool> log_scaling;
151  const RBParameters & mu = get_parameters();
152  unsigned int i=0;
153  for (const auto & pr : mu)
154  {
155  const std::string & param_name = pr.first;
156  log_scaling[param_name] = static_cast<bool>(infile("log_scaling", 0, i++));
157  }
158 
160  this->get_parameters_max(),
161  n_training_samples,
162  log_scaling,
163  deterministic_training); // use deterministic parameters
164 }
165 
167 {
168  // Print out info that describes the current setup
169  libMesh::out << std::endl << "RBSCMConstruction parameters:" << std::endl;
170  libMesh::out << "system name: " << this->name() << std::endl;
171  libMesh::out << "SCM Greedy tolerance: " << get_SCM_training_tolerance() << std::endl;
172  if (rb_scm_eval)
173  {
174  libMesh::out << "A_q operators attached: " << get_rb_theta_expansion().get_n_A_terms() << std::endl;
175  }
176  else
177  {
178  libMesh::out << "RBThetaExpansion member is not set yet" << std::endl;
179  }
180  libMesh::out << "Number of parameters: " << get_n_params() << std::endl;
181  for (const auto & pr : get_parameters())
182  {
183  const std::string & param_name = pr.first;
184  libMesh::out << "Parameter " << param_name
185  << ": Min = " << get_parameter_min(param_name)
186  << ", Max = " << get_parameter_max(param_name) << std::endl;
187  }
189  libMesh::out << "n_training_samples: " << get_n_training_samples() << std::endl;
190  libMesh::out << std::endl;
191 }
192 
194 {
195  // Clear SCM data vectors
196  rb_scm_eval->B_min.clear();
197  rb_scm_eval->B_max.clear();
198  rb_scm_eval->C_J.clear();
200  for (auto & vec : rb_scm_eval->SCM_UB_vectors)
201  vec.clear();
202  rb_scm_eval->SCM_UB_vectors.clear();
203 
204  // Resize the bounding box vectors
205  rb_scm_eval->B_min.resize(get_rb_theta_expansion().get_n_A_terms());
206  rb_scm_eval->B_max.resize(get_rb_theta_expansion().get_n_A_terms());
207 }
208 
209 void RBSCMConstruction::add_scaled_symm_Aq(unsigned int q_a, Number scalar)
210 {
211  LOG_SCOPE("add_scaled_symm_Aq()", "RBSCMConstruction");
212  // Load the operators from the RBConstruction
213  EquationSystems & es = this->get_equation_systems();
215  rb_system.add_scaled_Aq(scalar, q_a, matrix_A.get(), true);
216 }
217 
219 {
220  // Load the operators from the RBConstruction
221  EquationSystems & es = this->get_equation_systems();
223 
224  matrix_B->zero();
225  matrix_B->close();
226  matrix_B->add(1.,*rb_system.get_inner_product_matrix());
227 }
228 
230 {
231  LOG_SCOPE("perform_SCM_greedy()", "RBSCMConstruction");
232 
233  // initialize rb_scm_eval's parameters
235 
236 #ifdef LIBMESH_ENABLE_CONSTRAINTS
237  // Get a list of constrained dofs from rb_system
238  std::set<dof_id_type> constrained_dofs_set;
239  EquationSystems & es = this->get_equation_systems();
241 
242  for (dof_id_type i=0; i<rb_system.n_dofs(); i++)
243  {
244  if (rb_system.get_dof_map().is_constrained_dof(i))
245  {
246  constrained_dofs_set.insert(i);
247  }
248  }
249 
250  // Use these constrained dofs to identify which dofs we want to "get rid of"
251  // (i.e. condense) in our eigenproblems.
252  this->initialize_condensed_dofs(constrained_dofs_set);
253 #endif // LIBMESH_ENABLE_CONSTRAINTS
254 
255  // Copy the inner product matrix over from rb_system to be used as matrix_B
256  load_matrix_B();
257 
259 
261  // This loads the new parameter into current_parameters
262  enrich_C_J(0);
263 
264  unsigned int SCM_iter=0;
265  while (true)
266  {
267  // matrix_A is reinitialized for the current parameters
268  // on each call to evaluate_stability_constant
270 
271  std::pair<unsigned int,Real> SCM_error_pair = compute_SCM_bounds_on_training_set();
272 
273  libMesh::out << "SCM iteration " << SCM_iter
274  << ", max_SCM_error = " << SCM_error_pair.second << std::endl;
275 
276  if (SCM_error_pair.second < SCM_training_tolerance)
277  {
278  libMesh::out << std::endl << "SCM tolerance of " << SCM_training_tolerance << " reached."
279  << std::endl << std::endl;
280  break;
281  }
282 
283  // If we need another SCM iteration, then enrich C_J
284  enrich_C_J(SCM_error_pair.first);
285 
286  libMesh::out << std::endl << "-----------------------------------" << std::endl << std::endl;
287 
288  SCM_iter++;
289  }
290 }
291 
293 {
294  LOG_SCOPE("compute_SCM_bounding_box()", "RBSCMConstruction");
295 
296  // Resize the bounding box vectors
297  rb_scm_eval->B_min.resize(get_rb_theta_expansion().get_n_A_terms());
298  rb_scm_eval->B_max.resize(get_rb_theta_expansion().get_n_A_terms());
299 
300  for (unsigned int q=0; q<get_rb_theta_expansion().get_n_A_terms(); q++)
301  {
302  matrix_A->zero();
303  add_scaled_symm_Aq(q, 1.);
304 
305  // Compute B_min(q)
306  eigen_solver->set_position_of_spectrum(SMALLEST_REAL);
308 
309  solve();
310  // TODO: Assert convergence for eigensolver
311 
312  unsigned int nconv = get_n_converged();
313  if (nconv != 0)
314  {
315  std::pair<Real, Real> eval = get_eigenpair(0);
316 
317  // ensure that the eigenvalue is real
318  libmesh_assert_less (eval.second, TOLERANCE);
319 
320  rb_scm_eval->set_B_min(q, eval.first);
321  libMesh::out << std::endl << "B_min("<<q<<") = " << rb_scm_eval->get_B_min(q) << std::endl;
322  }
323  else
324  libmesh_error_msg("Eigen solver for computing B_min did not converge");
325 
326  // Compute B_max(q)
327  eigen_solver->set_position_of_spectrum(LARGEST_REAL);
329 
330  solve();
331  // TODO: Assert convergence for eigensolver
332 
333  nconv = get_n_converged();
334  if (nconv != 0)
335  {
336  std::pair<Real, Real> eval = get_eigenpair(0);
337 
338  // ensure that the eigenvalue is real
339  libmesh_assert_less (eval.second, TOLERANCE);
340 
341  rb_scm_eval->set_B_max(q,eval.first);
342  libMesh::out << "B_max("<<q<<") = " << rb_scm_eval->get_B_max(q) << std::endl;
343  }
344  else
345  libmesh_error_msg("Eigen solver for computing B_max did not converge");
346  }
347 }
348 
350 {
351  LOG_SCOPE("evaluate_stability_constant()", "RBSCMConstruction");
352 
353  // Get current index of C_J
354  const unsigned int j = cast_int<unsigned int>(rb_scm_eval->C_J.size()-1);
355 
356  eigen_solver->set_position_of_spectrum(SMALLEST_REAL);
357 
358  // We assume B is set in system assembly
359  // For coercive problems, B is set to the inner product matrix
360  // For non-coercive time-dependent problems, B is set to the mass matrix
361 
362  // Set matrix A corresponding to mu_star
363  matrix_A->zero();
364  for (unsigned int q=0; q<get_rb_theta_expansion().get_n_A_terms(); q++)
365  {
367  }
368 
370  solve();
371  // TODO: Assert convergence for eigensolver
372 
373  unsigned int nconv = get_n_converged();
374  if (nconv != 0)
375  {
376  std::pair<Real, Real> eval = get_eigenpair(0);
377 
378  // ensure that the eigenvalue is real
379  libmesh_assert_less (eval.second, TOLERANCE);
380 
381  // Store the coercivity constant corresponding to mu_star
383  libMesh::out << std::endl << "Stability constant for C_J("<<j<<") = "
384  << rb_scm_eval->get_C_J_stability_constraint(j) << std::endl << std::endl;
385 
386  // Compute and store the vector y = (y_1, \ldots, y_Q) for the
387  // eigenvector currently stored in eigen_system.solution.
388  // We use this later to compute the SCM upper bounds.
390 
391  for (unsigned int q=0; q<get_rb_theta_expansion().get_n_A_terms(); q++)
392  {
393  Real norm_Aq2 = libmesh_real( Aq_inner_product(q, *solution, *solution) );
394 
395  rb_scm_eval->set_SCM_UB_vector(j,q,norm_Aq2/norm_B2);
396  }
397  }
398  else
399  libmesh_error_msg("Error: Eigensolver did not converge in evaluate_stability_constant");
400 }
401 
403  const NumericVector<Number> & w) const
404 {
405  matrix_B->vector_mult(*inner_product_storage_vector, w);
406 
408 }
409 
411  const NumericVector<Number> & v,
412  const NumericVector<Number> & w)
413 {
414  if (q >= get_rb_theta_expansion().get_n_A_terms())
415  libmesh_error_msg("Error: We must have q < Q_a in Aq_inner_product.");
416 
417  matrix_A->zero();
418  add_scaled_symm_Aq(q, 1.);
419  matrix_A->vector_mult(*inner_product_storage_vector, w);
420 
422 }
423 
425 {
426  LOG_SCOPE("compute_SCM_bounds_on_training_set()", "RBSCMConstruction");
427 
428  // Now compute the maximum bound error over training_parameters
429  unsigned int new_C_J_index = 0;
430  Real max_SCM_error = 0.;
431 
433  for (unsigned int i=0; i<get_local_n_training_samples(); i++)
434  {
435  set_params_from_training_set(first_index+i);
437  Real LB = rb_scm_eval->get_SCM_LB();
438  Real UB = rb_scm_eval->get_SCM_UB();
439 
440  Real error_i = SCM_greedy_error_indicator(LB, UB);
441 
442  if (error_i > max_SCM_error)
443  {
444  max_SCM_error = error_i;
445  new_C_J_index = i;
446  }
447  }
448 
449  numeric_index_type global_index = first_index + new_C_J_index;
450  std::pair<numeric_index_type,Real> error_pair(global_index, max_SCM_error);
451  get_global_max_error_pair(this->comm(),error_pair);
452 
453  return error_pair;
454 }
455 
456 void RBSCMConstruction::enrich_C_J(unsigned int new_C_J_index)
457 {
458  LOG_SCOPE("enrich_C_J()", "RBSCMConstruction");
459 
461 
462  rb_scm_eval->C_J.push_back(get_parameters());
463 
464  libMesh::out << std::endl << "SCM: Added mu = (";
465 
466  bool first = true;
467  for (const auto & pr : get_parameters())
468  {
469  if (!first)
470  libMesh::out << ",";
471  const std::string & param_name = pr.first;
472  RBParameters C_J_params = rb_scm_eval->C_J[rb_scm_eval->C_J.size()-1];
473  libMesh::out << C_J_params.get_value(param_name);
474  first = false;
475  }
476  libMesh::out << ")" << std::endl;
477 
478  // Finally, resize C_J_stability_vector and SCM_UB_vectors
479  rb_scm_eval->C_J_stability_vector.push_back(0.);
480 
481  std::vector<Real> zero_vector(get_rb_theta_expansion().get_n_A_terms());
482  rb_scm_eval->SCM_UB_vectors.push_back(zero_vector);
483 }
484 
485 
486 } // namespace libMesh
487 
488 #endif // LIBMESH_HAVE_SLEPC && LIBMESH_HAVE_GLPK
libMesh::RBSCMConstruction::perform_SCM_greedy
virtual void perform_SCM_greedy()
Perform the SCM greedy algorithm to develop a lower bound over the training set.
Definition: rb_scm_construction.C:229
libMesh::RBParametrized::print_discrete_parameter_values
void print_discrete_parameter_values() const
Print out all the discrete parameter values.
Definition: rb_parametrized.C:389
libMesh::RBSCMConstruction::compute_SCM_bounds_on_training_set
virtual std::pair< unsigned int, Real > compute_SCM_bounds_on_training_set()
Compute upper and lower bounds for each SCM training point.
Definition: rb_scm_construction.C:424
libMesh::RBParametrized::get_parameters_min
const RBParameters & get_parameters_min() const
Get an RBParameters object that specifies the minimum allowable value for each parameter.
Definition: rb_parametrized.C:174
libMesh::RBConstructionBase< CondensedEigenSystem >::clear
virtual void clear()
Clear all the data structures associated with the system.
Definition: rb_construction_base.C:65
libMesh::dof_id_type
uint8_t dof_id_type
Definition: id_types.h:67
libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::System::get_equation_systems
const EquationSystems & get_equation_systems() const
Definition: system.h:720
libMesh::RBSCMConstruction::resize_SCM_vectors
virtual void resize_SCM_vectors()
Clear and resize the SCM data vectors.
Definition: rb_scm_construction.C:193
libMesh::RBSCMConstruction::load_matrix_B
virtual void load_matrix_B()
Copy over the matrix to store in matrix_B, usually this is the mass or inner-product matrix,...
Definition: rb_scm_construction.C:218
libMesh::EigenSystem::matrix_A
std::unique_ptr< SparseMatrix< Number > > matrix_A
The system matrix for standard eigenvalue problems.
Definition: eigen_system.h:165
libMesh::RBSCMEvaluation::C_J
std::vector< RBParameters > C_J
Vector storing the greedily selected parameters during SCM training.
Definition: rb_scm_evaluation.h:192
libMesh::RBSCMEvaluation::set_B_min
void set_B_min(unsigned int i, Real B_min_val)
Set B_min and B_max.
Definition: rb_scm_evaluation.C:150
libMesh::RBSCMConstruction::B_inner_product
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.
Definition: rb_scm_construction.C:402
libMesh::RBConstructionBase< CondensedEigenSystem >::get_first_local_training_index
numeric_index_type get_first_local_training_index() const
Get the first local index of the training parameters.
Definition: rb_construction_base.C:116
libMesh::RBConstructionBase< CondensedEigenSystem >::inner_product_storage_vector
std::unique_ptr< NumericVector< Number > > inner_product_storage_vector
We keep an extra temporary vector that is useful for performing inner products (avoids unnecessary me...
Definition: rb_construction_base.h:254
libMesh::RBSCMEvaluation::B_max
std::vector< Real > B_max
Definition: rb_scm_evaluation.h:186
libMesh::CondensedEigenSystem::solve
virtual void solve() override
Override to solve the condensed eigenproblem with the dofs in local_non_condensed_dofs_vector strippe...
Definition: condensed_eigen_system.C:91
libMesh::SMALLEST_REAL
Definition: enum_eigen_solver_type.h:79
libMesh::RBSCMConstruction::set_rb_scm_evaluation
void set_rb_scm_evaluation(RBSCMEvaluation &rb_scm_eval_in)
Set the RBSCMEvaluation object.
Definition: rb_scm_construction.C:76
libMesh::RBSCMConstruction::SCM_training_tolerance
Real SCM_training_tolerance
Tolerance which controls when to terminate the SCM Greedy.
Definition: rb_scm_construction.h:225
libMesh::libmesh_real
T libmesh_real(T a)
Definition: libmesh_common.h:166
libMesh::LARGEST_REAL
Definition: enum_eigen_solver_type.h:78
libMesh::RBConstructionBase< CondensedEigenSystem >::set_training_random_seed
void set_training_random_seed(unsigned int seed)
Set the seed that is used to randomly generate training parameters.
Definition: rb_construction_base.C:630
libMesh::RBParametrized::get_n_params
unsigned int get_n_params() const
Get the number of parameters.
Definition: rb_parametrized.C:115
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::EquationSystems::get_system
const T_sys & get_system(const std::string &name) const
Definition: equation_systems.h:757
libMesh::RBSCMConstruction::clear
virtual void clear() override
Clear all the data structures associated with the system.
Definition: rb_scm_construction.C:71
libMesh::RBConstructionBase< CondensedEigenSystem >::set_params_from_training_set_and_broadcast
virtual void set_params_from_training_set_and_broadcast(unsigned int index)
Load the specified training parameter and then broadcast to all processors.
Definition: rb_construction_base.C:160
libMesh::TOLERANCE
static const Real TOLERANCE
Definition: libmesh_common.h:128
libMesh::ParallelObject::comm
const Parallel::Communicator & comm() const
Definition: parallel_object.h:94
libMesh::RBSCMConstruction::set_eigensolver_properties
virtual void set_eigensolver_properties(int)
This function is called before truth eigensolves in compute_SCM_bounding_box and evaluate_stability_c...
Definition: rb_scm_construction.h:127
libMesh::RBParameters
This class is part of the rbOOmit framework.
Definition: rb_parameters.h:42
libMesh::EigenSystem::eigen_solver
std::unique_ptr< EigenSolver< Number > > eigen_solver
The EigenSolver, defining which interface, i.e solver package to use.
Definition: eigen_system.h:191
libMesh::RBSCMEvaluation::SCM_UB_vectors
std::vector< std::vector< Real > > SCM_UB_vectors
This matrix stores the infimizing vectors y_1( ),...,y_Q_a( ), for each in C_J, which are used in co...
Definition: rb_scm_evaluation.h:206
libMesh::CondensedEigenSystem::get_eigenpair
virtual std::pair< Real, Real > get_eigenpair(dof_id_type i) override
Override get_eigenpair() to retrieve the eigenpair for the condensed eigensolve.
Definition: condensed_eigen_system.C:177
libMesh::EigenSystem::set_eigenproblem_type
void set_eigenproblem_type(EigenProblemType ept)
Sets the type of the current eigen problem.
Definition: eigen_system.C:84
libMesh::RBSCMConstruction::RBSCMConstruction
RBSCMConstruction(EquationSystems &es, const std::string &name_in, const unsigned int number_in)
Constructor.
Definition: rb_scm_construction.C:47
libMesh::DofMap::is_constrained_dof
bool is_constrained_dof(const dof_id_type dof) const
Definition: dof_map.h:1962
libMesh::RBParametrized::set_parameters
void set_parameters(const RBParameters &params)
Set the current parameters to params.
Definition: rb_parametrized.C:155
libMesh::RBSCMConstruction::get_rb_theta_expansion
RBThetaExpansion & get_rb_theta_expansion()
Get a reference to the RBThetaExpansion object.
Definition: rb_scm_construction.C:89
libMesh::RBConstructionBase< CondensedEigenSystem >::get_global_max_error_pair
static void get_global_max_error_pair(const Parallel::Communicator &communicator, std::pair< numeric_index_type, Real > &error_pair)
Static function to return the error pair (index,error) that is corresponds to the largest error on al...
Definition: rb_construction_base.C:85
libMesh::RBSCMEvaluation::get_SCM_UB
virtual Real get_SCM_UB()
Evaluate single SCM upper bound.
Definition: rb_scm_evaluation.C:282
libMesh::RBParametrized::get_parameters_max
const RBParameters & get_parameters_max() const
Get an RBParameters object that specifies the maximum allowable value for each parameter.
Definition: rb_parametrized.C:182
libMesh::RBParametrized::get_parameters
const RBParameters & get_parameters() const
Get the current parameters.
Definition: rb_parametrized.C:166
libMesh::RBSCMConstruction::set_SCM_training_tolerance
void set_SCM_training_tolerance(Real SCM_training_tolerance_in)
Definition: rb_scm_construction.h:140
libMesh::RBSCMConstruction::Aq_inner_product
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.
Definition: rb_scm_construction.C:410
libMesh::NumericVector< Number >
libMesh::EigenSystem::get_n_converged
unsigned int get_n_converged() const
Definition: eigen_system.h:129
libMesh::EigenSystem::matrix_B
std::unique_ptr< SparseMatrix< Number > > matrix_B
A second system matrix for generalized eigenvalue problems.
Definition: eigen_system.h:170
libMesh::System::assemble_before_solve
bool assemble_before_solve
Flag which tells the system to whether or not to call the user assembly function during each call to ...
Definition: system.h:1493
libMesh::RBSCMEvaluation::C_J_stability_vector
std::vector< Real > C_J_stability_vector
Vector storing the (truth) stability values at the parameters in C_J.
Definition: rb_scm_evaluation.h:198
libMesh::RBSCMConstruction::process_parameters_file
virtual void process_parameters_file(const std::string &parameters_filename)
Read in the parameters from file specified by parameters_filename and set the this system's member va...
Definition: rb_scm_construction.C:94
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::RBSCMConstruction::get_rb_scm_evaluation
RBSCMEvaluation & get_rb_scm_evaluation()
Get a reference to the RBSCMEvaluation object.
Definition: rb_scm_construction.C:81
libMesh::RBParametrized::get_parameter_max
Real get_parameter_max(const std::string &param_name) const
Get maximum allowable value of parameter param_name.
Definition: rb_parametrized.C:198
libMesh::RBSCMConstruction::compute_SCM_bounding_box
virtual void compute_SCM_bounding_box()
Compute the SCM bounding box.
Definition: rb_scm_construction.C:292
libMesh::RBConstructionBase< CondensedEigenSystem >::get_n_training_samples
numeric_index_type get_n_training_samples() const
Get the total number of training samples.
Definition: rb_construction_base.C:98
libMesh::RBConstructionBase< CondensedEigenSystem >::set_params_from_training_set
void set_params_from_training_set(unsigned int index)
Set parameters to the RBParameters stored in index index of the training set.
Definition: rb_construction_base.C:130
libMesh::RBSCMEvaluation::set_SCM_UB_vector
void set_SCM_UB_vector(unsigned int j, unsigned int q, Real y_q)
Set entries of SCM_UB_vector, which stores the vector y, corresponding to the minimizing eigenvectors...
Definition: rb_scm_evaluation.C:100
libMesh::RBSCMConstruction::get_SCM_training_tolerance
Real get_SCM_training_tolerance() const
Get/set SCM_training_tolerance: tolerance for SCM greedy.
Definition: rb_scm_construction.h:139
libMesh::RBConstructionBase< CondensedEigenSystem >
libMesh::RBSCMEvaluation::B_min
std::vector< Real > B_min
B_min, B_max define the bounding box.
Definition: rb_scm_evaluation.h:185
libMesh::RBParameters::set_value
void set_value(const std::string &param_name, Real value)
Set the value of the specified parameter.
Definition: rb_parameters.C:47
libMesh::RBConstruction
This class is part of the rbOOmit framework.
Definition: rb_construction.h:53
libMesh::RBSCMConstruction::enrich_C_J
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...
Definition: rb_scm_construction.C:456
libMesh::RBSCMEvaluation::get_B_max
Real get_B_max(unsigned int i) const
Definition: rb_scm_evaluation.C:142
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::RBSCMConstruction::~RBSCMConstruction
virtual ~RBSCMConstruction()
Destructor.
Definition: rb_scm_construction.C:65
libMesh::CondensedEigenSystem::initialize_condensed_dofs
void initialize_condensed_dofs(const std::set< dof_id_type > &global_condensed_dofs_set=std::set< dof_id_type >())
Loop over the dofs on each processor to initialize the list of non-condensed dofs.
Definition: condensed_eigen_system.C:47
libMesh::numeric_index_type
dof_id_type numeric_index_type
Definition: id_types.h:99
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::RBSCMConstruction::evaluate_stability_constant
virtual void evaluate_stability_constant()
Compute the stability constant for current_parameters by solving a generalized eigenvalue problem ove...
Definition: rb_scm_construction.C:349
libMesh::RBParametrized::initialize_parameters
void initialize_parameters(const RBParameters &mu_min_in, const RBParameters &mu_max_in, const std::map< std::string, std::vector< Real >> &discrete_parameter_values)
Initialize the parameter ranges and set current_parameters.
Definition: rb_parametrized.C:60
libMesh::System::solution
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Definition: system.h:1539
libMesh::RBSCMConstruction::SCM_greedy_error_indicator
virtual Real SCM_greedy_error_indicator(Real LB, Real UB)
Helper function which provides an error indicator to be used in the SCM greedy.
Definition: rb_scm_construction.h:218
libMesh::RBSCMConstruction::rb_scm_eval
RBSCMEvaluation * rb_scm_eval
The current RBSCMEvaluation object we are using to perform the Evaluation stage of the SCM.
Definition: rb_scm_construction.h:238
libMesh::System::name
const std::string & name() const
Definition: system.h:2067
libMesh::RBSCMConstruction::RB_system_name
std::string RB_system_name
The name of the associated RB system.
Definition: rb_scm_construction.h:230
libMesh::RBSCMConstruction::print_info
virtual void print_info()
Print out info that describes the current setup of this RBSCMConstruction.
Definition: rb_scm_construction.C:166
libMesh::GHEP
Definition: enum_eigen_solver_type.h:58
libMesh::RBSCMEvaluation
This class is part of the rbOOmit framework.
Definition: rb_scm_evaluation.h:52
libMesh::RBThetaExpansion::get_n_A_terms
unsigned int get_n_A_terms() const
Get Q_a, the number of terms in the affine expansion for the bilinear form.
Definition: rb_theta_expansion.C:35
libMesh::RBSCMConstruction::attach_deflation_space
virtual void attach_deflation_space()
Attach the deflation space defined by the specified vector, can be useful in solving constrained eige...
Definition: rb_scm_construction.h:156
libMesh::RBSCMEvaluation::get_SCM_LB
virtual Real get_SCM_LB()
Evaluate single SCM lower bound.
Definition: rb_scm_evaluation.C:166
libMesh::RBSCMEvaluation::get_B_min
Real get_B_min(unsigned int i) const
Get B_min and B_max.
Definition: rb_scm_evaluation.C:133
libMesh::RBSCMEvaluation::set_C_J_stability_constraint
void set_C_J_stability_constraint(unsigned int j, Real stability_constraint_in)
Set stability constraints (i.e.
Definition: rb_scm_evaluation.C:80
libMesh::RBSCMEvaluation::set_B_max
void set_B_max(unsigned int i, Real B_max_val)
Definition: rb_scm_evaluation.C:158
libMesh::System::get_dof_map
const DofMap & get_dof_map() const
Definition: system.h:2099
libMesh::System::n_dofs
dof_id_type n_dofs() const
Definition: system.C:150
libMesh::RBConstructionBase< CondensedEigenSystem >::get_local_n_training_samples
numeric_index_type get_local_n_training_samples() const
Get the total number of training samples local to this processor.
Definition: rb_construction_base.C:109
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::RBThetaExpansion
This class stores the set of RBTheta functor objects that define the "parameter-dependent expansion" ...
Definition: rb_theta_expansion.h:44
libMesh::out
OStreamProxy out
libMesh::RBSCMConstruction::add_scaled_symm_Aq
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.
Definition: rb_scm_construction.C:209
libMesh::RBConstructionBase< CondensedEigenSystem >::initialize_training_parameters
virtual void initialize_training_parameters(const RBParameters &mu_min, const RBParameters &mu_max, unsigned int n_training_parameters, std::map< std::string, bool > log_param_scale, bool deterministic=true)
Initialize the parameter ranges and indicate whether deterministic or random training parameters shou...
Definition: rb_construction_base.C:181
libMesh::NumericVector::dot
virtual T dot(const NumericVector< T > &v) const =0
libMesh::RBSCMEvaluation::get_C_J_stability_constraint
Real get_C_J_stability_constraint(unsigned int j) const
Get stability constraints (i.e.
Definition: rb_scm_evaluation.C:92
libMesh::RBParameters::get_value
Real get_value(const std::string &param_name) const
Get the value of the specific parameter.
Definition: rb_parameters.C:41
libMesh::RBSCMEvaluation::get_rb_theta_expansion
RBThetaExpansion & get_rb_theta_expansion()
Get a reference to the rb_theta_expansion.
Definition: rb_scm_evaluation.C:72
libMesh::RBParametrized::get_parameter_min
Real get_parameter_min(const std::string &param_name) const
Get minimum allowable value of parameter param_name.
Definition: rb_parametrized.C:190