G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsStaticNewton.hpp
Go to the documentation of this file.
1
14#pragma once
15
16namespace gismo
17{
18
19template <class T>
21{
22 Base::defaultOptions();
23 m_options.setString("Solver","CGDiagonal"); // The CG solver is robust for membrane models, where zero-blocks in the matrix might occur.
24 m_options.addReal("Relaxation","Relaxation parameter",1);
25};
26
27template <class T>
29{
30 Base::getOptions();
31
32 // if (m_solverType!=solver::LDLT && m_stabilityMethod==stabmethod::Determinant)
33 // {
34 // gsWarn<<"Determinant method cannot be used with solvers other than LDLT. Bifurcation method will be set to 'Eigenvalue'.\n";
35 // m_stabilityMethod = stabmethod::Eigenvalue;
36 // }
37
38 m_relax = m_options.getReal("Relaxation");
39};
40
41template <class T>
43{
44 gsInfo<<"\t";
45 gsInfo<<std::setw(4)<<std::left<<"It.";
46 gsInfo<<std::setw(17)<<std::left<<"|R|";
47 gsInfo<<std::setw(17)<<std::left<<"|R|/|R0|";
48 gsInfo<<std::setw(17)<<std::left<<"|DU|";
49 gsInfo<<std::setw(17)<<std::left<<"|dU|";
50 gsInfo<<std::setw(17)<<std::left<<"|dU|/|DU|";
51 gsInfo<<std::setw(17)<<std::left<<"|dU|/|U+DU|";
52 gsInfo<<std::setw(17)<<std::left<<"log(Ri/R0):";
53 gsInfo<<std::setw(17)<<std::left<<"log(Ri+1/R0)";
54 gsInfo<<"\n";
55}
56
57template <class T>
59{
60 gsInfo<<"\t";
61 gsInfo<<std::setw(4)<<std::left<<k;
62 gsInfo<<std::setw(17)<<std::left<<m_residual;
63 gsInfo<<std::setw(17)<<std::left<<m_residual/m_residualIni;
64 gsInfo<<std::setw(17)<<std::left<<m_DeltaU.norm();
65 gsInfo<<std::setw(17)<<std::left<<m_relax * m_deltaU.norm();
66 gsInfo<<std::setw(17)<<std::left<<m_relax * m_deltaU.norm()/m_DeltaU.norm();
67 gsInfo<<std::setw(17)<<std::left<<m_relax * m_deltaU.norm()/(m_U+m_DeltaU).norm();
68 gsInfo<<std::setw(17)<<std::left<<math::log10(m_residualOld/m_residualIni);
69 gsInfo<<std::setw(17)<<std::left<<math::log10(m_residual/m_residualIni);
70 gsInfo<<"\n";
71
72}
73
74template <class T>
76{
77 try
78 {
79 _solveLinear();
80 m_status = gsStatus::Success;
81 }
82 catch (int errorCode)
83 {
84 if (errorCode==1)
85 m_status = gsStatus::NotConverged;
86 else if (errorCode==2)
87 m_status = gsStatus::AssemblyError;
88 else if (errorCode==3)
89 m_status = gsStatus::SolverError;
90 else
91 m_status = gsStatus::OtherError;
92 }
93 catch (...)
94 {
95 m_status = gsStatus::OtherError;
96 }
97 return m_status;
98};
99
100template <class T>
102{
103 this->getOptions();
104 if (m_verbose==2)
105 {
106 gsInfo<<"Matrix: \n"<<m_linear.toDense()<<"\n";
107 gsInfo<<"Vector: \n"<<m_force<<"\n";
108 }
109 _factorizeMatrix( m_linear );
110 m_U = _solveSystem(m_force);
111 return m_U;
112};
113
114template <class T>
116{
117 try
118 {
119 _solveNonlinear();
120 m_status = gsStatus::Success;
121 }
122 catch (int errorCode)
123 {
124 if (errorCode==1)
125 m_status = gsStatus::NotConverged;
126 else if (errorCode==2)
127 m_status = gsStatus::AssemblyError;
128 else if (errorCode==3)
129 m_status = gsStatus::SolverError;
130 else
131 m_status = gsStatus::OtherError;
132 }
133 catch (...)
134 {
135 m_status = gsStatus::OtherError;
136 }
137 return m_status;
138}
139
140
141template <class T>
143{
144 this->getOptions();
145 // m_start: true -> m_U given
146 // m_headstart: true -> m_DeltaU given
147
148 if (m_DeltaU.norm()==0 && m_DeltaU.rows()==0)
149 {
150 m_deltaU = m_DeltaU = this->_solveLinear();
151 m_U.setZero(); // Needed because linear solve modifies m_U.
152 m_headstart = true; // due to this, the relative residual is based on the solution of the linear solve
153 }
154 _start();
155
156 gsSparseMatrix<T> jacMat;
157
158 if (m_verbose>0) { initOutput(); }
159
160 for (m_numIterations = 0; m_numIterations != m_maxIterations; ++m_numIterations)
161 {
162 jacMat = this->_computeJacobian(m_U+m_DeltaU,m_deltaU);
163 if (m_verbose==2)
164 {
165 gsInfo<<"Matrix: \n"<<jacMat.toDense()<<"\n";
166 gsInfo<<"Vector: \n"<<m_R<<"\n";
167 }
168
169 _factorizeMatrix(jacMat);
170 m_deltaU = _solveSystem(m_R); // this is the UPDATE
171 m_DeltaU += m_relax * m_deltaU;
172
173 m_R = this->_computeResidual(m_U+m_DeltaU);
174 m_residual = m_R.norm();
175
176 if (m_verbose>0) { stepOutput(m_numIterations); }
177
178 m_residualOld = m_residual;
179
180 if (m_relax * m_deltaU.norm()/m_DeltaU.norm() < m_tolU && m_residual/m_residualIni < m_tolF)
181 {
182 m_U+=m_DeltaU;
183 break;
184 }
185 else if (m_numIterations+1 == m_maxIterations)
186 {
187 m_U+=m_DeltaU;
188 gsInfo<<"Maximum iterations reached. Solution did not converge\n";
189 throw 1;
190 }
191 }
192 return m_U;
193};
194
195template <class T>
197{
198 gsVector<T> resVec;
199 if (!m_residualFun(U, resVec))
200 throw 2;
201 return resVec;
202}
203
204template <class T>
205gsSparseMatrix<T> gsStaticNewton<T>::_computeJacobian(const gsVector<T> & U, const gsVector<T> & deltaU)
206{
207 // Compute Jacobian
208 gsSparseMatrix<T> m;
209 if (!m_dnonlinear(U,deltaU,m))
210 throw 2;
211 return m;
212}
213
214template <class T>
216{
217 m_solver->compute(jacMat);
218 if (m_solver->info()!=gsEigen::ComputationInfo::Success)
219 {
220 gsInfo<<"Solver error with code "<<m_solver->info()<<". See Eigen documentation on ComputationInfo \n"
221 <<gsEigen::ComputationInfo::Success<<": Success"<<"\n"
222 <<gsEigen::ComputationInfo::NumericalIssue<<": NumericalIssue"<<"\n"
223 <<gsEigen::ComputationInfo::NoConvergence<<": NoConvergence"<<"\n"
224 <<gsEigen::ComputationInfo::InvalidInput<<": InvalidInput"<<"\n";
225 throw 3;
226 }
227}
228
229template <class T>
231{
232 try
233 {
234 return m_solver->solve(F);
235 }
236 catch (...)
237 {
238 throw 3;
239 }
240}
241
242template <class T>
244{
245 m_dofs = m_force.rows();
246 // resets m_U, m_DeltaU, m_deltaU, m_R, m_L, m_DeltaL, m_deltaL and m_headstart
247 Base::reset();
248}
249
250template <class T>
252{
253 this->reset();
254 if( m_dnonlinear==nullptr || m_residualFun==nullptr)
255 m_NL=false;
256 else
257 m_NL = true;
258
259 m_stabilityMethod = 0;
260 m_start = false;
261
262 if (m_dofs==0)
263 gsWarn<<"The number of degrees of freedom is equal to zero. This can lead to bad initialization.\n";
264
265 m_residual = m_residualIni = m_residualOld = 0;
266
267 defaultOptions();
268
269 m_status = gsStatus::NotStarted;
270}
271
272template <class T>
274{
275 // Define residual measures:
276 // residual = current residual
277 // residual0 = residual on m_U
278 // residualOld = residual in previous step
279
280 if (!m_headstart) // no headstart
281 {
282 // We can reset the update to ensure we properly restart
283 if (m_dofs==0) gsWarn<<"The number of degrees of freedom is equal to zero. This can lead to bad initialization.\n";
284 m_DeltaU.setZero(m_dofs);
285 // Compute current residual and its norm
286 m_R = this->_computeResidual(m_U);
287 m_residual = m_R.norm();
288 // If the residual is 0 (e.g. with purely displacment loading), we set it to 1 to allow divisions
289 if (m_residual==0) m_residual=1;
290 // All residual norms are equal
291 m_residualIni = m_residualOld = m_residual;
292 }
293 else
294 {
295 // If we have a headstart, we need to compute Residual0 on the solution m_U
296 // Compute current residual and its norm
297 m_R = this->_computeResidual(m_U + m_DeltaU);
298 m_residual = m_R.norm();
299 // If the residual is 0 (e.g. with purely displacment loading), we set it to 1 to allow divisions
300 if (m_residual==0) m_residual=1;
301 // The previous step residual is the same as the residual
302 m_residualOld = m_residual;
303 // Residual0 is the residual without m_DeltaU
304 m_residualIni = this->_computeResidual(m_U).norm();
305 // If the residual is 0 (e.g. with purely displacment loading), we set it to 1 to allow divisions
306 if (m_residualIni==0) m_residualIni=1;
307
308 // Now we can reset the headstart
309 m_headstart = false;
310 }
311
312 // Reset incremental update
313 m_deltaU.setZero(m_dofs);
314
315}
316
317} // namespace gismo
Sparse matrix class, based on gsEigen::SparseMatrix.
Definition gsSparseMatrix.h:139
Static solver using a newton method.
Definition gsStaticNewton.h:30
gsVector< T > _solveLinear()
Perform a linear solve.
Definition gsStaticNewton.hpp:101
void _start()
Starts the method.
Definition gsStaticNewton.hpp:273
void stepOutput(index_t k) override
See gsStaticBase.
Definition gsStaticNewton.hpp:58
gsVector< T > _solveNonlinear()
Perform a nonlinear solve.
Definition gsStaticNewton.hpp:142
void _init()
Initializes the method.
Definition gsStaticNewton.hpp:251
void initOutput() override
See gsStaticBase.
Definition gsStaticNewton.hpp:42
gsVector< T > _solveSystem(const gsVector< T > &F)
Solves the system with RHS F and LHS the Jacobian.
Definition gsStaticNewton.hpp:230
void defaultOptions() override
See gsStaticBase.
Definition gsStaticNewton.hpp:20
void getOptions() override
See gsStaticBase.
Definition gsStaticNewton.hpp:28
void _factorizeMatrix(const gsSparseMatrix< T > &jacMat) const
Factorizes the jacMat.
Definition gsStaticNewton.hpp:215
void reset() override
See gsStaticBase.
Definition gsStaticNewton.hpp:243
gsStatus solveLinear()
Perform a linear solve.
Definition gsStaticNewton.hpp:75
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
#define index_t
Definition gsConfig.h:32
#define gsWarn
Definition gsDebug.h:50
#define gsInfo
Definition gsDebug.h:43
The G+Smo namespace, containing all definitions for the library.
gsStatus
Definition gsStructuralAnalysisTypes.h:21
@ Success
Successful.
@ OtherError
Other error.
@ SolverError
Assembly problem in step.
@ AssemblyError
Assembly problem in step.
@ NotConverged
Step did not converge.
@ NotStarted
ALM has not started yet.