G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsDynamicNewmark.h
Go to the documentation of this file.
1 
17 #pragma once
18 #include <gsCore/gsLinearAlgebra.h>
19 
21 #include <gsIO/gsOptionList.h>
22 
23 namespace gismo
24 {
25 
33 template <class T, bool _NL>
35 {
36  typedef gsDynamicBase<T> Base;
37 
38 protected:
39 
40  typedef typename gsStructuralAnalysisOps<T>::Force_t Force_t;
41  typedef typename gsStructuralAnalysisOps<T>::TForce_t TForce_t;
42  typedef typename gsStructuralAnalysisOps<T>::Residual_t Residual_t;
43  typedef typename gsStructuralAnalysisOps<T>::TResidual_t TResidual_t;
44  typedef typename gsStructuralAnalysisOps<T>::Mass_t Mass_t;
45  typedef typename gsStructuralAnalysisOps<T>::TMass_t TMass_t;
46  typedef typename gsStructuralAnalysisOps<T>::Damping_t Damping_t;
47  typedef typename gsStructuralAnalysisOps<T>::TDamping_t TDamping_t;
48  typedef typename gsStructuralAnalysisOps<T>::Stiffness_t Stiffness_t;
49  typedef typename gsStructuralAnalysisOps<T>::Jacobian_t Jacobian_t;
50  typedef typename gsStructuralAnalysisOps<T>::TJacobian_t TJacobian_t;
51 
52 public:
53 
54  virtual ~gsDynamicNewmark() {};
55 
58  const Mass_t & Mass,
59  const Damping_t & Damping,
60  const Stiffness_t & Stiffness,
61  const Force_t & Force
62  )
63  :
64  Base(Mass,Damping,Stiffness,Force)
65  {
66  this->defaultOptions();
67  }
68 
71  const Mass_t & Mass,
72  const Damping_t & Damping,
73  const Stiffness_t & Stiffness,
74  const TForce_t & TForce
75  )
76  :
77  Base(Mass,Damping,Stiffness,TForce)
78  {
79  this->defaultOptions();
80  }
81 
84  const Mass_t & Mass,
85  const Damping_t & Damping,
86  const Jacobian_t & Jacobian,
87  const Residual_t & Residual
88  )
89  :
90  Base(Mass,Damping,Jacobian,Residual)
91  {
92  this->defaultOptions();
93  }
94 
97  const Mass_t & Mass,
98  const Damping_t & Damping,
99  const Jacobian_t & Jacobian,
100  const TResidual_t & TResidual
101  )
102  :
103  Base(Mass,Damping,Jacobian,TResidual)
104  {
105  this->defaultOptions();
106  }
107 
110  const Mass_t & Mass,
111  const Damping_t & Damping,
112  const TJacobian_t & TJacobian,
113  const TResidual_t & TResidual
114  )
115  :
116  Base(Mass,Damping,TJacobian,TResidual)
117  {
118  this->defaultOptions();
119  }
120 
123  const TMass_t & TMass,
124  const TDamping_t & TDamping,
125  const TJacobian_t & TJacobian,
126  const TResidual_t & TResidual
127  )
128  :
129  Base(TMass,TDamping,TJacobian,TResidual)
130  {
131  this->defaultOptions();
132  }
133 
135  virtual void defaultOptions() override;
136 
137 // General functions
138 protected:
139 
140  gsStatus _step(const T t, const T dt, gsVector<T> & U, gsVector<T> & V, gsVector<T> & A) const override;
141 
142  void _initOutput() const;
143  void _stepOutput(const index_t it, const T resnorm, const T updatenorm) const;
144 
145  gsSparseMatrix<T> m_sysmat;
146 
147  using Base::_computeForce;
149  using Base::_computeMass;
151  using Base::_computeDamping;
153 
154 protected:
155 
156  using Base::m_solver;
157 
158  using Base::m_numIterations;
159 
160  using Base::m_options;
161 
162 
163 private:
164  template <bool _nonlinear>
165  typename std::enable_if<(_nonlinear==false), gsStatus>::type
166  _step_impl(const T t, const T dt, gsVector<T> & U, gsVector<T> & V, gsVector<T> & A) const;
167 
168  template <bool _nonlinear>
169  typename std::enable_if<(_nonlinear==true), gsStatus>::type
170  _step_impl(const T t, const T dt, gsVector<T> & U, gsVector<T> & V, gsVector<T> & A) const;
171 };
172 
173 } // namespace gismo
174 
175 #ifndef GISMO_BUILD_LIB
176 #include GISMO_HPP_HEADER(gsDynamicNewmark.hpp)
177 #endif
std::function< bool(gsVector< T > const &, const T, gsSparseMatrix< T > &) > TJacobian_t
Jacobian.
Definition: gsStructuralAnalysisTypes.h:85
virtual void _computeForce(const T time, gsVector< T > &F) const
Compute the residual.
Definition: gsDynamicBase.h:261
virtual void _computeMassInverse(const gsSparseMatrix< T > &M, gsSparseMatrix< T > &Minv) const
Compute the mass matrix.
Definition: gsDynamicBase.h:282
gsDynamicNewmark(const Mass_t &Mass, const Damping_t &Damping, const Jacobian_t &Jacobian, const Residual_t &Residual)
Constructor.
Definition: gsDynamicNewmark.h:83
std::function< bool(gsVector< T > const &, const T, gsVector< T > &)> TResidual_t
Time-dependent Residual Fint(t)-Fext(t)
Definition: gsStructuralAnalysisTypes.h:69
Base class to perform time integration of second-order structural dynamics systems.
gsStatus _step(const T t, const T dt, gsVector< T > &U, gsVector< T > &V, gsVector< T > &A) const override
Initialize the ALM.
Definition: gsDynamicNewmark.hpp:153
Performs the arc length method to solve a nonlinear system of equations.
Definition: gsDynamicBase.h:34
virtual void _computeMass(const T time, gsSparseMatrix< T > &M) const
Compute the mass matrix.
Definition: gsDynamicBase.h:275
gsDynamicNewmark(const Mass_t &Mass, const Damping_t &Damping, const TJacobian_t &TJacobian, const TResidual_t &TResidual)
Constructor.
Definition: gsDynamicNewmark.h:109
#define index_t
Definition: gsConfig.h:32
gsStatus
Definition: gsStructuralAnalysisTypes.h:20
std::function< bool(gsVector< T > const &, gsVector< T > &)> Residual_t
Residual, Fint-Fext.
Definition: gsStructuralAnalysisTypes.h:65
std::function< bool(gsSparseMatrix< T > &) > Mass_t
Mass matrix.
Definition: gsStructuralAnalysisTypes.h:72
std::function< bool(gsVector< T > const &, gsSparseMatrix< T > &) > Jacobian_t
Jacobian.
Definition: gsStructuralAnalysisTypes.h:83
std::function< bool(gsVector< T > const &, gsSparseMatrix< T > &) > Damping_t
Damping matrix.
Definition: gsStructuralAnalysisTypes.h:76
Provides a list of labeled parameters/options that can be set and accessed easily.
virtual void _computeDamping(const gsVector< T > &U, const T time, gsSparseMatrix< T > &C) const
Compute the damping matrix.
Definition: gsDynamicBase.h:297
gsDynamicNewmark(const TMass_t &TMass, const TDamping_t &TDamping, const TJacobian_t &TJacobian, const TResidual_t &TResidual)
Constructor.
Definition: gsDynamicNewmark.h:122
std::function< bool(gsVector< T > &)> Force_t
Force.
Definition: gsStructuralAnalysisTypes.h:60
std::function< bool(const T, gsSparseMatrix< T > &) > TMass_t
Time-dependent mass matrix.
Definition: gsStructuralAnalysisTypes.h:74
std::function< bool(gsSparseMatrix< T > &) > Stiffness_t
Stiffness matrix.
Definition: gsStructuralAnalysisTypes.h:81
virtual void defaultOptions() override
Set default options.
Definition: gsDynamicNewmark.hpp:20
gsDynamicNewmark(const Mass_t &Mass, const Damping_t &Damping, const Jacobian_t &Jacobian, const TResidual_t &TResidual)
Constructor.
Definition: gsDynamicNewmark.h:96
std::function< bool(const T, gsVector< T > &)> TForce_t
Time-dependent force.
Definition: gsStructuralAnalysisTypes.h:62
virtual void _computeJacobian(const gsVector< T > &U, const T time, gsSparseMatrix< T > &K) const
Compute the Jacobian matrix.
Definition: gsDynamicBase.h:304
index_t m_numIterations
Number of iterations performed.
Definition: gsDynamicBase.h:348
This is the main header file that collects wrappers of Eigen for linear algebra.
gsDynamicNewmark(const Mass_t &Mass, const Damping_t &Damping, const Stiffness_t &Stiffness, const Force_t &Force)
Constructor.
Definition: gsDynamicNewmark.h:57
Performs the arc length method to solve a nonlinear system of equations.
Definition: gsDynamicNewmark.h:34
std::function< bool(gsVector< T > const &, const T, gsSparseMatrix< T > &) > TDamping_t
Time-dependent Damping matrix.
Definition: gsStructuralAnalysisTypes.h:78
gsDynamicNewmark(const Mass_t &Mass, const Damping_t &Damping, const Stiffness_t &Stiffness, const TForce_t &TForce)
Constructor.
Definition: gsDynamicNewmark.h:70
virtual void _computeResidual(const gsVector< T > &U, const T time, gsVector< T > &R) const
Compute the residual.
Definition: gsDynamicBase.h:268