Go to the documentation of this file.
23 #include "libmesh/newmark_system.h"
24 #include "libmesh/equation_systems.h"
25 #include "libmesh/sparse_matrix.h"
26 #include "libmesh/libmesh_logging.h"
27 #include "libmesh/numeric_vector.h"
46 const std::string & name_in,
47 const unsigned int number_in) :
49 _a_0 (1./(_default_alpha*_default_timestep*_default_timestep)),
50 _a_1 (_default_delta/(_default_alpha*_default_timestep)),
51 _a_2 (1./(_default_alpha*_default_timestep)),
52 _a_3 (1./(2.*_default_alpha)-1.),
53 _a_4 (_default_delta/_default_alpha -1.),
54 _a_5 (_default_timestep/2.*(_default_delta/_default_alpha-2.)),
55 _a_6 (_default_timestep*(1.-_default_delta)),
56 _a_7 (_default_delta*_default_timestep),
57 _finished_assemble (false)
130 libmesh_not_implemented();
159 LOG_SCOPE(
"initial_conditions ()",
"NewmarkSystem");
192 LOG_SCOPE(
"update_rhs ()",
"NewmarkSystem");
225 LOG_SCOPE(
"update_u_v_a ()",
"NewmarkSystem");
258 libmesh_assert_not_equal_to (delta_T, 0.);
273 _a_0 = 1./(alpha*delta_T*delta_T);
274 _a_1 = delta/(alpha*delta_T);
275 _a_2 = 1./(alpha*delta_T);
276 _a_3 = 1./(2.*alpha)-1.;
277 _a_4 = delta/alpha -1.;
278 _a_5 = delta_T/2.*(delta/alpha-2.);
279 _a_6 = delta_T*(1.-delta);
280 _a_7 = delta*delta_T;
virtual void close()=0
Calls the SparseMatrix's internal assembly routines, ensuring that the values are consistent across p...
void initial_conditions()
Apply initial conditions.
virtual void zero()=0
Set all entries to zero.
const EquationSystems & get_equation_systems() const
virtual void add(const numeric_index_type i, const T value)=0
Adds value to each entry of the vector.
static const Real _default_alpha
Default Newmark alpha.
NumericVector< Number > * rhs
The system matrix.
virtual void clear() override
Clear all the data structures associated with the system.
void update_u_v_a()
Update displacement, velocity and acceleration.
The libMesh namespace provides an interface to certain functionality in the library.
void set_newmark_parameters(const Real delta_T=_default_timestep, const Real alpha=_default_alpha, const Real delta=_default_delta)
Set the time step size and the newmark parameter alpha and delta and calculate the constant parameter...
SparseMatrix< Number > & add_matrix(const std::string &mat_name)
Adds the additional matrix mat_name to this system.
~NewmarkSystem()
Destructor.
static const Real _default_timestep
Default Newmark time step.
Real _a_0
Constants used for the time integration.
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices)
Computes , where v is a pointer and each dof_indices[i] specifies where to add value v[i].
void compute_matrix()
Compute the global matrix by adding up scaled mass damping and stiffness matrix.
NewmarkSystem(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
virtual void reinit() override
Reinitializes the member data fields associated with the system, so that, e.g., assemble() may be use...
SparseMatrix< Number > * matrix
The system matrix.
virtual void assemble() override
Prepares matrix and _dof_map for matrix assembly.
virtual void add(const numeric_index_type i, const numeric_index_type j, const T value)=0
Add value to the element (i,j).
This is the EquationSystems class.
T & set(const std::string &)
void update_rhs()
Update the rhs.
virtual void clear() override
Clear all the data structures associated with the system.
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
bool _finished_assemble
true if the matrix assembly is finished.
virtual void assemble() override
Assemble the linear system.
virtual void scale(const T factor)=0
Scale each element of the vector by the given factor.
NumericVector< Number > & add_vector(const std::string &vec_name, const bool projections=true, const ParallelType type=PARALLEL)
Adds the additional vector vec_name to this system.
virtual void user_initialization()
Calls user's attached initialization function, or is overridden by the user in derived classes.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
const SparseMatrix< Number > & get_matrix(const std::string &mat_name) const
const NumericVector< Number > & get_vector(const std::string &vec_name) const
virtual void zero()=0
Set all entries to 0.
static const Real _default_delta
Default Newmark delta.
Parameters parameters
Data structure holding arbitrary parameters.