libMesh
Functions
reduced_basis_ex7.C File Reference

Go to the source code of this file.

Functions

int main (int argc, char **argv)
 

Function Documentation

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 64 of file reduced_basis_ex7.C.

References libMesh::EquationSystems::add_system(), libMesh::ParallelObject::comm(), libMesh::command_line_next(), libMesh::default_solver_package(), dim, libMesh::err, libMesh::RBConstruction::get_rb_evaluation(), libMesh::TriangleWrapper::init(), libMesh::EquationSystems::init(), libMesh::RBConstruction::initialize_rb_construction(), libMesh::RBEvaluation::legacy_write_offline_data_to_files(), libMesh::RBConstruction::load_rb_solution(), mesh, libMesh::out, libMesh::PETSC_SOLVERS, libMesh::RBConstruction::print_basis_function_orthogonality(), libMesh::RBConstruction::print_info(), libMesh::EquationSystems::print_info(), libMesh::MeshBase::print_info(), libMesh::RBConstruction::process_parameters_file(), libMesh::MeshBase::read(), libMesh::RBDataDeserialization::RBEvaluationDeserialization::read_from_file(), libMesh::Real, libMesh::RBConstruction::set_rb_evaluation(), libMesh::RBParameters::set_value(), libMesh::RBConstruction::train_reduced_basis(), libMesh::ExodusII_IO::write_equation_systems(), libMesh::RBEvaluation::write_out_basis_functions(), and libMesh::RBDataSerialization::RBEvaluationSerialization::write_to_file().

65 {
66  // Initialize libMesh.
67  LibMeshInit init (argc, argv);
68 
69 #ifndef LIBMESH_USE_COMPLEX_NUMBERS
70  libmesh_example_requires(false, "--enable-complex");
71 #else
72 
73 #if !defined(LIBMESH_HAVE_XDR)
74  // We need XDR support to write out reduced bases
75  libmesh_example_requires(false, "--enable-xdr");
76 #elif defined(LIBMESH_DEFAULT_SINGLE_PRECISION)
77  // XDR binary support requires double precision
78  libmesh_example_requires(false, "--disable-singleprecision");
79 #elif defined(LIBMESH_DEFAULT_TRIPLE_PRECISION)
80  // long double doesn't work on the previous examples, don't try it
81  // here
82  libmesh_example_requires(false, "double precision");
83 #endif
84  // FIXME: This example currently segfaults with Trilinos?
85  libmesh_example_requires(libMesh::default_solver_package() == PETSC_SOLVERS, "--enable-petsc");
86 
87 #if LIBMESH_HAVE_PETSC
88  // lets parse the input to check which set of arguments was provided.
89  // Skip the runs with wrong arguments:
90  GetPot input(argc, argv);
91 
92  // skip if the specification of the solver-package is given in wrong syntax:
93 #if PETSC_VERSION_LESS_THAN(3,9,0)
94  if (input.search("-pc_factor_mat_solver_type"))
95  {
96  libMesh::out << "LibMesh was configured with PETSc < 3.9, but command-line options "
97  << "use syntax understood by later versions only. Skipping this example."
98  << std::endl;
99  return 77;
100  }
101  input.search("-pc_factor_mat_solver_package");
102 #else
103  if (input.search("-pc_factor_mat_solver_package"))
104  {
105  libMesh::out << "LibMesh was configured with PETSc >= 3.9, but command-line options "
106  << "use deprecated syntax. Skipping now."
107  << std::endl;
108  return 77;
109  }
110  input.search("-pc_factor_mat_solver_type");
111 #endif
112 
113  // check that we have the required solver:
114  std::string solver;
115  solver = input.next("invalid_solver");
116  if (solver == "mumps")
117  {
118 #ifndef LIBMESH_PETSC_HAVE_MUMPS
119  libmesh_example_requires(false, "PETSc compiled with MUMPS support");
120 #endif
121  }
122  else if (solver == "superlu")
123  {
124 #ifndef LIBMESH_PETSC_HAVE_SUPERLU_DIST
125  libmesh_example_requires(false, "PETSc compiled with SuperLU support");
126 #endif
127  }
128  else
129  {
130  libMesh::err << "Error: Solver " << solver << " is unknown."
131  << std::endl;
132  return 1;
133  }
134 
135 #endif //LIBMESH_HAVE_PETSC
136 
137  // Skip this 2D example if libMesh was compiled as 1D-only.
138  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
139 
140  // Parse the input file (reduced_basis_ex7.in) using GetPot
141  std::string parameters_filename = "reduced_basis_ex7.in";
142  GetPot infile(parameters_filename);
143 
144  const unsigned int dim = 2; // The number of spatial dimensions
145 
146  bool store_basis_functions = infile("store_basis_functions", true); // Do we write the RB basis functions to disk?
147 
148  // Read the "online_mode" flag from the command line
149  const int online_mode =
150  libMesh::command_line_next("-online_mode", 0);
151 
152  // Create a mesh on the default MPI communicator.
153  Mesh mesh(init.comm(), dim);
154  mesh.read("horn.msh");
155 
156  // Create an equation systems object.
157  EquationSystems equation_systems (mesh);
158 
159  // We override RBConstruction with SimpleRBConstruction in order to
160  // specialize a few functions for this particular problem.
161  SimpleRBConstruction & rb_con =
162  equation_systems.add_system<SimpleRBConstruction> ("Acoustics");
163 
164  // Initialize the data structures for the equation system.
165  equation_systems.init ();
166 
167  // Print out some information about the "truth" discretization
168  equation_systems.print_info();
169  mesh.print_info();
170 
171  // Build a new RBEvaluation object which will be used to perform
172  // Reduced Basis calculations. This is required in both the
173  // "Offline" and "Online" stages.
174  SimpleRBEvaluation rb_eval(mesh.comm());
175 
176  // We need to give the RBConstruction object a pointer to
177  // our RBEvaluation object
178  rb_con.set_rb_evaluation(rb_eval);
179 
180  if (!online_mode) // Perform the Offline stage of the RB method
181  {
182  // Read in the data that defines this problem from the specified text file
183  rb_con.process_parameters_file(parameters_filename);
184 
185  // Print out info that describes the current setup of rb_con
186  rb_con.print_info();
187 
188  // Prepare rb_con for the Construction stage of the RB method.
189  // This sets up the necessary data structures and performs
190  // initial assembly of the "truth" affine expansion of the PDE.
192 
193  // Compute the reduced basis space by computing "snapshots", i.e.
194  // "truth" solves, at well-chosen parameter values and employing
195  // these snapshots as basis functions.
196  libmesh_try
197  {
198  rb_con.train_reduced_basis();
199  }
200  libmesh_catch (LogicError & e)
201  {
202  libMesh::err << "\n\n"
203  << "********************************************************************************\n"
204  << "Training reduced basis failed, this example requires a direct solver.\n"
205  << "Try running with -ksp_type preonly -pc_type lu instead.\n"
206  << "********************************************************************************"
207  << std::endl;
208  return 77;
209  }
210 
211  // Write out the data that will subsequently be required for the Evaluation stage
212 #if defined(LIBMESH_HAVE_CAPNPROTO)
214  rb_eval_writer.write_to_file("rb_eval.bin");
215 #else
217 #endif
218 
219  // If requested, write out the RB basis functions for visualization purposes
220  if (store_basis_functions)
221  {
222  // Write out the basis functions
224  }
225 
226  // Basis functions should be orthonormal, so
227  // print out the inner products to check this
229  }
230  else // Perform the Online stage of the RB method
231  {
232  // Read in the reduced basis data
233 #if defined(LIBMESH_HAVE_CAPNPROTO)
235  rb_eval_reader.read_from_file("rb_eval.bin", /*read_error_bound_data*/ true);
236 #else
237  rb_eval.legacy_read_offline_data_from_files();
238 #endif
239 
240  // Read in online_N and initialize online parameters
241  Real online_frequency = infile("online_frequency", 0.);
242  RBParameters online_mu;
243  online_mu.set_value("frequency", online_frequency);
244  rb_eval.set_parameters(online_mu);
245  rb_eval.print_parameters();
246 
247  // Now do the Online solve using the precomputed reduced basis
248  rb_eval.rb_solve(rb_eval.get_n_basis_functions());
249 
250  if (store_basis_functions)
251  {
252  // Read in the basis functions
253  rb_eval.read_in_basis_functions(rb_con);
254 
255  // Plot the solution
256  rb_con.load_rb_solution();
257 
258  //Output the variable in ExodusII format (single precision)
259  ExodusII_IO(mesh, /*single_precision=*/true).write_equation_systems ("RB_sol_float.exo", equation_systems);
260  }
261 
262  // Now do a sweep over frequencies and write out the reflection coefficient for each frequency
263  std::ofstream reflection_coeffs_out("reflection_coefficients.dat");
264 
265  Real n_frequencies = infile("n_frequencies", 0.);
266  Real delta_f = (rb_eval.get_parameter_max("frequency") - rb_eval.get_parameter_min("frequency")) / (n_frequencies-1);
267  for (unsigned int freq_i=0; freq_i<n_frequencies; freq_i++)
268  {
269  Real frequency = rb_eval.get_parameter_min("frequency") + freq_i * delta_f;
270  online_mu.set_value("frequency", frequency);
271  rb_eval.set_parameters(online_mu);
272  rb_eval.rb_solve(rb_eval.get_n_basis_functions());
273 
274  Number complex_one(1., 0.);
275  reflection_coeffs_out << frequency << " " << std::abs(rb_eval.RB_outputs[0] - complex_one) << std::endl;
276  }
277  reflection_coeffs_out.close();
278 
279  }
280 
281 #endif
282 
283  return 0;
284 }
OStreamProxy err
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.
T command_line_next(std::string name, T default_value)
Use GetPot&#39;s search()/next() functions to get following arguments from the command line...
Definition: libmesh.C:1078
This is the EquationSystems class.
virtual void read(const std::string &name, void *mesh_data=nullptr, bool skip_renumber_nodes_and_elements=false, bool skip_find_neighbors=false)=0
Interfaces for reading/writing a mesh to/from a file.
virtual void write_out_basis_functions(System &sys, const std::string &directory_name="offline_data", const bool write_binary_basis_functions=true)
Write out all the basis functions to file.
virtual Real train_reduced_basis(const bool resize_rb_eval_data=true)
Train the reduced basis.
This class serializes an RBEvaluation object using the Cap&#39;n Proto library.
unsigned int dim
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:52
MeshBase & mesh
const Parallel::Communicator & comm() const
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:90
void init()
Initializes degrees of freedom on the current mesh.
Definition: system.C:197
virtual void print_info() const
Print out info that describes the current setup of this RBConstruction.
void print_basis_function_orthogonality() const
Print out a matrix that shows the orthogonality of the RB basis functions.
SolverPackage default_solver_package()
Definition: libmesh.C:1117
virtual void legacy_write_offline_data_to_files(const std::string &directory_name="offline_data", const bool write_binary_data=true)
Write out all the data to text files in order to segregate the Offline stage from the Online stage...
virtual void write_equation_systems(const std::string &fname, const EquationSystems &es, const std::set< std::string > *system_names=nullptr) override
Writes out the solution for no specific time or timestep.
Definition: exodusII_io.C:2033
void print_info(std::ostream &os=libMesh::out, const unsigned int verbosity=0, const bool global=true) const
Prints relevant information about the mesh.
Definition: mesh_base.C:1562
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
This class is part of the rbOOmit framework.
Definition: rb_parameters.h:52
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void write_to_file(const std::string &path, bool use_packing=false)
Write the Cap&#39;n&#39;Proto buffer to disk.
A class to represent the internal "this should never happen" errors, to be thrown by "libmesh_error()...
void set_value(const std::string &param_name, Real value)
Set the value of the specified parameter.
OStreamProxy out
virtual void process_parameters_file(const std::string &parameters_filename)
Read in from the file specified by parameters_filename and set the this system&#39;s member variables acc...
This class de-serializes an RBEvaluation object using the Cap&#39;n Proto library.
RBEvaluation & get_rb_evaluation()
Get a reference to the RBEvaluation object.
virtual void load_rb_solution()
Load the RB solution from the most recent solve with rb_eval into this system&#39;s solution vector...
void set_rb_evaluation(RBEvaluation &rb_eval_in)
Set the RBEvaluation object.
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50