G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsDynamicBathe.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 Bathe's method, such that 0 =< 2 beta =< 1",0.25);
24 m_options.addReal("delta","Beta parameter for Bathe's method, such that 0 =< delta =< 1",0.50);
25 m_options.addReal("gamma","Beta parameter for Bathe's method, such that 0 =< gamma =< 1",0.50);
26}
27
28
29template <class T, bool _NL>
30template <bool _nonlinear>
31typename std::enable_if<(_nonlinear==false), gsStatus>::type
32gsDynamicBathe<T,_NL>::_step_impl(const T t, const T dt, gsVector<T> & U, gsVector<T> & V, gsVector<T> & A) const
33{
34 T gamma = m_options.getReal("gamma");
35
36 gsVector<T> Uold = U;
37 gsVector<T> Vold = V;
38 gsVector<T> Aold = A;
39
40 gsVector<T> Ustep = U;
41 gsVector<T> Vstep = V;
42 gsVector<T> Astep = A;
43
45 this->_stageOutput(1);
46 Base::_step(t,gamma*dt,Ustep,Vstep,Astep);
47
50 gsSparseMatrix<T> M, C, K;
51
52 // Computed at t=t0+dt
53 this->_computeMass(t+dt,M);
54 this->_computeForce(t+dt,F);
55 this->_computeDamping(U,t+dt,C);
56 this->_computeJacobian(U,t+dt,K);
57
58 T c1 = (1-gamma)/(gamma*dt);
59 T c2 = (-1)/((1-gamma)*gamma*dt);
60 T c3 = (2-gamma)/((1-gamma)*dt);
61
62 gsSparseMatrix<T> lhs = K + c3*c3*M + c3*C;
63 gsMatrix<T> rhs = F - M*(c1*Vold+c2*Vstep+c1*c3*Uold+c3*c2*Ustep) - C*(c2*Ustep+c1*Uold);
64
65 this->_stageOutput(2);
66 this->_initOutput();
67 typename gsSparseSolver<T>::uPtr solver;
68 solver = gsSparseSolver<T>::get( m_options.getString("Solver") );
69 solver->compute(lhs);
70 U = solver->solve(rhs);
71 V = c1*Uold + c2*Ustep + c3*U;
72 A = c1*Vold + c2*Vstep + c3*V;
73
74 this->_stepOutput(0,(F - K*U - M * A - C * V).norm(),U.norm());
75 if (math::isinf(U.norm()) || math::isnan(U.norm()))
77 else
78 return gsStatus::Success;
79}
80
81
82template <class T, bool _NL>
83template <bool _nonlinear>
84typename std::enable_if<(_nonlinear==true), gsStatus>::type
85gsDynamicBathe<T,_NL>::_step_impl(const T t, const T dt, gsVector<T> & U, gsVector<T> & V, gsVector<T> & A) const
86{
87 T gamma = m_options.getReal("gamma");
88
89 T c1 = (1-gamma)/(gamma*dt);
90 T c2 = (-1)/((1-gamma)*gamma*dt);
91 T c3 = (2-gamma)/((1-gamma)*dt);
92
93
94 gsVector<T> Uold = U;
95 gsVector<T> Vold = V;
96 gsVector<T> Aold = A;
97
98 gsVector<T> Ustep = U;
99 gsVector<T> Vstep = V;
100 gsVector<T> Astep = A;
101
103 this->_stageOutput(1);
104 Base::_step(t,gamma*dt,Ustep,Vstep,Astep);
105
107 gsVector<T> R;
108 gsSparseMatrix<T> M, C, K;
110 gsMatrix<T> rhs;
111 gsMatrix<T> dU;
112
113 // Computed at t=t0+dt
114 this->_computeMass(t+dt,M);
115 this->_computeDamping(U,t+dt,C);
116 this->_computeResidual(U,t+dt,R);
117
118 U = Ustep;
119 V = Vstep;
120 A = Astep;
121
122 rhs = R - M*(c1*Vold+c2*Vstep+c1*c3*Uold+c3*c2*Ustep+c3*c3*U) - C*(c1*Uold+c2*Ustep+c3*U);
123
124 T tolU = m_options.getReal("TolU");
125 T tolF = m_options.getReal("TolF");
126 T updateNorm = 10.0*tolU;
127 T residualNorm = rhs.norm();
128 T residualNorm0 = (residualNorm!=0) ? residualNorm : 1;
129 this->_initOutput();
130 T Unorm, dUnorm;
131 for (index_t numIterations = 0; numIterations < m_options.getInt("MaxIter"); ++numIterations)
132 {
133 if ((!m_options.getSwitch("Quasi")) || ((numIterations==0) || (numIterations % m_options.getInt("QuasiIterations") == 0)))
134 {
135 // Computed at t=t0+dt
136 this->_computeDamping(U,t+dt,C);
137 this->_computeJacobian(U,t+dt,K);
138 }
139
140 lhs = K + c3*c3*M + c3*C;
141
142 typename gsSparseSolver<T>::uPtr solver;
143 solver = gsSparseSolver<T>::get( m_options.getString("Solver") );
144 solver->compute(lhs);
145 // U = solver->solve(rhs);
146 // V = delta/( dt*alpha ) * (U-Uold) - Vold;
147 // A = 1/( dt*dt*alpha ) * (U-Uold-Vold*dt) - Aold;
148
149 dU = solver->solve(rhs);
150 U += dU;
151
152 Unorm = U.norm();
153 dUnorm = dU.norm();
154 updateNorm = (Unorm!=0) ? dUnorm/Unorm : dUnorm;
155
156 this->_computeResidual(U,t+dt,R);
157 rhs = R - M*(c1*Vold+c2*Vstep+c1*c3*Uold+c3*c2*Ustep+c3*c3*U) - C*(c1*Uold+c2*Ustep+c3*U);
158 residualNorm = rhs.norm() / residualNorm0;
159 updateNorm = dU.norm() / U.norm();
160
161 this->_stepOutput(numIterations,residualNorm,updateNorm);
162
163 if ( (updateNorm<tolU && residualNorm<tolF) )
164 {
165 V = c1*Uold + c2*Ustep + c3*U;
166 A = c1*Vold + c2*Vstep + c3*V;
167 return gsStatus::Success;
168 }
169 }
170
171 V = c1*Uold + c2*Ustep + c3*U;
172 A = c1*Vold + c2*Vstep + c3*V;
173 gsInfo<<"maximum iterations reached. Solution did not converge\n";
175}
176
177template <class T, bool _NL>
179 gsVector<T> & U, gsVector<T> & V,
180 gsVector<T> & A) const
181{
183 status = _step_impl<_NL>(t,dt,U,V,A);
184 return status;
185}
186
187template <class T, bool _NL>
189{
190 if (m_options.getSwitch("Verbose"))
191 {
192 gsInfo<<"Stage "<<stage<<":";
193 gsInfo<<"\n";
194 }
195}
196
197template <class T, bool _NL>
198void gsDynamicBathe<T,_NL>::_initOutput() const
199{
200 if (m_options.getSwitch("Verbose"))
201 {
202 gsInfo<<"\t";
203 gsInfo<<std::setw(4)<<std::left<<"It.";
204 gsInfo<<std::setw(17)<<std::left<<"|R|/|R0|";
205 gsInfo<<std::setw(17)<<std::left<<"|dU|/|U0|"<<"\n";
206 }
207}
208
209template <class T, bool _NL>
210void gsDynamicBathe<T,_NL>::_stepOutput(const index_t it, const T resnorm, const T updatenorm) const
211{
212 if (m_options.getSwitch("Verbose"))
213 {
214 gsInfo<<"\t";
215 gsInfo<<std::setw(4)<<std::left<<it;
216 gsInfo<<std::setw(17)<<std::left<<resnorm;
217 gsInfo<<std::setw(17)<<std::left<<updatenorm<<"\n";
218 }
219}
220
221} // namespace gismo
Performs the arc length method to solve a nonlinear system of equations.
Definition gsDynamicBathe.h:36
gsStatus _step(const T t, const T dt, gsVector< T > &U, gsVector< T > &V, gsVector< T > &A) const override
Initialize the ALM.
Definition gsDynamicBathe.hpp:178
virtual void defaultOptions() override
Set default options.
Definition gsDynamicBathe.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
c1
See gismo/filedata/surfaces/simple.xml for the geometry.
Definition example_shell3D.py:40
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.