libMesh
reduced_basis_ex2.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 // rbOOmit: An implementation of the Certified Reduced Basis method.
19 // Copyright (C) 2009, 2010 David J. Knezevic
20 // This file is part of rbOOmit.
21 
22 
23 
24 // <h1>Reduced Basis Example 2 - Successive Constraint Method</h1>
25 // \author David Knezevic
26 // \date 2011
27 //
28 // In this example we extend reduced_basis_ex1 to solve a steady
29 // convection-diffusion problem on the unit square via the Reduced
30 // Basis Method. In this case, we modify the PDE so that it no longer
31 // has a parameter-independent coercivity constant. Therefore, in
32 // order to obtain an error bound, we need to employ the Successive
33 // Constraint Method (SCM) implemented in
34 // RBSCMConstruction/RBSCMEvaluation to obtain a parameter-dependent
35 // lower bound for the coercivity constant.
36 //
37 // The PDE being solved is div(k*grad(u)) + Beta*grad(u) = f
38 // k is the diffusion coefficient :
39 // - constant in the domain 0<=x<0.5 , its value is given by the first parameter mu[0]
40 // - constant in the domain 0.5<=x<=1 , its value is given by the second parameter mu[1]
41 // Beta is the convection velocity :
42 // - constant in the whole domain
43 // - equal to zero in the y-direction
44 // - its value in the x-direction is given by the third (and last) parameter mu[2]
45 // Boundary conditions :
46 // - dyu=0 on top and bottom
47 // - u=0 on the left side
48 // - dxu + Beta*u = 0 on the right side
49 
50 // C++ include files that we need
51 #include <iostream>
52 #include <algorithm>
53 #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC
54 #include <cmath>
55 #include <set>
56 
57 // Basic include file needed for the mesh functionality.
58 #include "libmesh/libmesh.h"
59 #include "libmesh/mesh.h"
60 #include "libmesh/mesh_generation.h"
61 #include "libmesh/exodusII_io.h"
62 #include "libmesh/equation_systems.h"
63 #include "libmesh/dof_map.h"
64 #include "libmesh/getpot.h"
65 #include "libmesh/elem.h"
66 #include "libmesh/rb_data_serialization.h"
67 #include "libmesh/rb_data_deserialization.h"
68 #include "libmesh/enum_solver_package.h"
69 
70 // local includes
71 #include "rb_classes.h"
72 #include "assembly.h"
73 
74 // Bring in everything from the libMesh namespace
75 using namespace libMesh;
76 
77 // The main program.
78 int main (int argc, char ** argv)
79 {
80  // Initialize libMesh.
81  LibMeshInit init (argc, argv);
82 
83  // This example requires SLEPc and GLPK
84 #if !defined(LIBMESH_HAVE_SLEPC) || !defined(LIBMESH_HAVE_GLPK)
85  libmesh_example_requires(false, "--enable-slepc --enable-glpk");
86 #else
87 
88 #if !defined(LIBMESH_HAVE_XDR)
89  // We need XDR support to write out reduced bases
90  libmesh_example_requires(false, "--enable-xdr");
91 #elif defined(LIBMESH_DEFAULT_SINGLE_PRECISION)
92  // XDR binary support requires double precision
93  libmesh_example_requires(false, "--disable-singleprecision");
94 #endif
95  // FIXME: This example currently segfaults with Trilinos?
96  libmesh_example_requires(libMesh::default_solver_package() == PETSC_SOLVERS, "--enable-petsc");
97 
98  // FIXME: This example needs lapack, which is serial only
99  libmesh_example_requires(init.comm().size() == 1, "mpirun -np 1");
100 
101  // Skip this 2D example if libMesh was compiled as 1D-only.
102  libmesh_example_requires(2 <= LIBMESH_DIM, "2D support");
103 
104 #ifndef LIBMESH_ENABLE_DIRICHLET
105  libmesh_example_requires(false, "--enable-dirichlet");
106 #else
107 
108  // Parse the input file (reduced_basis_ex2.in) using GetPot
109  std::string parameters_filename = "reduced_basis_ex2.in";
110  GetPot infile(parameters_filename);
111 
112  unsigned int n_elem = infile("n_elem", 1); // Determines the number of elements in the "truth" mesh
113  const unsigned int dim = 2; // The number of spatial dimensions
114 
115  bool store_basis_functions = infile("store_basis_functions", true); // Do we write the RB basis functions to disk?
116 
117  // Read the "online_mode" flag from the command line
118  const int online_mode = libMesh::command_line_next("-online_mode", 0);
119 
120  // Build a mesh on the default MPI communicator.
121  Mesh mesh (init.comm(), dim);
123  n_elem, n_elem,
124  0., 1.,
125  0., 1.,
126  QUAD4);
127 
128  // Create an equation systems object.
129  EquationSystems equation_systems (mesh);
130 
131  // We override RBConstruction with SimpleRBConstruction in order to
132  // specialize a few functions for this particular problem.
133  SimpleRBConstruction & rb_con =
134  equation_systems.add_system<SimpleRBConstruction> ("RBConvectionDiffusion");
135 
136  // Initialize the SCM Construction object
137  RBSCMConstruction & rb_scm_con =
138  equation_systems.add_system<RBSCMConstruction> ("RBSCMConvectionDiffusion");
139  rb_scm_con.set_RB_system_name("RBConvectionDiffusion");
140  rb_scm_con.add_variable("p", FIRST);
141 
142  // Initialize the data structures for the equation system.
143  equation_systems.init ();
144 
145  // Print out some information about the "truth" discretization
146  equation_systems.print_info();
147  mesh.print_info();
148 
149  // Set parameters for the eigenvalue problems that will be solved by rb_scm_con
150  equation_systems.parameters.set<unsigned int>("eigenpairs") = 1;
151  equation_systems.parameters.set<unsigned int>("basis vectors") = 3;
152  equation_systems.parameters.set<unsigned int>
153  ("linear solver maximum iterations") = 1000;
154 
155  // Build a new RBEvaluation object which will be used to perform
156  // Reduced Basis calculations. This is required in both the
157  // "Offline" and "Online" stages.
158  SimpleRBEvaluation rb_eval(mesh.comm());
159 
160  // We need to give the RBConstruction object a pointer to
161  // our RBEvaluation object
162  rb_con.set_rb_evaluation(rb_eval);
163 
164  // We also need a SCM evaluation object to perform SCM calculations
165  RBSCMEvaluation rb_scm_eval(mesh.comm());
166  rb_scm_eval.set_rb_theta_expansion(rb_eval.get_rb_theta_expansion());
167 
168  // Tell rb_eval about rb_scm_eval
169  rb_eval.rb_scm_eval = &rb_scm_eval;
170 
171  // Finally, need to give rb_scm_con and rb_eval a pointer to the
172  // SCM evaluation object, rb_scm_eval
173  rb_scm_con.set_rb_scm_evaluation(rb_scm_eval);
174 
175  if (!online_mode) // Perform the Offline stage of the RB method
176  {
177  // Read in the data that defines this problem from the specified text file
178  rb_con.process_parameters_file(parameters_filename);
179  rb_scm_con.process_parameters_file(parameters_filename);
180 
181  // Print out info that describes the current setup of rb_con
182  rb_con.print_info();
183  rb_scm_con.print_info();
184 
185  // Prepare rb_con for the Construction stage of the RB method.
186  // This sets up the necessary data structures and performs
187  // initial assembly of the "truth" affine expansion of the PDE.
189 
190  // Perform the SCM Greedy algorithm to derive the data required
191  // for rb_scm_eval to provide a coercivity lower bound.
192  rb_scm_con.perform_SCM_greedy();
193 
194  // Compute the reduced basis space by computing "snapshots", i.e.
195  // "truth" solves, at well-chosen parameter values and employing
196  // these snapshots as basis functions.
197  rb_con.train_reduced_basis();
198 
199  // Write out the data that will subsequently be required for the Evaluation stage
200 #if defined(LIBMESH_HAVE_CAPNPROTO)
202  rb_eval_writer.write_to_file("rb_eval.bin");
203  RBDataSerialization::RBSCMEvaluationSerialization rb_scm_eval_writer(rb_scm_con.get_rb_scm_evaluation());
204  rb_scm_eval_writer.write_to_file("rb_scm_eval.bin");
205 #else
207  rb_scm_con.get_rb_scm_evaluation().legacy_write_offline_data_to_files("scm_data");
208 #endif
209 
210  // If requested, write out the RB basis functions for visualization purposes
211  if (store_basis_functions)
212  {
213  // Write out the basis functions
214  rb_con.get_rb_evaluation().write_out_basis_functions(rb_con, "rb_data");
215  }
216  }
217  else // Perform the Online stage of the RB method
218  {
219 
220  // Read in the reduced basis data
221 #if defined(LIBMESH_HAVE_CAPNPROTO)
223  rb_eval_reader.read_from_file("rb_eval.bin", /*read_error_bound_data*/ true);
224  RBDataDeserialization::RBSCMEvaluationDeserialization rb_scm_eval_reader(rb_scm_eval);
225  rb_scm_eval_reader.read_from_file("rb_scm_eval.bin");
226 #else
227  rb_eval.legacy_read_offline_data_from_files("rb_data");
228  rb_scm_eval.legacy_read_offline_data_from_files("scm_data");
229 #endif
230 
231  // Read in online_N and initialize online parameters
232  unsigned int online_N = infile("online_N", 1);
233  Real online_mu_0 = infile("online_mu_0", 0.);
234  Real online_mu_1 = infile("online_mu_1", 0.);
235  Real online_mu_2 = infile("online_mu_2", 0.);
236  RBParameters online_mu;
237  online_mu.set_value("mu_0", online_mu_0);
238  online_mu.set_value("mu_1", online_mu_1);
239  online_mu.set_value("mu_2", online_mu_2);
240  rb_eval.set_parameters(online_mu);
241  rb_eval.print_parameters();
242 
243  // Now do the Online solve using the precomputed reduced basis
244  rb_eval.rb_solve(online_N);
245 
246  // Print out outputs as well as the corresponding output error bounds.
247  libMesh::out << "output 1, value = " << rb_eval.RB_outputs[0]
248  << ", bound = " << rb_eval.RB_output_error_bounds[0]
249  << std::endl;
250  libMesh::out << "output 2, value = " << rb_eval.RB_outputs[1]
251  << ", bound = " << rb_eval.RB_output_error_bounds[1]
252  << std::endl;
253  libMesh::out << "output 3, value = " << rb_eval.RB_outputs[2]
254  << ", bound = " << rb_eval.RB_output_error_bounds[2]
255  << std::endl;
256  libMesh::out << "output 4, value = " << rb_eval.RB_outputs[3]
257  << ", bound = " << rb_eval.RB_output_error_bounds[3]
258  << std::endl << std::endl;
259 
260  if (store_basis_functions)
261  {
262  // Read in the basis functions
263  rb_eval.read_in_basis_functions(rb_con, "rb_data");
264 
265  // Plot the solution
266  rb_con.load_rb_solution();
267 #ifdef LIBMESH_HAVE_EXODUS_API
268  ExodusII_IO(mesh).write_equation_systems ("RB_sol.e", equation_systems);
269 #endif
270 
271  // Plot the first basis function that was generated from the train_reduced_basis
272  // call in the Offline stage
273  rb_con.load_basis_function(0);
274 #ifdef LIBMESH_HAVE_EXODUS_API
275  ExodusII_IO(mesh).write_equation_systems ("bf0.e", equation_systems);
276 #endif
277  }
278  }
279 
280 #endif // LIBMESH_ENABLE_DIRICHLET
281 
282  return 0;
283 
284 #endif // LIBMESH_HAVE_SLEPC && LIBMESH_HAVE_GLPK
285 }
int main(int argc, char **argv)
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 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.
void set_rb_theta_expansion(RBThetaExpansion &rb_theta_expansion_in)
Set the RBThetaExpansion object.
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.
void read_from_file(const std::string &path, bool use_packing=false)
Read the Cap&#39;n&#39;Proto buffer from disk.
dof_id_type n_elem(const MeshBase::const_element_iterator &begin, const MeshBase::const_element_iterator &end)
Count up the number of elements of a specific type (as defined by an iterator range).
Definition: mesh_tools.C:969
unsigned int dim
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:52
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
MeshBase & mesh
const Parallel::Communicator & comm() const
void build_square(UnstructuredMesh &mesh, const unsigned int nx, const unsigned int ny, const Real xmin=0., const Real xmax=1., const Real ymin=0., const Real ymax=1., const ElemType type=INVALID_ELEM, const bool gauss_lobatto_grid=false)
A specialized build_cube() for 2D meshes.
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:90
The libMesh namespace provides an interface to certain functionality in the library.
virtual void print_info() const
Print out info that describes the current setup of this RBConstruction.
SolverPackage default_solver_package()
Definition: libmesh.C:1117
virtual void load_basis_function(unsigned int i)
Load the i^th RB function into the RBConstruction solution vector.
This class is part of the rbOOmit framework.
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 read_from_file(const std::string &path, bool read_error_bound_data, bool use_packing=false)
Read the Cap&#39;n&#39;Proto buffer from disk.
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
void write_to_file(const std::string &path, bool use_packing=false)
Write the Cap&#39;n&#39;Proto buffer to disk.
This class is part of the rbOOmit framework.
Definition: rb_parameters.h:52
This class de-serializes a RBSCMEvaluation object using the Cap&#39;n Proto library.
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.
void set_value(const std::string &param_name, Real value)
Set the value of the specified parameter.
T & set(const std::string &)
Definition: parameters.h:469
OStreamProxy out
Parameters parameters
Data structure holding arbitrary parameters.
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 init()
Initialize all the systems.
virtual void load_rb_solution()
Load the RB solution from the most recent solve with rb_eval into this system&#39;s solution vector...
This class serializes an RBSCMEvaluation object using the Cap&#39;n Proto library.
void set_rb_evaluation(RBEvaluation &rb_eval_in)
Set the RBEvaluation object.
void set_RB_system_name(const std::string &new_name)
Set the name of the associated RB system — we need this to load the (symmetrized) affine operators...
virtual System & add_system(std::string_view system_type, std::string_view name)
Add the system of type system_type named name to the systems array.
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
This class is part of the rbOOmit framework.