G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsALMRiks.h
Go to the documentation of this file.
1 
14 #pragma once
15 
17 
18 namespace gismo
19 {
20 
28 template <class T>
29 class gsALMRiks : public gsALMBase<T>
30 {
31 
32  typedef gsALMBase<T> Base;
33 
34  typedef typename Base::ALResidual_t ALResidual_t;
35  typedef typename Base::Jacobian_t Jacobian_t;
36  typedef typename Base::dJacobian_t dJacobian_t;
37 
38 public:
39 
40  using Base::setLength;
41 
42 protected:
43 
45  using Base::getOptions;
46  using Base::computeJacobian;
47  using Base::computeResidual;
49  using Base::computeUt;
50  using Base::computeUbar;
52  using Base::computeLength;
53 
54 public:
55 
57  gsALMRiks( const Jacobian_t &Jacobian,
58  const ALResidual_t&ALResidual,
59  const gsVector<T> &Force )
60  : Base(Jacobian,ALResidual,Force)
61  {
63  getOptions();
64 
65  initMethods();
66  }
67 
69  gsALMRiks( const dJacobian_t &dJacobian,
70  const ALResidual_t&ALResidual,
71  const gsVector<T> &Force )
72  : Base(dJacobian,ALResidual,Force)
73  {
75  getOptions();
76 
77  initMethods();
78  }
79 
80 public:
81  T distance(const gsVector<T>& DeltaU, const T DeltaL) const
82  {
83  return math::pow(m_phi * math::pow(m_DeltaU.norm(),2.0) + (1.0-m_phi) * math::pow(m_DeltaL,2.0),0.5);
84  }
85 
86 protected:
87 
89  void initMethods();
91  void initiateStep();
93  void iterationFinish();
94 
96  void quasiNewtonPredictor();
98  void quasiNewtonIteration();
99 
101  void predictor();
102  void predictorGuess();
104  void iteration();
105 
107  void initOutput();
109  void stepOutput();
110 
111 protected:
112 
113  // Number of degrees of freedom
114  using Base::m_numDof;
115 
116  using Base::m_jacobian;
117  using Base::m_djacobian;
118  using Base::m_residualFun;
119  using Base::m_forcing;
120 
122  using Base::m_options;
123 
125  using Base::m_numIterations;
126 
128  using Base::m_maxIterations;
129 
131  using Base::m_arcLength;
132  using Base::m_arcLength_prev;
133 
135  using Base::m_verbose;
136 
138  using Base::m_note;
139 
141  using Base::m_converged;
142 
144  using Base::m_residueF;
145 
147  using Base::m_residueU;
148 
150  using Base::m_residueL;
151 
153  using Base::m_indicator;
154  using Base::m_negatives;
155 
157  using Base::m_relax;
158 
159  // Previous update
160  using Base::m_DeltaUold;
161  using Base::m_DeltaLold;
163  using Base::m_U;
164  using Base::m_Uprev;
165  using Base::m_Uguess;
167  using Base::m_DeltaU;
169  using Base::m_deltaUbar;
171  using Base::m_deltaUt;
173  using Base::m_deltaU;
174 
176  using Base::m_L;
177  using Base::m_Lprev;
178  using Base::m_Lguess;
180  using Base::m_DeltaL;
182  using Base::m_deltaL;
184  using Base::m_deltaLs;
185 
187  using Base::m_jacMat;
188  using Base::m_detKT;
189 
191  T m_phi;
192 };
193 
194 } // namespace gismo
195 
196 #ifndef GISMO_BUILD_LIB
197 #include GISMO_HPP_HEADER(gsALMRiks.hpp)
198 #endif
virtual void defaultOptions()
Set default options.
Definition: gsALMBase.hpp:23
gsALMRiks(const dJacobian_t &dJacobian, const ALResidual_t &ALResidual, const gsVector< T > &Force)
Constructor using the jacobian that takes the solution and the solution step.
Definition: gsALMRiks.h:69
Performs the Riks arc length method to solve a nonlinear equation system.
Definition: gsALMRiks.h:29
gsVector< T > m_deltaLs
Vector with lambda updates.
Definition: gsALMBase.h:504
void iterationFinish()
See gsALMBase.
Definition: gsALMRiks.hpp:205
gsALMRiks(const Jacobian_t &Jacobian, const ALResidual_t &ALResidual, const gsVector< T > &Force)
Constructor.
Definition: gsALMRiks.h:57
virtual gsStatus computeStability(bool jacobian=true, T shift=-1e2)
Calculates the stability of the solution x.
Definition: gsALMBase.hpp:521
virtual void computeUbar()
Compute .
Definition: gsALMBase.hpp:249
void initMethods()
See gsALMBase.
Definition: gsALMRiks.hpp:23
void quasiNewtonIteration()
See gsALMBase.
Definition: gsALMRiks.hpp:48
gsVector< T > m_deltaUbar
u_bar
Definition: gsALMBase.h:491
bool m_converged
Convergence result.
Definition: gsALMBase.h:454
void initiateStep()
See gsALMBase.
Definition: gsALMRiks.hpp:87
gsVector< T > m_deltaUt
u_t
Definition: gsALMBase.h:493
virtual void computeUt()
Compute .
Definition: gsALMBase.hpp:255
T m_residueU
Displacement residuum.
Definition: gsALMBase.h:464
T m_residueF
Force residuum.
Definition: gsALMBase.h:460
T m_relax
Relaxation factor.
Definition: gsALMBase.h:479
T m_deltaL
Update of update of lambda.
Definition: gsALMBase.h:502
virtual void computeResidualNorms()
Compute the residual error norms.
Definition: gsALMBase.hpp:151
index_t m_maxIterations
Maximum number of Arc Length iterations allowed.
Definition: gsALMBase.h:421
gsVector< T > m_DeltaU
Update of displacement vector.
Definition: gsALMBase.h:489
gsSparseMatrix< T > m_jacMat
Jacobian matrix.
Definition: gsALMBase.h:510
T m_DeltaL
Update of lambdaGeneralizedSelfAdjointEigenSolver.
Definition: gsALMBase.h:500
virtual void setLength(T length)
Set arc length to length.
Definition: gsALMBase.h:118
void iteration()
See gsALMBase.
Definition: gsALMRiks.hpp:69
T m_residueL
Load residuum.
Definition: gsALMBase.h:468
void stepOutput()
See gsALMBase.
Definition: gsALMRiks.hpp:243
Base class to perform the arc length method to solve a nonlinear equation system. ...
virtual void computeLength()
Compute the adaptive arc-length.
Definition: gsALMBase.hpp:106
gsVector< T > m_U
Displacement vector (present, at previously converged point)
Definition: gsALMBase.h:487
T m_L
Lambda (present, at previously converged point)
Definition: gsALMBase.h:498
index_t m_numIterations
Number of Arc Length iterations performed.
Definition: gsALMBase.h:418
T m_arcLength
Length of the step in the u,f plane.
Definition: gsALMBase.h:427
void initOutput()
See gsALMBase.
Definition: gsALMRiks.hpp:219
void predictor()
See gsALMBase.
Definition: gsALMRiks.hpp:98
void quasiNewtonPredictor()
See gsALMBase.
Definition: gsALMRiks.hpp:38
Performs the arc length method to solve a nonlinear system of equations.
Definition: gsALMBase.h:37
T m_phi
Scaling parameter.
Definition: gsALMRiks.h:191
T m_indicator
Indicator for bifurcation.
Definition: gsALMBase.h:475
virtual void getOptions()
Apply options.
Definition: gsALMBase.hpp:57
gsVector< T > m_deltaU
Update of update of displacement vector.
Definition: gsALMBase.h:495