SlepcEigenSolver< T > Class Template Reference
#include <slepc_eigen_solver.h>

Public Member Functions | |
| SlepcEigenSolver () | |
| ~SlepcEigenSolver () | |
| void | clear () |
| void | init () |
| std::pair< unsigned int, unsigned int > | solve_standard (SparseMatrix< T > &matrix_A, int nev, int ncv, const double tol, const unsigned int m_its) |
| std::pair< unsigned int, unsigned int > | solve_standard (ShellMatrix< T > &shell_matrix, int nev, int ncv, const double tol, const unsigned int m_its) |
| std::pair< unsigned int, unsigned int > | solve_generalized (SparseMatrix< T > &matrix_A, SparseMatrix< T > &matrix_B, int nev, int ncv, const double tol, const unsigned int m_its) |
| std::pair< unsigned int, unsigned int > | solve_generalized (ShellMatrix< T > &matrix_A, SparseMatrix< T > &matrix_B, int nev, int ncv, const double tol, const unsigned int m_its) |
| std::pair< unsigned int, unsigned int > | solve_generalized (SparseMatrix< T > &matrix_A, ShellMatrix< T > &matrix_B, int nev, int ncv, const double tol, const unsigned int m_its) |
| std::pair< unsigned int, unsigned int > | solve_generalized (ShellMatrix< T > &matrix_A, ShellMatrix< T > &matrix_B, int nev, int ncv, const double tol, const unsigned int m_its) |
| std::pair< Real, Real > | get_eigenpair (unsigned int i, NumericVector< T > &solution_in) |
| Real | get_relative_error (unsigned int i) |
| void | attach_deflation_space (NumericVector< T > &deflation_vector) |
| bool | initialized () const |
| EigenSolverType | eigen_solver_type () const |
| EigenProblemType | eigen_problem_type () const |
| PositionOfSpectrum | position_of_spectrum () const |
| void | set_eigensolver_type (const EigenSolverType est) |
| void | set_eigenproblem_type (EigenProblemType ept) |
| void | set_position_of_spectrum (PositionOfSpectrum pos) |
Static Public Member Functions | |
| static AutoPtr< EigenSolver< T > > | build (const SolverPackage solver_package=SLEPC_SOLVERS) |
| static std::string | get_info () |
| static void | print_info () |
| static unsigned int | n_objects () |
Protected Types | |
| typedef std::map< std::string, std::pair< unsigned int, unsigned int > > | Counts |
Protected Member Functions | |
| void | increment_constructor_count (const std::string &name) |
| void | increment_destructor_count (const std::string &name) |
Protected Attributes | |
| EigenSolverType | _eigen_solver_type |
| EigenProblemType | _eigen_problem_type |
| PositionOfSpectrum | _position_of_spectrum |
| bool | _is_initialized |
Static Protected Attributes | |
| static Counts | _counts |
| static Threads::atomic < unsigned int > | _n_objects |
| static Threads::spin_mutex | _mutex |
Private Member Functions | |
| std::pair< unsigned int, unsigned int > | _solve_standard_helper (Mat mat, int nev, int ncv, const double tol, const unsigned int m_its) |
| std::pair< unsigned int, unsigned int > | _solve_generalized_helper (Mat mat_A, Mat mat_B, int nev, int ncv, const double tol, const unsigned int m_its) |
| void | set_slepc_solver_type () |
| void | set_slepc_problem_type () |
| void | set_slepc_position_of_spectrum () |
Static Private Member Functions | |
| static PetscErrorCode | _petsc_shell_matrix_mult (Mat mat, Vec arg, Vec dest) |
| static PetscErrorCode | _petsc_shell_matrix_get_diagonal (Mat mat, Vec dest) |
Private Attributes | |
| EPS | _eps |
Detailed Description
template<typename T>
class SlepcEigenSolver< T >
This class provides an interface to the SLEPc eigenvalue solver library www.grycap.upv.es/slepc/.
Definition at line 68 of file slepc_eigen_solver.h.
Member Typedef Documentation
typedef std::map<std::string, std::pair<unsigned int, unsigned int> > ReferenceCounter::Counts [protected, inherited] |
Data structure to log the information. The log is identified by the class name.
Definition at line 105 of file reference_counter.h.
Constructor & Destructor Documentation
| SlepcEigenSolver< T >::SlepcEigenSolver | ( | ) | [inline] |
Constructor. Initializes Petsc data structures
Definition at line 274 of file slepc_eigen_solver.h.
References EigenSolver< T >::_eigen_problem_type, EigenSolver< T >::_eigen_solver_type, libMeshEnums::ARNOLDI, and libMeshEnums::NHEP.
00275 { 00276 this->_eigen_solver_type = ARNOLDI; 00277 this->_eigen_problem_type = NHEP; 00278 }
| SlepcEigenSolver< T >::~SlepcEigenSolver | ( | ) | [inline] |
Destructor.
Definition at line 284 of file slepc_eigen_solver.h.
References SlepcEigenSolver< T >::clear().
00285 { 00286 this->clear (); 00287 }
Member Function Documentation
| PetscErrorCode SlepcEigenSolver< T >::_petsc_shell_matrix_get_diagonal | ( | Mat | mat, | |
| Vec | dest | |||
| ) | [inline, static, private] |
Internal function if shell matrix mode is used, this just calls the shell matrix's get_diagonal function. Required in order to use Jacobi preconditioning.
Definition at line 727 of file slepc_eigen_solver.C.
References libMesh::COMM_WORLD, and ShellMatrix< T >::get_diagonal().
Referenced by SlepcEigenSolver< T >::solve_generalized(), and SlepcEigenSolver< T >::solve_standard().
00728 { 00729 /* Get the matrix context. */ 00730 int ierr=0; 00731 void* ctx; 00732 ierr = MatShellGetContext(mat,&ctx); 00733 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00734 00735 /* Get user shell matrix object. */ 00736 const ShellMatrix<T>& shell_matrix = *static_cast<const ShellMatrix<T>*>(ctx); 00737 00738 /* Make \p NumericVector instances around the vector. */ 00739 PetscVector<T> dest_global(dest); 00740 00741 /* Call the user function. */ 00742 shell_matrix.get_diagonal(dest_global); 00743 00744 return ierr; 00745 }
| PetscErrorCode SlepcEigenSolver< T >::_petsc_shell_matrix_mult | ( | Mat | mat, | |
| Vec | arg, | |||
| Vec | dest | |||
| ) | [inline, static, private] |
Internal function if shell matrix mode is used, this just calls the shell matrix's matrix multiplication function. See PetscLinearSolver for a similar implementation.
Definition at line 705 of file slepc_eigen_solver.C.
References libMesh::COMM_WORLD, and ShellMatrix< T >::vector_mult().
Referenced by SlepcEigenSolver< T >::solve_generalized(), and SlepcEigenSolver< T >::solve_standard().
00706 { 00707 /* Get the matrix context. */ 00708 int ierr=0; 00709 void* ctx; 00710 ierr = MatShellGetContext(mat,&ctx); 00711 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00712 00713 /* Get user shell matrix object. */ 00714 const ShellMatrix<T>& shell_matrix = *static_cast<const ShellMatrix<T>*>(ctx); 00715 00716 /* Make \p NumericVector instances around the vectors. */ 00717 PetscVector<T> arg_global(arg); 00718 PetscVector<T> dest_global(dest); 00719 00720 /* Call the user function. */ 00721 shell_matrix.vector_mult(dest_global,arg_global); 00722 00723 return ierr; 00724 }
| std::pair< unsigned int, unsigned int > SlepcEigenSolver< T >::_solve_generalized_helper | ( | Mat | mat_A, | |
| Mat | mat_B, | |||
| int | nev, | |||
| int | ncv, | |||
| const double | tol, | |||
| const unsigned int | m_its | |||
| ) | [inline, private] |
Helper function that actually performs the generalized eigensolve.
Definition at line 429 of file slepc_eigen_solver.C.
References SlepcEigenSolver< T >::_eps, libMesh::COMM_WORLD, SlepcEigenSolver< T >::set_slepc_position_of_spectrum(), and SlepcEigenSolver< T >::set_slepc_problem_type().
Referenced by SlepcEigenSolver< T >::solve_generalized().
00435 { 00436 START_LOG("solve_generalized()", "SlepcEigenSolver"); 00437 00438 int ierr=0; 00439 00440 // converged eigen pairs and number of iterations 00441 int nconv=0; 00442 int its=0; 00443 00444 #ifdef DEBUG 00445 // The relative error. 00446 PetscReal error, re, im; 00447 00448 // Pointer to vectors of the real parts, imaginary parts. 00449 PetscScalar kr, ki; 00450 #endif 00451 00452 // Set operators. 00453 ierr = EPSSetOperators (_eps, mat_A, mat_B); 00454 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00455 00456 //set the problem type and the position of the spectrum 00457 set_slepc_problem_type(); 00458 set_slepc_position_of_spectrum(); 00459 00460 // Set eigenvalues to be computed. 00461 #if SLEPC_VERSION_LESS_THAN(3,0,0) 00462 ierr = EPSSetDimensions (_eps, nev, ncv); 00463 #else 00464 ierr = EPSSetDimensions (_eps, nev, ncv, PETSC_DECIDE); 00465 #endif 00466 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00467 00468 00469 // Set the tolerance and maximum iterations. 00470 ierr = EPSSetTolerances (_eps, tol, m_its); 00471 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00472 00473 // Set runtime options, e.g., 00474 // -eps_type <type>, -eps_nev <nev>, -eps_ncv <ncv> 00475 // Similar to PETSc, these options will override those specified 00476 // above as long as EPSSetFromOptions() is called _after_ any 00477 // other customization routines. 00478 ierr = EPSSetFromOptions (_eps); 00479 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00480 00481 // Solve the eigenproblem. 00482 ierr = EPSSolve (_eps); 00483 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00484 00485 // Get the number of iterations. 00486 ierr = EPSGetIterationNumber (_eps, &its); 00487 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00488 00489 // Get number of converged eigenpairs. 00490 ierr = EPSGetConverged(_eps,&nconv); 00491 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00492 00493 00494 #ifdef DEBUG 00495 // ierr = PetscPrintf(libMesh::COMM_WORLD, 00496 // "\n Number of iterations: %d\n" 00497 // " Number of converged eigenpairs: %d\n\n", its, nconv); 00498 00499 // Display eigenvalues and relative errors. 00500 ierr = PetscPrintf(libMesh::COMM_WORLD, 00501 " k ||Ax-kx||/|kx|\n" 00502 " ----------------- -----------------\n" ); 00503 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00504 00505 for(int i=0; i<nconv; i++ ) 00506 { 00507 ierr = EPSGetEigenpair(_eps, i, &kr, &ki, PETSC_NULL, PETSC_NULL); 00508 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00509 00510 ierr = EPSComputeRelativeError(_eps, i, &error); 00511 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00512 00513 #ifdef LIBMESH_USE_COMPLEX_NUMBERS 00514 re = PetscRealPart(kr); 00515 im = PetscImaginaryPart(kr); 00516 #else 00517 re = kr; 00518 im = ki; 00519 #endif 00520 00521 if (im != .0) 00522 { 00523 ierr = PetscPrintf(libMesh::COMM_WORLD," %9f%+9f i %12f\n", re, im, error); 00524 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00525 } 00526 else 00527 { 00528 ierr = PetscPrintf(libMesh::COMM_WORLD," %12f %12f\n", re, error); 00529 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00530 } 00531 } 00532 00533 ierr = PetscPrintf(libMesh::COMM_WORLD,"\n" ); 00534 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00535 #endif // DEBUG 00536 00537 STOP_LOG("solve_generalized()", "SlepcEigenSolver"); 00538 00539 // return the number of converged eigenpairs 00540 // and the number of iterations 00541 return std::make_pair(nconv, its); 00542 00543 }
| std::pair< unsigned int, unsigned int > SlepcEigenSolver< T >::_solve_standard_helper | ( | Mat | mat, | |
| int | nev, | |||
| int | ncv, | |||
| const double | tol, | |||
| const unsigned int | m_its | |||
| ) | [inline, private] |
Helper function that actually performs the standard eigensolve.
Definition at line 151 of file slepc_eigen_solver.C.
References libMesh::COMM_WORLD.
Referenced by SlepcEigenSolver< T >::solve_standard().
00156 { 00157 START_LOG("solve_standard()", "SlepcEigenSolver"); 00158 00159 int ierr=0; 00160 00161 // converged eigen pairs and number of iterations 00162 int nconv=0; 00163 int its=0; 00164 00165 #ifdef DEBUG 00166 // The relative error. 00167 PetscReal error, re, im; 00168 00169 // Pointer to vectors of the real parts, imaginary parts. 00170 PetscScalar kr, ki; 00171 #endif 00172 00173 // Set operators. 00174 ierr = EPSSetOperators (_eps, mat, PETSC_NULL); 00175 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00176 00177 //set the problem type and the position of the spectrum 00178 set_slepc_problem_type(); 00179 set_slepc_position_of_spectrum(); 00180 00181 // Set eigenvalues to be computed. 00182 #if SLEPC_VERSION_LESS_THAN(3,0,0) 00183 ierr = EPSSetDimensions (_eps, nev, ncv); 00184 #else 00185 ierr = EPSSetDimensions (_eps, nev, ncv, PETSC_DECIDE); 00186 #endif 00187 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00188 // Set the tolerance and maximum iterations. 00189 ierr = EPSSetTolerances (_eps, tol, m_its); 00190 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00191 00192 // Set runtime options, e.g., 00193 // -eps_type <type>, -eps_nev <nev>, -eps_ncv <ncv> 00194 // Similar to PETSc, these options will override those specified 00195 // above as long as EPSSetFromOptions() is called _after_ any 00196 // other customization routines. 00197 ierr = EPSSetFromOptions (_eps); 00198 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00199 00200 // Solve the eigenproblem. 00201 ierr = EPSSolve (_eps); 00202 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00203 00204 // Get the number of iterations. 00205 ierr = EPSGetIterationNumber (_eps, &its); 00206 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00207 00208 // Get number of converged eigenpairs. 00209 ierr = EPSGetConverged(_eps,&nconv); 00210 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00211 00212 00213 #ifdef DEBUG 00214 // ierr = PetscPrintf(libMesh::COMM_WORLD, 00215 // "\n Number of iterations: %d\n" 00216 // " Number of converged eigenpairs: %d\n\n", its, nconv); 00217 00218 // Display eigenvalues and relative errors. 00219 ierr = PetscPrintf(libMesh::COMM_WORLD, 00220 " k ||Ax-kx||/|kx|\n" 00221 " ----------------- -----------------\n" ); 00222 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00223 00224 for(int i=0; i<nconv; i++ ) 00225 { 00226 ierr = EPSGetEigenpair(_eps, i, &kr, &ki, PETSC_NULL, PETSC_NULL); 00227 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00228 00229 ierr = EPSComputeRelativeError(_eps, i, &error); 00230 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00231 00232 #ifdef LIBMESH_USE_COMPLEX_NUMBERS 00233 re = PetscRealPart(kr); 00234 im = PetscImaginaryPart(kr); 00235 #else 00236 re = kr; 00237 im = ki; 00238 #endif 00239 00240 if (im != .0) 00241 { 00242 ierr = PetscPrintf(libMesh::COMM_WORLD," %9f%+9f i %12f\n", re, im, error); 00243 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00244 } 00245 else 00246 { 00247 ierr = PetscPrintf(libMesh::COMM_WORLD," %12f %12f\n", re, error); 00248 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00249 } 00250 } 00251 00252 ierr = PetscPrintf(libMesh::COMM_WORLD,"\n" ); 00253 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00254 #endif // DEBUG 00255 00256 00257 STOP_LOG("solve_standard()", "SlepcEigenSolver"); 00258 00259 // return the number of converged eigenpairs 00260 // and the number of iterations 00261 return std::make_pair(nconv, its); 00262 00263 }
| void SlepcEigenSolver< T >::attach_deflation_space | ( | NumericVector< T > & | deflation_vector | ) | [inline, virtual] |
Attach a deflation space defined by a single vector.
Implements EigenSolver< T >.
Definition at line 693 of file slepc_eigen_solver.C.
References SlepcEigenSolver< T >::_eps, libMesh::COMM_WORLD, and SlepcEigenSolver< T >::init().
00694 { 00695 this->init(); 00696 00697 int ierr = 0; 00698 Vec deflation_vector = (libmesh_cast_ptr<PetscVector<T>*>(&deflation_vector_in))->vec(); 00699 Vec* deflation_space = &deflation_vector; 00700 ierr = EPSAttachDeflationSpace(_eps, 1, deflation_space, PETSC_FALSE); 00701 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00702 }
| AutoPtr< EigenSolver< T > > EigenSolver< T >::build | ( | const SolverPackage | solver_package = SLEPC_SOLVERS |
) | [inline, static, inherited] |
Builds an EigenSolver using the linear solver package specified by solver_package
Definition at line 36 of file eigen_solver.C.
References libMeshEnums::SLEPC_SOLVERS.
00037 { 00038 // Build the appropriate solver 00039 switch (solver_package) 00040 { 00041 00042 00043 00044 #ifdef LIBMESH_HAVE_SLEPC 00045 case SLEPC_SOLVERS: 00046 { 00047 AutoPtr<EigenSolver<T> > ap(new SlepcEigenSolver<T>); 00048 return ap; 00049 } 00050 #endif 00051 00052 00053 default: 00054 std::cerr << "ERROR: Unrecognized eigen solver package: " 00055 << solver_package 00056 << std::endl; 00057 libmesh_error(); 00058 } 00059 00060 AutoPtr<EigenSolver<T> > ap(NULL); 00061 return ap; 00062 }
| void SlepcEigenSolver< T >::clear | ( | ) | [inline, virtual] |
Release all memory and clear data structures.
Reimplemented from EigenSolver< T >.
Definition at line 37 of file slepc_eigen_solver.C.
References EigenSolver< T >::_eigen_solver_type, SlepcEigenSolver< T >::_eps, EigenSolver< T >::_is_initialized, libMeshEnums::ARNOLDI, libMesh::COMM_WORLD, EigenSolver< T >::initialized(), and libMeshEnums::KRYLOVSCHUR.
Referenced by SlepcEigenSolver< T >::~SlepcEigenSolver().
00038 { 00039 if (this->initialized()) 00040 { 00041 this->_is_initialized = false; 00042 00043 int ierr=0; 00044 00045 ierr = EPSDestroy(_eps); 00046 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00047 00048 // SLEPc default eigenproblem solver 00049 #if SLEPC_VERSION_LESS_THAN(2,3,2) 00050 this->_eigen_solver_type = ARNOLDI; 00051 #else 00052 // Krylov-Schur showed up as of Slepc 2.3.2 00053 this->_eigen_solver_type = KRYLOVSCHUR; 00054 #endif 00055 } 00056 }
| EigenProblemType EigenSolver< T >::eigen_problem_type | ( | ) | const [inline, inherited] |
Returns the type of the eigen problem.
Definition at line 97 of file eigen_solver.h.
References EigenSolver< T >::_eigen_problem_type.
00097 { return _eigen_problem_type;}
| EigenSolverType EigenSolver< T >::eigen_solver_type | ( | ) | const [inline, inherited] |
Returns the type of eigensolver to use.
Definition at line 92 of file eigen_solver.h.
References EigenSolver< T >::_eigen_solver_type.
00092 { return _eigen_solver_type; }
| std::pair< Real, Real > SlepcEigenSolver< T >::get_eigenpair | ( | unsigned int | i, | |
| NumericVector< T > & | solution_in | |||
| ) | [inline, virtual] |
This function returns the real and imaginary part of the ith eigenvalue and copies the respective eigenvector to the solution vector. Note that also in case of purely real matrix entries the eigenpair may be complex values.
Implements EigenSolver< T >.
Definition at line 649 of file slepc_eigen_solver.C.
References SlepcEigenSolver< T >::_eps, PetscVector< T >::close(), libMesh::COMM_WORLD, and PetscVector< T >::vec().
00651 { 00652 int ierr=0; 00653 00654 PetscReal re, im; 00655 00656 // Make sure the NumericVector passed in is really a PetscVector 00657 PetscVector<T>* solution = libmesh_cast_ptr<PetscVector<T>*>(&solution_in); 00658 00659 // real and imaginary part of the ith eigenvalue. 00660 PetscScalar kr, ki; 00661 00662 solution->close(); 00663 00664 ierr = EPSGetEigenpair(_eps, i, &kr, &ki, solution->vec(), PETSC_NULL); 00665 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00666 00667 #ifdef LIBMESH_USE_COMPLEX_NUMBERS 00668 re = PetscRealPart(kr); 00669 im = PetscImaginaryPart(kr); 00670 #else 00671 re = kr; 00672 im = ki; 00673 #endif 00674 00675 return std::make_pair(re, im); 00676 }
| std::string ReferenceCounter::get_info | ( | ) | [static, inherited] |
Gets a string containing the reference information.
Definition at line 45 of file reference_counter.C.
References ReferenceCounter::_counts, and QuadratureRules::name().
Referenced by ReferenceCounter::print_info().
00046 { 00047 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG) 00048 00049 std::ostringstream out; 00050 00051 out << '\n' 00052 << " ---------------------------------------------------------------------------- \n" 00053 << "| Reference count information |\n" 00054 << " ---------------------------------------------------------------------------- \n"; 00055 00056 for (Counts::iterator it = _counts.begin(); 00057 it != _counts.end(); ++it) 00058 { 00059 const std::string name(it->first); 00060 const unsigned int creations = it->second.first; 00061 const unsigned int destructions = it->second.second; 00062 00063 out << "| " << name << " reference count information:\n" 00064 << "| Creations: " << creations << '\n' 00065 << "| Destructions: " << destructions << '\n'; 00066 } 00067 00068 out << " ---------------------------------------------------------------------------- \n"; 00069 00070 return out.str(); 00071 00072 #else 00073 00074 return ""; 00075 00076 #endif 00077 }
| Real SlepcEigenSolver< T >::get_relative_error | ( | unsigned int | i | ) | [inline] |
- Returns:
- the relative error ||A*x-lambda*x||/|lambda*x| of the ith eigenpair. (or the equivalent for a general eigenvalue problem)
Definition at line 680 of file slepc_eigen_solver.C.
References SlepcEigenSolver< T >::_eps, and libMesh::COMM_WORLD.
00681 { 00682 int ierr=0; 00683 PetscReal error; 00684 00685 ierr = EPSComputeRelativeError(_eps, i, &error); 00686 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00687 00688 return error; 00689 }
| void ReferenceCounter::increment_constructor_count | ( | const std::string & | name | ) | [inline, protected, inherited] |
Increments the construction counter. Should be called in the constructor of any derived class that will be reference counted.
Definition at line 149 of file reference_counter.h.
References ReferenceCounter::_counts, and Threads::spin_mtx.
Referenced by ReferenceCountedObject< SparseMatrix< T > >::ReferenceCountedObject().
00150 { 00151 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 00152 std::pair<unsigned int, unsigned int>& p = _counts[name]; 00153 00154 p.first++; 00155 }
| void ReferenceCounter::increment_destructor_count | ( | const std::string & | name | ) | [inline, protected, inherited] |
Increments the destruction counter. Should be called in the destructor of any derived class that will be reference counted.
Definition at line 167 of file reference_counter.h.
References ReferenceCounter::_counts, and Threads::spin_mtx.
Referenced by ReferenceCountedObject< SparseMatrix< T > >::~ReferenceCountedObject().
00168 { 00169 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 00170 std::pair<unsigned int, unsigned int>& p = _counts[name]; 00171 00172 p.second++; 00173 }
| void SlepcEigenSolver< T >::init | ( | ) | [inline, virtual] |
Initialize data structures if not done so already.
Implements EigenSolver< T >.
Definition at line 61 of file slepc_eigen_solver.C.
References SlepcEigenSolver< T >::_eps, EigenSolver< T >::_is_initialized, libMesh::COMM_WORLD, EigenSolver< T >::initialized(), and SlepcEigenSolver< T >::set_slepc_solver_type().
Referenced by SlepcEigenSolver< T >::attach_deflation_space(), SlepcEigenSolver< T >::solve_generalized(), and SlepcEigenSolver< T >::solve_standard().
00062 { 00063 00064 int ierr=0; 00065 00066 // Initialize the data structures if not done so already. 00067 if (!this->initialized()) 00068 { 00069 this->_is_initialized = true; 00070 00071 // Create the eigenproblem solver context 00072 ierr = EPSCreate (libMesh::COMM_WORLD, &_eps); 00073 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00074 00075 // Set user-specified solver 00076 set_slepc_solver_type(); 00077 } 00078 }
| bool EigenSolver< T >::initialized | ( | ) | const [inline, inherited] |
- Returns:
- true if the data structures are initialized, false otherwise.
Definition at line 76 of file eigen_solver.h.
References EigenSolver< T >::_is_initialized.
Referenced by SlepcEigenSolver< T >::clear(), and SlepcEigenSolver< T >::init().
00076 { return _is_initialized; }
| static unsigned int ReferenceCounter::n_objects | ( | ) | [inline, static, inherited] |
Prints the number of outstanding (created, but not yet destroyed) objects.
Definition at line 76 of file reference_counter.h.
References ReferenceCounter::_n_objects.
00077 { return _n_objects; }
| PositionOfSpectrum EigenSolver< T >::position_of_spectrum | ( | ) | const [inline, inherited] |
Returns the position of the spectrum to compute.
Definition at line 102 of file eigen_solver.h.
References EigenSolver< T >::_position_of_spectrum.
00103 { return _position_of_spectrum;}
| void ReferenceCounter::print_info | ( | ) | [static, inherited] |
Prints the reference information to std::cout.
Definition at line 83 of file reference_counter.C.
References ReferenceCounter::get_info().
00084 { 00085 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG) 00086 00087 std::cout << ReferenceCounter::get_info(); 00088 00089 #endif 00090 }
| void EigenSolver< T >::set_eigenproblem_type | ( | EigenProblemType | ept | ) | [inline, inherited] |
Sets the type of the eigenproblem.
Definition at line 114 of file eigen_solver.h.
References EigenSolver< T >::_eigen_problem_type.
00115 {_eigen_problem_type = ept;}
| void EigenSolver< T >::set_eigensolver_type | ( | const EigenSolverType | est | ) | [inline, inherited] |
Sets the type of eigensolver to use.
Definition at line 108 of file eigen_solver.h.
References EigenSolver< T >::_eigen_solver_type.
00109 { _eigen_solver_type = est; }
| void EigenSolver< T >::set_position_of_spectrum | ( | PositionOfSpectrum | pos | ) | [inline, inherited] |
Sets the position of the spectrum.
Definition at line 120 of file eigen_solver.h.
References EigenSolver< T >::_position_of_spectrum.
00121 {_position_of_spectrum= pos;}
| void SlepcEigenSolver< T >::set_slepc_position_of_spectrum | ( | ) | [inline, private] |
Tells Slepc to compute the spectrum at the position stored in _position_of_spectrum
Definition at line 616 of file slepc_eigen_solver.C.
References SlepcEigenSolver< T >::_eps, EigenSolver< T >::_position_of_spectrum, libMesh::COMM_WORLD, libMeshEnums::LARGEST_IMAGINARY, libMeshEnums::LARGEST_MAGNITUDE, libMeshEnums::LARGEST_REAL, libMeshEnums::SMALLEST_IMAGINARY, libMeshEnums::SMALLEST_MAGNITUDE, and libMeshEnums::SMALLEST_REAL.
Referenced by SlepcEigenSolver< T >::_solve_generalized_helper().
00617 { 00618 int ierr = 0; 00619 00620 switch (this->_position_of_spectrum) 00621 { 00622 case LARGEST_MAGNITUDE: 00623 ierr = EPSSetWhichEigenpairs (_eps, EPS_LARGEST_MAGNITUDE); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00624 case SMALLEST_MAGNITUDE: 00625 ierr = EPSSetWhichEigenpairs (_eps, EPS_SMALLEST_MAGNITUDE); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00626 case LARGEST_REAL: 00627 ierr = EPSSetWhichEigenpairs (_eps, EPS_LARGEST_REAL); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00628 case SMALLEST_REAL: 00629 ierr = EPSSetWhichEigenpairs (_eps, EPS_SMALLEST_REAL); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00630 case LARGEST_IMAGINARY: 00631 ierr = EPSSetWhichEigenpairs (_eps, EPS_LARGEST_IMAGINARY); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00632 case SMALLEST_IMAGINARY: 00633 ierr = EPSSetWhichEigenpairs (_eps, EPS_SMALLEST_IMAGINARY); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00634 00635 00636 default: 00637 std::cerr << "ERROR: Unsupported SLEPc position of spectrum: " 00638 << this->_position_of_spectrum << std::endl; 00639 libmesh_error(); 00640 } 00641 }
| void SlepcEigenSolver< T >::set_slepc_problem_type | ( | ) | [inline, private] |
Tells Slepc to deal with the type of problem stored in _eigen_problem_type
Definition at line 591 of file slepc_eigen_solver.C.
References EigenSolver< T >::_eigen_problem_type, SlepcEigenSolver< T >::_eps, libMesh::COMM_WORLD, libMeshEnums::GHEP, libMeshEnums::GNHEP, libMeshEnums::HEP, and libMeshEnums::NHEP.
Referenced by SlepcEigenSolver< T >::_solve_generalized_helper().
00592 { 00593 int ierr = 0; 00594 00595 switch (this->_eigen_problem_type) 00596 { 00597 case NHEP: 00598 ierr = EPSSetProblemType (_eps, EPS_NHEP); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00599 case GNHEP: 00600 ierr = EPSSetProblemType (_eps, EPS_GNHEP); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00601 case HEP: 00602 ierr = EPSSetProblemType (_eps, EPS_HEP); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00603 case GHEP: 00604 ierr = EPSSetProblemType (_eps, EPS_GHEP); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00605 00606 default: 00607 std::cerr << "ERROR: Unsupported SLEPc Eigen Problem: " 00608 << this->_eigen_problem_type << std::endl 00609 << "Continuing with SLEPc defaults" << std::endl; 00610 } 00611 }
| void SlepcEigenSolver< T >::set_slepc_solver_type | ( | ) | [inline, private] |
Tells Slepc to use the user-specified solver stored in _eigen_solver_type
Definition at line 556 of file slepc_eigen_solver.C.
References EigenSolver< T >::_eigen_solver_type, SlepcEigenSolver< T >::_eps, libMeshEnums::ARNOLDI, libMesh::COMM_WORLD, libMeshEnums::KRYLOVSCHUR, libMeshEnums::LANCZOS, libMeshEnums::LAPACK, libMeshEnums::POWER, and libMeshEnums::SUBSPACE.
Referenced by SlepcEigenSolver< T >::init().
00557 { 00558 int ierr = 0; 00559 00560 switch (this->_eigen_solver_type) 00561 { 00562 case POWER: 00563 ierr = EPSSetType (_eps, (char*) EPSPOWER); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00564 case SUBSPACE: 00565 ierr = EPSSetType (_eps, (char*) EPSSUBSPACE); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00566 case LAPACK: 00567 ierr = EPSSetType (_eps, (char*) EPSLAPACK); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00568 case ARNOLDI: 00569 ierr = EPSSetType (_eps, (char*) EPSARNOLDI); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00570 case LANCZOS: 00571 ierr = EPSSetType (_eps, (char*) EPSLANCZOS); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00572 #if !SLEPC_VERSION_LESS_THAN(2,3,2) 00573 // EPSKRYLOVSCHUR added in 2.3.2 00574 case KRYLOVSCHUR: 00575 ierr = EPSSetType (_eps, (char*) EPSKRYLOVSCHUR); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00576 #endif 00577 // case ARPACK: 00578 // ierr = EPSSetType (_eps, (char*) EPSARPACK); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00579 00580 default: 00581 std::cerr << "ERROR: Unsupported SLEPc Eigen Solver: " 00582 << this->_eigen_solver_type << std::endl 00583 << "Continuing with SLEPc defaults" << std::endl; 00584 } 00585 }
| std::pair< unsigned int, unsigned int > SlepcEigenSolver< T >::solve_generalized | ( | ShellMatrix< T > & | matrix_A, | |
| ShellMatrix< T > & | matrix_B, | |||
| int | nev, | |||
| int | ncv, | |||
| const double | tol, | |||
| const unsigned int | m_its | |||
| ) | [inline, virtual] |
Solve generalized eigenproblem when both matrix_A and matrix_B are of type ShellMatrix. When using this function, one should use the command line options: -st_ksp_type gmres -st_pc_type none or -st_ksp_type gmres -st_pc_type jacobi or similar.
Implements EigenSolver< T >.
Definition at line 377 of file slepc_eigen_solver.C.
References SlepcEigenSolver< T >::_petsc_shell_matrix_get_diagonal(), SlepcEigenSolver< T >::_petsc_shell_matrix_mult(), SlepcEigenSolver< T >::_solve_generalized_helper(), libMesh::COMM_WORLD, SlepcEigenSolver< T >::init(), ShellMatrix< T >::m(), and ShellMatrix< T >::n().
00383 { 00384 this->init (); 00385 00386 int ierr=0; 00387 00388 // Prepare the matrix. 00389 Mat mat_A; 00390 ierr = MatCreateShell(libMesh::COMM_WORLD, 00391 shell_matrix_A.m(), // Specify the number of local rows 00392 shell_matrix_A.n(), // Specify the number of local columns 00393 PETSC_DETERMINE, 00394 PETSC_DETERMINE, 00395 const_cast<void*>(static_cast<const void*>(&shell_matrix_A)), 00396 &mat_A); 00397 00398 Mat mat_B; 00399 ierr = MatCreateShell(libMesh::COMM_WORLD, 00400 shell_matrix_B.m(), // Specify the number of local rows 00401 shell_matrix_B.n(), // Specify the number of local columns 00402 PETSC_DETERMINE, 00403 PETSC_DETERMINE, 00404 const_cast<void*>(static_cast<const void*>(&shell_matrix_B)), 00405 &mat_B); 00406 00407 /* Note that the const_cast above is only necessary because PETSc 00408 does not accept a const void*. Inside the member function 00409 _petsc_shell_matrix() below, the pointer is casted back to a 00410 const ShellMatrix<T>*. */ 00411 00412 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00413 ierr = MatShellSetOperation(mat_A,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult)); 00414 ierr = MatShellSetOperation(mat_A,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal)); 00415 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00416 00417 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00418 ierr = MatShellSetOperation(mat_B,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult)); 00419 ierr = MatShellSetOperation(mat_B,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal)); 00420 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00421 00422 return _solve_generalized_helper (mat_A, mat_B, nev, ncv, tol, m_its); 00423 }
| std::pair< unsigned int, unsigned int > SlepcEigenSolver< T >::solve_generalized | ( | SparseMatrix< T > & | matrix_A, | |
| ShellMatrix< T > & | matrix_B, | |||
| int | nev, | |||
| int | ncv, | |||
| const double | tol, | |||
| const unsigned int | m_its | |||
| ) | [inline, virtual] |
Solve generalized eigenproblem when matrix_A is of type SparseMatrix, matrix_B is of type ShellMatrix. When using this function, one should use the command line options: -st_ksp_type gmres -st_pc_type none or -st_ksp_type gmres -st_pc_type jacobi or similar.
Implements EigenSolver< T >.
Definition at line 335 of file slepc_eigen_solver.C.
References SlepcEigenSolver< T >::_petsc_shell_matrix_get_diagonal(), SlepcEigenSolver< T >::_petsc_shell_matrix_mult(), SlepcEigenSolver< T >::_solve_generalized_helper(), PetscMatrix< T >::close(), libMesh::COMM_WORLD, SlepcEigenSolver< T >::init(), ShellMatrix< T >::m(), PetscMatrix< T >::mat(), and ShellMatrix< T >::n().
00341 { 00342 this->init (); 00343 00344 int ierr=0; 00345 00346 PetscMatrix<T>* matrix_A = libmesh_cast_ptr<PetscMatrix<T>*>(&matrix_A_in); 00347 00348 // Close the matrix and vectors in case this wasn't already done. 00349 matrix_A->close (); 00350 00351 // Prepare the matrix. 00352 Mat mat_B; 00353 ierr = MatCreateShell(libMesh::COMM_WORLD, 00354 shell_matrix_B.m(), // Specify the number of local rows 00355 shell_matrix_B.n(), // Specify the number of local columns 00356 PETSC_DETERMINE, 00357 PETSC_DETERMINE, 00358 const_cast<void*>(static_cast<const void*>(&shell_matrix_B)), 00359 &mat_B); 00360 00361 00362 /* Note that the const_cast above is only necessary because PETSc 00363 does not accept a const void*. Inside the member function 00364 _petsc_shell_matrix() below, the pointer is casted back to a 00365 const ShellMatrix<T>*. */ 00366 00367 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00368 ierr = MatShellSetOperation(mat_B,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult)); 00369 ierr = MatShellSetOperation(mat_B,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal)); 00370 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00371 00372 return _solve_generalized_helper (matrix_A->mat(), mat_B, nev, ncv, tol, m_its); 00373 }
| std::pair< unsigned int, unsigned int > SlepcEigenSolver< T >::solve_generalized | ( | ShellMatrix< T > & | matrix_A, | |
| SparseMatrix< T > & | matrix_B, | |||
| int | nev, | |||
| int | ncv, | |||
| const double | tol, | |||
| const unsigned int | m_its | |||
| ) | [inline, virtual] |
Solve generalized eigenproblem when matrix_A is of type ShellMatrix, matrix_B is of type SparseMatrix.
Implements EigenSolver< T >.
Definition at line 294 of file slepc_eigen_solver.C.
References SlepcEigenSolver< T >::_petsc_shell_matrix_get_diagonal(), SlepcEigenSolver< T >::_petsc_shell_matrix_mult(), SlepcEigenSolver< T >::_solve_generalized_helper(), libMesh::COMM_WORLD, SlepcEigenSolver< T >::init(), ShellMatrix< T >::m(), and ShellMatrix< T >::n().
00300 { 00301 this->init (); 00302 00303 int ierr=0; 00304 00305 // Prepare the matrix. 00306 Mat mat_A; 00307 ierr = MatCreateShell(libMesh::COMM_WORLD, 00308 shell_matrix_A.m(), // Specify the number of local rows 00309 shell_matrix_A.n(), // Specify the number of local columns 00310 PETSC_DETERMINE, 00311 PETSC_DETERMINE, 00312 const_cast<void*>(static_cast<const void*>(&shell_matrix_A)), 00313 &mat_A); 00314 00315 PetscMatrix<T>* matrix_B = libmesh_cast_ptr<PetscMatrix<T>*>(&matrix_B_in); 00316 00317 // Close the matrix and vectors in case this wasn't already done. 00318 matrix_B->close (); 00319 00320 /* Note that the const_cast above is only necessary because PETSc 00321 does not accept a const void*. Inside the member function 00322 _petsc_shell_matrix() below, the pointer is casted back to a 00323 const ShellMatrix<T>*. */ 00324 00325 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00326 ierr = MatShellSetOperation(mat_A,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult)); 00327 ierr = MatShellSetOperation(mat_A,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal)); 00328 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00329 00330 return _solve_generalized_helper (mat_A, matrix_B->mat(), nev, ncv, tol, m_its); 00331 }
| std::pair< unsigned int, unsigned int > SlepcEigenSolver< T >::solve_generalized | ( | SparseMatrix< T > & | matrix_A, | |
| SparseMatrix< T > & | matrix_B, | |||
| int | nev, | |||
| int | ncv, | |||
| const double | tol, | |||
| const unsigned int | m_its | |||
| ) | [inline, virtual] |
This function calls the SLEPc solver to compute the eigenpairs for the generalized eigenproblem defined by the matrix_A and matrix_B, which are of type SparseMatrix. The argument nev is the number of eigenpairs to be computed and ncv is the number of basis vectors to be used in the solution procedure. Return values are the number of converged eigen values and the number of the iterations carried out by the eigen solver.
Implements EigenSolver< T >.
Definition at line 271 of file slepc_eigen_solver.C.
References SlepcEigenSolver< T >::_solve_generalized_helper(), PetscMatrix< T >::close(), SlepcEigenSolver< T >::init(), and PetscMatrix< T >::mat().
00277 { 00278 this->init (); 00279 00280 // Make sure the data passed in are really of Petsc types 00281 PetscMatrix<T>* matrix_A = libmesh_cast_ptr<PetscMatrix<T>*>(&matrix_A_in); 00282 PetscMatrix<T>* matrix_B = libmesh_cast_ptr<PetscMatrix<T>*>(&matrix_B_in); 00283 00284 // Close the matrix and vectors in case this wasn't already done. 00285 matrix_A->close (); 00286 matrix_B->close (); 00287 00288 00289 return _solve_generalized_helper (matrix_A->mat(), matrix_B->mat(), nev, ncv, tol, m_its); 00290 }
| std::pair< unsigned int, unsigned int > SlepcEigenSolver< T >::solve_standard | ( | ShellMatrix< T > & | shell_matrix, | |
| int | nev, | |||
| int | ncv, | |||
| const double | tol, | |||
| const unsigned int | m_its | |||
| ) | [inline, virtual] |
Same as above except that matrix_A is a ShellMatrix in this case.
Implements EigenSolver< T >.
Definition at line 114 of file slepc_eigen_solver.C.
References SlepcEigenSolver< T >::_petsc_shell_matrix_get_diagonal(), SlepcEigenSolver< T >::_petsc_shell_matrix_mult(), SlepcEigenSolver< T >::_solve_standard_helper(), libMesh::COMM_WORLD, SlepcEigenSolver< T >::init(), ShellMatrix< T >::m(), and ShellMatrix< T >::n().
00119 { 00120 this->init (); 00121 00122 int ierr=0; 00123 00124 // Prepare the matrix. 00125 Mat mat; 00126 ierr = MatCreateShell(libMesh::COMM_WORLD, 00127 shell_matrix.m(), // Specify the number of local rows 00128 shell_matrix.n(), // Specify the number of local columns 00129 PETSC_DETERMINE, 00130 PETSC_DETERMINE, 00131 const_cast<void*>(static_cast<const void*>(&shell_matrix)), 00132 &mat); 00133 00134 /* Note that the const_cast above is only necessary because PETSc 00135 does not accept a const void*. Inside the member function 00136 _petsc_shell_matrix() below, the pointer is casted back to a 00137 const ShellMatrix<T>*. */ 00138 00139 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00140 ierr = MatShellSetOperation(mat,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult)); 00141 ierr = MatShellSetOperation(mat,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal)); 00142 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00143 00144 00145 return _solve_standard_helper(mat, nev, ncv, tol, m_its); 00146 }
| std::pair< unsigned int, unsigned int > SlepcEigenSolver< T >::solve_standard | ( | SparseMatrix< T > & | matrix_A, | |
| int | nev, | |||
| int | ncv, | |||
| const double | tol, | |||
| const unsigned int | m_its | |||
| ) | [inline, virtual] |
This function calls the SLEPc solver to compute the eigenpairs of the SparseMatrix matrix_A. nev is the number of eigenpairs to be computed and ncv is the number of basis vectors to be used in the solution procedure. Return values are the number of converged eigen values and the number of the iterations carried out by the eigen solver.
Implements EigenSolver< T >.
Definition at line 84 of file slepc_eigen_solver.C.
References SlepcEigenSolver< T >::_solve_standard_helper(), PetscMatrix< T >::close(), SlepcEigenSolver< T >::init(), and PetscMatrix< T >::mat().
00089 { 00090 // START_LOG("solve_standard()", "SlepcEigenSolver"); 00091 00092 this->init (); 00093 00094 // Make sure the SparseMatrix passed in is really a PetscMatrix 00095 PetscMatrix<T>* matrix_A = libmesh_cast_ptr<PetscMatrix<T>*>(&matrix_A_in); 00096 00097 // Close the matrix and vectors in case this wasn't already done. 00098 matrix_A->close (); 00099 00100 // just for debugging, remove this 00101 // char mat_file[] = "matA.petsc"; 00102 // PetscViewer petsc_viewer; 00103 // ierr = PetscViewerBinaryOpen(libMesh::COMM_WORLD, mat_file, PETSC_FILE_CREATE, &petsc_viewer); 00104 // CHKERRABORT(libMesh::COMM_WORLD,ierr); 00105 // ierr = MatView(matrix_A->mat(),petsc_viewer); 00106 // CHKERRABORT(libMesh::COMM_WORLD,ierr); 00107 00108 return _solve_standard_helper(matrix_A->mat(), nev, ncv, tol, m_its); 00109 }
Member Data Documentation
ReferenceCounter::Counts ReferenceCounter::_counts [static, protected, inherited] |
Actually holds the data.
Definition at line 110 of file reference_counter.h.
Referenced by ReferenceCounter::get_info(), ReferenceCounter::increment_constructor_count(), and ReferenceCounter::increment_destructor_count().
EigenProblemType EigenSolver< T >::_eigen_problem_type [protected, inherited] |
Enum stating which type of eigen problem we deal with.
Definition at line 215 of file eigen_solver.h.
Referenced by EigenSolver< T >::eigen_problem_type(), EigenSolver< T >::set_eigenproblem_type(), SlepcEigenSolver< T >::set_slepc_problem_type(), and SlepcEigenSolver< T >::SlepcEigenSolver().
EigenSolverType EigenSolver< T >::_eigen_solver_type [protected, inherited] |
Enum stating which type of eigensolver to use.
Definition at line 210 of file eigen_solver.h.
Referenced by SlepcEigenSolver< T >::clear(), EigenSolver< T >::eigen_solver_type(), EigenSolver< T >::set_eigensolver_type(), SlepcEigenSolver< T >::set_slepc_solver_type(), and SlepcEigenSolver< T >::SlepcEigenSolver().
EPS SlepcEigenSolver< T >::_eps [private] |
Eigenproblem solver context
Definition at line 266 of file slepc_eigen_solver.h.
Referenced by SlepcEigenSolver< T >::_solve_generalized_helper(), SlepcEigenSolver< T >::attach_deflation_space(), SlepcEigenSolver< T >::clear(), SlepcEigenSolver< T >::get_eigenpair(), SlepcEigenSolver< T >::get_relative_error(), SlepcEigenSolver< T >::init(), SlepcEigenSolver< T >::set_slepc_position_of_spectrum(), SlepcEigenSolver< T >::set_slepc_problem_type(), and SlepcEigenSolver< T >::set_slepc_solver_type().
bool EigenSolver< T >::_is_initialized [protected, inherited] |
Flag indicating if the data structures have been initialized.
Definition at line 225 of file eigen_solver.h.
Referenced by SlepcEigenSolver< T >::clear(), SlepcEigenSolver< T >::init(), and EigenSolver< T >::initialized().
Threads::spin_mutex ReferenceCounter::_mutex [static, protected, inherited] |
Mutual exclusion object to enable thread-safe reference counting.
Definition at line 123 of file reference_counter.h.
Threads::atomic< unsigned int > ReferenceCounter::_n_objects [static, protected, inherited] |
The number of objects. Print the reference count information when the number returns to 0.
Definition at line 118 of file reference_counter.h.
Referenced by ReferenceCounter::n_objects(), ReferenceCounter::ReferenceCounter(), and ReferenceCounter::~ReferenceCounter().
PositionOfSpectrum EigenSolver< T >::_position_of_spectrum [protected, inherited] |
Enum stating where to evaluate the spectrum.
Definition at line 220 of file eigen_solver.h.
Referenced by EigenSolver< T >::position_of_spectrum(), EigenSolver< T >::set_position_of_spectrum(), and SlepcEigenSolver< T >::set_slepc_position_of_spectrum().
The documentation for this class was generated from the following files: