G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsStaticOpt.h
Go to the documentation of this file.
1
16
17#pragma once
18namespace gismo
19{
20
38template <typename T>
40{
41protected:
42
43 typedef gsOptProblem<T> Base;
44 typedef typename gsStructuralAnalysisOps<T>::Residual_t Residual_t;
45 typedef typename gsStructuralAnalysisOps<T>::ALResidual_t ALResidual_t;
46
47public:
56 :
57 m_residualFun(residualFun),
58 m_L(L) // unused, in fact, but needs to be assigned
59 {
61 m_curDesign.resize(numDesignVars,1);
62 m_curDesign.setZero();
63 }
64
73 :
74 m_ALresidualFun(ALresidualFun),
75 m_L(L)
76 {
78 m_residualFun = [this](gsVector<T> const & x, gsVector<T> & result) -> bool
79 {
80 return m_ALresidualFun(x,m_L,result);
81 };
82 m_curDesign.resize(numDesignVars,1);
83 m_curDesign.setZero();
84 }
85
86public:
93 T evalObj( const gsAsConstVector<T> & u ) const
94 {
95 gsVector<T> result;
96 m_residualFun(u,result);
97 return result.norm();
98 }
99
106 void gradObj_into( const gsAsConstVector<T> & u, gsAsVector<T> & result) const
107 {
108 gsVector<T> tmp;
109 result.resize(u.rows());
110 m_residualFun(u,tmp);
111 result = -tmp;
112 }
113
114private:
115
116 Residual_t m_residualFun;
117 ALResidual_t m_ALresidualFun;
118 const T & m_L;
119 // Lastly, we forward the memebers of the base clase gsOptProblem
123
126
129
130 using Base::m_conJacRows;
131 using Base::m_conJacCols;
132
133 using Base::m_curDesign;
134};
135
148template <class T, class Optimizer = gsGradientDescent<T>>
149class gsStaticOpt : public gsStaticBase<T>
150{
151protected:
152
153 typedef gsStaticBase<T> Base;
154
155 typedef typename Base::Residual_t Residual_t;
156 typedef typename Base::ALResidual_t ALResidual_t;
157
158public:
159
166 gsStaticOpt( const Residual_t &Residual, const index_t numDofs )
167 :
168 m_optimizer(&m_optProblem),
169 m_optProblem(Residual,m_L,numDofs)
170 {
171 m_dofs = numDofs;
172 this->_init();
173 }
174
181 gsStaticOpt( const ALResidual_t & ALResidual, const index_t numDofs )
182 :
183 m_optimizer(&m_optProblem),
184 m_optProblem(ALResidual,m_L,numDofs)
185 {
186 m_L = 1.0;
187 m_dofs = numDofs;
188 this->_init();
189 }
190
191public:
197 gsOptionList & optimizerOptions() { return m_optimizer.options(); }
198
200public:
202 gsStatus solve() override;
203
204 // /// See \ref gsStaticBase
205 // void initialize() override;
206
207 // /// See \ref gsStaticBase
208 // void initOutput() override;
209
210 // /// See \ref gsStaticBase
211 // void stepOutput(index_t k) override;
212
214 void defaultOptions() override;
215
216 // /// See \ref gsStaticBase
217 // void reset() override;
218
220 void getOptions() override;
221
227
233
234protected:
236 void _solve();
238 void _init();
239
240
241protected:
242 Optimizer m_optimizer;
243 gsOptProblemStatic<T> m_optProblem;
244
245 // // Solution
246 using Base::m_U;
247 using Base::m_DeltaU;
248 // using Base::m_deltaU;
249
250 using Base::m_L;
251 // using Base::m_DeltaL;
252 // using Base::m_deltaL;
253
254 // Iterations
255 using Base::m_numIterations;
256 using Base::m_maxIterations;
257
258 // // Residuals
259 // using Base::m_R;
260
261 // // Tolerances
262 using Base::m_tolF;
263 using Base::m_tolU;
264
265 // Options
266 using Base::m_options;
267 using Base::m_verbose;
268
269 // DoFs
270 using Base::m_dofs;
271
272 // Headstart
273 using Base::m_headstart;
274
275 // Solver status
276 using Base::m_status;
277};
278
279
280
281} //namespace
282
283#ifndef GISMO_BUILD_LIB
284#include GISMO_HPP_HEADER(gsStaticOpt.hpp)
285#endif
Creates a mapped object or data pointer to a const vector without copying data.
Definition gsAsMatrix.h:285
Creates a mapped object or data pointer to a vector without copying data.
Definition gsAsMatrix.h:239
A class representing a static optimization problem.
Definition gsStaticOpt.h:40
T evalObj(const gsAsConstVector< T > &u) const
Evaluates the objective function.
Definition gsStaticOpt.h:93
void gradObj_into(const gsAsConstVector< T > &u, gsAsVector< T > &result) const
Computes the gradient of the objective function.
Definition gsStaticOpt.h:106
gsOptProblemStatic(const typename gsStructuralAnalysisOps< T >::Residual_t &residualFun, const T &L, index_t numDesignVars)
Constructor for gsOptProblemStatic.
Definition gsStaticOpt.h:55
gsOptProblemStatic(const typename gsStructuralAnalysisOps< T >::ALResidual_t &ALresidualFun, const T &L, index_t numDesignVars)
Constructor for gsOptProblemStatic with arc-length residual function.
Definition gsStaticOpt.h:72
int m_numDesignVars
Number of design variables.
Definition gsOptProblem.h:171
gsMatrix< T > m_curDesign
Current design variables (and starting point )
Definition gsOptProblem.h:198
Class defining an optimization problem.
Definition gsOptProblem.h:25
int m_numConstraints
Number of constraints.
Definition gsOptProblem.h:174
int numDesignVars() const
Callback function is executed after every iteration. Returning false causes premature termination of ...
Definition gsOptProblem.h:119
gsVector< T > m_conUpperBounds
Upper bounds for the constraints.
Definition gsOptProblem.h:189
gsVector< T > m_conLowerBounds
Lower bounds for the constraints.
Definition gsOptProblem.h:186
gsVector< T > m_desLowerBounds
Lower bounds for the design variables.
Definition gsOptProblem.h:180
int m_numConJacNonZero
Number of nonzero entries in the Constraint Jacobian.
Definition gsOptProblem.h:177
std::vector< index_t > m_conJacRows
Constraint Jacobian non-zero entries rows.
Definition gsOptProblem.h:192
int m_numDesignVars
Number of design variables.
Definition gsOptProblem.h:171
std::vector< index_t > m_conJacCols
Constraint Jacobian non-zero entries columns.
Definition gsOptProblem.h:195
gsMatrix< T > m_curDesign
Current design variables (and starting point )
Definition gsOptProblem.h:198
gsVector< T > m_desUpperBounds
Upper bounds for the design variables.
Definition gsOptProblem.h:183
Class which holds a list of parameters/options, and provides easy access to them.
Definition gsOptionList.h:33
Sparse matrix class, based on gsEigen::SparseMatrix.
Definition gsSparseMatrix.h:139
Base class for static solvers.
Definition gsStaticBase.h:38
virtual index_t numDofs()
Returns the number of DoFs of the system.
Definition gsStaticBase.h:132
Static solver using the Dynamic Relaxation method.
Definition gsStaticOpt.h:150
gsStaticOpt(const ALResidual_t &ALResidual, const index_t numDofs)
Constructor for gsStaticOpt with arc-length residual function.
Definition gsStaticOpt.h:181
void _init()
Initializes the method.
Definition gsStaticOpt.hpp:40
gsStaticOpt(const Residual_t &Residual, const index_t numDofs)
Constructor for gsStaticOpt.
Definition gsStaticOpt.h:166
void defaultOptions() override
See gsStaticBase.
Definition gsStaticOpt.hpp:34
void getOptions() override
See gsStaticBase.
Definition gsStaticOpt.hpp:55
gsStatus solve() override
gsStaticBase base functions
Definition gsStaticOpt.hpp:92
T indicator(const gsSparseMatrix< T > &, T)
Returns the stability indicator.
Definition gsStaticOpt.h:223
gsOptionList & optimizerOptions()
Returns the optimizer options.
Definition gsStaticOpt.h:197
void _solve()
See solve()
gsVector< T > stabilityVec(const gsSparseMatrix< T > &, T)
Returns the stability vector.
Definition gsStaticOpt.h:229
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
#define index_t
Definition gsConfig.h:32
#define GISMO_NO_IMPLEMENTATION
Definition gsDebug.h:129
Provides declaration of the gradient descent method.
Base class for static solvers.
The G+Smo namespace, containing all definitions for the library.
gsStatus
Definition gsStructuralAnalysisTypes.h:21
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 &, gsVector< T > &)> Residual_t
Residual, Fint-Fext.
Definition gsStructuralAnalysisTypes.h:68