Go to the documentation of this file.
   44 #include "libmesh/libmesh.h" 
   45 #include "libmesh/mesh.h" 
   46 #include "libmesh/mesh_generation.h" 
   47 #include "libmesh/exodusII_io.h" 
   48 #include "libmesh/equation_systems.h" 
   49 #include "libmesh/dof_map.h" 
   50 #include "libmesh/getpot.h" 
   51 #include "libmesh/elem.h" 
   52 #include "libmesh/quadrature_gauss.h" 
   53 #include "libmesh/libmesh_logging.h" 
   54 #include "libmesh/rb_data_serialization.h" 
   55 #include "libmesh/rb_data_deserialization.h" 
   56 #include "libmesh/enum_solver_package.h" 
   57 #include "libmesh/enum_solver_type.h" 
   60 #include "rb_classes.h" 
   64 #define BOUNDARY_ID_MIN_Z 0 
   65 #define BOUNDARY_ID_MIN_Y 1 
   66 #define BOUNDARY_ID_MAX_X 2 
   67 #define BOUNDARY_ID_MAX_Y 3 
   68 #define BOUNDARY_ID_MIN_X 4 
   69 #define BOUNDARY_ID_MAX_Z 5 
   70 #define NODE_BOUNDARY_ID 10 
   78                          const std::string & filename);
 
   84 int main(
int argc, 
char ** argv)
 
   91                            "--enable-petsc, --enable-trilinos, or --enable-eigen");
 
   93 #if !defined(LIBMESH_HAVE_XDR) 
   95   libmesh_example_requires(
false, 
"--enable-xdr");
 
   96 #elif defined(LIBMESH_DEFAULT_SINGLE_PRECISION) 
   98   libmesh_example_requires(
false, 
"--disable-singleprecision");
 
   99 #elif defined(LIBMESH_DEFAULT_TRIPLE_PRECISION) 
  101   libmesh_example_requires(
false, 
"double precision");
 
  113   const unsigned int dim = 3;
 
  114   libmesh_example_requires(
dim == LIBMESH_DIM, 
"3D support");
 
  116 #ifndef LIBMESH_ENABLE_DIRICHLET 
  117   libmesh_example_requires(
false, 
"--enable-dirichlet");
 
  120   std::string parameters_filename = 
"reduced_basis_ex5.in";
 
  121   GetPot infile(parameters_filename);
 
  123   unsigned int n_elem_x = infile(
"n_elem_x", 0);
 
  124   unsigned int n_elem_y = infile(
"n_elem_y", 0);
 
  125   unsigned int n_elem_z = infile(
"n_elem_z", 0);
 
  126   Real x_size           = infile(
"x_size", 0.);
 
  127   Real y_size           = infile(
"y_size", 0.);
 
  128   Real z_size           = infile(
"z_size", 0.);
 
  130   bool store_basis_functions = infile(
"store_basis_functions", 
true);
 
  133   GetPot command_line(argc, argv);
 
  135   if (command_line.search(1, 
"-online_mode"))
 
  136     online_mode = command_line.next(online_mode);
 
  155       unsigned int side_max_x = 0, side_max_y = 0, side_max_z = 0;
 
  156       bool found_side_max_x = 
false, found_side_max_y = 
false, found_side_max_z = 
false;
 
  157       for (
auto side : elem->side_index_range())
 
  162               found_side_max_x = 
true;
 
  168               found_side_max_y = 
true;
 
  174               found_side_max_z = 
true;
 
  181       if (found_side_max_x && found_side_max_y && found_side_max_z)
 
  183           for (
auto n : elem->node_index_range())
 
  185               if (elem->is_node_on_side(n, side_max_x) &&
 
  186                   elem->is_node_on_side(n, side_max_y) &&
 
  187                   elem->is_node_on_side(n, side_max_z))
 
  220   equation_systems.
init ();
 
  255 #if defined(LIBMESH_HAVE_CAPNPROTO) 
  263       if (store_basis_functions)
 
  272 #if defined(LIBMESH_HAVE_CAPNPROTO) 
  276       rb_eval.legacy_read_offline_data_from_files();
 
  280       Real online_x_scaling = infile(
"online_x_scaling", 0.);
 
  281       Real online_load_Fx   = infile(
"online_load_Fx",   0.);
 
  282       Real online_load_Fy   = infile(
"online_load_Fy",   0.);
 
  283       Real online_load_Fz   = infile(
"online_load_Fz",   0.);
 
  284       Real online_point_load_Fx   = infile(
"online_point_load_Fx",   0.);
 
  285       Real online_point_load_Fy   = infile(
"online_point_load_Fy",   0.);
 
  286       Real online_point_load_Fz   = infile(
"online_point_load_Fz",   0.);
 
  288       online_mu.
set_value(
"x_scaling", online_x_scaling);
 
  289       online_mu.
set_value(
"load_Fx",   online_load_Fx);
 
  290       online_mu.
set_value(
"load_Fy",   online_load_Fy);
 
  291       online_mu.
set_value(
"load_Fz",   online_load_Fz);
 
  292       online_mu.
set_value(
"point_load_Fx",   online_point_load_Fx);
 
  293       online_mu.
set_value(
"point_load_Fy",   online_point_load_Fy);
 
  294       online_mu.
set_value(
"point_load_Fz",   online_point_load_Fz);
 
  295       rb_eval.set_parameters(online_mu);
 
  296       rb_eval.print_parameters();
 
  299       rb_eval.rb_solve(rb_eval.get_n_basis_functions());
 
  301       if (store_basis_functions)
 
  304           rb_eval.read_in_basis_functions(rb_con);
 
  309           const RBParameters & rb_eval_params = rb_eval.get_parameters();
 
  314 #endif // LIBMESH_ENABLE_DIRICHLET 
  319 #ifdef LIBMESH_ENABLE_DIRICHLET 
  323                          const std::string & filename)
 
  334 #ifdef LIBMESH_HAVE_EXODUS_API 
  345   LOG_SCOPE(
"compute_stresses()", 
"main");
 
  353   unsigned int displacement_vars[3];
 
  363   fe->attach_quadrature_rule (&qrule);
 
  365   const std::vector<Real> & JxW = fe->get_JxW();
 
  366   const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
 
  371   unsigned int sigma_vars[3][3];
 
  381   unsigned int vonMises_var = stress_system.
variable_number (
"vonMises");
 
  384   std::vector<std::vector<dof_id_type>> dof_indices_var(system.
n_vars());
 
  385   std::vector<dof_id_type> stress_dof_indices_var;
 
  392       for (
unsigned int var=0; var<3; var++)
 
  393         dof_map.
dof_indices (elem, dof_indices_var[var], displacement_vars[var]);
 
  399       for (
unsigned int qp=0; qp<qrule.n_points(); qp++)
 
  400         for (
unsigned int C_i=0; C_i<3; C_i++)
 
  401           for (
unsigned int C_j=0; C_j<3; C_j++)
 
  402             for (
unsigned int C_k=0; C_k<3; C_k++)
 
  404                 const unsigned int n_var_dofs = dof_indices_var[C_k].size();
 
  408                 for (
unsigned int l=0; l<n_var_dofs; l++)
 
  411                 for (
unsigned int C_l=0; C_l<3; C_l++)
 
  412                   elem_sigma(C_i,C_j) += JxW[qp] * (
elasticity_tensor(C_i, C_j, C_k, C_l) * displacement_gradient(C_l));
 
  416       elem_sigma.
scale(1./elem->volume());
 
  419       for (
unsigned int i=0; i<3; i++)
 
  420         for (
unsigned int j=0; j<3; j++)
 
  422             stress_dof_map.dof_indices (elem, stress_dof_indices_var, sigma_vars[i][j]);
 
  428             if ((stress_system.
solution->first_local_index() <= dof_index) &&
 
  429                 (dof_index < stress_system.solution->last_local_index()))
 
  431                 stress_system.
solution->set(dof_index, elem_sigma(i,j));
 
  438                                              pow(elem_sigma(1,1) - elem_sigma(2,2),2.) +
 
  439                                              pow(elem_sigma(2,2) - elem_sigma(0,0),2.) +
 
  440                                              6.*(
pow(elem_sigma(0,1),2.) + 
pow(elem_sigma(1,2),2.) + 
pow(elem_sigma(2,0),2.))
 
  442       stress_dof_map.dof_indices (elem, stress_dof_indices_var, vonMises_var);
 
  444       if ((stress_system.
solution->first_local_index() <= dof_index) &&
 
  445           (dof_index < stress_system.solution->last_local_index()))
 
  447           stress_system.
solution->set(dof_index, vonMises_value);
 
  457 #endif // LIBMESH_ENABLE_DIRICHLET 
  
Real elasticity_tensor(unsigned int i, unsigned int j, unsigned int k, unsigned int l)
 
unsigned int n_vars() const
 
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
 
virtual System & add_system(const std::string &system_type, const std::string &name)
Add the system of type system_type named name to the systems array.
 
const MeshBase & get_mesh() const
 
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
 
This class implements specific orders of Gauss quadrature.
 
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Fills the vector di with the global degree of freedom indices for the element.
 
virtual void process_parameters_file(const std::string ¶meters_filename)
Read in from the file specified by parameters_filename and set the this system's member variables acc...
 
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 add_node(const Node *node, const boundary_id_type id)
Add Node node with boundary id id to the boundary information data structures.
 
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
 
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.
 
void set_rb_evaluation(RBEvaluation &rb_eval_in)
Set the RBEvaluation object.
 
The libMesh namespace provides an interface to certain functionality in the library.
 
const T_sys & get_system(const std::string &name) const
 
MetaPhysicL::DualNumber< T, D > sqrt(const MetaPhysicL::DualNumber< T, D > &in)
 
const Parallel::Communicator & comm() const
 
This class is part of the rbOOmit framework.
 
This class serializes an RBEvaluation object using the Cap'n Proto library.
 
SolverPackage default_solver_package()
 
unsigned int mesh_dimension() const
 
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
 
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
 
Number current_solution(const dof_id_type global_dof_number) const
 
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
 
virtual SimpleRange< element_iterator > element_ptr_range()=0
 
This is the MeshBase class.
 
void set_solver_type(const SolverType st)
Sets the type of solver to use.
 
void write_to_file(const std::string &path)
Write the Cap'n'Proto buffer to disk.
 
virtual LinearSolver< Number > * get_linear_solver() const override
 
virtual SimpleRange< node_iterator > node_ptr_range()=0
 
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
 
unsigned int add_variable(const std::string &var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variable var to the list of variables for this system.
 
virtual void init()
Initialize all the systems.
 
void add_scaled(const TypeVector< T2 > &, const T &)
Add a scaled value to this vector without creating a temporary.
 
int main(int argc, char **argv)
 
RBEvaluation & get_rb_evaluation()
Get a reference to the RBEvaluation object.
 
void set_value(const std::string ¶m_name, Real value)
Set the value of the specified parameter.
 
double pow(double a, int b)
 
const FEType & variable_type(const unsigned int c) const
 
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
 
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
 
This is the EquationSystems class.
 
virtual void write_equation_systems(const std::string &, const EquationSystems &, const std::set< std::string > *system_names=nullptr)
This method implements writing a mesh with data to a specified file where the data is taken from the ...
 
virtual void print_info()
Print out info that describes the current setup of this RBConstruction.
 
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.
 
Manages consistently variables, degrees of freedom, and coefficient vectors for explicit systems.
 
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
 
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
 
This class handles the numbering of degrees of freedom on a mesh.
 
void compute_stresses(EquationSystems &es)
 
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
 
virtual void load_rb_solution()
Load the RB solution from the most recent solve with rb_eval into this system's solution vector.
 
This class de-serializes an RBEvaluation object using the Cap'n Proto library.
 
void scale_mesh_and_plot(EquationSystems &es, const RBParameters &mu, const std::string &filename)
 
unsigned short int variable_number(const std::string &var) const
 
void read_from_file(const std::string &path, bool read_error_bound_data)
Write the Cap'n'Proto buffer to disk.
 
const DofMap & get_dof_map() const
 
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
 
Order default_quadrature_order() const
 
bool has_boundary_id(const Node *const node, const boundary_id_type id) const
 
virtual void update()
Update the local values to reflect the solution on neighboring processors.
 
void scale(const T factor)
Multiplies every element in the matrix by factor.
 
virtual Real train_reduced_basis(const bool resize_rb_eval_data=true)
Train the reduced basis.
 
Real get_value(const std::string ¶m_name) const
Get the value of the specific parameter.