34 #include "libmesh/equation_systems.h" 35 #include "libmesh/error_vector.h" 36 #include "libmesh/getpot.h" 37 #include "libmesh/exodusII_io.h" 38 #include "libmesh/kelly_error_estimator.h" 39 #include "libmesh/mesh.h" 40 #include "libmesh/mesh_generation.h" 41 #include "libmesh/enum_solver_package.h" 42 #include "libmesh/enum_solver_type.h" 43 #include "libmesh/numeric_vector.h" 46 #include "elasticity_system.h" 47 #include "libmesh/diff_solver.h" 48 #include "libmesh/newmark_solver.h" 49 #include "libmesh/steady_solver.h" 50 #include "libmesh/euler_solver.h" 51 #include "libmesh/euler2_solver.h" 52 #include "libmesh/elem.h" 53 #include "libmesh/newton_solver.h" 54 #include "libmesh/eigen_sparse_linear_solver.h" 62 int main (
int argc,
char ** argv)
69 "--enable-petsc, --enable-trilinos, or --enable-eigen");
72 libmesh_example_requires(LIBMESH_DIM > 2,
"3D support");
75 #ifndef LIBMESH_ENABLE_DIRICHLET 76 libmesh_example_requires(
false,
"--enable-dirichlet");
80 GetPot infile(
"fem_system_ex3.in");
83 infile.parse_command_line(argc, argv);
86 const Real deltat = infile(
"deltat", 0.25);
87 unsigned int n_timesteps = infile(
"n_timesteps", 1);
89 #ifdef LIBMESH_HAVE_EXODUS_API 90 const unsigned int write_interval = infile(
"write_interval", 1);
94 const unsigned int dim = 3;
97 libmesh_example_requires(
dim == LIBMESH_DIM,
"3D support");
99 const unsigned int nx = infile(
"nx", 32);
100 const unsigned int ny = infile(
"ny", 8);
101 const unsigned int nz = infile(
"nz", 4);
121 for (
const auto & elem :
mesh.element_ptr_range())
124 side_max_x = 0, side_min_y = 0,
125 side_max_y = 0, side_max_z = 0;
127 found_side_max_x =
false, found_side_max_y =
false,
128 found_side_min_y =
false, found_side_max_z =
false;
129 for (
auto side : elem->side_index_range())
134 found_side_max_x =
true;
140 found_side_min_y =
true;
146 found_side_max_y =
true;
152 found_side_max_z =
true;
159 if (found_side_max_x && found_side_max_y && found_side_max_z)
160 for (
auto n : elem->node_index_range())
161 if (elem->is_node_on_side(n, side_max_x) &&
162 elem->is_node_on_side(n, side_max_y) &&
163 elem->is_node_on_side(n, side_max_z))
169 if (found_side_max_x && found_side_min_y)
170 for (
auto e : elem->edge_index_range())
171 if (elem->is_edge_on_side(e, side_max_x) &&
172 elem->is_edge_on_side(e, side_min_y))
184 const std::string time_solver = infile(
"time_solver",
"DIE!");
189 if (time_solver ==
"newmark")
203 system.
time_solver = std::make_unique<NewmarkSolver>(system);
206 else if (time_solver ==
"euler")
208 system.
time_solver = std::make_unique<EulerSolver>(system);
210 euler_solver.theta = infile(
"theta", 1.0);
213 else if (time_solver ==
"euler2")
215 system.
time_solver = std::make_unique<Euler2Solver>(system);
217 euler_solver.theta = infile(
"theta", 1.0);
220 else if (time_solver ==
"steady")
222 system.
time_solver = std::make_unique<SteadySolver>(system);
223 libmesh_assert_equal_to (n_timesteps, 1);
226 libmesh_error_msg(std::string(
"ERROR: invalid time_solver ")+time_solver);
229 equation_systems.
init ();
236 solver.
quiet = infile(
"solver_quiet",
true);
255 (time_solver ==
"euler" || time_solver ==
"euler2"))
257 #ifdef LIBMESH_HAVE_EIGEN_SPARSE 262 if (eigen_linear_solver)
267 if (time_solver ==
"newmark")
278 #ifdef LIBMESH_HAVE_EXODUS_API 281 std::ostringstream file_name;
284 file_name << std::string(
"out.")+time_solver+std::string(
".e-s.")
296 #endif // #ifdef LIBMESH_HAVE_EXODUS_API 300 for (
unsigned int t_step=0; t_step != n_timesteps; ++t_step)
316 if (time_solver ==
"newmark")
322 #ifdef LIBMESH_HAVE_EXODUS_API 324 if ((t_step+1)%write_interval == 0)
326 std::ostringstream file_name;
329 file_name << std::string(
"out.")+time_solver+std::string(
".e-s.")
341 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
Real time
For time-dependent problems, this is the time t at the beginning of the current timestep.
LinearSolver< Number > & get_linear_solver()
This is the EquationSystems class.
const boundary_id_type edge_boundary_id
bool has_boundary_id(const Node *const node, const boundary_id_type id) const
int main(int argc, char **argv)
bool quiet
The DiffSolver should not print anything to libMesh::out unless quiet is set to false; default is tru...
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
unsigned int max_nonlinear_iterations
The DiffSolver should exit in failure if max_nonlinear_iterations is exceeded and continue_after_max_...
Real absolute_residual_tolerance
The DiffSolver should exit after the residual is reduced to either less than absolute_residual_tolera...
std::unique_ptr< TimeSolver > time_solver
A pointer to the solver object we're going to use.
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
This class provides an interface to Eigen iterative solvers that is compatible with the libMesh Linea...
This class defines a Newmark time integrator for second order (in time) DifferentiableSystems.
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
The libMesh namespace provides an interface to certain functionality in the library.
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
boundary_id_type boundary_id_max_x
This class defines a solver which uses the default libMesh linear solver in a quasiNewton method to h...
virtual void compute_initial_accel()
This method uses the specified initial displacement and velocity to compute the initial acceleration ...
SolverPackage default_solver_package()
This is a generic class that defines a solver to handle ImplicitSystem classes, including NonlinearIm...
void add_node(const Node *node, const boundary_id_type id)
Add Node node with boundary id id to the boundary information data structures.
This class defines a theta-method (defaulting to Backward Euler with theta = 1.0) solver to handle ti...
boundary_id_type boundary_id_max_z
const boundary_id_type node_boundary_id
Real deltat
For time-dependent problems, this is the amount delta t to advance the solution in time...
unsigned int max_linear_iterations
Each linear solver step should exit after max_linear_iterations is exceeded.
void print_info(std::ostream &os=libMesh::out, const unsigned int verbosity=0, const bool global=true) const
Prints relevant information about the mesh.
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
boundary_id_type boundary_id_max_y
unsigned int add_variable(std::string_view 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.
void set_solver_type(const SolverType st)
Sets the type of solver to use.
This class defines a theta-method Euler (defaulting to Backward Euler with theta = 1...
bool verbose
The DiffSolver may print a lot more to libMesh::out if verbose is set to true; default is false...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
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.
virtual void solve() override
Invokes the solver associated with the system.
virtual void init()
Initialize all the systems.
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.
boundary_id_type boundary_id_min_y
double initial_linear_tolerance
Any required linear solves will at first be done with this tolerance; the DiffSolver may tighten the ...
Manages consistently variables, degrees of freedom, and coefficient vectors for explicit systems...
Real relative_residual_tolerance
Real relative_step_tolerance
const NumericVector< Number > & get_vector(std::string_view vec_name) const
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...