G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsALMConsistentCrisfield.hpp
Go to the documentation of this file.
1
14#pragma once
15
16#include <typeinfo>
17#include <gsStructuralAnalysis/src/gsALMSolvers/gsALMHelper.h>
18
19namespace gismo
20{
21
22template <class T>
24{
25 Base::defaultOptions();
26 m_options.addReal("Scaling","Set Scaling factor Phi",-1);
27}
28
29template <class T>
31{
32 Base::getOptions();
33 m_phi = m_options.getReal("Scaling");
34 m_phi_user = m_phi == -1 ? false : true;
35}
36
38template <class T>
40{
41 m_numDof = m_forcing.size();
42 m_DeltaU = m_U = gsVector<T>::Zero(m_numDof);
43 m_DeltaL = m_L = 0.0;
44
45 m_Uprev = gsVector<T>::Zero(m_numDof);
46 m_Lprev = 0.0;
47}
48
49// ------------------------------------------------------------------------------------------------------------
50// ---------------------------------------Consistent Crisfield's method----------------------------------------
51// ------------------------------------------------------------------------------------------------------------
52
54template <class T>
56{
57 m_jacMat = computeJacobian();
58 this->factorizeMatrix(m_jacMat);
59 computeUt(); // rhs does not depend on solution
60 computeUbar(); // rhs contains residual and should be computed every time
61
62}
63
65template <class T>
67{
68 m_jacMat = computeJacobian();
69 this->factorizeMatrix(m_jacMat);
70 computeUt(); // rhs does not depend on solution
71}
72
73
74template <class T>
75void gsALMConsistentCrisfield<T>::iteration() // see Carrera1994 eqs 20 and 21
76{
77 computeUbar(); // rhs contains residual and should be computed every time
78
79 T cres = m_DeltaU.dot(m_DeltaU) + m_phi*m_phi*m_DeltaL*m_DeltaL * m_forcing.dot(m_forcing) - m_arcLength*m_arcLength;
80 T num = cres + (2*m_DeltaU).dot(m_deltaUbar);
81 T denum = (2*m_DeltaU).dot(m_deltaUt) + m_phi*m_phi*2*m_DeltaL*m_forcing.dot(m_forcing);
82 m_deltaL = - num / denum;
83 m_deltaU = m_deltaL * m_deltaUt + m_deltaUbar;
84
85 m_DeltaU += m_deltaU;
86 m_DeltaL += m_deltaL;
87}
88
89template <class T>
91{
92 // (m_U,m_L) is the present solution (iteratively updated)
93 // (m_Uprev,m_Lprev) is the previously converged solution before (m_Lprev,m_Uprev)
94
95 // Reset step
96 m_DeltaU = m_deltaU = gsVector<T>::Zero(m_numDof);
97 m_DeltaL = m_deltaL = 0.0;
98
99}
100
101template <class T>
103{
104 m_jacMat = computeJacobian();
105 this->factorizeMatrix(m_jacMat);
106
107 // Check if the solution on start and prev are similar.
108 // Then compute predictor of the method
109 T tol = 1e-10;
110
111 if ( ((m_U-m_Uprev).norm() < tol) && ((m_L - m_Lprev) * (m_L - m_Lprev) < tol ) )
112 {
113 m_note+= "predictor\t";
114 T DL = 1.;
115 m_deltaUt = this->solveSystem(m_forcing);
116 m_deltaU = m_deltaUt / math::sqrt( m_deltaUt.dot(m_deltaUt) + m_DeltaL*DL );
117 m_deltaL = DL / math::sqrt( m_deltaUt.dot(m_deltaUt) + m_DeltaL*DL );
118
119 if (!m_phi_user)
120 m_phi = math::pow( m_deltaUt.dot(m_deltaUt) / m_forcing.dot(m_forcing),0.5);
121 }
122 else
123 {
124 m_deltaL = 1./m_arcLength_prev*(m_L - m_Lprev);
125 m_deltaU = 1./m_arcLength_prev*(m_U - m_Uprev);
126
127 if (!m_phi_user)
128 m_phi = math::pow( m_deltaUt.dot(m_deltaUt) / m_forcing.dot(m_forcing),0.5);
129 }
130
131 // Update iterative step
132 m_deltaL *= m_arcLength;
133 m_deltaU *= m_arcLength;
134
135 // Update load step
136 m_DeltaU += m_deltaU;
137 m_DeltaL += m_deltaL;
138}
139
140template <class T>
142{
143 m_jacMat = computeJacobian();
144 this->factorizeMatrix(m_jacMat);
145
146 // Check if the solution on start and prev are similar.
147 // Then compute predictor of the method
148 T tol = 1e-10;
149
150 if ( ((m_Uguess-m_U).norm() < tol) && ((m_Lguess - m_L) * (m_Lguess - m_L) < tol ) )
151 {
152 m_note+= "predictor\t";
153 T DL = 1.;
154 m_deltaUt = this->solveSystem(m_forcing);
155 m_deltaU = m_deltaUt / math::sqrt( m_deltaUt.dot(m_deltaUt) + m_DeltaL*DL );
156 m_deltaL = DL / math::sqrt( m_deltaUt.dot(m_deltaUt) + m_DeltaL*DL );
157
158 if (!m_phi_user)
159 m_phi = math::pow( m_deltaUt.dot(m_deltaUt) / m_forcing.dot(m_forcing),0.5);
160 }
161
162 else
163 {
164 m_deltaL = 1./m_arcLength_prev*(m_Lguess - m_L);
165 m_deltaU = 1./m_arcLength_prev*(m_Uguess - m_U);
166
167 if (!m_phi_user)
168 m_phi = math::pow( m_deltaUt.dot(m_deltaUt) / m_forcing.dot(m_forcing),0.5);
169 }
170
171 // Update iterative step
172 m_deltaL *= m_arcLength;
173 m_deltaU *= m_arcLength;
174
175 // Update load step
176 m_DeltaU += m_deltaU;
177 m_DeltaL += m_deltaL;
178}
179
180template <class T>
182{
183 m_converged = true;
184 m_Uprev = m_U;
185 m_Lprev = m_L;
186 m_U += m_DeltaU;
187 m_L += m_DeltaL;
188 m_DeltaUold = m_DeltaU;
189 m_DeltaLold = m_DeltaL;
190}
191
192
193// ------------------------------------------------------------------------------------------------------------
194// ---------------------------------------Output functions-----------------------------------------------------
195// ------------------------------------------------------------------------------------------------------------
196
197template <class T>
199{
200 gsInfo<<"\t";
201 gsInfo<<std::setw(4)<<std::left<<"It.";
202 gsInfo<<std::setw(17)<<std::left<<"Res. F";
203 gsInfo<<std::setw(17)<<std::left<<"|dU|/|Du|";
204 gsInfo<<std::setw(17)<<std::left<<"dL/DL";
205 gsInfo<<std::setw(17)<<std::left<<"|U|";
206 gsInfo<<std::setw(17)<<std::left<<"L";
207 gsInfo<<std::setw(17)<<std::left<<"|DU|";
208 gsInfo<<std::setw(17)<<std::left<<"DL";
209 gsInfo<<std::setw(17)<<std::left<<"|dU|";
210 gsInfo<<std::setw(17)<<std::left<<"dL";
211 gsInfo<<std::setw(17)<<std::left<<"ds²";
212 gsInfo<<std::setw(17)<<std::left<<"|dU|²";
213 gsInfo<<std::setw(17)<<std::left<<"dL²";
214 gsInfo<<std::setw(17)<<std::left<<"Dmin";
215 gsInfo<<std::setw(17)<<std::left<<"m_note";
216 gsInfo<<"\n";
217
218 m_note = "";
219}
220
221template <class T>
223{
224 computeStability(false);
225
226 gsInfo<<"\t";
227 gsInfo<<std::setw(4)<<std::left<<m_numIterations;
228 gsInfo<<std::setw(17)<<std::left<<m_residueF;
229 gsInfo<<std::setw(17)<<std::left<<m_residueU;
230 gsInfo<<std::setw(17)<<std::left<<m_residueL;
231 gsInfo<<std::setw(17)<<std::left<<(m_U+m_DeltaU).norm();
232 gsInfo<<std::setw(17)<<std::left<<(m_L + m_DeltaL);
233 gsInfo<<std::setw(17)<<std::left<<m_DeltaU.norm();
234 gsInfo<<std::setw(17)<<std::left<<m_DeltaL;
235 gsInfo<<std::setw(17)<<std::left<<m_deltaU.norm();
236 gsInfo<<std::setw(17)<<std::left<<m_deltaL;
237 gsInfo<<std::setw(17)<<std::left<<m_phi * math::pow(m_DeltaU.norm(),2.0) + (1.0-m_phi) * math::pow(m_DeltaL,2.0);
238 gsInfo<<std::setw(17)<<std::left<<m_phi * math::pow(m_DeltaU.norm(),2.0);
239 gsInfo<<std::setw(17)<<std::left<<(1-m_phi) * math::pow(m_DeltaL,2.0);
240 gsInfo<<std::setw(17)<<std::left<<m_indicator;
241 gsInfo<<std::setw(17)<<std::left<<m_note;
242 gsInfo<<"\n";
243
244 m_note = "";
245}
246
247} // namespace gismo
Performs the Consistent Crisfield arc length method to solve a nonlinear equation system.
Definition gsALMConsistentCrisfield.h:35
void defaultOptions()
See gsALMBase.
Definition gsALMConsistentCrisfield.hpp:23
void predictor()
See gsALMBase.
Definition gsALMConsistentCrisfield.hpp:102
void initOutput()
See gsALMBase.
Definition gsALMConsistentCrisfield.hpp:198
void stepOutput()
See gsALMBase.
Definition gsALMConsistentCrisfield.hpp:222
void iteration()
See gsALMBase.
Definition gsALMConsistentCrisfield.hpp:75
void quasiNewtonPredictor()
See gsALMBase.
Definition gsALMConsistentCrisfield.hpp:55
void initiateStep()
See gsALMBase.
Definition gsALMConsistentCrisfield.hpp:90
void quasiNewtonIteration()
See gsALMBase.
Definition gsALMConsistentCrisfield.hpp:66
void initMethods()
See gsALMBase.
Definition gsALMConsistentCrisfield.hpp:39
void iterationFinish()
See gsALMBase.
Definition gsALMConsistentCrisfield.hpp:181
void getOptions()
See gsALMBase.
Definition gsALMConsistentCrisfield.hpp:30
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
#define gsInfo
Definition gsDebug.h:43
The G+Smo namespace, containing all definitions for the library.