G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsStaticDR.hpp
1
14#pragma once
15
16namespace gismo
17{
18
19template <class T>
21{
22 Base::defaultOptions();
23 m_options.addReal("damping","damping factor",1.0);
24 m_options.addReal("alpha","mass coefficient",2.0);
25 m_options.addReal("tolE","Kinetic energy tolerance",1e-6);
26 m_options.addInt("ResetIt","Reset rate of velocities if damping is zero",-1);
27}
28
29template <class T>
31{
32 Base::getOptions();
33 m_c = m_options.getReal("damping");
34 m_alpha = m_options.getReal("alpha");
35 m_tolE = m_options.getReal("tolE");
36 m_resetIterations = m_options.getInt("ResetIt");
37}
38
39template <class T>
41{
42 gsInfo<<"\t";
43 gsInfo<<std::setw(4)<<std::left<<"It.";
44 gsInfo<<std::setw(17)<<std::left<<"|R|";
45 gsInfo<<std::setw(17)<<std::left<<"|R|/|R0|";
46 gsInfo<<std::setw(17)<<std::left<<" Ek";
47 gsInfo<<std::setw(17)<<std::left<<" Ek/Ek0";
48 gsInfo<<std::setw(17)<<std::left<<"|dU|";
49 gsInfo<<std::setw(17)<<std::left<<"|dU|/|DU|";
50 gsInfo<<std::setw(17)<<std::left<<"|dU|/|U+DU|";
51 gsInfo<<std::setw(17)<<std::left<<"|dV|";
52 gsInfo<<"\n";
53}
54
55template <class T>
57{
58 gsInfo<<"\t";
59 gsInfo<<std::setw(4)<<std::left<<k;
60 gsInfo<<std::setw(17)<<std::left<<m_residual;
61 gsInfo<<std::setw(17)<<std::left<<m_residual/m_residualIni;
62 gsInfo<<std::setw(17)<<std::left<<m_Ek;
63 gsInfo<<std::setw(17)<<std::left<<m_Ek/m_Ek0;
64 gsInfo<<std::setw(17)<<std::left<<m_deltaU.norm();
65 gsInfo<<std::setw(17)<<std::left<<m_deltaU.norm()/m_DeltaU.norm();
66 gsInfo<<std::setw(17)<<std::left<<m_deltaU.norm()/(m_U+m_DeltaU).norm();
67 gsInfo<<std::setw(17)<<std::left<<m_v.norm();
68 gsInfo<<"\n";
69}
70
71template <class T>
73{
74 // try
75 // {
76 _solve();
77 m_status = gsStatus::Success;
78 // }
79 // catch (int errorCode)
80 // {
81 // if (errorCode==1)
82 // m_status = gsStatus::NotConverged;
83 // else if (errorCode==2)
84 // m_status = gsStatus::AssemblyError;
85 // else if (errorCode==3)
86 // m_status = gsStatus::SolverError;
87 // else
88 // m_status = gsStatus::OtherError;
89 // }
90 // catch (...)
91 // {
92 // m_status = gsStatus::OtherError;
93 // }
94 return m_status;
95}
96
97template <class T>
99{
100 m_Eks.clear();
101 m_Eks.reserve(m_maxIterations);
102
103 if (m_verbose) initOutput();
104
105 _start(); // initializes m_Ek, etc
106
107 m_Ek0 = m_Ek;
108 m_Eks.push_back(m_Ek);
109
110 if (m_verbose != 0) stepOutput(0);
111 index_t resetIt = 0;
112 for (m_numIterations=1; m_numIterations!=m_maxIterations; m_numIterations++, resetIt++)
113 {
114 _iteration();
115 if ((m_c==0 && m_Ek_prev > m_Ek) || resetIt==m_resetIterations)// || (m_Ek/m_Ek_prev > 1/m_tolE && m_Ek_prev!=0))
116 {
117 resetIt = 0;
118 _peak();
119 }
120
121 if (m_verbose!=0)
122 if (m_numIterations % m_verbose == 0 || m_verbose==-1 ) stepOutput(m_numIterations);
123
124 m_Eks.push_back(m_Ek);
125
126 m_residualOld = m_residual;
127
128 if (m_residual/m_residualIni < m_tolF && m_Ek/m_Ek0 < m_tolE && m_deltaU.norm()/m_DeltaU.norm() < m_tolU)
129 {
130 m_U += m_DeltaU;
131 gsDebug <<"Converged: \n";
132 gsDebug <<"\t |R|/|R0| = "<<m_residual/m_residualIni<<" < tolF = "<<m_tolF<<"\n";
133 gsDebug <<"\t |E|/|E0| = "<<m_Ek/m_Ek0 <<" < tolE = "<<m_tolE<<"\n";
134 gsDebug <<"\t |U|/|U0| = "<<m_deltaU.norm()/m_DeltaU.norm()<<" < tolF = "<<m_tolU<<"\n";
135 break;
136 }
137 if (m_numIterations==m_maxIterations-1)
138 {
139 m_U += m_DeltaU;
140 gsInfo<<"Maximum iterations reached. Solution did not converge\n";
141 throw 1;
142 }
143 }
144 gsInfo<<"\n";
145};
146
147template <class T>
149{
150 this->reset();
151 getOptions();
152 m_massInv *= 1./m_alpha;
153 m_damp = m_c * m_mass;
154}
155
156template <class T>
158{
159 m_dofs = m_mass.rows();
160 m_massInv = m_mass.array().inverse();
161 m_Ek = m_Ek0 = 0.0;
162 // resets m_U, m_DeltaU, m_deltaU, m_R, m_L, m_DeltaL, m_deltaL and m_headstart
163 Base::reset();
164 m_v.setZero(m_dofs);
165 m_status = gsStatus::NotStarted;
166}
167
168template <class T>
170{
171 this->reset();
172 m_dt = 1.0;
173 defaultOptions();
174}
175
176template <class T>
178{
179 gsVector<T> resVec;
180 if (!m_residualFun(U, resVec))
181 throw 2;
182 return resVec;
183}
184
185template <class T>
187{
188 m_Ek_prev = m_Ek;
189 m_R = _computeResidual(m_U+m_DeltaU) - m_damp.cwiseProduct(m_v);
190 m_residual = m_R.norm();
191//----------------------------------------------------------------------------------
192 m_v += m_dt * m_massInv.cwiseProduct(m_R); // Velocities at t+dt/2
193 m_deltaU = m_dt * m_v; // Velocities at t+dt
194//----------------------------------------------------------------------------------
195 // m_deltaU = m_deltaU + m_dt*m_dt*m_massInv*m_R;
196 // m_v = m_deltaU/m_dt;
197//----------------------------------------------------------------------------------
198
199 m_DeltaU += m_deltaU; // Velocities at t+dt
200 m_Ek = m_v.transpose() * m_mass.cwiseProduct(m_v);
201}
202
203template <class T>
205{
206 m_R = _computeResidual(m_U+m_DeltaU)- m_damp.cwiseProduct(m_v);
207 m_deltaU = - 1.5 * m_dt * m_v + m_dt*m_dt / 2. * m_massInv.cwiseProduct(m_R);
208 m_DeltaU += m_deltaU;
209
210 m_v.setZero();
211 m_Ek = 0.0;
212 // m_v = 0.5 * m_dt * m_massInv * m_R; // Velocities at dt/2
213}
214
215template <class T>
217{
218 m_v.setZero();
219
220//----------------------------------------------------------------------------------
221 // Define residual measures:
222 // m_residual = current residual
223 // m_residualIni = residual on m_U
224 // m_residualOld = residual in previous step
225 if (!m_headstart) // no headstart
226 {
227 // We can reset the update to ensure we properly restart
228 m_DeltaU.setZero(m_dofs);
229 // Compute current residual and its norm
230 m_R = _computeResidual(m_U);
231 // m_residual = m_R.norm();
232 m_residual = m_forcing.norm();
233 // If the residual is 0 (e.g. with purely displacment loading), we set it to 1 to allow divisions
234 if (m_residual==0) m_residual=1;
235 // All residual norms are equal
236 m_residualIni = m_residualOld = m_residual;
237 }
238 else
239 {
240 // Compute current residual and its norm
241 m_R = _computeResidual(m_U + m_DeltaU);
242 // m_residual = m_R.norm();
243 m_residual = m_forcing.norm();
244 // If the residual is 0 (e.g. with purely displacment loading), we set it to 1 to allow divisions
245 if (m_residual==0) m_residual=1;
246 // The previous step residual is the same as the residual
247 m_residualOld = m_residual;
248 // Residual0 is the residual without m_DeltaU
249 m_residualIni = _computeResidual(m_U).norm();
250 // If the residual is 0 (e.g. with purely displacment loading), we set it to 1 to allow divisions
251 if (m_residualIni==0) m_residualIni=1;
252
253 // Now we can reset the headstart
254 m_headstart = false;
255 }
256
257//----------------------------------------------------------------------------------
258 // m_R = _computeResidual(m_U+m_DeltaU);
259 m_deltaU = 0.5*m_dt*m_dt*m_massInv.cwiseProduct(m_R);
260 m_v = m_deltaU/m_dt;
261//----------------------------------------------------------------------------------
262
263 m_DeltaU += m_deltaU;
264
265 m_Ek = m_v.transpose() * m_mass.cwiseProduct(m_v);
266}
267
268} // namespace gismo
Static solver using the Dynamic Relaxation method.
Definition gsStaticDR.h:28
void _start()
Starts the method.
Definition gsStaticDR.hpp:216
void stepOutput(index_t k) override
See gsStaticBase.
Definition gsStaticDR.hpp:56
void _peak()
Identifies a peak.
Definition gsStaticDR.hpp:204
void initialize() override
See gsStaticBase.
Definition gsStaticDR.hpp:148
void _init()
Initializes the method.
Definition gsStaticDR.hpp:169
void initOutput() override
See gsStaticBase.
Definition gsStaticDR.hpp:40
void defaultOptions() override
See gsStaticBase.
Definition gsStaticDR.hpp:20
void getOptions() override
See gsStaticBase.
Definition gsStaticDR.hpp:30
gsStatus solve() override
gsStaticBase base functions
Definition gsStaticDR.hpp:72
void _iteration()
Performs an iteration.
Definition gsStaticDR.hpp:186
void _solve()
See solve()
Definition gsStaticDR.hpp:98
void reset() override
See gsStaticBase.
Definition gsStaticDR.hpp:157
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
#define index_t
Definition gsConfig.h:32
#define gsDebug
Definition gsDebug.h:61
#define gsInfo
Definition gsDebug.h:43
The G+Smo namespace, containing all definitions for the library.
gsStatus
Definition gsStructuralAnalysisTypes.h:21
@ Success
Successful.
@ NotStarted
ALM has not started yet.