ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | Private Types | Private Attributes | List of all members
ModuleBase::Chebyshev< REAL, Device > Class Template Reference

A class to treat the Chebyshev expansion. More...

#include <math_chebyshev.h>

Collaboration diagram for ModuleBase::Chebyshev< REAL, Device >:

Public Member Functions

 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)
 

Public Attributes

int norder
 
int norder2
 
REAL * coef_real = nullptr
 
std::complex< REAL > * coef_complex = nullptr
 
REAL * coefr_cpu = nullptr
 
std::complex< REAL > * coefc_cpu = nullptr
 
FFTW< REAL > fftw
 
REAL * polytrace
 
bool getcoef_real
 
bool getcoef_complex
 

Private Types

using ct_Device = typename container::PsiToContainer< Device >::type
 
using resmem_complex_op = base_device::memory::resize_memory_op< std::complex< REAL >, Device >
 
using resmem_var_op = base_device::memory::resize_memory_op< REAL, Device >
 
using delmem_complex_op = base_device::memory::delete_memory_op< std::complex< REAL >, Device >
 
using delmem_var_op = base_device::memory::delete_memory_op< REAL, Device >
 
using syncmem_var_h2d_op = base_device::memory::synchronize_memory_op< REAL, Device, base_device::DEVICE_CPU >
 
using syncmem_var_d2h_op = base_device::memory::synchronize_memory_op< REAL, base_device::DEVICE_CPU, Device >
 
using syncmem_complex_h2d_op = base_device::memory::synchronize_memory_op< std::complex< REAL >, Device, base_device::DEVICE_CPU >
 
using syncmem_complex_d2h_op = base_device::memory::synchronize_memory_op< std::complex< REAL >, base_device::DEVICE_CPU, Device >
 
using memcpy_var_op = base_device::memory::synchronize_memory_op< REAL, Device, Device >
 
using memcpy_complex_op = base_device::memory::synchronize_memory_op< std::complex< REAL >, Device, Device >
 
using setmem_complex_op = base_device::memory::set_memory_op< std::complex< REAL >, Device >
 

Private Attributes

Device * ctx = {}
 
base_device::DEVICE_CPU * cpu_ctx = {}
 

Detailed Description

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)

  1. 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]

  2. 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>

  3. che.tracepolyA(hpsi, psi_in, npw, npwx, nbands) //calculate \sum_i^{nbands} <psi_i|T_n(H)|psi_i>
  4. 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) }
  5. che.recurs_complex(hpsi, vp1, v, vm1, npw) //calculate vp1: |vp1> = 2 H|v> - |vm1>;

Member Typedef Documentation

◆ ct_Device

template<typename REAL , typename Device = base_device::DEVICE_CPU>
using ModuleBase::Chebyshev< REAL, Device >::ct_Device = typename container::PsiToContainer<Device>::type
private

◆ delmem_complex_op

template<typename REAL , typename Device = base_device::DEVICE_CPU>
using ModuleBase::Chebyshev< REAL, Device >::delmem_complex_op = base_device::memory::delete_memory_op<std::complex<REAL>, Device>
private

◆ delmem_var_op

template<typename REAL , typename Device = base_device::DEVICE_CPU>
using ModuleBase::Chebyshev< REAL, Device >::delmem_var_op = base_device::memory::delete_memory_op<REAL, Device>
private

◆ memcpy_complex_op

template<typename REAL , typename Device = base_device::DEVICE_CPU>
using ModuleBase::Chebyshev< REAL, Device >::memcpy_complex_op = base_device::memory::synchronize_memory_op<std::complex<REAL>, Device, Device>
private

◆ memcpy_var_op

template<typename REAL , typename Device = base_device::DEVICE_CPU>
using ModuleBase::Chebyshev< REAL, Device >::memcpy_var_op = base_device::memory::synchronize_memory_op<REAL, Device, Device>
private

◆ resmem_complex_op

template<typename REAL , typename Device = base_device::DEVICE_CPU>
using ModuleBase::Chebyshev< REAL, Device >::resmem_complex_op = base_device::memory::resize_memory_op<std::complex<REAL>, Device>
private

◆ resmem_var_op

template<typename REAL , typename Device = base_device::DEVICE_CPU>
using ModuleBase::Chebyshev< REAL, Device >::resmem_var_op = base_device::memory::resize_memory_op<REAL, Device>
private

◆ setmem_complex_op

template<typename REAL , typename Device = base_device::DEVICE_CPU>
using ModuleBase::Chebyshev< REAL, Device >::setmem_complex_op = base_device::memory::set_memory_op<std::complex<REAL>, Device>
private

◆ syncmem_complex_d2h_op

template<typename REAL , typename Device = base_device::DEVICE_CPU>
using ModuleBase::Chebyshev< REAL, Device >::syncmem_complex_d2h_op = base_device::memory::synchronize_memory_op<std::complex<REAL>, base_device::DEVICE_CPU, Device>
private

◆ syncmem_complex_h2d_op

template<typename REAL , typename Device = base_device::DEVICE_CPU>
using ModuleBase::Chebyshev< REAL, Device >::syncmem_complex_h2d_op = base_device::memory::synchronize_memory_op<std::complex<REAL>, Device, base_device::DEVICE_CPU>
private

◆ syncmem_var_d2h_op

template<typename REAL , typename Device = base_device::DEVICE_CPU>
using ModuleBase::Chebyshev< REAL, Device >::syncmem_var_d2h_op = base_device::memory::synchronize_memory_op<REAL, base_device::DEVICE_CPU, Device>
private

◆ syncmem_var_h2d_op

template<typename REAL , typename Device = base_device::DEVICE_CPU>
using ModuleBase::Chebyshev< REAL, Device >::syncmem_var_h2d_op = base_device::memory::synchronize_memory_op<REAL, Device, base_device::DEVICE_CPU>
private

Constructor & Destructor Documentation

◆ Chebyshev()

template<typename REAL , typename Device >
ModuleBase::Chebyshev< REAL, Device >::Chebyshev ( const int  norder)
Here is the call graph for this function:

◆ ~Chebyshev()

template<typename REAL , typename Device >
ModuleBase::Chebyshev< REAL, Device >::~Chebyshev ( )

Member Function Documentation

◆ calcoef_complex()

template<typename REAL , typename Device >
void ModuleBase::Chebyshev< REAL, Device >::calcoef_complex ( std::function< std::complex< REAL >(std::complex< REAL >)>  fun)

◆ calcoef_pair()

template<typename REAL , typename Device >
void ModuleBase::Chebyshev< REAL, Device >::calcoef_pair ( std::function< REAL(REAL)>  fun1,
std::function< REAL(REAL)>  fun2 
)
Here is the caller graph for this function:

◆ calcoef_real()

template<typename REAL , typename Device >
void ModuleBase::Chebyshev< REAL, Device >::calcoef_real ( std::function< REAL(REAL)>  fun)
Here is the caller graph for this function:

◆ calfinalvec_complex()

template<typename REAL , typename Device >
void ModuleBase::Chebyshev< REAL, Device >::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 
)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ calfinalvec_real() [1/2]

template<typename REAL , typename Device = base_device::DEVICE_CPU>
void ModuleBase::Chebyshev< REAL, Device >::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 
)
Here is the caller graph for this function:

◆ calfinalvec_real() [2/2]

template<typename REAL , typename Device >
void ModuleBase::Chebyshev< REAL, Device >::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 
)
Here is the call graph for this function:

◆ calpolyvec_complex()

template<typename REAL , typename Device = base_device::DEVICE_CPU>
void ModuleBase::Chebyshev< REAL, Device >::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 
)
Here is the caller graph for this function:

◆ calpolyvec_real()

template<typename REAL , typename Device = base_device::DEVICE_CPU>
void ModuleBase::Chebyshev< REAL, Device >::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 
)

◆ checkconverge()

template<typename REAL , typename Device >
bool ModuleBase::Chebyshev< REAL, Device >::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 
)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ddot_real()

template<typename REAL , typename Device >
REAL ModuleBase::Chebyshev< REAL, Device >::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 
)

◆ getpolyval()

template<typename REAL , typename Device >
void ModuleBase::Chebyshev< REAL, Device >::getpolyval ( REAL  x,
REAL *  polyval,
const int  N 
)

◆ recurs()

template<typename REAL , typename Device >
REAL ModuleBase::Chebyshev< REAL, Device >::recurs ( const REAL  x,
const REAL  Tn,
const REAL  Tn_1 
)
inline

◆ recurs_complex()

template<typename REAL , typename Device >
void ModuleBase::Chebyshev< REAL, Device >::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 
)

◆ recurs_real()

template<typename REAL , typename Device = base_device::DEVICE_CPU>
void ModuleBase::Chebyshev< REAL, Device >::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 
)

◆ tracepolyA()

template<typename REAL , typename Device >
void ModuleBase::Chebyshev< REAL, Device >::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 
)
Here is the caller graph for this function:

Member Data Documentation

◆ coef_complex

template<typename REAL , typename Device = base_device::DEVICE_CPU>
std::complex<REAL>* ModuleBase::Chebyshev< REAL, Device >::coef_complex = nullptr

◆ coef_real

template<typename REAL , typename Device = base_device::DEVICE_CPU>
REAL* ModuleBase::Chebyshev< REAL, Device >::coef_real = nullptr

◆ coefc_cpu

template<typename REAL , typename Device = base_device::DEVICE_CPU>
std::complex<REAL>* ModuleBase::Chebyshev< REAL, Device >::coefc_cpu = nullptr

◆ coefr_cpu

template<typename REAL , typename Device = base_device::DEVICE_CPU>
REAL* ModuleBase::Chebyshev< REAL, Device >::coefr_cpu = nullptr

◆ cpu_ctx

template<typename REAL , typename Device = base_device::DEVICE_CPU>
base_device::DEVICE_CPU* ModuleBase::Chebyshev< REAL, Device >::cpu_ctx = {}
private

◆ ctx

template<typename REAL , typename Device = base_device::DEVICE_CPU>
Device* ModuleBase::Chebyshev< REAL, Device >::ctx = {}
private

◆ fftw

template<typename REAL , typename Device = base_device::DEVICE_CPU>
FFTW<REAL> ModuleBase::Chebyshev< REAL, Device >::fftw

◆ getcoef_complex

template<typename REAL , typename Device = base_device::DEVICE_CPU>
bool ModuleBase::Chebyshev< REAL, Device >::getcoef_complex

◆ getcoef_real

template<typename REAL , typename Device = base_device::DEVICE_CPU>
bool ModuleBase::Chebyshev< REAL, Device >::getcoef_real

◆ norder

template<typename REAL , typename Device = base_device::DEVICE_CPU>
int ModuleBase::Chebyshev< REAL, Device >::norder

◆ norder2

template<typename REAL , typename Device = base_device::DEVICE_CPU>
int ModuleBase::Chebyshev< REAL, Device >::norder2

◆ polytrace

template<typename REAL , typename Device = base_device::DEVICE_CPU>
REAL* ModuleBase::Chebyshev< REAL, Device >::polytrace

The documentation for this class was generated from the following files: