Go to the documentation of this file.
   33 #include "libmesh/equation_systems.h" 
   34 #include "libmesh/error_vector.h" 
   35 #include "libmesh/getpot.h" 
   36 #include "libmesh/exodusII_io.h" 
   37 #include "libmesh/kelly_error_estimator.h" 
   38 #include "libmesh/mesh.h" 
   39 #include "libmesh/mesh_generation.h" 
   40 #include "libmesh/enum_solver_package.h" 
   41 #include "libmesh/enum_solver_type.h" 
   42 #include "libmesh/auto_ptr.h"  
   46 #include "libmesh/diff_solver.h" 
   47 #include "libmesh/newmark_solver.h" 
   48 #include "libmesh/steady_solver.h" 
   49 #include "libmesh/euler_solver.h" 
   50 #include "libmesh/euler2_solver.h" 
   51 #include "libmesh/elem.h" 
   52 #include "libmesh/newton_solver.h" 
   53 #include "libmesh/eigen_sparse_linear_solver.h" 
   61 int main (
int argc, 
char ** argv)
 
   68                            "--enable-petsc, --enable-trilinos, or --enable-eigen");
 
   71   libmesh_example_requires(LIBMESH_DIM > 2, 
"3D support");
 
   74 #ifndef LIBMESH_ENABLE_DIRICHLET 
   75   libmesh_example_requires(
false, 
"--enable-dirichlet");
 
   79   GetPot infile(
"fem_system_ex3.in");
 
   82   infile.parse_command_line(argc, argv);
 
   85   const Real deltat           = infile(
"deltat", 0.25);
 
   86   unsigned int n_timesteps    = infile(
"n_timesteps", 1);
 
   88 #ifdef LIBMESH_HAVE_EXODUS_API 
   89   const unsigned int write_interval    = infile(
"write_interval", 1);
 
   93   const unsigned int dim = 3;
 
   96   libmesh_example_requires(
dim == LIBMESH_DIM, 
"3D support");
 
  119         side_max_x = 0, side_min_y = 0,
 
  120         side_max_y = 0, side_max_z = 0;
 
  122         found_side_max_x = 
false, found_side_max_y = 
false,
 
  123         found_side_min_y = 
false, found_side_max_z = 
false;
 
  124       for (
auto side : elem->side_index_range())
 
  129               found_side_max_x = 
true;
 
  135               found_side_min_y = 
true;
 
  141               found_side_max_y = 
true;
 
  147               found_side_max_z = 
true;
 
  154       if (found_side_max_x && found_side_max_y && found_side_max_z)
 
  155         for (
auto n : elem->node_index_range())
 
  156           if (elem->is_node_on_side(n, side_max_x) &&
 
  157               elem->is_node_on_side(n, side_max_y) &&
 
  158               elem->is_node_on_side(n, side_max_z))
 
  164       if (found_side_max_x && found_side_min_y)
 
  165         for (
auto e : elem->edge_index_range())
 
  166           if (elem->is_edge_on_side(e, side_max_x) &&
 
  167               elem->is_edge_on_side(e, side_min_y))
 
  179   std::string time_solver = infile(
"time_solver",
"DIE!");
 
  184   if( time_solver == std::string(
"newmark") )
 
  199   if (time_solver == std::string(
"newmark"))
 
  200     system.
time_solver = libmesh_make_unique<NewmarkSolver>(system);
 
  202   else if( time_solver == std::string(
"euler") )
 
  204       system.
time_solver = libmesh_make_unique<EulerSolver>(system);
 
  206       euler_solver.
theta = infile(
"theta", 1.0);
 
  209   else if( time_solver == std::string(
"euler2") )
 
  211       system.
time_solver = libmesh_make_unique<Euler2Solver>(system);
 
  213       euler_solver.
theta = infile(
"theta", 1.0);
 
  216   else if( time_solver == std::string(
"steady"))
 
  218       system.
time_solver = libmesh_make_unique<SteadySolver>(system);
 
  219       libmesh_assert_equal_to (n_timesteps, 1);
 
  222     libmesh_error_msg(std::string(
"ERROR: invalid time_solver ")+time_solver);
 
  225   equation_systems.
init ();
 
  232   solver.
quiet = infile(
"solver_quiet", 
true);
 
  249   NewtonSolver * newton_solver = dynamic_cast<NewtonSolver *>( &solver );
 
  251       (time_solver == std::string(
"euler") || time_solver == std::string(
"euler2") ) )
 
  253 #ifdef LIBMESH_HAVE_EIGEN_SPARSE 
  258       if( eigen_linear_solver )
 
  263   if( time_solver == std::string(
"newmark") )
 
  274 #ifdef LIBMESH_HAVE_EXODUS_API 
  277     std::ostringstream file_name;
 
  280     file_name << std::string(
"out.")+time_solver+std::string(
".e-s.")
 
  292 #endif // #ifdef LIBMESH_HAVE_EXODUS_API 
  296   for (
unsigned int t_step=0; t_step != n_timesteps; ++t_step)
 
  312       if( time_solver == std::string(
"newmark") )
 
  318 #ifdef LIBMESH_HAVE_EXODUS_API 
  320       if ((t_step+1)%write_interval == 0)
 
  322           std::ostringstream file_name;
 
  325           file_name << std::string(
"out.")+time_solver+std::string(
".e-s.")
 
  337 #endif // #ifdef LIBMESH_HAVE_EXODUS_API 
  
This class defines a theta-method Euler (defaulting to Backward Euler with theta = 1....
 
Real deltat
For time-dependent problems, this is the amount delta t to advance the solution in time.
 
virtual void solve() override
Invokes the solver associated with the system.
 
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 BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
 
void add_node(const Node *node, const boundary_id_type id)
Add Node node with boundary id id to the boundary information data structures.
 
The libMesh namespace provides an interface to certain functionality in the library.
 
Real theta
The value for the theta method to employ: 1.0 corresponds to backwards Euler, 0.0 corresponds to forw...
 
This is a generic class that defines a solver to handle ImplicitSystem classes, including NonlinearIm...
 
SolverPackage default_solver_package()
 
void write_timestep(const std::string &fname, const EquationSystems &es, const int timestep, const Real time, const std::set< std::string > *system_names=nullptr)
Writes out the solution at a specific timestep.
 
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
 
std::unique_ptr< TimeSolver > time_solver
A pointer to the solver object we're going to use.
 
This class defines a theta-method (defaulting to Backward Euler with theta = 1.0) solver to handle ti...
 
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
 
virtual SimpleRange< element_iterator > element_ptr_range()=0
 
This class defines a Newmark time integrator for second order (in time) DifferentiableSystems.
 
void set_solver_type(const SolverType st)
Sets the type of solver to use.
 
Real relative_step_tolerance
 
bool verbose
The DiffSolver may print a lot more to libMesh::out if verbose is set to true; default is false.
 
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.
 
int main(int argc, char **argv)
 
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
 
This class defines a solver which uses the default libMesh linear solver in a quasiNewton method to h...
 
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.
 
Real time
For time-dependent problems, this is the time t at the beginning of the current timestep.
 
bool quiet
The DiffSolver should not print anything to libMesh::out unless quiet is set to false; default is tru...
 
Manages consistently variables, degrees of freedom, and coefficient vectors for explicit systems.
 
virtual void compute_initial_accel()
This method uses the specified initial displacement and velocity to compute the initial acceleration ...
 
Real relative_residual_tolerance
 
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
 
Real absolute_residual_tolerance
The DiffSolver should exit after the residual is reduced to either less than absolute_residual_tolera...
 
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
 
void add_edge(const dof_id_type elem, const unsigned short int edge, const boundary_id_type id)
Add edge edge of element number elem with boundary id id to the boundary information data structure.
 
This class provides an interface to Eigen iterative solvers that is compatible with the libMesh Linea...
 
LinearSolver< Number > & get_linear_solver()
 
unsigned int max_linear_iterations
Each linear solver step should exit after max_linear_iterations is exceeded.
 
Real theta
The value for the theta method to employ: 1.0 corresponds to backwards Euler, 0.0 corresponds to forw...
 
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
 
unsigned int max_nonlinear_iterations
The DiffSolver should exit in failure if max_nonlinear_iterations is exceeded and continue_after_max_...
 
bool has_boundary_id(const Node *const node, const boundary_id_type id) const
 
const NumericVector< Number > & get_vector(const std::string &vec_name) const
 
double initial_linear_tolerance
Any required linear solves will at first be done with this tolerance; the DiffSolver may tighten the ...