libMesh::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) |
| std::pair< Real, Real > | get_eigenvalue (unsigned int i) |
| 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 (std::ostream &out=libMesh::out) |
| static unsigned int | n_objects () |
| static void | enable_print_counter_info () |
| static void | disable_print_counter_info () |
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 |
| static bool | _enable_print_counter = true |
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 libMesh::SlepcEigenSolver< T >
This class provides an interface to the SLEPc eigenvalue solver library www.grycap.upv.es/slepc/.
Definition at line 51 of file slepc_eigen_solver.h.
Member Typedef Documentation
typedef std::map<std::string, std::pair<unsigned int, unsigned int> > libMesh::ReferenceCounter::Counts [protected, inherited] |
Data structure to log the information. The log is identified by the class name.
Definition at line 113 of file reference_counter.h.
Constructor & Destructor Documentation
| libMesh::SlepcEigenSolver< T >::SlepcEigenSolver | ( | ) | [inline] |
Constructor. Initializes Petsc data structures
Definition at line 262 of file slepc_eigen_solver.h.
References libMesh::EigenSolver< T >::_eigen_problem_type, libMesh::EigenSolver< T >::_eigen_solver_type, libMeshEnums::ARNOLDI, and libMeshEnums::NHEP.
00263 { 00264 this->_eigen_solver_type = ARNOLDI; 00265 this->_eigen_problem_type = NHEP; 00266 }
| libMesh::SlepcEigenSolver< T >::~SlepcEigenSolver | ( | ) | [inline] |
Destructor.
Definition at line 272 of file slepc_eigen_solver.h.
References libMesh::SlepcEigenSolver< T >::clear().
00273 { 00274 this->clear (); 00275 }
Member Function Documentation
| PetscErrorCode libMesh::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 759 of file slepc_eigen_solver.C.
References libMesh::COMM_WORLD, and libMesh::ShellMatrix< T >::get_diagonal().
Referenced by libMesh::SlepcEigenSolver< T >::solve_generalized(), and libMesh::SlepcEigenSolver< T >::solve_standard().
00760 { 00761 /* Get the matrix context. */ 00762 int ierr=0; 00763 void* ctx; 00764 ierr = MatShellGetContext(mat,&ctx); 00765 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00766 00767 /* Get user shell matrix object. */ 00768 const ShellMatrix<T>& shell_matrix = *static_cast<const ShellMatrix<T>*>(ctx); 00769 00770 /* Make \p NumericVector instances around the vector. */ 00771 PetscVector<T> dest_global(dest); 00772 00773 /* Call the user function. */ 00774 shell_matrix.get_diagonal(dest_global); 00775 00776 return ierr; 00777 }
| PetscErrorCode libMesh::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 737 of file slepc_eigen_solver.C.
References libMesh::COMM_WORLD, and libMesh::ShellMatrix< T >::vector_mult().
Referenced by libMesh::SlepcEigenSolver< T >::solve_generalized(), and libMesh::SlepcEigenSolver< T >::solve_standard().
00738 { 00739 /* Get the matrix context. */ 00740 int ierr=0; 00741 void* ctx; 00742 ierr = MatShellGetContext(mat,&ctx); 00743 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00744 00745 /* Get user shell matrix object. */ 00746 const ShellMatrix<T>& shell_matrix = *static_cast<const ShellMatrix<T>*>(ctx); 00747 00748 /* Make \p NumericVector instances around the vectors. */ 00749 PetscVector<T> arg_global(arg); 00750 PetscVector<T> dest_global(dest); 00751 00752 /* Call the user function. */ 00753 shell_matrix.vector_mult(dest_global,arg_global); 00754 00755 return ierr; 00756 }
| std::pair< unsigned int, unsigned int > libMesh::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 432 of file slepc_eigen_solver.C.
References libMesh::SlepcEigenSolver< T >::_eps, libMesh::COMM_WORLD, libMesh::SlepcEigenSolver< T >::set_slepc_position_of_spectrum(), and libMesh::SlepcEigenSolver< T >::set_slepc_problem_type().
Referenced by libMesh::SlepcEigenSolver< T >::solve_generalized().
00438 { 00439 START_LOG("solve_generalized()", "SlepcEigenSolver"); 00440 00441 int ierr=0; 00442 00443 // converged eigen pairs and number of iterations 00444 int nconv=0; 00445 int its=0; 00446 00447 #ifdef DEBUG 00448 // The relative error. 00449 PetscReal error, re, im; 00450 00451 // Pointer to vectors of the real parts, imaginary parts. 00452 PetscScalar kr, ki; 00453 #endif 00454 00455 // Set operators. 00456 ierr = EPSSetOperators (_eps, mat_A, mat_B); 00457 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00458 00459 //set the problem type and the position of the spectrum 00460 set_slepc_problem_type(); 00461 set_slepc_position_of_spectrum(); 00462 00463 // Set eigenvalues to be computed. 00464 #if SLEPC_VERSION_LESS_THAN(3,0,0) 00465 ierr = EPSSetDimensions (_eps, nev, ncv); 00466 #else 00467 ierr = EPSSetDimensions (_eps, nev, ncv, PETSC_DECIDE); 00468 #endif 00469 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00470 00471 00472 // Set the tolerance and maximum iterations. 00473 ierr = EPSSetTolerances (_eps, tol, m_its); 00474 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00475 00476 // Set runtime options, e.g., 00477 // -eps_type <type>, -eps_nev <nev>, -eps_ncv <ncv> 00478 // Similar to PETSc, these options will override those specified 00479 // above as long as EPSSetFromOptions() is called _after_ any 00480 // other customization routines. 00481 ierr = EPSSetFromOptions (_eps); 00482 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00483 00484 // Solve the eigenproblem. 00485 ierr = EPSSolve (_eps); 00486 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00487 00488 // Get the number of iterations. 00489 ierr = EPSGetIterationNumber (_eps, &its); 00490 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00491 00492 // Get number of converged eigenpairs. 00493 ierr = EPSGetConverged(_eps,&nconv); 00494 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00495 00496 00497 #ifdef DEBUG 00498 // ierr = PetscPrintf(libMesh::COMM_WORLD, 00499 // "\n Number of iterations: %d\n" 00500 // " Number of converged eigenpairs: %d\n\n", its, nconv); 00501 00502 // Display eigenvalues and relative errors. 00503 ierr = PetscPrintf(libMesh::COMM_WORLD, 00504 " k ||Ax-kx||/|kx|\n" 00505 " ----------------- -----------------\n" ); 00506 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00507 00508 for(int i=0; i<nconv; i++ ) 00509 { 00510 ierr = EPSGetEigenpair(_eps, i, &kr, &ki, PETSC_NULL, PETSC_NULL); 00511 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00512 00513 ierr = EPSComputeRelativeError(_eps, i, &error); 00514 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00515 00516 #ifdef LIBMESH_USE_COMPLEX_NUMBERS 00517 re = PetscRealPart(kr); 00518 im = PetscImaginaryPart(kr); 00519 #else 00520 re = kr; 00521 im = ki; 00522 #endif 00523 00524 if (im != .0) 00525 { 00526 ierr = PetscPrintf(libMesh::COMM_WORLD," %9f%+9f i %12f\n", re, im, error); 00527 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00528 } 00529 else 00530 { 00531 ierr = PetscPrintf(libMesh::COMM_WORLD," %12f %12f\n", re, error); 00532 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00533 } 00534 } 00535 00536 ierr = PetscPrintf(libMesh::COMM_WORLD,"\n" ); 00537 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00538 #endif // DEBUG 00539 00540 STOP_LOG("solve_generalized()", "SlepcEigenSolver"); 00541 00542 // return the number of converged eigenpairs 00543 // and the number of iterations 00544 return std::make_pair(nconv, its); 00545 00546 }
| std::pair< unsigned int, unsigned int > libMesh::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 154 of file slepc_eigen_solver.C.
References libMesh::COMM_WORLD.
Referenced by libMesh::SlepcEigenSolver< T >::solve_standard().
00159 { 00160 START_LOG("solve_standard()", "SlepcEigenSolver"); 00161 00162 int ierr=0; 00163 00164 // converged eigen pairs and number of iterations 00165 int nconv=0; 00166 int its=0; 00167 00168 #ifdef DEBUG 00169 // The relative error. 00170 PetscReal error, re, im; 00171 00172 // Pointer to vectors of the real parts, imaginary parts. 00173 PetscScalar kr, ki; 00174 #endif 00175 00176 // Set operators. 00177 ierr = EPSSetOperators (_eps, mat, PETSC_NULL); 00178 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00179 00180 //set the problem type and the position of the spectrum 00181 set_slepc_problem_type(); 00182 set_slepc_position_of_spectrum(); 00183 00184 // Set eigenvalues to be computed. 00185 #if SLEPC_VERSION_LESS_THAN(3,0,0) 00186 ierr = EPSSetDimensions (_eps, nev, ncv); 00187 #else 00188 ierr = EPSSetDimensions (_eps, nev, ncv, PETSC_DECIDE); 00189 #endif 00190 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00191 // Set the tolerance and maximum iterations. 00192 ierr = EPSSetTolerances (_eps, tol, m_its); 00193 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00194 00195 // Set runtime options, e.g., 00196 // -eps_type <type>, -eps_nev <nev>, -eps_ncv <ncv> 00197 // Similar to PETSc, these options will override those specified 00198 // above as long as EPSSetFromOptions() is called _after_ any 00199 // other customization routines. 00200 ierr = EPSSetFromOptions (_eps); 00201 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00202 00203 // Solve the eigenproblem. 00204 ierr = EPSSolve (_eps); 00205 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00206 00207 // Get the number of iterations. 00208 ierr = EPSGetIterationNumber (_eps, &its); 00209 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00210 00211 // Get number of converged eigenpairs. 00212 ierr = EPSGetConverged(_eps,&nconv); 00213 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00214 00215 00216 #ifdef DEBUG 00217 // ierr = PetscPrintf(libMesh::COMM_WORLD, 00218 // "\n Number of iterations: %d\n" 00219 // " Number of converged eigenpairs: %d\n\n", its, nconv); 00220 00221 // Display eigenvalues and relative errors. 00222 ierr = PetscPrintf(libMesh::COMM_WORLD, 00223 " k ||Ax-kx||/|kx|\n" 00224 " ----------------- -----------------\n" ); 00225 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00226 00227 for(int i=0; i<nconv; i++ ) 00228 { 00229 ierr = EPSGetEigenpair(_eps, i, &kr, &ki, PETSC_NULL, PETSC_NULL); 00230 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00231 00232 ierr = EPSComputeRelativeError(_eps, i, &error); 00233 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00234 00235 #ifdef LIBMESH_USE_COMPLEX_NUMBERS 00236 re = PetscRealPart(kr); 00237 im = PetscImaginaryPart(kr); 00238 #else 00239 re = kr; 00240 im = ki; 00241 #endif 00242 00243 if (im != .0) 00244 { 00245 ierr = PetscPrintf(libMesh::COMM_WORLD," %9f%+9f i %12f\n", re, im, error); 00246 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00247 } 00248 else 00249 { 00250 ierr = PetscPrintf(libMesh::COMM_WORLD," %12f %12f\n", re, error); 00251 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00252 } 00253 } 00254 00255 ierr = PetscPrintf(libMesh::COMM_WORLD,"\n" ); 00256 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00257 #endif // DEBUG 00258 00259 00260 STOP_LOG("solve_standard()", "SlepcEigenSolver"); 00261 00262 // return the number of converged eigenpairs 00263 // and the number of iterations 00264 return std::make_pair(nconv, its); 00265 00266 }
| void libMesh::SlepcEigenSolver< T >::attach_deflation_space | ( | NumericVector< T > & | deflation_vector | ) | [inline, virtual] |
Attach a deflation space defined by a single vector.
Implements libMesh::EigenSolver< T >.
Definition at line 721 of file slepc_eigen_solver.C.
References libMesh::SlepcEigenSolver< T >::_eps, libMesh::COMM_WORLD, and libMesh::SlepcEigenSolver< T >::init().
00722 { 00723 this->init(); 00724 00725 int ierr = 0; 00726 Vec deflation_vector = (libmesh_cast_ptr<PetscVector<T>*>(&deflation_vector_in))->vec(); 00727 Vec* deflation_space = &deflation_vector; 00728 #if SLEPC_VERSION_LESS_THAN(3,1,0) 00729 ierr = EPSAttachDeflationSpace(_eps, 1, deflation_space, PETSC_FALSE); 00730 #else 00731 ierr = EPSSetDeflationSpace(_eps, 1, deflation_space); 00732 #endif 00733 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00734 }
| AutoPtr< EigenSolver< T > > libMesh::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 37 of file eigen_solver.C.
References libMesh::err, and libMeshEnums::SLEPC_SOLVERS.
00038 { 00039 // Build the appropriate solver 00040 switch (solver_package) 00041 { 00042 00043 00044 00045 #ifdef LIBMESH_HAVE_SLEPC 00046 case SLEPC_SOLVERS: 00047 { 00048 AutoPtr<EigenSolver<T> > ap(new SlepcEigenSolver<T>); 00049 return ap; 00050 } 00051 #endif 00052 00053 00054 default: 00055 libMesh::err << "ERROR: Unrecognized eigen solver package: " 00056 << solver_package 00057 << std::endl; 00058 libmesh_error(); 00059 } 00060 00061 AutoPtr<EigenSolver<T> > ap(NULL); 00062 return ap; 00063 }
| void libMesh::SlepcEigenSolver< T >::clear | ( | ) | [inline, virtual] |
Release all memory and clear data structures.
Reimplemented from libMesh::EigenSolver< T >.
Definition at line 40 of file slepc_eigen_solver.C.
References libMesh::EigenSolver< T >::_eigen_solver_type, libMesh::SlepcEigenSolver< T >::_eps, libMesh::EigenSolver< T >::_is_initialized, libMeshEnums::ARNOLDI, libMesh::COMM_WORLD, libMesh::EigenSolver< T >::initialized(), and libMeshEnums::KRYLOVSCHUR.
Referenced by libMesh::SlepcEigenSolver< T >::~SlepcEigenSolver().
00041 { 00042 if (this->initialized()) 00043 { 00044 this->_is_initialized = false; 00045 00046 int ierr=0; 00047 00048 ierr = LibMeshEPSDestroy(&_eps); 00049 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00050 00051 // SLEPc default eigenproblem solver 00052 #if SLEPC_VERSION_LESS_THAN(2,3,2) 00053 this->_eigen_solver_type = ARNOLDI; 00054 #else 00055 // Krylov-Schur showed up as of Slepc 2.3.2 00056 this->_eigen_solver_type = KRYLOVSCHUR; 00057 #endif 00058 } 00059 }
| void libMesh::ReferenceCounter::disable_print_counter_info | ( | ) | [static, inherited] |
Definition at line 106 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter.
00107 { 00108 _enable_print_counter = false; 00109 return; 00110 }
| EigenProblemType libMesh::EigenSolver< T >::eigen_problem_type | ( | ) | const [inline, inherited] |
Returns the type of the eigen problem.
Definition at line 98 of file eigen_solver.h.
00098 { return _eigen_problem_type;}
| EigenSolverType libMesh::EigenSolver< T >::eigen_solver_type | ( | ) | const [inline, inherited] |
Returns the type of eigensolver to use.
Definition at line 93 of file eigen_solver.h.
00093 { return _eigen_solver_type; }
| void libMesh::ReferenceCounter::enable_print_counter_info | ( | ) | [static, inherited] |
Methods to enable/disable the reference counter output from print_info()
Definition at line 100 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter.
00101 { 00102 _enable_print_counter = true; 00103 return; 00104 }
| std::pair< Real, Real > libMesh::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 libMesh::EigenSolver< T >.
Definition at line 652 of file slepc_eigen_solver.C.
References libMesh::SlepcEigenSolver< T >::_eps, libMesh::PetscVector< T >::close(), libMesh::COMM_WORLD, and libMesh::PetscVector< T >::vec().
00654 { 00655 int ierr=0; 00656 00657 PetscReal re, im; 00658 00659 // Make sure the NumericVector passed in is really a PetscVector 00660 PetscVector<T>* solution = libmesh_cast_ptr<PetscVector<T>*>(&solution_in); 00661 00662 // real and imaginary part of the ith eigenvalue. 00663 PetscScalar kr, ki; 00664 00665 solution->close(); 00666 00667 ierr = EPSGetEigenpair(_eps, i, &kr, &ki, solution->vec(), PETSC_NULL); 00668 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00669 00670 #ifdef LIBMESH_USE_COMPLEX_NUMBERS 00671 re = PetscRealPart(kr); 00672 im = PetscImaginaryPart(kr); 00673 #else 00674 re = kr; 00675 im = ki; 00676 #endif 00677 00678 return std::make_pair(re, im); 00679 }
| std::pair< Real, Real > libMesh::SlepcEigenSolver< T >::get_eigenvalue | ( | unsigned int | i | ) | [inline, virtual] |
Same as above, but does not copy the eigenvector.
Implements libMesh::EigenSolver< T >.
Definition at line 683 of file slepc_eigen_solver.C.
References libMesh::SlepcEigenSolver< T >::_eps, and libMesh::COMM_WORLD.
00684 { 00685 int ierr=0; 00686 00687 PetscReal re, im; 00688 00689 // real and imaginary part of the ith eigenvalue. 00690 PetscScalar kr, ki; 00691 00692 ierr = EPSGetEigenvalue(_eps, i, &kr, &ki); 00693 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00694 00695 #ifdef LIBMESH_USE_COMPLEX_NUMBERS 00696 re = PetscRealPart(kr); 00697 im = PetscImaginaryPart(kr); 00698 #else 00699 re = kr; 00700 im = ki; 00701 #endif 00702 00703 return std::make_pair(re, im); 00704 }
| std::string libMesh::ReferenceCounter::get_info | ( | ) | [static, inherited] |
Gets a string containing the reference information.
Definition at line 47 of file reference_counter.C.
References libMesh::ReferenceCounter::_counts, and libMesh::Quality::name().
Referenced by libMesh::ReferenceCounter::print_info().
00048 { 00049 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG) 00050 00051 std::ostringstream oss; 00052 00053 oss << '\n' 00054 << " ---------------------------------------------------------------------------- \n" 00055 << "| Reference count information |\n" 00056 << " ---------------------------------------------------------------------------- \n"; 00057 00058 for (Counts::iterator it = _counts.begin(); 00059 it != _counts.end(); ++it) 00060 { 00061 const std::string name(it->first); 00062 const unsigned int creations = it->second.first; 00063 const unsigned int destructions = it->second.second; 00064 00065 oss << "| " << name << " reference count information:\n" 00066 << "| Creations: " << creations << '\n' 00067 << "| Destructions: " << destructions << '\n'; 00068 } 00069 00070 oss << " ---------------------------------------------------------------------------- \n"; 00071 00072 return oss.str(); 00073 00074 #else 00075 00076 return ""; 00077 00078 #endif 00079 }
| Real libMesh::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 708 of file slepc_eigen_solver.C.
References libMesh::SlepcEigenSolver< T >::_eps, and libMesh::COMM_WORLD.
00709 { 00710 int ierr=0; 00711 PetscReal error; 00712 00713 ierr = EPSComputeRelativeError(_eps, i, &error); 00714 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00715 00716 return error; 00717 }
| void libMesh::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 163 of file reference_counter.h.
References libMesh::ReferenceCounter::_counts, and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().
00164 { 00165 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 00166 std::pair<unsigned int, unsigned int>& p = _counts[name]; 00167 00168 p.first++; 00169 }
| void libMesh::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 176 of file reference_counter.h.
References libMesh::ReferenceCounter::_counts, and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().
00177 { 00178 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 00179 std::pair<unsigned int, unsigned int>& p = _counts[name]; 00180 00181 p.second++; 00182 }
| void libMesh::SlepcEigenSolver< T >::init | ( | ) | [inline, virtual] |
Initialize data structures if not done so already.
Implements libMesh::EigenSolver< T >.
Definition at line 64 of file slepc_eigen_solver.C.
References libMesh::SlepcEigenSolver< T >::_eps, libMesh::EigenSolver< T >::_is_initialized, libMesh::COMM_WORLD, libMesh::EigenSolver< T >::initialized(), and libMesh::SlepcEigenSolver< T >::set_slepc_solver_type().
Referenced by libMesh::SlepcEigenSolver< T >::attach_deflation_space(), libMesh::SlepcEigenSolver< T >::solve_generalized(), and libMesh::SlepcEigenSolver< T >::solve_standard().
00065 { 00066 00067 int ierr=0; 00068 00069 // Initialize the data structures if not done so already. 00070 if (!this->initialized()) 00071 { 00072 this->_is_initialized = true; 00073 00074 // Create the eigenproblem solver context 00075 ierr = EPSCreate (libMesh::COMM_WORLD, &_eps); 00076 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00077 00078 // Set user-specified solver 00079 set_slepc_solver_type(); 00080 } 00081 }
| bool libMesh::EigenSolver< T >::initialized | ( | ) | const [inline, inherited] |
- Returns:
- true if the data structures are initialized, false otherwise.
Definition at line 77 of file eigen_solver.h.
Referenced by libMesh::SlepcEigenSolver< T >::clear(), and libMesh::SlepcEigenSolver< T >::init().
00077 { return _is_initialized; }
| static unsigned int libMesh::ReferenceCounter::n_objects | ( | ) | [inline, static, inherited] |
Prints the number of outstanding (created, but not yet destroyed) objects.
Definition at line 79 of file reference_counter.h.
References libMesh::ReferenceCounter::_n_objects.
00080 { return _n_objects; }
| PositionOfSpectrum libMesh::EigenSolver< T >::position_of_spectrum | ( | ) | const [inline, inherited] |
Returns the position of the spectrum to compute.
Definition at line 103 of file eigen_solver.h.
00104 { return _position_of_spectrum;}
| void libMesh::ReferenceCounter::print_info | ( | std::ostream & | out = libMesh::out |
) | [static, inherited] |
Prints the reference information, by default to libMesh::out.
Definition at line 88 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter, and libMesh::ReferenceCounter::get_info().
00089 { 00090 if( _enable_print_counter ) out_stream << ReferenceCounter::get_info(); 00091 }
| void libMesh::EigenSolver< T >::set_eigenproblem_type | ( | EigenProblemType | ept | ) | [inline, inherited] |
Sets the type of the eigenproblem.
Definition at line 115 of file eigen_solver.h.
00116 {_eigen_problem_type = ept;}
| void libMesh::EigenSolver< T >::set_eigensolver_type | ( | const EigenSolverType | est | ) | [inline, inherited] |
Sets the type of eigensolver to use.
Definition at line 109 of file eigen_solver.h.
00110 { _eigen_solver_type = est; }
| void libMesh::EigenSolver< T >::set_position_of_spectrum | ( | PositionOfSpectrum | pos | ) | [inline, inherited] |
Sets the position of the spectrum.
Definition at line 121 of file eigen_solver.h.
00122 {_position_of_spectrum= pos;}
| void libMesh::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 619 of file slepc_eigen_solver.C.
References libMesh::SlepcEigenSolver< T >::_eps, libMesh::EigenSolver< T >::_position_of_spectrum, libMesh::COMM_WORLD, libMesh::err, libMeshEnums::LARGEST_IMAGINARY, libMeshEnums::LARGEST_MAGNITUDE, libMeshEnums::LARGEST_REAL, libMeshEnums::SMALLEST_IMAGINARY, libMeshEnums::SMALLEST_MAGNITUDE, and libMeshEnums::SMALLEST_REAL.
Referenced by libMesh::SlepcEigenSolver< T >::_solve_generalized_helper().
00620 { 00621 int ierr = 0; 00622 00623 switch (this->_position_of_spectrum) 00624 { 00625 case LARGEST_MAGNITUDE: 00626 ierr = EPSSetWhichEigenpairs (_eps, EPS_LARGEST_MAGNITUDE); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00627 case SMALLEST_MAGNITUDE: 00628 ierr = EPSSetWhichEigenpairs (_eps, EPS_SMALLEST_MAGNITUDE); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00629 case LARGEST_REAL: 00630 ierr = EPSSetWhichEigenpairs (_eps, EPS_LARGEST_REAL); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00631 case SMALLEST_REAL: 00632 ierr = EPSSetWhichEigenpairs (_eps, EPS_SMALLEST_REAL); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00633 case LARGEST_IMAGINARY: 00634 ierr = EPSSetWhichEigenpairs (_eps, EPS_LARGEST_IMAGINARY); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00635 case SMALLEST_IMAGINARY: 00636 ierr = EPSSetWhichEigenpairs (_eps, EPS_SMALLEST_IMAGINARY); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00637 00638 00639 default: 00640 libMesh::err << "ERROR: Unsupported SLEPc position of spectrum: " 00641 << this->_position_of_spectrum << std::endl; 00642 libmesh_error(); 00643 } 00644 }
| void libMesh::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 594 of file slepc_eigen_solver.C.
References libMesh::EigenSolver< T >::_eigen_problem_type, libMesh::SlepcEigenSolver< T >::_eps, libMesh::COMM_WORLD, libMesh::err, libMeshEnums::GHEP, libMeshEnums::GNHEP, libMeshEnums::HEP, and libMeshEnums::NHEP.
Referenced by libMesh::SlepcEigenSolver< T >::_solve_generalized_helper().
00595 { 00596 int ierr = 0; 00597 00598 switch (this->_eigen_problem_type) 00599 { 00600 case NHEP: 00601 ierr = EPSSetProblemType (_eps, EPS_NHEP); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00602 case GNHEP: 00603 ierr = EPSSetProblemType (_eps, EPS_GNHEP); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00604 case HEP: 00605 ierr = EPSSetProblemType (_eps, EPS_HEP); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00606 case GHEP: 00607 ierr = EPSSetProblemType (_eps, EPS_GHEP); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00608 00609 default: 00610 libMesh::err << "ERROR: Unsupported SLEPc Eigen Problem: " 00611 << this->_eigen_problem_type << std::endl 00612 << "Continuing with SLEPc defaults" << std::endl; 00613 } 00614 }
| void libMesh::SlepcEigenSolver< T >::set_slepc_solver_type | ( | ) | [inline, private] |
Tells Slepc to use the user-specified solver stored in _eigen_solver_type
Definition at line 559 of file slepc_eigen_solver.C.
References libMesh::EigenSolver< T >::_eigen_solver_type, libMesh::SlepcEigenSolver< T >::_eps, libMeshEnums::ARNOLDI, libMesh::COMM_WORLD, libMesh::err, libMeshEnums::KRYLOVSCHUR, libMeshEnums::LANCZOS, libMeshEnums::LAPACK, libMeshEnums::POWER, and libMeshEnums::SUBSPACE.
Referenced by libMesh::SlepcEigenSolver< T >::init().
00560 { 00561 int ierr = 0; 00562 00563 switch (this->_eigen_solver_type) 00564 { 00565 case POWER: 00566 ierr = EPSSetType (_eps, (char*) EPSPOWER); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00567 case SUBSPACE: 00568 ierr = EPSSetType (_eps, (char*) EPSSUBSPACE); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00569 case LAPACK: 00570 ierr = EPSSetType (_eps, (char*) EPSLAPACK); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00571 case ARNOLDI: 00572 ierr = EPSSetType (_eps, (char*) EPSARNOLDI); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00573 case LANCZOS: 00574 ierr = EPSSetType (_eps, (char*) EPSLANCZOS); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00575 #if !SLEPC_VERSION_LESS_THAN(2,3,2) 00576 // EPSKRYLOVSCHUR added in 2.3.2 00577 case KRYLOVSCHUR: 00578 ierr = EPSSetType (_eps, (char*) EPSKRYLOVSCHUR); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00579 #endif 00580 // case ARPACK: 00581 // ierr = EPSSetType (_eps, (char*) EPSARPACK); CHKERRABORT(libMesh::COMM_WORLD,ierr); return; 00582 00583 default: 00584 libMesh::err << "ERROR: Unsupported SLEPc Eigen Solver: " 00585 << this->_eigen_solver_type << std::endl 00586 << "Continuing with SLEPc defaults" << std::endl; 00587 } 00588 }
| std::pair< unsigned int, unsigned int > libMesh::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 libMesh::EigenSolver< T >.
Definition at line 380 of file slepc_eigen_solver.C.
References libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_mult(), libMesh::SlepcEigenSolver< T >::_solve_generalized_helper(), libMesh::COMM_WORLD, libMesh::SlepcEigenSolver< T >::init(), libMesh::ShellMatrix< T >::m(), and libMesh::ShellMatrix< T >::n().
00386 { 00387 this->init (); 00388 00389 int ierr=0; 00390 00391 // Prepare the matrix. 00392 Mat mat_A; 00393 ierr = MatCreateShell(libMesh::COMM_WORLD, 00394 shell_matrix_A.m(), // Specify the number of local rows 00395 shell_matrix_A.n(), // Specify the number of local columns 00396 PETSC_DETERMINE, 00397 PETSC_DETERMINE, 00398 const_cast<void*>(static_cast<const void*>(&shell_matrix_A)), 00399 &mat_A); 00400 00401 Mat mat_B; 00402 ierr = MatCreateShell(libMesh::COMM_WORLD, 00403 shell_matrix_B.m(), // Specify the number of local rows 00404 shell_matrix_B.n(), // Specify the number of local columns 00405 PETSC_DETERMINE, 00406 PETSC_DETERMINE, 00407 const_cast<void*>(static_cast<const void*>(&shell_matrix_B)), 00408 &mat_B); 00409 00410 /* Note that the const_cast above is only necessary because PETSc 00411 does not accept a const void*. Inside the member function 00412 _petsc_shell_matrix() below, the pointer is casted back to a 00413 const ShellMatrix<T>*. */ 00414 00415 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00416 ierr = MatShellSetOperation(mat_A,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult)); 00417 ierr = MatShellSetOperation(mat_A,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal)); 00418 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00419 00420 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00421 ierr = MatShellSetOperation(mat_B,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult)); 00422 ierr = MatShellSetOperation(mat_B,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal)); 00423 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00424 00425 return _solve_generalized_helper (mat_A, mat_B, nev, ncv, tol, m_its); 00426 }
| std::pair< unsigned int, unsigned int > libMesh::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 libMesh::EigenSolver< T >.
Definition at line 338 of file slepc_eigen_solver.C.
References libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_mult(), libMesh::SlepcEigenSolver< T >::_solve_generalized_helper(), libMesh::PetscMatrix< T >::close(), libMesh::COMM_WORLD, libMesh::SlepcEigenSolver< T >::init(), libMesh::ShellMatrix< T >::m(), libMesh::PetscMatrix< T >::mat(), and libMesh::ShellMatrix< T >::n().
00344 { 00345 this->init (); 00346 00347 int ierr=0; 00348 00349 PetscMatrix<T>* matrix_A = libmesh_cast_ptr<PetscMatrix<T>*>(&matrix_A_in); 00350 00351 // Close the matrix and vectors in case this wasn't already done. 00352 matrix_A->close (); 00353 00354 // Prepare the matrix. 00355 Mat mat_B; 00356 ierr = MatCreateShell(libMesh::COMM_WORLD, 00357 shell_matrix_B.m(), // Specify the number of local rows 00358 shell_matrix_B.n(), // Specify the number of local columns 00359 PETSC_DETERMINE, 00360 PETSC_DETERMINE, 00361 const_cast<void*>(static_cast<const void*>(&shell_matrix_B)), 00362 &mat_B); 00363 00364 00365 /* Note that the const_cast above is only necessary because PETSc 00366 does not accept a const void*. Inside the member function 00367 _petsc_shell_matrix() below, the pointer is casted back to a 00368 const ShellMatrix<T>*. */ 00369 00370 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00371 ierr = MatShellSetOperation(mat_B,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult)); 00372 ierr = MatShellSetOperation(mat_B,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal)); 00373 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00374 00375 return _solve_generalized_helper (matrix_A->mat(), mat_B, nev, ncv, tol, m_its); 00376 }
| std::pair< unsigned int, unsigned int > libMesh::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 libMesh::EigenSolver< T >.
Definition at line 297 of file slepc_eigen_solver.C.
References libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_mult(), libMesh::SlepcEigenSolver< T >::_solve_generalized_helper(), libMesh::COMM_WORLD, libMesh::SlepcEigenSolver< T >::init(), libMesh::ShellMatrix< T >::m(), and libMesh::ShellMatrix< T >::n().
00303 { 00304 this->init (); 00305 00306 int ierr=0; 00307 00308 // Prepare the matrix. 00309 Mat mat_A; 00310 ierr = MatCreateShell(libMesh::COMM_WORLD, 00311 shell_matrix_A.m(), // Specify the number of local rows 00312 shell_matrix_A.n(), // Specify the number of local columns 00313 PETSC_DETERMINE, 00314 PETSC_DETERMINE, 00315 const_cast<void*>(static_cast<const void*>(&shell_matrix_A)), 00316 &mat_A); 00317 00318 PetscMatrix<T>* matrix_B = libmesh_cast_ptr<PetscMatrix<T>*>(&matrix_B_in); 00319 00320 // Close the matrix and vectors in case this wasn't already done. 00321 matrix_B->close (); 00322 00323 /* Note that the const_cast above is only necessary because PETSc 00324 does not accept a const void*. Inside the member function 00325 _petsc_shell_matrix() below, the pointer is casted back to a 00326 const ShellMatrix<T>*. */ 00327 00328 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00329 ierr = MatShellSetOperation(mat_A,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult)); 00330 ierr = MatShellSetOperation(mat_A,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal)); 00331 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00332 00333 return _solve_generalized_helper (mat_A, matrix_B->mat(), nev, ncv, tol, m_its); 00334 }
| std::pair< unsigned int, unsigned int > libMesh::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 libMesh::EigenSolver< T >.
Definition at line 274 of file slepc_eigen_solver.C.
References libMesh::SlepcEigenSolver< T >::_solve_generalized_helper(), libMesh::PetscMatrix< T >::close(), libMesh::SlepcEigenSolver< T >::init(), and libMesh::PetscMatrix< T >::mat().
00280 { 00281 this->init (); 00282 00283 // Make sure the data passed in are really of Petsc types 00284 PetscMatrix<T>* matrix_A = libmesh_cast_ptr<PetscMatrix<T>*>(&matrix_A_in); 00285 PetscMatrix<T>* matrix_B = libmesh_cast_ptr<PetscMatrix<T>*>(&matrix_B_in); 00286 00287 // Close the matrix and vectors in case this wasn't already done. 00288 matrix_A->close (); 00289 matrix_B->close (); 00290 00291 00292 return _solve_generalized_helper (matrix_A->mat(), matrix_B->mat(), nev, ncv, tol, m_its); 00293 }
| std::pair< unsigned int, unsigned int > libMesh::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 libMesh::EigenSolver< T >.
Definition at line 117 of file slepc_eigen_solver.C.
References libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_mult(), libMesh::SlepcEigenSolver< T >::_solve_standard_helper(), libMesh::COMM_WORLD, libMesh::SlepcEigenSolver< T >::init(), libMesh::ShellMatrix< T >::m(), and libMesh::ShellMatrix< T >::n().
00122 { 00123 this->init (); 00124 00125 int ierr=0; 00126 00127 // Prepare the matrix. 00128 Mat mat; 00129 ierr = MatCreateShell(libMesh::COMM_WORLD, 00130 shell_matrix.m(), // Specify the number of local rows 00131 shell_matrix.n(), // Specify the number of local columns 00132 PETSC_DETERMINE, 00133 PETSC_DETERMINE, 00134 const_cast<void*>(static_cast<const void*>(&shell_matrix)), 00135 &mat); 00136 00137 /* Note that the const_cast above is only necessary because PETSc 00138 does not accept a const void*. Inside the member function 00139 _petsc_shell_matrix() below, the pointer is casted back to a 00140 const ShellMatrix<T>*. */ 00141 00142 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00143 ierr = MatShellSetOperation(mat,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult)); 00144 ierr = MatShellSetOperation(mat,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal)); 00145 CHKERRABORT(libMesh::COMM_WORLD,ierr); 00146 00147 00148 return _solve_standard_helper(mat, nev, ncv, tol, m_its); 00149 }
| std::pair< unsigned int, unsigned int > libMesh::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 libMesh::EigenSolver< T >.
Definition at line 87 of file slepc_eigen_solver.C.
References libMesh::SlepcEigenSolver< T >::_solve_standard_helper(), libMesh::PetscMatrix< T >::close(), libMesh::SlepcEigenSolver< T >::init(), and libMesh::PetscMatrix< T >::mat().
00092 { 00093 // START_LOG("solve_standard()", "SlepcEigenSolver"); 00094 00095 this->init (); 00096 00097 // Make sure the SparseMatrix passed in is really a PetscMatrix 00098 PetscMatrix<T>* matrix_A = libmesh_cast_ptr<PetscMatrix<T>*>(&matrix_A_in); 00099 00100 // Close the matrix and vectors in case this wasn't already done. 00101 matrix_A->close (); 00102 00103 // just for debugging, remove this 00104 // char mat_file[] = "matA.petsc"; 00105 // PetscViewer petsc_viewer; 00106 // ierr = PetscViewerBinaryOpen(libMesh::COMM_WORLD, mat_file, PETSC_FILE_CREATE, &petsc_viewer); 00107 // CHKERRABORT(libMesh::COMM_WORLD,ierr); 00108 // ierr = MatView(matrix_A->mat(),petsc_viewer); 00109 // CHKERRABORT(libMesh::COMM_WORLD,ierr); 00110 00111 return _solve_standard_helper(matrix_A->mat(), nev, ncv, tol, m_its); 00112 }
Member Data Documentation
ReferenceCounter::Counts libMesh::ReferenceCounter::_counts [static, protected, inherited] |
Actually holds the data.
Definition at line 118 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::get_info(), libMesh::ReferenceCounter::increment_constructor_count(), and libMesh::ReferenceCounter::increment_destructor_count().
EigenProblemType libMesh::EigenSolver< T >::_eigen_problem_type [protected, inherited] |
Enum stating which type of eigen problem we deal with.
Definition at line 222 of file eigen_solver.h.
Referenced by libMesh::EigenSolver< Number >::eigen_problem_type(), libMesh::EigenSolver< Number >::set_eigenproblem_type(), libMesh::SlepcEigenSolver< T >::set_slepc_problem_type(), and libMesh::SlepcEigenSolver< T >::SlepcEigenSolver().
EigenSolverType libMesh::EigenSolver< T >::_eigen_solver_type [protected, inherited] |
Enum stating which type of eigensolver to use.
Definition at line 217 of file eigen_solver.h.
Referenced by libMesh::SlepcEigenSolver< T >::clear(), libMesh::EigenSolver< Number >::eigen_solver_type(), libMesh::EigenSolver< Number >::set_eigensolver_type(), libMesh::SlepcEigenSolver< T >::set_slepc_solver_type(), and libMesh::SlepcEigenSolver< T >::SlepcEigenSolver().
bool libMesh::ReferenceCounter::_enable_print_counter = true [static, protected, inherited] |
Flag to control whether reference count information is printed when print_info is called.
Definition at line 137 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::disable_print_counter_info(), libMesh::ReferenceCounter::enable_print_counter_info(), and libMesh::ReferenceCounter::print_info().
EPS libMesh::SlepcEigenSolver< T >::_eps [private] |
Eigenproblem solver context
Definition at line 254 of file slepc_eigen_solver.h.
Referenced by libMesh::SlepcEigenSolver< T >::_solve_generalized_helper(), libMesh::SlepcEigenSolver< T >::attach_deflation_space(), libMesh::SlepcEigenSolver< T >::clear(), libMesh::SlepcEigenSolver< T >::get_eigenpair(), libMesh::SlepcEigenSolver< T >::get_eigenvalue(), libMesh::SlepcEigenSolver< T >::get_relative_error(), libMesh::SlepcEigenSolver< T >::init(), libMesh::SlepcEigenSolver< T >::set_slepc_position_of_spectrum(), libMesh::SlepcEigenSolver< T >::set_slepc_problem_type(), and libMesh::SlepcEigenSolver< T >::set_slepc_solver_type().
bool libMesh::EigenSolver< T >::_is_initialized [protected, inherited] |
Flag indicating if the data structures have been initialized.
Definition at line 232 of file eigen_solver.h.
Referenced by libMesh::SlepcEigenSolver< T >::clear(), libMesh::SlepcEigenSolver< T >::init(), and libMesh::EigenSolver< Number >::initialized().
Threads::spin_mutex libMesh::ReferenceCounter::_mutex [static, protected, inherited] |
Mutual exclusion object to enable thread-safe reference counting.
Definition at line 131 of file reference_counter.h.
Threads::atomic< unsigned int > libMesh::ReferenceCounter::_n_objects [static, protected, inherited] |
The number of objects. Print the reference count information when the number returns to 0.
Definition at line 126 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::n_objects(), libMesh::ReferenceCounter::ReferenceCounter(), and libMesh::ReferenceCounter::~ReferenceCounter().
PositionOfSpectrum libMesh::EigenSolver< T >::_position_of_spectrum [protected, inherited] |
Enum stating where to evaluate the spectrum.
Definition at line 227 of file eigen_solver.h.
Referenced by libMesh::EigenSolver< Number >::position_of_spectrum(), libMesh::EigenSolver< Number >::set_position_of_spectrum(), and libMesh::SlepcEigenSolver< T >::set_slepc_position_of_spectrum().
The documentation for this class was generated from the following files:
Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:39 UTC
Hosted By: