|
| Chebyshev (const int norder) |
|
| ~Chebyshev () |
|
void | calcoef_real (std::function< REAL(REAL)> fun) |
|
void | calcoef_complex (std::function< std::complex< REAL >(std::complex< REAL >)> fun) |
|
void | calcoef_pair (std::function< REAL(REAL)> fun1, std::function< REAL(REAL)> fun2) |
|
void | calfinalvec_real (std::function< void(REAL *, REAL *, const int)> funA, REAL *wavein, REAL *waveout, const int N, const int LDA=1, const int m=1) |
|
void | calfinalvec_real (std::function< void(std::complex< REAL > *, std::complex< REAL > *, const int)> funA, std::complex< REAL > *wavein, std::complex< REAL > *waveout, const int N, const int LDA=1, const int m=1) |
|
void | calfinalvec_complex (std::function< void(std::complex< REAL > *, std::complex< REAL > *, const int)> funA, std::complex< REAL > *wavein, std::complex< REAL > *waveout, const int N, const int LDA=1, const int m=1) |
|
void | tracepolyA (std::function< void(std::complex< REAL > *in, std::complex< REAL > *out, const int)> funA, std::complex< REAL > *wavein, const int N, const int LDA=1, const int m=1) |
|
void | getpolyval (REAL x, REAL *polyval, const int N) |
|
void | calpolyvec_real (std::function< void(REAL *in, REAL *out, const int)> funA, REAL *wavein, REAL *waveout, const int N, const int LDA=1, const int m=1) |
|
void | calpolyvec_complex (std::function< void(std::complex< REAL > *in, std::complex< REAL > *out, const int)> funA, std::complex< REAL > *wavein, std::complex< REAL > *waveout, const int N, const int LDA=1, const int m=1) |
|
void | recurs_real (std::function< void(REAL *in, REAL *out, const int)> funA, REAL *arraynp1, REAL *arrayn, REAL *arrayn_1, const int N, const int LDA=1, const int m=1) |
|
void | recurs_complex (std::function< void(std::complex< REAL > *in, std::complex< REAL > *out, const int)> funA, std::complex< REAL > *arraynp1, std::complex< REAL > *arrayn, std::complex< REAL > *arrayn_1, const int N, const int LDA=1, const int m=1) |
|
REAL | recurs (const REAL x, const REAL Tn, const REAL Tn_1) |
|
bool | checkconverge (std::function< void(std::complex< REAL > *in, std::complex< REAL > *out, const int)> funA, std::complex< REAL > *wavein, const int N, const int LDA, REAL &tmax, REAL &tmin, REAL stept) |
|
REAL | ddot_real (const std::complex< REAL > *psi_L, const std::complex< REAL > *psi_R, const int N, const int LDA=1, const int m=1) |
|
template<typename REAL, typename Device = base_device::DEVICE_CPU>
class ModuleBase::Chebyshev< REAL, Device >
A class to treat the Chebyshev expansion.
- Author
- qianrui on 2022-05-18
Math: I. Chebyshev polynomial: T_0(x) = 1; T_1(x) = x; T_2(x) = 2x^2 -1; T_3(x) = 4x^3-3x; T_{n+2}(x) = 2xT_{n+1}(x) - T_n(x)
II. Any analytical function f(x) can be expanded by Chebyshev polynomial: f(x) = \sum_{n=0}^{norder-1} C_n[f]*T_n(x) (|x| < 1), where C_n[f] = \frac{2-\delta_{n0}}{\pi} \int_0^\pi f(cos(\theta))cos(n\theta) d\theta Here C_n can be calculate with FFT.
III. Any functions of linear Operator or matrix f(\hat{A}) or f(A) can also be expanded as well: f(A) = \sum_{n=0}^{norder-1} C_n[f]*T_n(A) (|all eigenvalues of A| < 1). f(A)v = \sum_{n=0}^{norder-1} C_n[f]*T_n(A)v, where v is column vector = \sum_{n=0}^{norder-1} C_n[f]*v_n, where v_n = T_n(A)v, v_0 = v v_{n+2} = 2Av_{n+1} - v_n
IV. v^+f(A)v = \sum_{n=0}^{norder-1} C_n[f]*v^+v_n = \sum_{n=0}^{norder-1} C_n[f] * w_n, where w_n = v^+ * v_n = v^+ * T_n(A) * v
USAGE: Chebyshev che(10); // constructe a chebyshev expansion of 10 orders (n=0,1,...,9)
che.calcoef_real(cos) // calculate C_n[f], where f is a.cos for(int i=0;i<10;++i) cout<<che.coef_real[i]<<endl; //Then we print C_n[f]
che.calcoef_complex(expi) // calculate C_n[g], where g is b.expi for(int i=0;i<10;++i) cout<<che.coef_complex[i]<<endl; //Then we print C_n[g]
che.calcoef_pair(cos, sin) // calculate C_n[g], where g is (c.cos, c.sin) for(int i=0;i<10;++i) cout<<che.coef_complex[i]<<endl; //Then we print C_n[g]
che.calcoef_real(Occupy::fd) che.calfinalvec_real(hpsi, psi_in, psi_out, npw); //calculate f(H)|psi>, where f is occ.fd and H is hamilt.hpsi
che.calcoef_complex(expi) che.calfinalvec_complex(hpsi, psi_in, psi_out, npw, npwx, nbands); //calculate exp(iH)|psi_i>
- che.tracepolyA(hpsi, psi_in, npw, npwx, nbands) //calculate \sum_i^{nbands} <psi_i|T_n(H)|psi_i>
- che.calcoef_complex(expi); //calculate C_n[exp(ix)] che.getpolyval(PI/4, T, norder); //get T_n(pi/4) std::complex<double> sum(0,0); for(int i = 0; i < norder ; ++i) { sum += che.coef_complex[i]*T[i]; //sum = exp(i*pi/4) = \sum_n C_n[exp(ix)]*T_n(pi/4) }
- che.recurs_complex(hpsi, vp1, v, vm1, npw) //calculate vp1: |vp1> = 2 H|v> - |vm1>;