G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsDynamicWilson.h
Go to the documentation of this file.
1 
17 #pragma once
18 #include <gsCore/gsLinearAlgebra.h>
19 
22 #include <gsIO/gsOptionList.h>
23 
24 namespace gismo
25 {
26 
34 template <class T, bool _NL>
35 class gsDynamicWilson : public gsDynamicNewmark<T,_NL>
36 {
38 
39 protected:
40 
41  typedef typename gsStructuralAnalysisOps<T>::Force_t Force_t;
42  typedef typename gsStructuralAnalysisOps<T>::TForce_t TForce_t;
43  typedef typename gsStructuralAnalysisOps<T>::Residual_t Residual_t;
44  typedef typename gsStructuralAnalysisOps<T>::TResidual_t TResidual_t;
45  typedef typename gsStructuralAnalysisOps<T>::Mass_t Mass_t;
46  typedef typename gsStructuralAnalysisOps<T>::TMass_t TMass_t;
47  typedef typename gsStructuralAnalysisOps<T>::Damping_t Damping_t;
48  typedef typename gsStructuralAnalysisOps<T>::TDamping_t TDamping_t;
49  typedef typename gsStructuralAnalysisOps<T>::Stiffness_t Stiffness_t;
50  typedef typename gsStructuralAnalysisOps<T>::Jacobian_t Jacobian_t;
51  typedef typename gsStructuralAnalysisOps<T>::TJacobian_t TJacobian_t;
52 
53 public:
54 
55  virtual ~gsDynamicWilson() {};
56 
59  const Mass_t & Mass,
60  const Damping_t & Damping,
61  const Stiffness_t & Stiffness,
62  const Force_t & Force
63  )
64  :
65  Base(Mass,Damping,Stiffness,Force)
66  {
67  this->defaultOptions();
68  }
69 
72  const Mass_t & Mass,
73  const Damping_t & Damping,
74  const Stiffness_t & Stiffness,
75  const TForce_t & TForce
76  )
77  :
78  Base(Mass,Damping,Stiffness,TForce)
79  {
80  this->defaultOptions();
81  }
82 
85  const Mass_t & Mass,
86  const Damping_t & Damping,
87  const Jacobian_t & Jacobian,
88  const Residual_t & Residual
89  )
90  :
91  Base(Mass,Damping,Jacobian,Residual)
92  {
93  this->defaultOptions();
94  }
95 
98  const Mass_t & Mass,
99  const Damping_t & Damping,
100  const Jacobian_t & Jacobian,
101  const TResidual_t & TResidual
102  )
103  :
104  Base(Mass,Damping,Jacobian,TResidual)
105  {
106  this->defaultOptions();
107  }
108 
111  const Mass_t & Mass,
112  const Damping_t & Damping,
113  const TJacobian_t & TJacobian,
114  const TResidual_t & TResidual
115  )
116  :
117  Base(Mass,Damping,TJacobian,TResidual)
118  {
119  this->defaultOptions();
120  }
121 
124  const TMass_t & TMass,
125  const TDamping_t & TDamping,
126  const TJacobian_t & TJacobian,
127  const TResidual_t & TResidual
128  )
129  :
130  Base(TMass,TDamping,TJacobian,TResidual)
131  {
132  this->defaultOptions();
133  }
134 
136  virtual void defaultOptions() override;
137 
138 // General functions
139 protected:
140 
141  gsStatus _step(const T t, const T dt, gsVector<T> & U, gsVector<T> & V, gsVector<T> & A) const override;
142 
143  void _initOutput() const;
144  void _stageOutput(index_t stage) const;
145  void _stepOutput(const index_t it, const T resnorm, const T updatenorm) const;
146 
147  gsSparseMatrix<T> m_sysmat;
148 
149  using Base::_computeForce;
151  using Base::_computeMass;
153  using Base::_computeDamping;
155 
156 protected:
157 
158  using Base::m_solver;
159 
160  using Base::m_numIterations;
161 
162  using Base::m_options;
163 
164 
165 private:
166  template <bool _nonlinear>
167  typename std::enable_if<(_nonlinear==false), gsStatus>::type
168  _step_impl(const T t, const T dt, gsVector<T> & U, gsVector<T> & V, gsVector<T> & A) const;
169 
170  template <bool _nonlinear>
171  typename std::enable_if<(_nonlinear==true), gsStatus>::type
172  _step_impl(const T t, const T dt, gsVector<T> & U, gsVector<T> & V, gsVector<T> & A) const;
173 };
174 
175 } // namespace gismo
176 
177 #ifndef GISMO_BUILD_LIB
178 #include GISMO_HPP_HEADER(gsDynamicWilson.hpp)
179 #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
Class to perform time integration of second-order structural dynamics systems using the Explicit Eule...
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.
virtual void _computeMass(const T time, gsSparseMatrix< T > &M) const
Compute the mass matrix.
Definition: gsDynamicBase.h:275
gsDynamicWilson(const Mass_t &Mass, const Damping_t &Damping, const TJacobian_t &TJacobian, const TResidual_t &TResidual)
Constructor.
Definition: gsDynamicWilson.h:110
#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
gsDynamicWilson(const Mass_t &Mass, const Damping_t &Damping, const Jacobian_t &Jacobian, const TResidual_t &TResidual)
Constructor.
Definition: gsDynamicWilson.h:97
gsStatus _step(const T t, const T dt, gsVector< T > &U, gsVector< T > &V, gsVector< T > &A) const override
Initialize the ALM.
Definition: gsDynamicWilson.hpp:87
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
gsDynamicWilson(const TMass_t &TMass, const TDamping_t &TDamping, const TJacobian_t &TJacobian, const TResidual_t &TResidual)
Constructor.
Definition: gsDynamicWilson.h:123
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
gsDynamicWilson(const Mass_t &Mass, const Damping_t &Damping, const Stiffness_t &Stiffness, const TForce_t &TForce)
Constructor.
Definition: gsDynamicWilson.h:71
gsDynamicWilson(const Mass_t &Mass, const Damping_t &Damping, const Jacobian_t &Jacobian, const Residual_t &Residual)
Constructor.
Definition: gsDynamicWilson.h:84
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: gsDynamicWilson.hpp:20
gsDynamicWilson(const Mass_t &Mass, const Damping_t &Damping, const Stiffness_t &Stiffness, const Force_t &Force)
Constructor.
Definition: gsDynamicWilson.h:58
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.
Performs the arc length method to solve a nonlinear system of equations.
Definition: gsDynamicWilson.h:35
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
virtual void _computeResidual(const gsVector< T > &U, const T time, gsVector< T > &R) const
Compute the residual.
Definition: gsDynamicBase.h:268