newmark_system.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #ifndef LIBMESH_NEWMARK_SYSTEM_H
21 #define LIBMESH_NEWMARK_SYSTEM_H
22 
23 // Local Includes
25 
26 // C++ includes
27 
28 namespace libMesh
29 {
30 
52 {
53 public:
54 
60  const std::string & name,
61  const unsigned int number);
62 
66  ~NewmarkSystem ();
67 
72 
77  virtual void clear () override;
78 
83  virtual void reinit () override;
84 
89  virtual void assemble () override;
90 
95  virtual std::string system_type () const override { return "Newmark"; }
96 
97 
98  //---------------------------------------------------------
99  // These members are specific to the Newmark system
100  //
101 
105  void initial_conditions ();
106 
111  void compute_matrix ();
112 
116  void update_rhs ();
117 
121  void update_u_v_a ();
122 
128  void set_newmark_parameters (const Real delta_T = _default_timestep,
129  const Real alpha = _default_alpha,
130  const Real delta = _default_delta);
131 
132 private:
133 
145 
150 
154  static const Real _default_alpha;
155 
159  static const Real _default_delta;
160 
164  static const Real _default_timestep;
165 };
166 
167 } // namespace libMesh
168 
169 #endif // LIBMESH_NEWMARK_SYSTEM_H
Manages multiples systems of equations.
void set_newmark_parameters(const Real delta_T=_default_timestep, const Real alpha=_default_alpha, const Real delta=_default_delta)
virtual std::string system_type() const override
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and linear solvers ...
NewmarkSystem sys_type
virtual void assemble() override
virtual void clear() override
unsigned int number() const
Definition: system.h:2025
virtual void reinit() override
static const Real _default_delta
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static const Real _default_timestep
NewmarkSystem(EquationSystems &es, const std::string &name, const unsigned int number)
static const Real _default_alpha
Implements the Newmark time integration scheme.
const std::string & name() const
Definition: system.h:2017