G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsStaticBase.h
Go to the documentation of this file.
1
14#include <typeinfo>
15
17#ifdef gsSpectra_ENABLED
18#include <gsSpectra/gsSpectra.h>
19#endif
20#include <gsIO/gsOptionList.h>
22
23#pragma once
24
25
26namespace gismo
27{
28
36template <class T>
38{
39protected:
40
41 typedef typename gsStructuralAnalysisOps<T>::Residual_t Residual_t;
42 typedef typename gsStructuralAnalysisOps<T>::ALResidual_t ALResidual_t;
43 typedef typename gsStructuralAnalysisOps<T>::Jacobian_t Jacobian_t;
44 typedef typename gsStructuralAnalysisOps<T>::dJacobian_t dJacobian_t;
45
46public:
47
48 virtual ~gsStaticBase() {};
49
51 virtual gsStatus solve() = 0;
52
54 virtual void initialize()
55 {
56 this->reset();
57 this->getOptions();
58 }
59
61 virtual void initOutput() {};
63 virtual void stepOutput(index_t /*k*/) {};
64
66 virtual void defaultOptions()
67 {
68 m_options.addReal("tol","Relative Tolerance",1e-6);
69 m_options.addReal("tolF","Residual relative tolerance",-1);
70 m_options.addReal("tolU","Solution step relative tolerance",-1);
71 m_options.addInt("maxIt","Maximum number of iterations",25);
72 m_options.addInt("verbose","Verbose output",0);
73 m_options.addInt ("BifurcationMethod","Bifurcation Identification based on: 0: Determinant; 1: Eigenvalue",stabmethod::Eigenvalue);
74 m_options.addString("Solver","Sparse linear solver", "SimplicialLDLT");
75 }
76
78 virtual void getOptions()
79 {
80 m_tolF = m_options.getReal("tolF")!=-1 ? m_options.getReal("tolF") : m_options.getReal("tol");
81 m_tolU = m_options.getReal("tolU")!=-1 ? m_options.getReal("tolU") : m_options.getReal("tol");
82 m_maxIterations = m_options.getInt("maxIt");
83 m_verbose = m_options.getInt("verbose");
84 m_stabilityMethod = m_options.getInt ("BifurcationMethod");
85 m_solver = gsSparseSolver<T>::get( m_options.askString("Solver","SimplicialLDLT") );
86 }
87
89 virtual void setOptions(gsOptionList & options) {m_options.update(options,gsOptionList::addIfUnknown); }
90
92 virtual const gsOptionList & options() const {return m_options;}
93 virtual gsOptionList & options() {return m_options;}
94
96 virtual gsVector<T> solution() const {return m_U;}
97
99 virtual gsVector<T> update() const {return m_DeltaU;}
100
102 virtual void setDisplacement(const gsVector<T> & displacement)
103 {
104 m_start = true;
105 m_U = displacement;
106 }
107
109 virtual void setLoad(const T L) { m_L = L;}
110
112 virtual void setSolution(const gsVector<T> & displacement, const T L)
113 {
114 this->setDisplacement(displacement);
115 this->setLoad(L);
116 }
117 virtual void setUpdate(const gsVector<T> & update)
118 {
119 m_headstart = true;
120 m_DeltaU = update;
121 }
122
124 virtual index_t iterations() const { return m_numIterations; }
125
127 virtual bool converged() const { return m_status==gsStatus::Success; }
128
130 virtual gsStatus status() const { return m_status; }
131
133 virtual index_t numDofs() { return m_dofs; }
134
136 virtual void reset()
137 {
138 m_U.setZero(m_dofs);
139 m_DeltaU.setZero(m_dofs);
140 m_deltaU.setZero(m_dofs);
141 m_R.setZero(m_dofs);
142 m_L = m_DeltaL = m_deltaL = 0.0;
143 m_headstart = false;
144 }
145
147 virtual T indicator(const gsSparseMatrix<T> & jacMat, T shift = -1e-2)
148 {
149 _computeStability(jacMat, shift);
150 return m_indicator;
151 }
152
154 virtual gsVector<T> stabilityVec(const gsSparseMatrix<T> & jacMat, T shift = -1e-2)
155 {
156 _computeStability(jacMat, shift);
157 return m_stabilityVec;
158 }
159
160protected:
162 virtual bool _computeStabilityDet(const gsSparseMatrix<T> & jacMat)
163 {
164 m_solver->compute(jacMat);
165 // If 1: matrix is not SPD
166 if (m_solver->info()!=gsEigen::ComputationInfo::Success)
167 {
168 gsInfo<<"Solver error with code "<<m_solver->info()<<". See Eigen documentation on ComputationInfo \n"
169 <<gsEigen::ComputationInfo::Success<<": Success"<<"\n"
170 <<gsEigen::ComputationInfo::NumericalIssue<<": NumericalIssue"<<"\n"
171 <<gsEigen::ComputationInfo::NoConvergence<<": NoConvergence"<<"\n"
172 <<gsEigen::ComputationInfo::InvalidInput<<": InvalidInput"<<"\n";
173 return false;
174 }
175
176 if ( auto * s = dynamic_cast<typename gsSparseSolver<T>::SimplicialLDLT*>(m_solver.get()) )
177 m_stabilityVec = s->vectorD();
178 return true;
179 }
180
182 virtual bool _computeStabilityEig(const gsSparseMatrix<T> & jacMat, T shift)
183 {
184#ifdef gsSpectra_ENABLED
185 index_t number = std::min(static_cast<index_t>(std::floor(jacMat.cols()/5.)),10);
186 /*
187 // Without shift!
188 // This one can sometimes not converge, because spectra is better at finding large values.
189 gsSpectraSymSolver<gsSparseMatrix<T>> es(jacMat,number,5*number);
190 es.init();
191 es.compute(Spectra::SortRule::SmallestAlge,1000,1e-6,Spectra::SortRule::SmallestAlge);
192 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.
193 */
194
195 // With shift!
196 // This one converges easier. However, a shift must be provided!
197 gsSpectraSymShiftSolver<gsSparseMatrix<T>> es(jacMat,number,5*number,shift);
198 es.init();
199 es.compute(Spectra::SortRule::LargestAlge,1000,1e-6,Spectra::SortRule::SmallestAlge);
200 if (es.info()!=Spectra::CompInfo::Successful)
201 {
202 gsWarn<<"Spectra did not converge!\n";
203 return false;
204 }
205
206 // if (es.info()==Spectra::CompInfo::NotComputed)
207 // if (es.info()==Spectra::CompInfo::NotConverging)
208 // if (es.info()==Spectra::CompInfo::NumericalIssue)
209 // gsEigen::SelfAdjointEigenSolver< gsMatrix<T> > es(jacMat);
210 m_stabilityVec = es.eigenvalues();
211#else
212 GISMO_UNUSED(shift);
213 gsEigen::SelfAdjointEigenSolver<gsMatrix<T>> es2(jacMat);
214 m_stabilityVec = es2.eigenvalues();
215#endif
216
217 m_indicator = m_stabilityVec.colwise().minCoeff()[0]; // This is required since D does not necessarily have one column.
218 return true;
219 }
220
222 virtual bool _computeStability (const gsSparseMatrix<T> & jacMat, T shift)
223 {
224 bool success = false;
225 if (m_stabilityMethod == stabmethod::Determinant)
226 success = _computeStabilityDet(jacMat);
227 else if (m_stabilityMethod == stabmethod::Eigenvalue)
228 success = _computeStabilityEig(jacMat, shift);
229 else
230 gsInfo<<"bifurcation method unknown!";
231
232 if (success)
233 m_indicator = m_stabilityVec.colwise().minCoeff()[0]; // This is required since D does not necessarily have one column
234 return success;
235 }
236
237protected:
238 mutable gsOptionList m_options;
239
240 T m_L, m_DeltaL, m_deltaL;
241 gsVector<T> m_U, m_DeltaU, m_deltaU;
242 gsVector<T> m_R;
243 T m_residual, m_residualIni, m_residualOld;
244
245 index_t m_numIterations;
246 index_t m_maxIterations;
247 bool m_start;
248 index_t m_headstart;
249
250 index_t m_verbose;
251
252 index_t m_dofs;
253
254 T m_tolF,m_tolU;
255
256 T m_indicator;
257
258 gsVector<T> m_stabilityVec;
259
260 mutable typename gsSparseSolver<T>::uPtr m_solver; // Cholesky by default
261
262 index_t m_stabilityMethod;
263
264 struct stabmethod
265 {
266 enum type
267 {
268 Determinant = 0,
269 Eigenvalue = 1,
270 };
271 };
272
273 gsStatus m_status;
274
275};
276
277
278} // namespace gismo
279
280
281#ifndef GISMO_BUILD_LIB
282#include GISMO_HPP_HEADER(gsStaticBase.hpp)
283#endif
Class which holds a list of parameters/options, and provides easy access to them.
Definition gsOptionList.h:33
void addInt(const std::string &label, const std::string &desc, const index_t &value)
Adds a option named label, with description desc and value value.
Definition gsOptionList.cpp:201
void update(const gsOptionList &other, updateType type=ignoreIfUnknown)
Updates the object using the data from other.
Definition gsOptionList.cpp:253
const index_t & getInt(const std::string &label) const
Reads value for option label from options.
Definition gsOptionList.cpp:37
Real getReal(const std::string &label) const
Reads value for option label from options.
Definition gsOptionList.cpp:44
void addReal(const std::string &label, const std::string &desc, const Real &value)
Adds a option named label, with description desc and value value.
Definition gsOptionList.cpp:211
std::string askString(const std::string &label, const std::string &value="") const
Reads value for option label from options.
Definition gsOptionList.cpp:106
void addString(const std::string &label, const std::string &desc, const std::string &value)
Adds a option named label, with description desc and value value.
Definition gsOptionList.cpp:190
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
Base class for static solvers.
Definition gsStaticBase.h:38
virtual void initOutput()
Initialize output.
Definition gsStaticBase.h:61
virtual gsStatus status() const
Returns the status.
Definition gsStaticBase.h:130
virtual const gsOptionList & options() const
Get options.
Definition gsStaticBase.h:92
virtual gsVector< T > update() const
Access the update.
Definition gsStaticBase.h:99
virtual void setDisplacement(const gsVector< T > &displacement)
Set the displacement.
Definition gsStaticBase.h:102
virtual index_t iterations() const
Returns the number of iterations.
Definition gsStaticBase.h:124
virtual bool _computeStabilityDet(const gsSparseMatrix< T > &jacMat)
Computes the stability vector using the determinant of the Jacobian.
Definition gsStaticBase.h:162
virtual bool converged() const
Returns whether the solver converged or not.
Definition gsStaticBase.h:127
virtual gsStatus solve()=0
Solve.
virtual void setLoad(const T L)
Set the load.
Definition gsStaticBase.h:109
virtual bool _computeStability(const gsSparseMatrix< T > &jacMat, T shift)
Computes the stability of the Jacobian, optionally applying a shift (if provided)
Definition gsStaticBase.h:222
virtual bool _computeStabilityEig(const gsSparseMatrix< T > &jacMat, T shift)
Computes the stability vector using the eigenvalues of the Jacobian, optionally applying a shift.
Definition gsStaticBase.h:182
virtual void initialize()
See gsStaticBase.
Definition gsStaticBase.h:54
virtual void reset()
Reset the stored solution.
Definition gsStaticBase.h:136
virtual void getOptions()
Apply the options.
Definition gsStaticBase.h:78
virtual T indicator(const gsSparseMatrix< T > &jacMat, T shift=-1e-2)
Returns the stability indicator.
Definition gsStaticBase.h:147
virtual void stepOutput(index_t)
Stepwise output.
Definition gsStaticBase.h:63
virtual index_t numDofs()
Returns the number of DoFs of the system.
Definition gsStaticBase.h:133
virtual gsVector< T > solution() const
Access the solution.
Definition gsStaticBase.h:96
virtual gsVector< T > stabilityVec(const gsSparseMatrix< T > &jacMat, T shift=-1e-2)
Returns the stability vector.
Definition gsStaticBase.h:154
virtual void setSolution(const gsVector< T > &displacement, const T L)
Set the displacement and the load.
Definition gsStaticBase.h:112
virtual void setOptions(gsOptionList &options)
Set the options from options.
Definition gsStaticBase.h:89
virtual void defaultOptions()
Get default options.
Definition gsStaticBase.h:66
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
This is the main header file that collects wrappers of Eigen for linear algebra.
Provides a list of labeled parameters/options that can be set and accessed easily.
Header file for using Spectra extension.
Provides a status object and typedefs.
The G+Smo namespace, containing all definitions for the library.
gsStatus
Definition gsStructuralAnalysisTypes.h:21
@ Success
Successful.
std::function< bool(gsVector< T > const &, gsVector< T > const &, gsSparseMatrix< T > &) > dJacobian_t
Jacobian with solution update as argument.
Definition gsStructuralAnalysisTypes.h:90
std::function< bool(gsVector< T > const &, const T, gsVector< T > &)> ALResidual_t
Arc-Length Residual, Fint-lambda*Fext.
Definition gsStructuralAnalysisTypes.h:70
std::function< bool(gsVector< T > const &, gsSparseMatrix< T > &) > Jacobian_t
Jacobian.
Definition gsStructuralAnalysisTypes.h:86
std::function< bool(gsVector< T > const &, gsVector< T > &)> Residual_t
Residual, Fint-Fext.
Definition gsStructuralAnalysisTypes.h:68