|
ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
|
Functions | |
| template<typename T > | |
| void | diago_hs_para (T *h, T *s, const int lda, const int nband, typename GetTypeReal< T >::type *const ekb, T *const wfc, const MPI_Comm &comm, const int diag_subspace, const int block_size=0) |
| Parallel do the generalized eigenvalue problem. | |
| template void | diago_hs_para< double > (double *h, double *s, const int lda, const int nband, typename GetTypeReal< double >::type *const ekb, double *const wfc, const MPI_Comm &comm, const int diag_subspace, const int block_size) |
| template void | diago_hs_para< std::complex< double > > (std::complex< double > *h, std::complex< double > *s, const int lda, const int nband, typename GetTypeReal< std::complex< double > >::type *const ekb, std::complex< double > *const wfc, const MPI_Comm &comm, const int diag_subspace, const int block_size) |
| template void | diago_hs_para< float > (float *h, float *s, const int lda, const int nband, typename GetTypeReal< float >::type *const ekb, float *const wfc, const MPI_Comm &comm, const int diag_subspace, const int block_size) |
| template void | diago_hs_para< std::complex< float > > (std::complex< float > *h, std::complex< float > *s, const int lda, const int nband, typename GetTypeReal< std::complex< float > >::type *const ekb, std::complex< float > *const wfc, const MPI_Comm &comm, const int diag_subspace, const int block_size) |
| template<typename T , typename Device > | |
| void | setup_diago_params_pw (const int istep, const int iter, const double ethr, const Input_para &inp) |
| Setup diagonalization parameters for PW method. | |
| template<typename T , typename Device > | |
| void | setup_diago_params_sdft (const int istep, const int iter, const double ethr, const Input_para &inp) |
| Setup diagonalization parameters for SDFT method. | |
| template void | setup_diago_params_pw< double, base_device::DEVICE_CPU > (const int istep, const int iter, const double ethr, const Input_para &inp) |
| template void | setup_diago_params_sdft< double, base_device::DEVICE_CPU > (const int istep, const int iter, const double ethr, const Input_para &inp) |
| void | pxxxgvx_ (const int *itype, const char *jobz, const char *range, const char *uplo, const int *n, double *A, const int *ia, const int *ja, const int *desca, double *B, const int *ib, const int *jb, const int *descb, const double *vl, const double *vu, const int *il, const int *iu, const double *abstol, int *m, int *nz, double *w, const double *orfac, double *Z, const int *iz, const int *jz, const int *descz, double *work, int *lwork, double *rwork, int *lrwork, int *iwork, int *liwork, int *ifail, int *iclustr, double *gap, int *info) |
| Wrapper function for Scalapack's generalized eigensolver routines. | |
| void | pxxxgvx_ (const int *itype, const char *jobz, const char *range, const char *uplo, const int *n, std::complex< double > *A, const int *ia, const int *ja, const int *desca, std::complex< double > *B, const int *ib, const int *jb, const int *descb, const double *vl, const double *vu, const int *il, const int *iu, const double *abstol, int *m, int *nz, double *w, const double *orfac, std::complex< double > *Z, const int *iz, const int *jz, const int *descz, std::complex< double > *work, int *lwork, double *rwork, int *lrwork, int *iwork, int *liwork, int *ifail, int *iclustr, double *gap, int *info) |
| void | pxxxgvx_ (const int *itype, const char *jobz, const char *range, const char *uplo, const int *n, float *A, const int *ia, const int *ja, const int *desca, float *B, const int *ib, const int *jb, const int *descb, const float *vl, const float *vu, const int *il, const int *iu, const float *abstol, int *m, int *nz, float *w, const float *orfac, float *Z, const int *iz, const int *jz, const int *descz, float *work, int *lwork, float *rwork, int *lrwork, int *iwork, int *liwork, int *ifail, int *iclustr, float *gap, int *info) |
| void | pxxxgvx_ (const int *itype, const char *jobz, const char *range, const char *uplo, const int *n, std::complex< float > *A, const int *ia, const int *ja, const int *desca, std::complex< float > *B, const int *ib, const int *jb, const int *descb, const float *vl, const float *vu, const int *il, const int *iu, const float *abstol, int *m, int *nz, float *w, const float *orfac, std::complex< float > *Z, const int *iz, const int *jz, const int *descz, std::complex< float > *work, int *lwork, float *rwork, int *lrwork, int *iwork, int *liwork, int *ifail, int *iclustr, float *gap, int *info) |
| void | pxxxgvx_post_processing (const int info, const std::vector< int > &ifail, const std::vector< int > &iclustr, const int M, const int NZ, const int nbands, int °eneracy_max) |
| void | get_lwork (int &lwork, std::vector< double > &work) |
| void | get_lwork (int &lwork, std::vector< float > &work) |
| void | get_lwork (int &lwork, std::vector< std::complex< double > > &work) |
| void | get_lwork (int &lwork, std::vector< std::complex< float > > &work) |
| template<typename T > | |
| void | pxxxgvx_diag (const int *const desc, const int ncol, const int nrow, const int nbands, const T *const h_mat, const T *const s_mat, typename GetTypeReal< T >::type *const ekb, T *const wfc_2d) |
| Wrapper function for Scalapack's generalized eigensolver routines: pdsygvx_, pzhegvx_, pssygvx_, pchegvx_. | |
| template void | pxxxgvx_diag (const int *const desc, const int ncol, const int nrow, const int nbands, const double *const h_mat, const double *const s_mat, double *const ekb, double *const wfc_2d) |
| template void | pxxxgvx_diag (const int *const desc, const int ncol, const int nrow, const int nbands, const std::complex< double > *const h_mat, const std::complex< double > *const s_mat, double *const ekb, std::complex< double > *const wfc_2d) |
| template void | pxxxgvx_diag (const int *const desc, const int ncol, const int nrow, const int nbands, const float *const h_mat, const float *const s_mat, float *const ekb, float *const wfc_2d) |
| template void | pxxxgvx_diag (const int *const desc, const int ncol, const int nrow, const int nbands, const std::complex< float > *const h_mat, const std::complex< float > *const s_mat, float *const ekb, std::complex< float > *const wfc_2d) |
| double | set_diagethr_ks (const std::string basis_type, const std::string esolver_type, const std::string calculation_in, const std::string init_chg_in, const std::string precision_flag_in, const int istep, const int iter, const double drho, const double pw_diag_thr_init, const double diag_ethr_in, const double nelec_in) |
| double | set_diagethr_sdft (const std::string basis_type, const std::string esolver_type, const std::string calculation_in, const std::string init_chg_in, const int istep, const int iter, const double drho, const double pw_diag_thr_init, const double diag_ethr_in, const int nband_in, const double stoiter_ks_ne_in) |
| double | reset_diag_ethr (std::ofstream &ofs_running, const std::string basis_type, const std::string esolver_type, const std::string precision_flag_in, const double hsover_error, const double drho_in, const double diag_ethr_in, const double nelec_in) |
| double | cal_hsolve_error (const std::string basis_type, const std::string esolver_type, const double diag_ethr_in, const double nelec_in) |
| double | get_real (const std::complex< double > &x) |
| float | get_real (const std::complex< float > &x) |
| double | get_real (const double &x) |
| float | get_real (const float &x) |
This is the module for wrapper of DeNse Generalized eigenValue (eXtended) HErmitian / SYmmetric
Tested function:
| double hsolver::cal_hsolve_error | ( | const std::string | basis_type, |
| const std::string | esolver_type, | ||
| const double | diag_ethr_in, | ||
| const double | nelec_in | ||
| ) |
| void hsolver::diago_hs_para | ( | T * | h, |
| T * | s, | ||
| const int | lda, | ||
| const int | nband, | ||
| typename GetTypeReal< T >::type *const | ekb, | ||
| T *const | wfc, | ||
| const MPI_Comm & | comm, | ||
| const int | diag_subspace, | ||
| const int | block_size = 0 |
||
| ) |
Parallel do the generalized eigenvalue problem.
| T | double or std::complex<double> or float or complex<float> |
| H | the hermitian matrix H. |
| S | the overlap matrix S. |
| lda | the leading dimension of H and S |
| nband | the number of bands to be calculated |
| ekb | to store the eigenvalues. |
| wfc | to store the eigenvectors |
| comm | the communicator |
| diag_subspace | the method to solve the generalized eigenvalue problem |
| block_size | the block size in 2d block cyclic distribution if use elpa or scalapack. |
| template void hsolver::diago_hs_para< double > | ( | double * | h, |
| double * | s, | ||
| const int | lda, | ||
| const int | nband, | ||
| typename GetTypeReal< double >::type *const | ekb, | ||
| double *const | wfc, | ||
| const MPI_Comm & | comm, | ||
| const int | diag_subspace, | ||
| const int | block_size | ||
| ) |
| template void hsolver::diago_hs_para< float > | ( | float * | h, |
| float * | s, | ||
| const int | lda, | ||
| const int | nband, | ||
| typename GetTypeReal< float >::type *const | ekb, | ||
| float *const | wfc, | ||
| const MPI_Comm & | comm, | ||
| const int | diag_subspace, | ||
| const int | block_size | ||
| ) |
| template void hsolver::diago_hs_para< std::complex< double > > | ( | std::complex< double > * | h, |
| std::complex< double > * | s, | ||
| const int | lda, | ||
| const int | nband, | ||
| typename GetTypeReal< std::complex< double > >::type *const | ekb, | ||
| std::complex< double > *const | wfc, | ||
| const MPI_Comm & | comm, | ||
| const int | diag_subspace, | ||
| const int | block_size | ||
| ) |
| template void hsolver::diago_hs_para< std::complex< float > > | ( | std::complex< float > * | h, |
| std::complex< float > * | s, | ||
| const int | lda, | ||
| const int | nband, | ||
| typename GetTypeReal< std::complex< float > >::type *const | ekb, | ||
| std::complex< float > *const | wfc, | ||
| const MPI_Comm & | comm, | ||
| const int | diag_subspace, | ||
| const int | block_size | ||
| ) |
| void hsolver::get_lwork | ( | int & | lwork, |
| std::vector< double > & | work | ||
| ) |
| void hsolver::get_lwork | ( | int & | lwork, |
| std::vector< float > & | work | ||
| ) |
| void hsolver::get_lwork | ( | int & | lwork, |
| std::vector< std::complex< double > > & | work | ||
| ) |
| void hsolver::get_lwork | ( | int & | lwork, |
| std::vector< std::complex< float > > & | work | ||
| ) |
|
inline |
|
inline |
|
inline |
|
inline |
| void hsolver::pxxxgvx_ | ( | const int * | itype, |
| const char * | jobz, | ||
| const char * | range, | ||
| const char * | uplo, | ||
| const int * | n, | ||
| double * | A, | ||
| const int * | ia, | ||
| const int * | ja, | ||
| const int * | desca, | ||
| double * | B, | ||
| const int * | ib, | ||
| const int * | jb, | ||
| const int * | descb, | ||
| const double * | vl, | ||
| const double * | vu, | ||
| const int * | il, | ||
| const int * | iu, | ||
| const double * | abstol, | ||
| int * | m, | ||
| int * | nz, | ||
| double * | w, | ||
| const double * | orfac, | ||
| double * | Z, | ||
| const int * | iz, | ||
| const int * | jz, | ||
| const int * | descz, | ||
| double * | work, | ||
| int * | lwork, | ||
| double * | rwork, | ||
| int * | lrwork, | ||
| int * | iwork, | ||
| int * | liwork, | ||
| int * | ifail, | ||
| int * | iclustr, | ||
| double * | gap, | ||
| int * | info | ||
| ) |
Wrapper function for Scalapack's generalized eigensolver routines.
| itype | Specifies the problem type to be solved. |
| jobz | Specifies whether to compute eigenvectors. |
| range | Specifies the range of eigenvalues to be found. |
| uplo | Specifies whether the upper or lower triangular part of the matrices is referenced. |
| n | The order of the matrices A and B. |
| A | The array containing the matrix A. |
| ia | The row index in the global array A. |
| ja | The column index in the global array A. |
| desca | The array descriptor for the distributed matrix A. |
| B | The array containing the matrix B. |
| ib | The row index in the global array B. |
| jb | The column index in the global array B. |
| descb | The array descriptor for the distributed matrix B. |
| vl | Lower bound of the interval to be searched for eigenvalues. |
| vu | Upper bound of the interval to be searched for eigenvalues. |
| il | Index of the smallest eigenvalue to be returned. |
| iu | Index of the largest eigenvalue to be returned. |
| abstol | The absolute error tolerance for the eigenvalues. |
| m | The total number of eigenvalues found. |
| nz | The total number of eigenvalues found in the interval (vl, vu]. |
| w | The array to store the eigenvalues. |
| orfac | The orthogonality factor. |
| Z | The array to store the eigenvectors. |
| iz | The row index in the global array Z. |
| jz | The column index in the global array Z. |
| descz | The array descriptor for the distributed matrix Z. |
| work | Workspace array. |
| lwork | The dimension of the array work. |
| rwork | Workspace array (not used in this function). |
| lrwork | The dimension of the array rwork (not used in this function). |
| iwork | Workspace array. |
| liwork | The dimension of the array iwork. |
| ifail | The array to store the indices of the eigenvectors that failed to converge. |
| iclustr | The array to store the indices of the eigenvalue clusters. |
| gap | The array to store the gaps between eigenvalue clusters. |
| info | Output status of the computation. |
| void hsolver::pxxxgvx_ | ( | const int * | itype, |
| const char * | jobz, | ||
| const char * | range, | ||
| const char * | uplo, | ||
| const int * | n, | ||
| float * | A, | ||
| const int * | ia, | ||
| const int * | ja, | ||
| const int * | desca, | ||
| float * | B, | ||
| const int * | ib, | ||
| const int * | jb, | ||
| const int * | descb, | ||
| const float * | vl, | ||
| const float * | vu, | ||
| const int * | il, | ||
| const int * | iu, | ||
| const float * | abstol, | ||
| int * | m, | ||
| int * | nz, | ||
| float * | w, | ||
| const float * | orfac, | ||
| float * | Z, | ||
| const int * | iz, | ||
| const int * | jz, | ||
| const int * | descz, | ||
| float * | work, | ||
| int * | lwork, | ||
| float * | rwork, | ||
| int * | lrwork, | ||
| int * | iwork, | ||
| int * | liwork, | ||
| int * | ifail, | ||
| int * | iclustr, | ||
| float * | gap, | ||
| int * | info | ||
| ) |
| void hsolver::pxxxgvx_ | ( | const int * | itype, |
| const char * | jobz, | ||
| const char * | range, | ||
| const char * | uplo, | ||
| const int * | n, | ||
| std::complex< double > * | A, | ||
| const int * | ia, | ||
| const int * | ja, | ||
| const int * | desca, | ||
| std::complex< double > * | B, | ||
| const int * | ib, | ||
| const int * | jb, | ||
| const int * | descb, | ||
| const double * | vl, | ||
| const double * | vu, | ||
| const int * | il, | ||
| const int * | iu, | ||
| const double * | abstol, | ||
| int * | m, | ||
| int * | nz, | ||
| double * | w, | ||
| const double * | orfac, | ||
| std::complex< double > * | Z, | ||
| const int * | iz, | ||
| const int * | jz, | ||
| const int * | descz, | ||
| std::complex< double > * | work, | ||
| int * | lwork, | ||
| double * | rwork, | ||
| int * | lrwork, | ||
| int * | iwork, | ||
| int * | liwork, | ||
| int * | ifail, | ||
| int * | iclustr, | ||
| double * | gap, | ||
| int * | info | ||
| ) |
| void hsolver::pxxxgvx_ | ( | const int * | itype, |
| const char * | jobz, | ||
| const char * | range, | ||
| const char * | uplo, | ||
| const int * | n, | ||
| std::complex< float > * | A, | ||
| const int * | ia, | ||
| const int * | ja, | ||
| const int * | desca, | ||
| std::complex< float > * | B, | ||
| const int * | ib, | ||
| const int * | jb, | ||
| const int * | descb, | ||
| const float * | vl, | ||
| const float * | vu, | ||
| const int * | il, | ||
| const int * | iu, | ||
| const float * | abstol, | ||
| int * | m, | ||
| int * | nz, | ||
| float * | w, | ||
| const float * | orfac, | ||
| std::complex< float > * | Z, | ||
| const int * | iz, | ||
| const int * | jz, | ||
| const int * | descz, | ||
| std::complex< float > * | work, | ||
| int * | lwork, | ||
| float * | rwork, | ||
| int * | lrwork, | ||
| int * | iwork, | ||
| int * | liwork, | ||
| int * | ifail, | ||
| int * | iclustr, | ||
| float * | gap, | ||
| int * | info | ||
| ) |
| template void hsolver::pxxxgvx_diag | ( | const int *const | desc, |
| const int | ncol, | ||
| const int | nrow, | ||
| const int | nbands, | ||
| const double *const | h_mat, | ||
| const double *const | s_mat, | ||
| double *const | ekb, | ||
| double *const | wfc_2d | ||
| ) |
| template void hsolver::pxxxgvx_diag | ( | const int *const | desc, |
| const int | ncol, | ||
| const int | nrow, | ||
| const int | nbands, | ||
| const float *const | h_mat, | ||
| const float *const | s_mat, | ||
| float *const | ekb, | ||
| float *const | wfc_2d | ||
| ) |
| template void hsolver::pxxxgvx_diag | ( | const int *const | desc, |
| const int | ncol, | ||
| const int | nrow, | ||
| const int | nbands, | ||
| const std::complex< double > *const | h_mat, | ||
| const std::complex< double > *const | s_mat, | ||
| double *const | ekb, | ||
| std::complex< double > *const | wfc_2d | ||
| ) |
| template void hsolver::pxxxgvx_diag | ( | const int *const | desc, |
| const int | ncol, | ||
| const int | nrow, | ||
| const int | nbands, | ||
| const std::complex< float > *const | h_mat, | ||
| const std::complex< float > *const | s_mat, | ||
| float *const | ekb, | ||
| std::complex< float > *const | wfc_2d | ||
| ) |
| void hsolver::pxxxgvx_diag | ( | const int *const | desc, |
| const int | ncol, | ||
| const int | nrow, | ||
| const int | nbands, | ||
| const T *const | h_mat, | ||
| const T *const | s_mat, | ||
| typename GetTypeReal< T >::type *const | ekb, | ||
| T *const | wfc_2d | ||
| ) |
Wrapper function for Scalapack's generalized eigensolver routines: pdsygvx_, pzhegvx_, pssygvx_, pchegvx_.
| desc | the descriptor of scalapack descriptor |
| ncol | the number of columns of the H/S matrix in current processor |
| nrow | the number of rows of the H/S matrix in current processor |
| nbands | the number of bands to be solved |
| h_mat | the Hamiltonian matrix |
| s_mat | the overlap matrix |
| ekb | the eigenvalues |
| wfc_2d | the eigenvectors in 2D block cyclic distribution |
| void hsolver::pxxxgvx_post_processing | ( | const int | info, |
| const std::vector< int > & | ifail, | ||
| const std::vector< int > & | iclustr, | ||
| const int | M, | ||
| const int | NZ, | ||
| const int | nbands, | ||
| int & | degeneracy_max | ||
| ) |
| double hsolver::reset_diag_ethr | ( | std::ofstream & | ofs_running, |
| const std::string | basis_type, | ||
| const std::string | esolver_type, | ||
| const std::string | precision_flag_in, | ||
| const double | hsover_error, | ||
| const double | drho_in, | ||
| const double | diag_ethr_in, | ||
| const double | nelec_in | ||
| ) |
| double hsolver::set_diagethr_ks | ( | const std::string | basis_type, |
| const std::string | esolver_type, | ||
| const std::string | calculation_in, | ||
| const std::string | init_chg_in, | ||
| const std::string | precision_flag_in, | ||
| const int | istep, | ||
| const int | iter, | ||
| const double | drho, | ||
| const double | pw_diag_thr_init, | ||
| const double | diag_ethr_in, | ||
| const double | nelec_in | ||
| ) |
| double hsolver::set_diagethr_sdft | ( | const std::string | basis_type, |
| const std::string | esolver_type, | ||
| const std::string | calculation_in, | ||
| const std::string | init_chg_in, | ||
| const int | istep, | ||
| const int | iter, | ||
| const double | drho, | ||
| const double | pw_diag_thr_init, | ||
| const double | diag_ethr_in, | ||
| const int | nband_in, | ||
| const double | stoiter_ks_ne_in | ||
| ) |
| void hsolver::setup_diago_params_pw | ( | const int | istep, |
| const int | iter, | ||
| const double | ethr, | ||
| const Input_para & | inp | ||
| ) |
Setup diagonalization parameters for PW method.
Template instantiation for CPU.
This function sets up the diagonalization parameters for plane wave method, including subspace diagonalization flag, SCF iteration number, diagonalization threshold, and maximum number of diagonalization steps.
| istep | Current ionic step |
| iter | Current SCF iteration |
| ethr | Diagonalization threshold |
| inp | Input parameters |
choose if psi should be diag in subspace be careful that istep start from 0 and iter start from 1
| template void hsolver::setup_diago_params_pw< double, base_device::DEVICE_CPU > | ( | const int | istep, |
| const int | iter, | ||
| const double | ethr, | ||
| const Input_para & | inp | ||
| ) |
| void hsolver::setup_diago_params_sdft | ( | const int | istep, |
| const int | iter, | ||
| const double | ethr, | ||
| const Input_para & | inp | ||
| ) |
Setup diagonalization parameters for SDFT method.
Template instantiation for GPU.
This function sets up the diagonalization parameters for stochastic DFT method, including subspace diagonalization flag, diagonalization threshold, and maximum number of diagonalization steps.
| istep | Current ionic step |
| iter | Current SCF iteration |
| ethr | Diagonalization threshold |
| inp | Input parameters |
Template instantiation for SDFT CPU
choose if psi should be diag in subspace be careful that istep start from 0 and iter start from 1
| template void hsolver::setup_diago_params_sdft< double, base_device::DEVICE_CPU > | ( | const int | istep, |
| const int | iter, | ||
| const double | ethr, | ||
| const Input_para & | inp | ||
| ) |