G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsALMBase.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 m_options.addInt ("MaxIter","Maximum iterations",100);
26 m_options.addReal("Tol","Tolerance",1e-6);
27 m_options.addReal("TolF","Tolerance",1e-3);
28 m_options.addReal("TolU","Tolerance",1e-6);
29 m_options.addReal("Perturbation","Set Perturbation factor Tau",1e3);
30
31 m_options.addReal("Length","Arclength",1e-2);
32
33 m_options.addSwitch("AdaptiveLength","Adaptive length",false);
34 m_options.addInt ("AdaptiveIterations","Desired iterations for adaptive length",10);
35
36 m_options.addSwitch("Quasi","Use Quasi Newton method",false);
37 m_options.addInt ("QuasiIterations","Number of iterations for quasi newton method",-1);
38
39 m_options.addInt ("BifurcationMethod","Bifurcation Identification based on: -1: nothing, 0: Determinant; 1: Eigenvalue",bifmethod::Eigenvalue);
40
41 m_options.addInt ("SingularPointFailure","What to do when a singular point determination fails?: 0 = Apply solution anyways; 1 = Proceed without singular point",SPfail::With);
42 m_options.addReal("SingularPointTestTol", "Detection tolerance for singular points. Product of the mode shape and the forcing should be below this tolerance", 1e-6);
43 m_options.addInt ("SingularPointTestIt" , "Number of iterations of the power method for the singular point test",5);
44
45 m_options.addReal("SingularPointComputeTolE", "Tolerance for the extended iterations to compute a bifurcation point", 1e-10);
46 m_options.addReal("SingularPointComputeTolB", "Tolerance for the bisection iterations to compute a bifurcation point. If tol = 0, no bi-section method is used.", 0);
47
48 m_options.addString("Solver","Sparse linear solver", "SimplicialLDLT");
49
50 m_options.addSwitch ("Verbose","Verbose output",false);
51
52 m_options.addReal("Relaxation","Set Relaxation factor alpha",1.0);
53
54}
55
56template <class T>
58{
59 m_maxIterations = m_options.getInt ("MaxIter");
60 m_tolerance = m_options.getReal("Tol");
61 m_toleranceF = m_options.getReal("TolF");
62 m_toleranceU = m_options.getReal("TolU");
63
64 m_tau = m_options.getReal("Perturbation");
65
66 m_adaptiveLength = m_options.getSwitch("AdaptiveLength");
67 m_desiredIterations = m_options.getInt ("AdaptiveIterations");
68
69 m_quasiNewton = m_options.getSwitch("Quasi");
70 m_quasiNewtonInterval = m_options.getInt ("QuasiIterations");
71
72 m_bifurcationMethod = m_options.getInt ("BifurcationMethod");
73 m_solver = gsSparseSolver<T>::get( m_options.getString("Solver") );
74 if (!dynamic_cast<typename gsSparseSolver<T>::SimplicialLDLT*>(m_solver.get()) && m_bifurcationMethod==bifmethod::Determinant)
75 {
76 gsWarn<<"Determinant method cannot be used with solvers other than LDLT. Bifurcation method will be set to 'Eigenvalue'.\n";
77 m_bifurcationMethod = bifmethod::Eigenvalue;
78 }
79
80 m_verbose = m_options.getSwitch ("Verbose");
81 m_relax = m_options.getReal("Relaxation");
82
83 m_arcLength = m_arcLength_prev = m_arcLength_ori = m_options.getReal("Length");
84
85 m_SPfail = m_options.getInt ("SingularPointFailure");
86 m_SPTestTol = m_options.getReal("SingularPointTestTol");
87 m_SPTestIt = m_options.getInt ("SingularPointTestIt");
88
89 m_SPCompTolE = m_options.getReal("SingularPointComputeTolE");
90 m_SPCompTolB = m_options.getReal("SingularPointComputeTolB");
91
92
93}
94
95template <class T>
96void gsALMBase<T>::init(bool stability)
97{
98 if (stability)
99 {
100 this->_computeStability(m_U);
101 m_stability = this->stability(); // requires Jabobian!!
102 }
103}
104
105template <class T>
107{
108 // m_desiredIterations := optimal number of iterations
109 T fac = m_desiredIterations / m_numIterations;
110 if (fac < 0.5)
111 fac = 0.5;
112 else if (fac > 2.0)
113 fac = 2.0;
114
115 m_arcLength_prev = m_arcLength;
116 m_arcLength = m_arcLength*fac;
117}
118
119template <class T>
121{
122 m_arcLength *= fac;
123 return m_arcLength;
124}
125
126template <class T>
128{
129 m_arcLength = m_arcLength_prev = m_arcLength_ori;
130 if (m_adaptiveLength)
131 computeLength();
132 return m_arcLength;
133}
134
135template <class T>
137{
138 gsVector<T> resVec;
139 if (!m_residualFun(U, L, resVec))
140 throw 2;
141 return resVec;
142}
143
144template <class T>
146{
147 m_resVec = this->computeResidual(m_U + m_DeltaU, m_L + m_DeltaL);
148}
149
150template <class T>
152{
153 if (m_numIterations ==0 ) // then define the residual
154 {
155 // m_basisResidualF = m_resVec.norm();
156 m_basisResidualF = ((m_L+m_DeltaL) * m_forcing).norm();
157 m_basisResidualU = m_DeltaU.norm();
158
159 // m_residueF = m_resVec.norm() / (((m_L+m_DeltaL) * m_forcing).norm());
160 m_residueF = m_resVec.norm() / (((m_L+m_DeltaL) * m_forcing).norm());
161 m_residueU = m_deltaU.norm()/m_DeltaU.norm();
162 // m_residueU = m_basisResidualU;
163 m_residueL = m_DeltaL;
164 }
165 else
166 {
167 // m_residueF = m_resVec.norm() / (((m_L+m_DeltaL) * m_forcing).norm());
168 m_residueF = m_resVec.norm() / m_basisResidualF;
169 m_residueU = m_deltaU.norm() / m_basisResidualU;
170 m_residueL = m_deltaL / m_DeltaL;
171 }
172 // gsInfo<<"m_DeltaU = "<<m_DeltaU.norm()
173 // <<"\t m_deltaU = "<<m_deltaU.norm()
174 // <<"\t m_deltaUbar = "<<m_deltaUbar.norm()
175 // <<"\t m_deltaUt = "<<m_deltaUt.norm()<<"\n";
176 // gsInfo<<"m_resVec.norm() = "<<m_resVec.norm()
177 // <<"\t ((m_L+m_DeltaL) * m_forcing).norm() = "<<((m_L+m_DeltaL) * m_forcing).norm()
178 // <<"\t m_deltaU.norm() = "<<m_deltaUbar.norm()
179 // <<"\t m_basisResidualU = "<<m_basisResidualU
180 // <<"\t m_deltaL = "<<m_deltaL
181 // <<"\t m_DeltaL = "<<m_DeltaL<<"\n";
182
183}
184
185template <class T>
187{
188 m_solver->compute(M);
189 if (m_solver->info()!=gsEigen::ComputationInfo::Success)
190 {
191 gsInfo<<"Solver error with code "<<m_solver->info()<<". See Eigen documentation on ComputationInfo \n"
192 <<gsEigen::ComputationInfo::Success<<": Success"<<"\n"
193 <<gsEigen::ComputationInfo::NumericalIssue<<": NumericalIssue"<<"\n"
194 <<gsEigen::ComputationInfo::NoConvergence<<": NoConvergence"<<"\n"
195 <<gsEigen::ComputationInfo::InvalidInput<<": InvalidInput"<<"\n";
196 throw 3;
197 }
198}
199
200template <class T>
202{
203 try
204 {
205 return m_solver->solve(F);
206 }
207 catch (...)
208 {
209 throw 3;
210 }
211}
212
213template <class T>
215{
216 // Compute Jacobian
218 m_note += "J";
219 if (!m_djacobian(U,deltaU,m))
220 throw 2;
221 this->factorizeMatrix(m);
222 return m;
223}
224
225template <class T>
227{
228 return this->_computeJacobian(U,deltaU);
229}
230
231template <class T>
233{
234 gsVector<T> DeltaU(U.rows());
235 DeltaU.setZero();
236 return this->computeJacobian(U,DeltaU);
237}
238
239template <class T>
241{
242 // Compute Jacobian
243 if (m_deltaU.rows() == 0)
244 m_deltaU = gsVector<T>::Zero(m_DeltaU.rows());
245 return this->computeJacobian(m_U + m_DeltaU, m_deltaU);
246}
247
248template <class T>
250{
251 m_deltaUbar = this->solveSystem(-m_resVec);
252}
253
254template <class T>
256{
257 m_deltaUt = this->solveSystem(m_forcing);
258}
259
260// template <class T>
261// void gsALMBase<T>::initiateStep()
262// {
263// m_converged = false;
264// m_numIterations = 0;
265// if (m_method == method::Riks)
266// initiateStepRiks();
267// else if (m_method == method::ConsistentCrisfield)
268// initiateStepConsistentCrisfield();
269// else if (m_method == method::ExplicitIterations)
270// initiateStepExplicitIterations();
271// else if (m_method == method::Crisfield)
272// initiateStepCrisfield();
273// else if (m_method == method::LoadControl)
274// initiateStepLC();
275// else
276// {
277// gsInfo<<"Error: Method unknown...\n Terminating process...\n";
278// std::terminate();
279// }
280// }
281
282// template <class T>
283// void gsALMBase<T>::predictor()
284// {
285// if (m_method == method::Riks)
286// predictorRiks();
287// else if (m_method == method::ConsistentCrisfield)
288// predictorConsistentCrisfield();
289// else if (m_method == method::ExplicitIterations)
290// predictorExplicitIterations();
291// else if (m_method == method::Crisfield)
292// predictorCrisfield();
293// else if (m_method == method::LoadControl)
294// predictorLC();
295// else
296// {
297// gsInfo<<"Error: Method unknown...\n Terminating process...\n";
298// std::terminate();
299// }
300// }
301
302// template <class T>
303// void gsALMBase<T>::iterationFinish()
304// {
305// if (m_method == method::Riks)
306// iterationFinishRiks();
307// else if (m_method == method::ConsistentCrisfield)
308// iterationFinishConsistentCrisfield();
309// else if (m_method == method::ExplicitIterations)
310// iterationFinishExplicitIterations();
311// else if (m_method == method::Crisfield)
312// iterationFinishCrisfield();
313// else if (m_method == method::LoadControl)
314// iterationFinishLC();
315// else
316// {
317// gsInfo<<"Error: Method unknown...\n Terminating process...\n";
318// std::terminate();
319// }
320// }
321
322
323// HV:
324// to do: make a stand-alone (static? const?) Unew,Lnew = step(Uold,Lold) function, which can be used in the bisection functions without problems
325// Some ideas: just make a (static) step function inside the classes and get rid of virtual sub-functions
326//
327template <class T>
330 try
331 {
332 _step();
333 m_status = gsStatus::Success;
334 }
335 catch (int errorCode)
336 {
337 if (errorCode==1)
339 else if (errorCode==2)
340 m_status = gsStatus::AssemblyError;
341 else if (errorCode==3)
342 m_status = gsStatus::SolverError;
343 else
344 m_status = gsStatus::OtherError;
345 }
346 catch (...)
347 {
348 m_status = gsStatus::OtherError;
349 }
350 return m_status;
351}
352
353template <class T>
355{
356 GISMO_ASSERT(m_initialized,"Arc-Length Method is not initialized! Call initialize()");
357
358 if (m_verbose)
359 initOutput();
360
361 m_converged = false;
362 m_numIterations = 0;
363 initiateStep();
364
365 if (m_Uguess.rows()!=0 && m_Uguess.cols()!=0 && (m_Uguess-m_U).norm()!=0 && (m_Lguess-m_L)!=0)
366 predictorGuess();
367 else
368 predictor();
369
370 computeResidual();
371 computeResidualNorms();
372
373 if (m_verbose)
374 stepOutput();
375
376 if (m_quasiNewton)
377 {
378 quasiNewtonPredictor();
379 }
380
381 m_stabilityPrev = m_stability;
382 for (m_numIterations = 1; m_numIterations < m_maxIterations; ++m_numIterations)
383 {
384 if ( (!m_quasiNewton) || ( ( m_quasiNewtonInterval>0 ) && ( m_numIterations % m_quasiNewtonInterval) < 1e-10 ) )
385 {
386 quasiNewtonIteration();
387 }
388
389 iteration();
390 computeStability(false);
391
392 computeResidual();
393 computeResidualNorms();
394 if (m_verbose)
395 stepOutput();
396
397 if ( m_residueF < m_toleranceF && m_residueU < m_toleranceU )
398 {
399 iterationFinish();
400 // Change arc length
401 if (m_adaptiveLength)
402 computeLength();
403 else
404 m_arcLength_prev = m_arcLength;
405
406 break;
407 }
408 else if (m_numIterations == m_maxIterations-1)
409 {
410 gsInfo<<"maximum iterations reached. Solution did not converge\n";
411 throw 1;
412 }
413
414 // GISMO_ERROR("maximum iterations reached. Solution did not converge");
415 }
416}
417
418// ------------------------------------------------------------------------------------------------------------
419// ---------------------------------------Singular point methods-----------------------------------------------
420// ------------------------------------------------------------------------------------------------------------
421template <class T>
422void gsALMBase<T>::_computeSingularPoint(const gsVector<T> & U, const T & L, bool switchBranch, bool jacobian, bool testPoint)
423{
424 // Determines if the point is a bifurcation point. If not, it assumes it is
425 bool test = (testPoint) ? this->_testSingularPoint(jacobian) : true;
426
427 m_U = U;
428 m_L = L;
429 if (test)
430 {
431 bool converged = false;
432 // First stage: bisection method
433 if (m_SPCompTolB != 0)
434 this->_bisectionSolve(m_U,m_L,m_SPCompTolB);
435
436 converged = this->_extendedSystemSolve(m_U, m_L, m_SPCompTolE);
437
438 if (switchBranch && (converged || m_SPfail==1))
439 this->switchBranch();
440
441 // here we assume that the stability of the singular point is
442 // equal to the one of the previous point...
443 // to avoid the algorithm to find a singular point again
444 m_stability = m_stabilityPrev;
445 }
446}
447
448// tolB and switchBranch will be defaulted
449template <class T>
450gsStatus gsALMBase<T>::computeSingularPoint(const gsVector<T> & U, const T & L, bool switchBranch, bool jacobian, bool testPoint)
451{
452 try
453 {
454 this->_computeSingularPoint(U,L,switchBranch,jacobian, testPoint);
455 }
456 catch (int errorCode)
457 {
458 if (errorCode==1)
459 m_status = gsStatus::NotConverged;
460 else if (errorCode==2)
461 m_status = gsStatus::AssemblyError;
462 else if (errorCode==3)
463 m_status = gsStatus::SolverError;
464 else
465 m_status = gsStatus::OtherError;
466 }
467 catch (...)
468 {
469 m_status = gsStatus::OtherError;
470 }
471 return m_status;
472}
473
474template <class T>
476{
477 // First, approximate the eigenvector of the Jacobian by a few arc length iterations
478 // Initiate m_V and m_DeltaVDET
479
480 if (jacobian)
481 {
482 m_jacMat = this->computeJacobian(m_U,m_deltaU);
483 this->factorizeMatrix(m_jacMat);
484 }
485
486 m_V = gsVector<T>::Ones(m_numDof);
487 m_V.normalize();
488 factorizeMatrix(m_jacMat);
489
490 for (index_t k = 0; k<m_SPTestIt; k++)
491 {
492 m_V = this->solveSystem(m_V);
493 m_V.normalize();
494 }
495
496 T dot = (abs(m_V.dot(m_forcing)));
497 if ( (abs(dot) / m_SPTestTol > 1e-1) && (abs(dot) / m_SPTestTol < 10) )
498 {
499 gsInfo<<"Warning: the singular point test is close to its tolerance. dot/tol = "<<abs(dot)/m_SPTestTol<<" dot = "<<dot<<"\t tolerance = "<<m_SPTestTol<<"\n";
500 }
501 if (dot < m_SPTestTol) // Bifurcation point
502 {
503 if (m_verbose) {gsInfo<<"\t Bifurcation point\n";}
504 return true;
505 }
506 else // Limit point
507 {
508 if (m_verbose) {gsInfo<<"\t Limit point\n";}
509 return false;
510 }
511}
512
513template <class T>
515{
516 // Controls the singular point test with the first two arguments
517 return this->_testSingularPoint(jacobian);
518}
519
520template <class T>
522{
523 try
524 {
525 _computeStability(m_U,jacobian,shift);
526 m_status = gsStatus::Success;
527 }
528 catch (int errorCode)
529 {
530 if (errorCode==1)
531 m_status = gsStatus::NotConverged;
532 else if (errorCode==2)
533 m_status = gsStatus::AssemblyError;
534 else if (errorCode==3)
535 m_status = gsStatus::SolverError;
536 else
537 m_status = gsStatus::OtherError;
538 }
539 catch (...)
540 {
541 m_status = gsStatus::OtherError;
542 }
543 return m_status;
544}
545
546template <class T>
547void gsALMBase<T>::_computeStability(const gsVector<T> & x, bool jacobian, T shift)
548{
549 if (jacobian)
550 {
551 gsVector<T> dx = gsVector<T>::Zero(x.size());
552 m_jacMat = this->computeJacobian(x,dx);
553 this->factorizeMatrix(m_jacMat);
554 } // otherwise the jacobian is already computed (on m_U+m_DeltaU)
555
556 // gsInfo<<"x = \n"<<x.transpose()<<"\n";
557 if (m_bifurcationMethod == bifmethod::Determinant)
558 {
559 if ( auto * s = dynamic_cast<typename gsSparseSolver<T>::SimplicialLDLT*>(m_solver.get()) )
560 {
561 factorizeMatrix(m_jacMat);
562 m_stabilityVec = s->vectorD();
563 }
564 else
565 {
566 gsWarn<<"Determinant stability method only works with SimplicialLDLT solver, current solver is "<<m_options.getString("Solver")<<"\n";
567 throw 3;
568 }
569 }
570 else if (m_bifurcationMethod == bifmethod::Eigenvalue)
571 {
572 #ifdef gsSpectra_ENABLED
573 index_t number = std::min(static_cast<index_t>(std::floor(m_jacMat.cols()/5.)),10);
574 /*
575 // Without shift!
576 // This one can sometimes not converge, because spectra is better at finding large values.
577 gsSpectraSymSolver<gsSparseMatrix<T>> es(m_jacMat,number,5*number);
578 es.init();
579 es.compute(Spectra::SortRule::SmallestAlge,1000,1e-6,Spectra::SortRule::SmallestAlge);
580 GISMO_ASSERT(es.info()==Spectra::CompInfo::Successful,"Spectra did not converge!"); // Reason for not converging can be due to the value of ncv (last input in the class member), which is too low.
581 */
582
583 // With shift!
584 // This one converges easier. However, a shift must be provided!
585 gsSpectraSymShiftSolver<gsSparseMatrix<T>> es(m_jacMat,number,5*number,shift);
586 es.init();
587 es.compute(Spectra::SortRule::LargestAlge,1000,1e-6,Spectra::SortRule::SmallestAlge);
588 if (es.info()!=Spectra::CompInfo::Successful)
589 {
590 gsWarn<<"Spectra did not converge!\n"; // Reason for not converging can be due to the value of ncv (last input in the class member), which is too low.
591 throw 3;
592 }
593
594 // if (es.info()==Spectra::CompInfo::NotComputed)
595 // if (es.info()==Spectra::CompInfo::NotConverging)
596 // if (es.info()==Spectra::CompInfo::NumericalIssue)
597 // gsEigen::SelfAdjointEigenSolver< gsMatrix<T> > es(m_jacMat);
598 m_stabilityVec = es.eigenvalues();
599 #else
600 GISMO_UNUSED(shift);
601 gsEigen::SelfAdjointEigenSolver<gsMatrix<T>> es2(m_jacMat);
602 m_stabilityVec = es2.eigenvalues();
603 #endif
604 }
605 else if (m_bifurcationMethod == bifmethod::Nothing)
606 {
607 m_stabilityVec = gsVector<T>::Zero(x.size());
608 }
609 else
610 gsInfo<<"bifurcation method unknown!";
611
612 m_negatives = countNegatives(m_stabilityVec);
613 m_indicator = m_stabilityVec.colwise().minCoeff()[0]; // This is required since D does not necessarily have one column.
614 m_stability = (m_indicator < 0) ? 1 : -1;
615}
616
617template <class T>
619{
620 return (m_indicator < 0) ? 1 : -1;
621}
622
623template <class T>
625{
626 return (m_stability*m_stabilityPrev < 0) ? true : false;
627}
628
629
630template <class T>
631bool gsALMBase<T>::_extendedSystemSolve(const gsVector<T> & U, const T L, const T tol)
632{
633 m_U = U;
634 m_L = L;
635 gsInfo<<"Extended iterations --- Starting with U.norm = "<<m_U.norm()<<" and L = "<<m_L<<"\n";
636
637 m_jacMat = this->computeJacobian(m_U,m_deltaU); // Jacobian evaluated on m_U
638 this->factorizeMatrix(m_jacMat);
639 m_basisResidualKTPhi = (m_jacMat*m_V).norm();
640
641 m_DeltaV = gsVector<T>::Zero(m_numDof);
642 m_DeltaU.setZero();
643 m_DeltaL = 0.0;
644
645 m_deltaV = gsVector<T>::Zero(m_numDof);
646 m_deltaU.setZero();
647 m_deltaL = 0.0;
648 if (m_verbose)
649 _initOutputExtended();
650 for (m_numIterations = 0; m_numIterations < m_maxIterations; ++m_numIterations)
651 {
652 _extendedSystemIteration();
653 m_DeltaU += m_deltaU;
654 m_DeltaL += m_deltaL;
655 m_DeltaV += m_deltaV;
656 // m_V.normalize();
657
658 // m_resVec = m_residualFun(m_U, m_L, m_forcing);
659 // m_residue = m_resVec.norm() / ( m_L * m_forcing.norm() );
660 // m_residue = (m_jacobian(m_U).toDense()*m_V).norm() / refError;
661 m_jacMat = this->computeJacobian(m_U+m_DeltaU,m_deltaU);
662 this->factorizeMatrix(m_jacMat);
663 m_residueKTPhi = (m_jacMat*(m_V+m_DeltaV)).norm(); // /m_basisResidualKTPhi;
664 m_resVec = this->computeResidual(m_U+m_DeltaU,m_L+m_DeltaL);
665 computeResidualNorms();
666 if (m_verbose)
667 _stepOutputExtended();
668
669 // termination criteria
670 // if ( m_residueF < m_toleranceF && m_residueU < m_toleranceU && m_residueKTPhi< tol )
671 if ( m_residueKTPhi< tol )
672 {
673 m_U += m_DeltaU;
674 m_L += m_DeltaL;
675 gsInfo<<"Iterations finished. U.norm() = "<<m_U.norm()<<"\t L = "<<m_L<<"\n";
676 break;
677 }
678 if (m_numIterations == m_maxIterations-1)
679 {
680 gsInfo<<"Warning: Extended iterations did not converge! \n";
681 if (m_SPfail==1)
682 {
683 m_U += m_DeltaU;
684 m_L += m_DeltaL;
685 gsInfo<<"Iterations finished; continuing with last solution U.norm() = "<<m_U.norm()<<"\t L = "<<m_L<<"\n";
686 }
687 else
688 {
689 m_DeltaV = gsVector<T>::Zero(m_numDof);
690 m_DeltaU.setZero();
691 m_DeltaL = 0.0;
692
693 m_deltaV = gsVector<T>::Zero(m_numDof);
694 m_deltaU.setZero();
695 m_deltaL = 0.0;
696 gsInfo<<"Iterations finished. continuing with original solution U.norm() = "<<m_U.norm()<<"\t L = "<<m_L<<"\n";
697 }
698 return false;
699 }
700 }
701 return true;
702}
703
704// TODO: Optimize memory
705template <class T>
707{
708 m_resVec = this->computeResidual(m_U+m_DeltaU, m_L+m_DeltaL);
709 m_jacMat = this->computeJacobian(m_U+m_DeltaU,m_deltaU);
710 this->factorizeMatrix(m_jacMat);
711
712 // m_jacMat = m_jacobian(m_U);
713
714 m_deltaUt = this->solveSystem(m_forcing); // DeltaV1
715 m_deltaUbar = this->solveSystem(-m_resVec); // DeltaV2
716
717 real_t eps = 1e-8;
718 gsSparseMatrix<T> jacMatEps = this->computeJacobian((m_U+m_DeltaU) + eps*(m_V+m_DeltaV));
719 m_note += "J"; // mark jacobian computation
720 gsVector<T> h1 = 1/eps * ( jacMatEps * m_deltaUt ) - 1/eps * m_forcing;
721 gsVector<T> h2 = m_jacMat * (m_V+m_DeltaV) + 1/eps * ( jacMatEps * m_deltaUbar + m_resVec );
722
723 this->factorizeMatrix(m_jacMat);
724
725 m_deltaVt = this->solveSystem(-h1); // DeltaV1
726 m_deltaVbar = this->solveSystem(-h2); // DeltaV2
727
728 m_deltaL = -( ( (m_V+m_DeltaV)/(m_V+m_DeltaV).norm() ).dot(m_deltaVbar) + (m_V+m_DeltaV).norm() - 1) / ( (m_V+m_DeltaV)/(m_V+m_DeltaV).norm() ).dot( m_deltaVt );
729
730 m_deltaU = m_deltaL * m_deltaUt + m_deltaUbar;
731 // gsInfo<<"m_DeltaU = \n"<<m_DeltaU<<"\n";
732
733 m_deltaV = m_deltaL * m_deltaVt + m_deltaVbar;
734 // gsInfo<<"m_DeltaV = \n"<<m_DeltaV<<"\n";
735}
736
737template <class T>
739{
740 this->_computeStability(x,jacobian);
741 m_note += "(" + std::to_string(m_negatives) + ")";
742 return m_negatives;
743}
744
745template <class T>
747{
748 this->_computeStability(x,jacobian);
749 // return m_stabilityVec.colwise().minCoeff()[0]; // This is required since D does not necessarily have one column.
750 return m_indicator;
751}
752
753template <class T>
754bool gsALMBase<T>::_bisectionSolve(const gsVector<T> & U, const T L, const T tol)
755{
756 m_U = U;
757 m_L = L;
758
759 gsVector<T> U_old = m_U;
760 T L_old = L;
761 // Store original arc length
762 real_t dL = m_arcLength;
763 bool adaptiveBool = m_adaptiveLength;
764 m_adaptiveLength = false;
765 bool converged = false;
766
767 T referenceError = _bisectionTerminationFunction(m_U, true);
768
769 gsInfo<<"Bisection iterations --- Starting with U.norm = "<<m_U.norm()<<" and L = "<<m_L<<"; Reference error = "<<referenceError<<"\n";
770
771 index_t fa, fb;
772 fa = _bisectionObjectiveFunction(m_U, false); // jacobian is already computed in the termination function
773
774 // m_arcLength = m_arcLength/2.0;
775
776 // Store termination function value of initial solution U
777 // T termU = _bisectionTerminationFunction(m_U);
778 for (index_t k = 1; k < m_maxIterations; ++k)
779 {
780 if (m_verbose) gsInfo<<"\t bisection iteration "<<k<<"\t arc length = "<<m_arcLength<<"; ";
781
782 // Reset start point
783 m_U = U_old;
784 m_L = L_old;
785
786 gsInfo<<"From U.norm = "<<m_U.norm()<<" and L = "<<m_L<<"\n";
787
788 // Make an arc length step; m_U and U_old and m_L and L_old are different
789 gsStatus status = step();
790
791 // Objective function on new point
792 fb = _bisectionObjectiveFunction(m_U, true); // m_u is the new solution
793
794 if (fa==fb && status==gsStatus::Success)
795 {
796 U_old = m_U;
797 L_old = m_L;
798 // Objective function on previous point
799 fa = fb;
800 }
801 else
802 m_arcLength = m_arcLength/2.0; // ONLY IF fa==fb??
803
804 // termination criteria
805 // relative error
806 T term = _bisectionTerminationFunction(m_U, false); // jacobian on the new point already computed
807
808 if (m_verbose) gsInfo<<"\t Finished. Relative error = "<< abs(term/referenceError)<<"\t"<<" obj.value = "<<term<<"\n";
809
810 // if ( m_arcLength < tol && abs(term) < tol )
811 if ( abs(term/referenceError) < tol )
812 {
813 converged = true;
814 break;
815 }
816
817 }
818 // Reset arc length
819 m_arcLength = dL;
820 m_adaptiveLength = adaptiveBool;
821
822 // Compute eigenvector
823 m_V = gsVector<T>::Ones(m_numDof);
824 m_V.normalize();
825 gsVector<T> Vold = gsVector<T>::Zero(m_numDof);
826 for (index_t k = 0; k<5; k++)
827 {
828 this->factorizeMatrix(m_jacMat);
829
830 m_V = this->solveSystem(m_V);
831 m_V.normalize();
832 if ( (m_V-Vold).norm() < m_tolerance )
833 {
834 converged = true;
835 break;
836 }
837 Vold = m_V;
838 }
839 return converged;
840}
841
842template <class T>
844{
845 m_V.normalize();
846 real_t lenPhi = m_V.norm();
847 real_t xi = lenPhi/m_tau;
848 // gsInfo<<xi<<"\n";
849 m_DeltaU = xi*m_V;
850 m_U = m_U + m_DeltaU;
851
852 m_DeltaLold = 0.0;
853 m_DeltaUold.setZero();
854 // m_DeltaLold = m_DeltaL;
855 // m_DeltaUold = m_DeltaU;
856}
857
858// ------------------------------------------------------------------------------------------------------------
859// ---------------------------------------Output functions-----------------------------------------------------
860// ------------------------------------------------------------------------------------------------------------
861
862template <class T>
864{
865 gsInfo<<"\t";
866 gsInfo<<std::setw(4)<<std::left<<"It.";
867 gsInfo<<std::setw(17)<<std::left<<"Res. F";
868 gsInfo<<std::setw(17)<<std::left<<"|dU|/|Du|";
869 gsInfo<<std::setw(17)<<std::left<<"dL/DL";
870 gsInfo<<std::setw(17)<<std::left<<"K_T * φ";
871 gsInfo<<std::setw(17)<<std::left<<"|U|";
872 gsInfo<<std::setw(17)<<std::left<<"|φ|";
873 gsInfo<<std::setw(17)<<std::left<<"L";
874 gsInfo<<std::setw(17)<<std::left<<"|DU|";
875 gsInfo<<std::setw(17)<<std::left<<"|Dφ|";
876 gsInfo<<std::setw(17)<<std::left<<"DL";
877 gsInfo<<std::setw(17)<<std::left<<"|dU|";
878 gsInfo<<std::setw(17)<<std::left<<"|dφ|";
879 gsInfo<<std::setw(17)<<std::left<<"dL";
880 gsInfo<<std::setw(17)<<std::left<<"Dmin";
881 gsInfo<<std::setw(17)<<std::left<<"note";
882 gsInfo<<"\n";
883
884 m_note = "";
885}
886
887template <class T>
889{
890 gsInfo<<"\t";
891 gsInfo<<std::setw(4)<<std::left<<m_numIterations;
892 gsInfo<<std::setw(17)<<std::left<<m_residueF;
893 gsInfo<<std::setw(17)<<std::left<<m_residueU;
894 gsInfo<<std::setw(17)<<std::left<<m_residueL;
895 gsInfo<<std::setw(17)<<std::left<<m_residueKTPhi;
896 gsInfo<<std::setw(17)<<std::left<<(m_U).norm();
897 gsInfo<<std::setw(17)<<std::left<<(m_V).norm();
898 gsInfo<<std::setw(17)<<std::left<<(m_L);
899 gsInfo<<std::setw(17)<<std::left<<m_DeltaU.norm();
900 gsInfo<<std::setw(17)<<std::left<<m_DeltaV.norm();
901 gsInfo<<std::setw(17)<<std::left<<m_DeltaL;
902 gsInfo<<std::setw(17)<<std::left<<m_deltaU.norm();
903 gsInfo<<std::setw(17)<<std::left<<m_deltaV.norm();
904 gsInfo<<std::setw(17)<<std::left<<m_deltaL;
905 gsInfo<<std::setw(17)<<std::left<<_bisectionTerminationFunction(m_U,false);
906 gsInfo<<std::setw(17)<<std::left<<m_note;
907 gsInfo<<"\n";
908
909 m_note = "";
910}
911
912} // namespace gismo
Performs the arc length method to solve a nonlinear system of equations.
Definition gsALMBase.h:38
virtual void computeLength()
Compute the adaptive arc-length.
Definition gsALMBase.hpp:106
virtual gsVector< T > solveSystem(const gsVector< T > &F)
Solve the system with right-hand side F.
Definition gsALMBase.hpp:201
virtual void _computeStability(const gsVector< T > &x, bool jacobian=true, T shift=-1e2)
See computeStability.
Definition gsALMBase.hpp:547
virtual void _initOutputExtended()
Initialize the output for extended iterations.
Definition gsALMBase.hpp:863
virtual void defaultOptions()
Set default options.
Definition gsALMBase.hpp:23
virtual bool stabilityChange() const
Checks if the stability of the system changed since the previously known solution.
Definition gsALMBase.hpp:624
virtual void _computeSingularPoint(const gsVector< T > &U, const T &L, bool switchBranch=false, bool jacobian=false, bool testPoint=true)
See computeSingularPoint.
Definition gsALMBase.hpp:422
virtual void _extendedSystemIteration()
Perform an extended system iteration.
Definition gsALMBase.hpp:706
virtual bool _extendedSystemSolve(const gsVector< T > &U, const T L, const T tol)
Perform an extended system solve to find a singular point.
Definition gsALMBase.hpp:631
virtual void init(bool stability)
Initialize the solver.
Definition gsALMBase.hpp:96
virtual gsSparseMatrix< T > _computeJacobian(const gsVector< T > &U, const gsVector< T > &dU)
Compute the jacobian matrix.
Definition gsALMBase.hpp:214
virtual gsStatus step()
Perform one arc-length step.
Definition gsALMBase.hpp:328
virtual T _bisectionTerminationFunction(const gsVector< T > &x, bool jacobian=true)
Returns the termination function for the bisection method given solution x.
Definition gsALMBase.hpp:746
virtual void computeResidualNorms()
Compute the residual error norms.
Definition gsALMBase.hpp:151
virtual void computeUt()
Compute .
Definition gsALMBase.hpp:255
virtual bool _testSingularPoint(bool jacobian=false)
Tests if a point is a bifurcation point.
Definition gsALMBase.hpp:475
virtual T resetLength()
Reset the length.
Definition gsALMBase.hpp:127
virtual void _step()
Implementation of step.
Definition gsALMBase.hpp:354
virtual gsStatus computeStability(bool jacobian=true, T shift=-1e2)
Calculates the stability of the solution x.
Definition gsALMBase.hpp:521
virtual index_t _bisectionObjectiveFunction(const gsVector< T > &x, bool jacobian=true)
Returns the objective function for the bisection method given solution x.
Definition gsALMBase.hpp:738
virtual gsStatus computeSingularPoint(const gsVector< T > &U, const T &L, bool switchBranch=false, bool jacobian=false, bool testPoint=true)
Calculates the singular point.
Definition gsALMBase.hpp:450
virtual T reduceLength(T fac=0.5)
Reduce the length by multiplication with a factor fac.
Definition gsALMBase.hpp:120
virtual void _stepOutputExtended()
Step output for extended iterations.
Definition gsALMBase.hpp:888
virtual void computeUbar()
Compute .
Definition gsALMBase.hpp:249
virtual void factorizeMatrix(const gsSparseMatrix< T > &M)
Factorize the matrix M.
Definition gsALMBase.hpp:186
virtual void switchBranch()
Switches branches.
Definition gsALMBase.hpp:843
virtual bool _bisectionSolve(const gsVector< T > &U, const T L, const T tol)
Perform a bisection system solve to find a singular point.
Definition gsALMBase.hpp:754
virtual index_t stability() const
Computes the stability: -1 if unstable, +1 if stable.
Definition gsALMBase.hpp:618
virtual void getOptions()
Apply options.
Definition gsALMBase.hpp:57
virtual bool isBifurcation(bool jacobian=false)
Returns true if the point is a bifurcation.
Definition gsALMBase.hpp:514
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
Shifted Eigenvalue solver for real symmetric matrices.
Definition gsSpectra.h:362
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 GISMO_UNUSED(x)
Definition gsDebug.h:112
#define gsInfo
Definition gsDebug.h:43
#define GISMO_ASSERT(cond, message)
Definition gsDebug.h:89
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.