Go to the documentation of this file.
23 #include "libmesh/equation_systems.h"
24 #include "libmesh/getpot.h"
25 #include "libmesh/libmesh.h"
26 #include "libmesh/libmesh_logging.h"
27 #include "libmesh/mesh.h"
28 #include "libmesh/mesh_generation.h"
29 #include "libmesh/numeric_vector.h"
30 #include "libmesh/string_to_enum.h"
31 #include "libmesh/time_solver.h"
32 #include "libmesh/transient_system.h"
33 #include "libmesh/vtk_io.h"
34 #include "libmesh/auto_ptr.h"
35 #include "libmesh/enum_solver_package.h"
59 ElemType eltype = Utility::string_to_enum<ElemType>(args(
"mesh/generation/element_type",
"hex8"));
60 int nx = args(
"mesh/generation/num_elem", 4, 0);
61 int ny = args(
"mesh/generation/num_elem", 4, 1);
62 int nz =
dim > 2 ? args(
"mesh/generation/num_elem", 4, 2) : 0;
63 double origx = args(
"mesh/generation/origin", -1.0, 0);
64 double origy = args(
"mesh/generation/origin", -1.0, 1);
65 double origz = args(
"mesh/generation/origin", 0.0, 2);
66 double sizex = args(
"mesh/generation/size", 2.0, 0);
67 double sizey = args(
"mesh/generation/size", 2.0, 1);
68 double sizez = args(
"mesh/generation/size", 2.0, 2);
70 origx, origx+sizex, origy, origy+sizey, origz, origz+sizez, eltype);
99 std::unique_ptr<VTKIO> io = libmesh_make_unique<VTKIO>(systems.
get_mesh());
101 Real duration = args(
"duration", 1.0);
103 for (
unsigned int t_step = 0; t_step < duration/solid_system.deltat; t_step++)
106 Real progress = t_step * solid_system.deltat / duration;
112 out <<
"===== Time Step " << std::setw(4) << t_step;
113 out <<
" (" << std::fixed << std::setprecision(2) << std::setw(6) << (progress*100.) <<
"%)";
114 out <<
", time = " << std::setw(7) << solid_system.time;
115 out <<
" =====" << std::endl;
118 aux_system.current_local_solution->close();
123 out <<
"Solving Solid" << std::endl;
124 solid_system.solve();
125 aux_system.current_local_solution->close();
126 *aux_system.solution = *aux_system.current_local_solution;
127 aux_system.solution->close();
130 out <<
"Doing a reinit of the equation systems" << std::endl;
133 if (t_step % args(
"output/frequency", 1) == 0)
135 std::stringstream file_name;
136 file_name << args(
"results_directory",
"./")
143 io->write_equation_systems(file_name.str(), systems);
146 out <<
"Advancing to next step" << std::endl;
147 solid_system.time_solver->advance_timestep();
153 int main(
int argc,
char ** argv)
160 "--enable-petsc, --enable-trilinos, or --enable-eigen");
163 #ifndef LIBMESH_HAVE_VTK
164 libmesh_example_requires(
false,
"--enable-vtk");
175 libMesh::out <<
"We skip fem_system_ex2 when using threads.\n"
181 GetPot args = GetPot(
"solid.in");
184 int dim = args(
"mesh/generation/dimension", 3);
185 libmesh_example_requires(
dim <= LIBMESH_DIM,
"3D support");
198 out <<
"Finished calculations" << std::endl;
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
Manages storage and variables for transient systems.
int main(int argc, char **argv)
The libMesh namespace provides an interface to certain functionality in the library.
const T_sys & get_system(const std::string &name) const
SolverPackage default_solver_package()
unsigned int mesh_dimension() const
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
void setup(EquationSystems &systems, Mesh &mesh, GetPot &args)
std::unique_ptr< NumericVector< Number > > older_local_solution
All the values I need to compute my contribution to the simulation at hand.
std::unique_ptr< NumericVector< Number > > current_local_solution
All the values I need to compute my contribution to the simulation at hand.
virtual void init()
Initialize all the systems.
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
This is the EquationSystems class.
T & set(const std::string &)
Manages consistently variables, degrees of freedom, and coefficient vectors for explicit systems.
virtual void reinit()
Handle any mesh changes and reinitialize all the systems on the updated mesh.
void run_timestepping(EquationSystems &systems, GetPot &args)
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::unique_ptr< NumericVector< Number > > old_local_solution
All the values I need to compute my contribution to the simulation at hand.
ElemType
Defines an enum for geometric element types.
Parameters parameters
Data structure holding arbitrary parameters.