G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsDynamicNewmark.hpp
Go to the documentation of this file.
1
14#pragma once
15
16namespace gismo
17{
18
19template <class T, bool _NL>
21{
22 Base::defaultOptions();
23 m_options.addReal("alpha","Beta parameter for Newmark's method, such that 0 =< 2 beta =< 1",0.25);
24 m_options.addReal("delta","Beta parameter for Newmark's method, such that 0 =< gamma =< 1",0.50);
25}
26
27
28template <class T, bool _NL>
29template <bool _nonlinear>
30typename std::enable_if<(_nonlinear==false), gsStatus>::type
31gsDynamicNewmark<T,_NL>::_step_impl(const T t, const T dt, gsVector<T> & U, gsVector<T> & V, gsVector<T> & A) const
32{
33 gsVector<T> Uold = U;
34 gsVector<T> Vold = V;
35 gsVector<T> Aold = A;
36
38 gsSparseMatrix<T> M, C, K;
39
40 T alpha = m_options.getReal("alpha");
41 T delta = m_options.getReal("delta");
42
43 // Computed at t=t0+dt
44 this->_computeMass(t+dt,M);
45 this->_computeForce(t+dt,F);
46
48 gsMatrix<T> rhs;
49 // predictors
50 U = Uold + dt*Vold + Aold*(0.5 - alpha)*dt*dt;
51 V = Vold + Aold*(1 - delta)*dt;
52
53 this->_computeDamping(U,t+dt,C);
54 this->_computeJacobian(U,t+dt,K);
55
56 // lhs and rhs
57 lhs = M + delta*dt*C + dt*dt*alpha*K;
58 rhs = F - K*U - C * V;
59
60 this->_initOutput();
61 typename gsSparseSolver<T>::uPtr solver;
62 solver = gsSparseSolver<T>::get( m_options.getString("Solver") );
63 solver->compute(lhs);
64
65 A = solver->solve(rhs);
66 V += A*delta*dt;
67 U += A*alpha*dt*dt;
68
69 this->_stepOutput(0,(F - K*U - M * A - C * V).norm(),U.norm());
70 if (math::isinf(U.norm()) || math::isnan(U.norm()))
72 else
73 return gsStatus::Success;
74}
75
76
77template <class T, bool _NL>
78template <bool _nonlinear>
79typename std::enable_if<(_nonlinear==true), gsStatus>::type
80gsDynamicNewmark<T,_NL>::_step_impl(const T t, const T dt, gsVector<T> & U, gsVector<T> & V, gsVector<T> & A) const
81{
82 //https://ethz.ch/content/dam/ethz/special-interest/baug/ibk/structural-mechanics-dam/education/femII/presentation_05_dynamics_v3.pdf
83 gsVector<T> Uold = U;
84 gsVector<T> Vold = V;
85 gsVector<T> Aold = A;
86
87 gsVector<T> R;
88 gsSparseMatrix<T> M, C, K;
89
90 T alpha = m_options.getReal("alpha");
91 T delta = m_options.getReal("delta");
92
93 gsSparseMatrix<T> lhs;
94 A.setZero();
95 U = Uold + dt*Vold + Aold*(0.5 - alpha)*dt*dt + alpha*dt*dt*A;
96 V = Vold + Aold*(1 - delta)*dt + delta*dt*A;
97 // Computed at t=t0+dt
98 this->_computeMass(t+dt,M);
99 this->_computeDamping(U,t+dt,C);
100 this->_computeResidual(U,t+dt,R);
101
102 gsMatrix<T> rhs = R - C * V - M*(A);
103
104 T tolU = m_options.getReal("TolU");
105 T tolF = m_options.getReal("TolF");
106 T updateNorm = 10.0*tolU;
107 T residualNorm = rhs.norm();
108 T residualNorm0 = (residualNorm!=0) ? residualNorm : 1;
109 gsVector<T> dA;
110 T Anorm, dAnorm;
111 this->_initOutput();
112 typename gsSparseSolver<T>::uPtr solver;
113 for (index_t numIterations = 0; numIterations < m_options.getInt("MaxIter"); ++numIterations)
114 {
115 if ((!m_options.getSwitch("Quasi")) || ((numIterations==0) || (numIterations % m_options.getInt("QuasiIterations") == 0)))
116 {
117 // Computed at t=t0+dt
118 this->_computeDamping(U,t+dt,C);
119 this->_computeJacobian(U,t+dt,K);
120 }
121
122 lhs = M + delta*dt*C + dt*dt*alpha*K;
123
124 solver = gsSparseSolver<T>::get( m_options.getString("Solver") );
125 solver->compute(lhs);
126
127 dA = solver->solve(rhs);
128 A += dA;
129 V += dA*delta*dt;
130 U += dA*alpha*dt*dt;
131
132 Anorm = A.norm();
133 dAnorm = dA.norm();
134 updateNorm = (Anorm!=0) ? dAnorm/Anorm : dAnorm;
135
136 this->_computeResidual(U,t+dt,R);
137 rhs = R - C*V - M*A;
138 residualNorm = rhs.norm() / residualNorm0;
139
140 this->_stepOutput(numIterations,residualNorm,updateNorm);
141
142 if ( (updateNorm<tolU && residualNorm<tolF) )
143 {
144 return gsStatus::Success;
145 }
146 }
147
148 gsInfo<<"maximum iterations reached. Solution did not converge\n";
150}
151
152template <class T, bool _NL>
154 gsVector<T> & U, gsVector<T> & V,
155 gsVector<T> & A) const
156{
158 status = _step_impl<_NL>(t,dt,U,V,A);
159 return status;
160}
161
162
163template <class T, bool _NL>
165{
166 if (m_options.getSwitch("Verbose"))
167 {
168 gsInfo<<"\t";
169 gsInfo<<std::setw(4)<<std::left<<"It.";
170 gsInfo<<std::setw(17)<<std::left<<"|R|/|R0|";
171 gsInfo<<std::setw(17)<<std::left<<"|dU|/|U0|"<<"\n";
172 }
173}
174
175template <class T, bool _NL>
176void gsDynamicNewmark<T,_NL>::_stepOutput(const index_t it, const T resnorm, const T updatenorm) const
177{
178 if (m_options.getSwitch("Verbose"))
179 {
180 gsInfo<<"\t";
181 gsInfo<<std::setw(4)<<std::left<<it;
182 gsInfo<<std::setw(17)<<std::left<<resnorm;
183 gsInfo<<std::setw(17)<<std::left<<updatenorm<<"\n";
184 }
185}
186
187} // namespace gismo
Performs the arc length method to solve a nonlinear system of equations.
Definition gsDynamicNewmark.h:35
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
virtual void defaultOptions() override
Set default options.
Definition gsDynamicNewmark.hpp:20
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
Sparse matrix class, based on gsEigen::SparseMatrix.
Definition gsSparseMatrix.h:139
Abstract class for solvers. The solver interface is base on 3 methods: -compute set the system matrix...
Definition gsSparseSolver.h:67
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
#define index_t
Definition gsConfig.h:32
#define gsInfo
Definition gsDebug.h:43
The G+Smo namespace, containing all definitions for the library.
gsStatus
Definition gsStructuralAnalysisTypes.h:21
@ Success
Successful.
@ NotConverged
Step did not converge.
@ NotStarted
ALM has not started yet.