G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsEigenProblemBase.hpp
Go to the documentation of this file.
1
14#include <typeinfo>
15#pragma once
16
17namespace gismo
18{
19
20template <class T>
21gsStatus gsEigenProblemBase<T>::compute()
22{
23 if (m_status==gsStatus::AssemblyError)
24 return m_status;
25
26 bool verbose = m_options.getSwitch("verbose");
27 if (verbose) { gsInfo<<"Solving eigenvalue problem" ; }
28 try
29 {
30 T shift = m_options.getReal("shift");
31 if (shift!=0.0)
32 m_eigSolver.compute(m_A-shift*m_B,m_B);
33 else
34 m_eigSolver.compute(m_A,m_B);
35
36 if (verbose) { gsInfo<<"." ; }
37 m_values = m_eigSolver.eigenvalues();
38 m_values.array() += shift;
39 if (verbose) { gsInfo<<"." ; }
40 m_vectors = m_eigSolver.eigenvectors();
41 if (verbose) { gsInfo<<"." ; }
42 if (verbose) { gsInfo<<"Finished\n" ; }
43 m_status = gsStatus::Success;
44 }
45 catch (...)
46 {
47 m_status = gsStatus::SolverError;
48 }
49 return m_status;
50};
51
52
53template <class T>
54gsStatus gsEigenProblemBase<T>::computeSparse(const index_t number)
55{
56 if (m_status==gsStatus::AssemblyError)
57 return m_status;
58
59 #ifdef gsSpectra_ENABLED
60 if (m_options.getInt("solver")==0)
61 return computeSparse_impl<Spectra::GEigsMode::Cholesky>(number);
62 else if (m_options.getInt("solver")==1)
63 return computeSparse_impl<Spectra::GEigsMode::RegularInverse>(number);
64 else if (m_options.getInt("solver")==2)
65 return computeSparse_impl<Spectra::GEigsMode::ShiftInvert>(number);
66 else if (m_options.getInt("solver")==3)
67 return computeSparse_impl<Spectra::GEigsMode::Buckling>(number);
68 else if (m_options.getInt("solver")==4)
69 return computeSparse_impl<Spectra::GEigsMode::Cayley>(number);
70 else
71 {
72 gsWarn<<"Solver index "<<m_options.getInt("solver")<<" unknown.\n";
73 m_status=gsStatus::NotStarted;
74 return m_status;
75 }
76 #else
77 GISMO_UNUSED(number);
78 gsWarn<<"Sparse solver is not implemented without gsSpectra. Please compile gismo with Spectra.\n";
79 m_status=gsStatus::NotStarted;
81 #endif
82};
83
84#ifdef gsSpectra_ENABLED
85template <class T>
86template< Spectra::GEigsMode _GEigsMode>
87typename std::enable_if<_GEigsMode==Spectra::GEigsMode::Cholesky ||
88 _GEigsMode==Spectra::GEigsMode::RegularInverse
89 ,
90 gsStatus>::type
91gsEigenProblemBase<T>::computeSparse_impl(index_t number)
92{
93 bool verbose = m_options.getSwitch("verbose");
94 T shift = m_options.getReal("shift");
95
96 Spectra::SortRule selectionRule = static_cast<Spectra::SortRule>(m_options.getInt("selectionRule"));
97 Spectra::SortRule sortRule = static_cast<Spectra::SortRule>(m_options.getInt("sortRule"));
98
99 index_t ncvFac = m_options.getInt("ncvFac");
100 T tol = m_options.getReal("tolerance");
101 if (verbose) { gsInfo<<"Solving eigenvalue problem" ; }
102
103 gsSparseMatrix<T> Atmp;
104 if (shift!=0.0)
105 Atmp = m_A-shift*m_B;
106 else
107 Atmp = m_A;
108
109 gsSpectraGenSymSolver<gsSparseMatrix<T>,_GEigsMode> solver(Atmp,m_B,number,ncvFac*number);
110
111 if (verbose) { gsInfo<<"." ; }
112 solver.init();
113 if (verbose) { gsInfo<<"." ; }
114 solver.compute(selectionRule,1000,tol,sortRule);
115
116 if (solver.info()==Spectra::CompInfo::Successful)
117 {
118 gsDebug<<"Spectra converged in "<<solver.num_iterations()<<" iterations and with "<<solver.num_operations()<<"operations. \n";
119 if (verbose) { gsInfo<<"." ; }
120 m_values = solver.eigenvalues();
121 m_values.array() += shift;
122 if (verbose) { gsInfo<<"." ; }
123 m_vectors = solver.eigenvectors();
124 if (verbose) { gsInfo<<"Finished\n" ; }
125
126 m_status = gsStatus::Success;
127 }
128 else if (solver.info()==Spectra::CompInfo::NotConverging)
129 {
130 gsWarn<<"Spectra did not converge! Error code: NotConverging\n";
131 m_status = gsStatus::NotConverged;
132 }
133 else if (solver.info()==Spectra::CompInfo::NumericalIssue)
134 {
135 gsWarn<<"Spectra did not converge! Error code: NumericalIssue\n";
136 m_status = gsStatus::SolverError;
137 }
138 else if (solver.info()==Spectra::CompInfo::NotComputed)
139 {
140 gsWarn<<"Spectra did not converge! Error code: NotComputed\n";
141 m_status = gsStatus::SolverError;
142 }
143 else
144 {
145 gsWarn<<"No error code known\n";
146 m_status = gsStatus::OtherError;
147 }
148 return m_status;
149}
150#endif
151
152#ifdef gsSpectra_ENABLED
153template <class T>
154template< Spectra::GEigsMode _GEigsMode>
155typename std::enable_if<_GEigsMode==Spectra::GEigsMode::ShiftInvert ||
156 _GEigsMode==Spectra::GEigsMode::Buckling ||
157 _GEigsMode==Spectra::GEigsMode::Cayley
158 ,
159 gsStatus>::type
160gsEigenProblemBase<T>::computeSparse_impl(index_t number)
161{
162 bool verbose = m_options.getSwitch("verbose");
163 T shift = m_options.getReal("shift");
164
165 Spectra::SortRule selectionRule = static_cast<Spectra::SortRule>(m_options.getInt("selectionRule"));
166 Spectra::SortRule sortRule = static_cast<Spectra::SortRule>(m_options.getInt("sortRule"));
167
168 if (selectionRule == Spectra::SortRule::SmallestMagn ) // =4
169 gsWarn<<"Selection Rule 'SmallestMagn' is selected, but for ShiftInvert, Buckling and Cayley it is advised to use 'LargestMagn'!\n";
170 if (selectionRule == Spectra::SortRule::SmallestAlge ) // =7
171 gsWarn<<"Selection Rule 'SmallestAlge' is selected, but for ShiftInvert, Buckling and Cayley it is advised to use 'LargestMagn'!\n";
172
173 index_t ncvFac = m_options.getInt("ncvFac");
174 T tol = m_options.getReal("tolerance");
175 if (verbose) { gsInfo<<"Solving eigenvalue problem" ; }
176 gsSpectraGenSymShiftSolver<gsSparseMatrix<T>,_GEigsMode> solver(m_A,m_B,number,ncvFac*number,shift);
177 if (verbose) { gsInfo<<"." ; }
178 solver.init();
179 if (verbose) { gsInfo<<"." ; }
180 solver.compute(selectionRule,1000,tol,sortRule);
181
182 if (solver.info()==Spectra::CompInfo::Successful)
183 {
184 gsDebug<<"Spectra converged in "<<solver.num_iterations()<<" iterations and with "<<solver.num_operations()<<"operations. \n";
185 if (verbose) { gsInfo<<"." ; }
186 m_values = solver.eigenvalues();
187 // m_values.array() += shift;
188 if (verbose) { gsInfo<<"." ; }
189 m_vectors = solver.eigenvectors();
190 if (verbose) { gsInfo<<"Finished\n" ; }
191
192 m_status = gsStatus::Success;
193 }
194 else if (solver.info()==Spectra::CompInfo::NotConverging)
195 {
196 gsWarn<<"Spectra did not converge! Error code: NotConverging\n";
197 m_status = gsStatus::NotConverged;
198 }
199 else if (solver.info()==Spectra::CompInfo::NumericalIssue)
200 {
201 gsWarn<<"Spectra did not converge! Error code: NumericalIssue\n";
202 m_status = gsStatus::SolverError;
203 }
204 else if (solver.info()==Spectra::CompInfo::NotComputed)
205 {
206 gsWarn<<"Spectra did not converge! Error code: NotComputed\n";
207 m_status = gsStatus::SolverError;
208 }
209 else
210 {
211 gsWarn<<"No error code known\n";
212 m_status = gsStatus::OtherError;
213 }
214 return m_status;
215};
216#endif
217
218template <class T>
219gsStatus gsEigenProblemBase<T>::computePower()
220{
221 if (m_status==gsStatus::AssemblyError)
222 return m_status;
223
224 bool verbose = m_options.getSwitch("verbose");
225 if (verbose) { gsInfo<<"Solving eigenvalue problem" ; }
226 gsMatrix<T> D = m_A.toDense().inverse() * (m_B);
227
228 gsVector<T> v(D.cols());
229 v.setOnes();
230 gsVector<T> v_old(D.cols());
231 v_old.setZero();
232
233 index_t kmax = 100;
234 real_t error,tol = m_options.getReal("tolerance");
235 for (index_t k=0; k!=kmax; k++)
236 {
237 v = D*v;
238 v.normalize();
239
240 error = (v-v_old).norm();
241
242 if ( error < tol )
243 {
244 m_status = gsStatus::Success;
245 break;
246 }
247 else if (k==kmax-1)
248 m_status = gsStatus::NotConverged;
249
250 v_old = v;
251 }
252
253 m_vectors = v;
254 m_values = (v.transpose() * v) / (v.transpose() * D * v);
255
256 if (verbose) { gsInfo<<"Finished\n" ; }
257 return m_status;
258};
259
260template <class T>
261std::vector<std::pair<T,gsMatrix<T>> > gsEigenProblemBase<T>::makeMode(int k) const
262{
263 std::vector<std::pair<T,gsMatrix<T>> > mode;
264 mode.push_back( std::make_pair( m_values.at(k), m_vectors.col(k) ) );
265 return mode;
266};
267
268} // namespace gismo
#define index_t
Definition gsConfig.h:32
#define gsDebug
Definition gsDebug.h:61
#define gsWarn
Definition gsDebug.h:50
#define GISMO_UNUSED(x)
Definition gsDebug.h:112
#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.