21gsStatus gsEigenProblemBase<T>::compute()
26 bool verbose = m_options.getSwitch(
"verbose");
27 if (verbose) {
gsInfo<<
"Solving eigenvalue problem" ; }
30 T shift = m_options.getReal(
"shift");
32 m_eigSolver.compute(m_A-shift*m_B,m_B);
34 m_eigSolver.compute(m_A,m_B);
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" ; }
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);
72 gsWarn<<
"Solver index "<<m_options.getInt(
"solver")<<
" unknown.\n";
78 gsWarn<<
"Sparse solver is not implemented without gsSpectra. Please compile gismo with Spectra.\n";
84#ifdef gsSpectra_ENABLED
86template< Spectra::GEigsMode _GEigsMode>
87typename std::enable_if<_GEigsMode==Spectra::GEigsMode::Cholesky ||
88 _GEigsMode==Spectra::GEigsMode::RegularInverse
91gsEigenProblemBase<T>::computeSparse_impl(
index_t number)
93 bool verbose = m_options.getSwitch(
"verbose");
94 T shift = m_options.getReal(
"shift");
96 Spectra::SortRule selectionRule =
static_cast<Spectra::SortRule
>(m_options.getInt(
"selectionRule"));
97 Spectra::SortRule sortRule =
static_cast<Spectra::SortRule
>(m_options.getInt(
"sortRule"));
99 index_t ncvFac = m_options.getInt(
"ncvFac");
100 T tol = m_options.getReal(
"tolerance");
101 if (verbose) {
gsInfo<<
"Solving eigenvalue problem" ; }
103 gsSparseMatrix<T> Atmp;
105 Atmp = m_A-shift*m_B;
109 gsSpectraGenSymSolver<gsSparseMatrix<T>,_GEigsMode> solver(Atmp,m_B,number,ncvFac*number);
111 if (verbose) {
gsInfo<<
"." ; }
113 if (verbose) {
gsInfo<<
"." ; }
114 solver.compute(selectionRule,1000,tol,sortRule);
116 if (solver.info()==Spectra::CompInfo::Successful)
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" ; }
128 else if (solver.info()==Spectra::CompInfo::NotConverging)
130 gsWarn<<
"Spectra did not converge! Error code: NotConverging\n";
133 else if (solver.info()==Spectra::CompInfo::NumericalIssue)
135 gsWarn<<
"Spectra did not converge! Error code: NumericalIssue\n";
138 else if (solver.info()==Spectra::CompInfo::NotComputed)
140 gsWarn<<
"Spectra did not converge! Error code: NotComputed\n";
145 gsWarn<<
"No error code known\n";
152#ifdef gsSpectra_ENABLED
154template< Spectra::GEigsMode _GEigsMode>
155typename std::enable_if<_GEigsMode==Spectra::GEigsMode::ShiftInvert ||
156 _GEigsMode==Spectra::GEigsMode::Buckling ||
157 _GEigsMode==Spectra::GEigsMode::Cayley
160gsEigenProblemBase<T>::computeSparse_impl(
index_t number)
162 bool verbose = m_options.getSwitch(
"verbose");
163 T shift = m_options.getReal(
"shift");
165 Spectra::SortRule selectionRule =
static_cast<Spectra::SortRule
>(m_options.getInt(
"selectionRule"));
166 Spectra::SortRule sortRule =
static_cast<Spectra::SortRule
>(m_options.getInt(
"sortRule"));
168 if (selectionRule == Spectra::SortRule::SmallestMagn )
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 )
171 gsWarn<<
"Selection Rule 'SmallestAlge' is selected, but for ShiftInvert, Buckling and Cayley it is advised to use 'LargestMagn'!\n";
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<<
"." ; }
179 if (verbose) {
gsInfo<<
"." ; }
180 solver.compute(selectionRule,1000,tol,sortRule);
182 if (solver.info()==Spectra::CompInfo::Successful)
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();
188 if (verbose) {
gsInfo<<
"." ; }
189 m_vectors = solver.eigenvectors();
190 if (verbose) {
gsInfo<<
"Finished\n" ; }
194 else if (solver.info()==Spectra::CompInfo::NotConverging)
196 gsWarn<<
"Spectra did not converge! Error code: NotConverging\n";
199 else if (solver.info()==Spectra::CompInfo::NumericalIssue)
201 gsWarn<<
"Spectra did not converge! Error code: NumericalIssue\n";
204 else if (solver.info()==Spectra::CompInfo::NotComputed)
206 gsWarn<<
"Spectra did not converge! Error code: NotComputed\n";
211 gsWarn<<
"No error code known\n";
219gsStatus gsEigenProblemBase<T>::computePower()
224 bool verbose = m_options.getSwitch(
"verbose");
225 if (verbose) {
gsInfo<<
"Solving eigenvalue problem" ; }
226 gsMatrix<T> D = m_A.toDense().inverse() * (m_B);
228 gsVector<T> v(D.cols());
230 gsVector<T> v_old(D.cols());
234 real_t error,tol = m_options.getReal(
"tolerance");
235 for (
index_t k=0; k!=kmax; k++)
240 error = (v-v_old).norm();
254 m_values = (v.transpose() * v) / (v.transpose() * D * v);
256 if (verbose) {
gsInfo<<
"Finished\n" ; }
261std::vector<std::pair<T,gsMatrix<T>> > gsEigenProblemBase<T>::makeMode(
int k)
const
263 std::vector<std::pair<T,gsMatrix<T>> > mode;
264 mode.push_back( std::make_pair( m_values.at(k), m_vectors.col(k) ) );
#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
@ SolverError
Assembly problem in step.
@ AssemblyError
Assembly problem in step.
@ NotConverged
Step did not converge.
@ NotStarted
ALM has not started yet.