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

#include <stress_func.h>

Inheritance diagram for Stress_Func< FPTYPE, Device >:
Collaboration diagram for Stress_Func< FPTYPE, Device >:

Public Member Functions

 Stress_Func ()
 
 ~Stress_Func ()
 
void stress_kin (ModuleBase::matrix &sigma, const ModuleBase::matrix &wg, ModuleSymmetry::Symmetry *p_symm, K_Vectors *p_kv, ModulePW::PW_Basis_K *wfc_basis, const UnitCell &ucell_in, const psi::Psi< std::complex< FPTYPE >, Device > *psi_in=nullptr)
 
void stress_har (const UnitCell &ucell, ModuleBase::matrix &sigma, ModulePW::PW_Basis *rho_basis, const bool is_pw, const Charge *const chr)
 
void stress_ewa (const UnitCell &ucell, ModuleBase::matrix &sigma, ModulePW::PW_Basis *rho_basis, const bool is_pw)
 
void stress_loc (const UnitCell &ucell, ModuleBase::matrix &sigma, ModulePW::PW_Basis *rho_basis, const ModuleBase::matrix &vloc, const Structure_Factor *p_sf, const bool is_pw, const Charge *const chr)
 
void dvloc_of_g (const int &msh, const FPTYPE *rab, const FPTYPE *r, const FPTYPE *vloc_at, const FPTYPE &zp, FPTYPE *dvloc, ModulePW::PW_Basis *rho_basis, const UnitCell &ucell_in)
 
void dvloc_coulomb (const UnitCell &ucell, const FPTYPE &zp, FPTYPE *dvloc, ModulePW::PW_Basis *rho_basis)
 compute the derivative of the coulomb potential in reciprocal space D V(g^2) / D g^2 = 4pi e^2/omegai /G^4
 
void stress_cc (ModuleBase::matrix &sigma, ModulePW::PW_Basis *rho_basis, UnitCell &ucell, const Structure_Factor *p_sf, const bool is_pw, const bool *numeric, const Charge *const chr)
 
void deriv_drhoc (const bool &numeric, const double &omega, const double &tpiba2, const int mesh, const FPTYPE *r, const FPTYPE *rab, const FPTYPE *rhoc, FPTYPE *drhocg, ModulePW::PW_Basis *rho_basis, int type)
 
void stress_gga (const UnitCell &ucell, ModuleBase::matrix &sigma, ModulePW::PW_Basis *rho_basis, const Charge *const chr)
 
void stress_mgga (const UnitCell &ucell, ModuleBase::matrix &sigma, const ModuleBase::matrix &wg, const ModuleBase::matrix &v_ofk, const Charge *const chr, K_Vectors *p_kv, ModulePW::PW_Basis_K *wfc_basis, const psi::Psi< std::complex< FPTYPE >, Device > *psi_in)
 
void stress_nl (ModuleBase::matrix &sigma, const ModuleBase::matrix &wg, const ModuleBase::matrix &ekb, Structure_Factor *p_sf, K_Vectors *p_kv, ModuleSymmetry::Symmetry *p_symm, ModulePW::PW_Basis_K *wfc_basis, const psi::Psi< std::complex< FPTYPE >, Device > *psi_in, const pseudopot_cell_vnl &nlpp_in, const UnitCell &ucell_in)
 This routine computes the atomic force of non-local pseudopotential Stress^{NL}_{ij} = -1/\Omega \sum_{n,k}f_{nk}\sum_I \sum_{lm,l'm'}D_{l,l'}^{I} [ \sum_G \langle c_{nk}(\mathbf{G+K})|\beta_{lm}^I(\mathbf{G+K})\rangle * \sum_{G'}\langle \partial \beta_{lm}^I(\mathbf{G+K})/\partial \varepsilon_{ij} |c_{nk}(\mathbf{G+K})\rangle ] there would be three parts in the above equation: (1) sum over becp and dbecp with D_{l,l'}^{I} --— first line in the above equation (2) calculate becp = <psi | beta> --— second line in the above equation (3) calculate dbecp = <psi | dbeta> --— third line in the above equation.
 
void stress_onsite (ModuleBase::matrix &sigma, const ModuleBase::matrix &wg, const ModulePW::PW_Basis_K *wfc_basis, const UnitCell &ucell_in, const psi::Psi< std::complex< FPTYPE >, Device > *psi_in, ModuleSymmetry::Symmetry *p_symm)
 This routine computes the stress contribution from the DFT+U and DeltaSpin calculations Stress^{NL}_{ij} = -1/\Omega \sum_{n,k}f_{nk}\sum_I \sum_{lm,l'm'}(V^U_{lmm'\sigma\sigma'} + f(\lambda,\sigma\sigma')) [ \sum_G \langle c_{nk}(\mathbf{G+K})|\alpha_{lm}^I(\mathbf{G+K})\rangle * \sum_{G'}\langle \partial \alpha_{lm}^I(\mathbf{G+K})/\partial \varepsilon_{ij} |c_{nk}(\mathbf{G+K})\rangle ] there would be three parts in the above equation: (1) sum over becp and dbecp with f(U+\lambda, \sigma\sigma', lmm')^{I} --— first line in the above equation (2) calculate becp = <psi | alpha> --— second line in the above equation (3) calculate dbecp = <psi | dalpha> --— third line in the above equation.
 
void get_dvnl1 (ModuleBase::ComplexMatrix &vkb, const int ik, const int ipol, Structure_Factor *p_sf, ModulePW::PW_Basis_K *wfc_basis)
 
void get_dvnl2 (ModuleBase::ComplexMatrix &vkb, const int ik, Structure_Factor *p_sf, ModulePW::PW_Basis_K *wfc_basis)
 
FPTYPE Polynomial_Interpolation_nl (const ModuleBase::realArray &table, const int &dim1, const int &dim2, const int &dim3, const FPTYPE &table_interval, const FPTYPE &x)
 
void dqvan2 (const pseudopot_cell_vnl &ppcell_in, const int ih, const int jh, const int itype, const int ipol, const int ng, const ModuleBase::Vector3< FPTYPE > *g, const FPTYPE *qnorm, const FPTYPE &tpiba, const ModuleBase::matrix &ylmk0, const ModuleBase::matrix &dylmk0, std::complex< FPTYPE > *dqg)
 Compute the derivatives of the radial Fourier transform of the Q functions.
 

Protected Attributes

Device * ctx = {}
 
base_device::DEVICE_CPU * cpu_ctx = {}
 
base_device::AbacusDevice_t device = {}
 
pseudopot_cell_vnlnlpp = nullptr
 
const UnitCellucell = nullptr
 

Private Types

using gemm_op = ModuleBase::gemm_op< std::complex< FPTYPE >, Device >
 
using cal_stress_nl_op = hamilt::cal_stress_nl_op< FPTYPE, Device >
 
using cal_dbecp_noevc_nl_op = hamilt::cal_dbecp_noevc_nl_op< FPTYPE, Device >
 
using resmem_complex_op = base_device::memory::resize_memory_op< std::complex< FPTYPE >, Device >
 
using resmem_complex_h_op = base_device::memory::resize_memory_op< std::complex< FPTYPE >, base_device::DEVICE_CPU >
 
using setmem_complex_op = base_device::memory::set_memory_op< std::complex< FPTYPE >, Device >
 
using delmem_complex_op = base_device::memory::delete_memory_op< std::complex< FPTYPE >, Device >
 
using delmem_complex_h_op = base_device::memory::delete_memory_op< std::complex< FPTYPE >, base_device::DEVICE_CPU >
 
using syncmem_complex_h2d_op = base_device::memory::synchronize_memory_op< std::complex< FPTYPE >, Device, base_device::DEVICE_CPU >
 
using syncmem_complex_d2h_op = base_device::memory::synchronize_memory_op< std::complex< FPTYPE >, base_device::DEVICE_CPU, Device >
 
using resmem_var_op = base_device::memory::resize_memory_op< FPTYPE, Device >
 
using resmem_var_h_op = base_device::memory::resize_memory_op< FPTYPE, base_device::DEVICE_CPU >
 
using setmem_var_op = base_device::memory::set_memory_op< FPTYPE, Device >
 
using delmem_var_op = base_device::memory::delete_memory_op< FPTYPE, Device >
 
using delmem_var_h_op = base_device::memory::delete_memory_op< FPTYPE, base_device::DEVICE_CPU >
 
using syncmem_var_h2d_op = base_device::memory::synchronize_memory_op< FPTYPE, Device, base_device::DEVICE_CPU >
 
using syncmem_var_d2h_op = base_device::memory::synchronize_memory_op< FPTYPE, base_device::DEVICE_CPU, Device >
 
using resmem_int_op = base_device::memory::resize_memory_op< int, Device >
 
using delmem_int_op = base_device::memory::delete_memory_op< int, Device >
 
using syncmem_int_h2d_op = base_device::memory::synchronize_memory_op< int, Device, base_device::DEVICE_CPU >
 
using cal_vq_op = hamilt::cal_vq_op< FPTYPE, Device >
 
using cal_vq_deri_op = hamilt::cal_vq_deri_op< FPTYPE, Device >
 
using cal_vkb_op = hamilt::cal_vkb_op< FPTYPE, Device >
 
using cal_vkb_deri_op = hamilt::cal_vkb_deri_op< FPTYPE, Device >
 

Member Typedef Documentation

◆ cal_dbecp_noevc_nl_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::cal_dbecp_noevc_nl_op = hamilt::cal_dbecp_noevc_nl_op<FPTYPE, Device>
private

◆ cal_stress_nl_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::cal_stress_nl_op = hamilt::cal_stress_nl_op<FPTYPE, Device>
private

◆ cal_vkb_deri_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::cal_vkb_deri_op = hamilt::cal_vkb_deri_op<FPTYPE, Device>
private

◆ cal_vkb_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::cal_vkb_op = hamilt::cal_vkb_op<FPTYPE, Device>
private

◆ cal_vq_deri_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::cal_vq_deri_op = hamilt::cal_vq_deri_op<FPTYPE, Device>
private

◆ cal_vq_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::cal_vq_op = hamilt::cal_vq_op<FPTYPE, Device>
private

◆ delmem_complex_h_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::delmem_complex_h_op = base_device::memory::delete_memory_op<std::complex<FPTYPE>, base_device::DEVICE_CPU>
private

◆ delmem_complex_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::delmem_complex_op = base_device::memory::delete_memory_op<std::complex<FPTYPE>, Device>
private

◆ delmem_int_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::delmem_int_op = base_device::memory::delete_memory_op<int, Device>
private

◆ delmem_var_h_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::delmem_var_h_op = base_device::memory::delete_memory_op<FPTYPE, base_device::DEVICE_CPU>
private

◆ delmem_var_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::delmem_var_op = base_device::memory::delete_memory_op<FPTYPE, Device>
private

◆ gemm_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::gemm_op = ModuleBase::gemm_op<std::complex<FPTYPE>, Device>
private

◆ resmem_complex_h_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::resmem_complex_h_op = base_device::memory::resize_memory_op<std::complex<FPTYPE>, base_device::DEVICE_CPU>
private

◆ resmem_complex_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::resmem_complex_op = base_device::memory::resize_memory_op<std::complex<FPTYPE>, Device>
private

◆ resmem_int_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::resmem_int_op = base_device::memory::resize_memory_op<int, Device>
private

◆ resmem_var_h_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::resmem_var_h_op = base_device::memory::resize_memory_op<FPTYPE, base_device::DEVICE_CPU>
private

◆ resmem_var_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::resmem_var_op = base_device::memory::resize_memory_op<FPTYPE, Device>
private

◆ setmem_complex_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::setmem_complex_op = base_device::memory::set_memory_op<std::complex<FPTYPE>, Device>
private

◆ setmem_var_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::setmem_var_op = base_device::memory::set_memory_op<FPTYPE, Device>
private

◆ syncmem_complex_d2h_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::syncmem_complex_d2h_op = base_device::memory::synchronize_memory_op<std::complex<FPTYPE>, base_device::DEVICE_CPU, Device>
private

◆ syncmem_complex_h2d_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::syncmem_complex_h2d_op = base_device::memory::synchronize_memory_op<std::complex<FPTYPE>, Device, base_device::DEVICE_CPU>
private

◆ syncmem_int_h2d_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::syncmem_int_h2d_op = base_device::memory::synchronize_memory_op<int, Device, base_device::DEVICE_CPU>
private

◆ syncmem_var_d2h_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::syncmem_var_d2h_op = base_device::memory::synchronize_memory_op<FPTYPE, base_device::DEVICE_CPU, Device>
private

◆ syncmem_var_h2d_op

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
using Stress_Func< FPTYPE, Device >::syncmem_var_h2d_op = base_device::memory::synchronize_memory_op<FPTYPE, Device, base_device::DEVICE_CPU>
private

Constructor & Destructor Documentation

◆ Stress_Func()

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
Stress_Func< FPTYPE, Device >::Stress_Func ( )
inline

◆ ~Stress_Func()

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
Stress_Func< FPTYPE, Device >::~Stress_Func ( )
inline

Member Function Documentation

◆ deriv_drhoc()

template<typename FPTYPE , typename Device >
void Stress_Func< FPTYPE, Device >::deriv_drhoc ( const bool &  numeric,
const double &  omega,
const double &  tpiba2,
const int  mesh,
const FPTYPE *  r,
const FPTYPE *  rab,
const FPTYPE *  rhoc,
FPTYPE *  drhocg,
ModulePW::PW_Basis rho_basis,
int  type 
)
Here is the call graph for this function:

◆ dqvan2()

template<typename FPTYPE , typename Device >
void Stress_Func< FPTYPE, Device >::dqvan2 ( const pseudopot_cell_vnl ppcell_in,
const int  ih,
const int  jh,
const int  itype,
const int  ipol,
const int  ng,
const ModuleBase::Vector3< FPTYPE > *  g,
const FPTYPE *  qnorm,
const FPTYPE &  tpiba,
const ModuleBase::matrix ylmk0,
const ModuleBase::matrix dylmk0,
std::complex< FPTYPE > *  dqg 
)

Compute the derivatives of the radial Fourier transform of the Q functions.

This routine computes the derivatives of the Fourier transform of the Q function needed in stress assuming that the radial fourier transform is already computed and stored in table qrad. The formula implemented here is:

dq(g,i,j) = sum_lm (-i)^l ap(lm,i,j) * ( yr_lm(g^) dqrad(g,l,i,j) + dyr_lm(g^) qrad(g,l,i,j))

Parameters
ih[in] the first index of Q
jh[in] the second index of Q
itype[in] the atomic type
ipol[in] the polarization of the derivative
ng[in] the number of G vectors
g[in] the G vectors
qnorm[in] the norm of q+g vectors
tpiba[in] 2pi/a factor, multiplies G vectors
ylmk0[in] the real spherical harmonics
dylmk0[in] derivetives of spherical harmonics
dqg[out] the Fourier transform of interest
Here is the call graph for this function:

◆ dvloc_coulomb()

template<typename FPTYPE , typename Device >
void Stress_Func< FPTYPE, Device >::dvloc_coulomb ( const UnitCell ucell,
const FPTYPE &  zp,
FPTYPE *  dvloc,
ModulePW::PW_Basis rho_basis 
)

compute the derivative of the coulomb potential in reciprocal space D V(g^2) / D g^2 = 4pi e^2/omegai /G^4

◆ dvloc_of_g()

template<typename FPTYPE , typename Device >
void Stress_Func< FPTYPE, Device >::dvloc_of_g ( const int &  msh,
const FPTYPE *  rab,
const FPTYPE *  r,
const FPTYPE *  vloc_at,
const FPTYPE &  zp,
FPTYPE *  dvloc,
ModulePW::PW_Basis rho_basis,
const UnitCell ucell_in 
)

◆ get_dvnl1()

template<typename FPTYPE , typename Device >
void Stress_Func< FPTYPE, Device >::get_dvnl1 ( ModuleBase::ComplexMatrix vkb,
const int  ik,
const int  ipol,
Structure_Factor p_sf,
ModulePW::PW_Basis_K wfc_basis 
)
Here is the call graph for this function:

◆ get_dvnl2()

template<typename FPTYPE , typename Device >
void Stress_Func< FPTYPE, Device >::get_dvnl2 ( ModuleBase::ComplexMatrix vkb,
const int  ik,
Structure_Factor p_sf,
ModulePW::PW_Basis_K wfc_basis 
)
Here is the call graph for this function:

◆ Polynomial_Interpolation_nl()

template<typename FPTYPE , typename Device >
FPTYPE Stress_Func< FPTYPE, Device >::Polynomial_Interpolation_nl ( const ModuleBase::realArray table,
const int &  dim1,
const int &  dim2,
const int &  dim3,
const FPTYPE &  table_interval,
const FPTYPE &  x 
)

◆ stress_cc()

template<typename FPTYPE , typename Device >
void Stress_Func< FPTYPE, Device >::stress_cc ( ModuleBase::matrix sigma,
ModulePW::PW_Basis rho_basis,
UnitCell ucell,
const Structure_Factor p_sf,
const bool  is_pw,
const bool *  numeric,
const Charge *const  chr 
)
Here is the call graph for this function:

◆ stress_ewa()

template<typename FPTYPE , typename Device >
void Stress_Func< FPTYPE, Device >::stress_ewa ( const UnitCell ucell,
ModuleBase::matrix sigma,
ModulePW::PW_Basis rho_basis,
const bool  is_pw 
)
Here is the call graph for this function:

◆ stress_gga()

template<typename FPTYPE , typename Device >
void Stress_Func< FPTYPE, Device >::stress_gga ( const UnitCell ucell,
ModuleBase::matrix sigma,
ModulePW::PW_Basis rho_basis,
const Charge *const  chr 
)
Here is the call graph for this function:

◆ stress_har()

template<typename FPTYPE , typename Device >
void Stress_Func< FPTYPE, Device >::stress_har ( const UnitCell ucell,
ModuleBase::matrix sigma,
ModulePW::PW_Basis rho_basis,
const bool  is_pw,
const Charge *const  chr 
)
Here is the call graph for this function:

◆ stress_kin()

template<typename FPTYPE , typename Device >
void Stress_Func< FPTYPE, Device >::stress_kin ( ModuleBase::matrix sigma,
const ModuleBase::matrix wg,
ModuleSymmetry::Symmetry p_symm,
K_Vectors p_kv,
ModulePW::PW_Basis_K wfc_basis,
const UnitCell ucell_in,
const psi::Psi< std::complex< FPTYPE >, Device > *  psi_in = nullptr 
)
Here is the call graph for this function:

◆ stress_loc()

template<typename FPTYPE , typename Device >
void Stress_Func< FPTYPE, Device >::stress_loc ( const UnitCell ucell,
ModuleBase::matrix sigma,
ModulePW::PW_Basis rho_basis,
const ModuleBase::matrix vloc,
const Structure_Factor p_sf,
const bool  is_pw,
const Charge *const  chr 
)
Here is the call graph for this function:

◆ stress_mgga()

template<typename FPTYPE , typename Device >
void Stress_Func< FPTYPE, Device >::stress_mgga ( const UnitCell ucell,
ModuleBase::matrix sigma,
const ModuleBase::matrix wg,
const ModuleBase::matrix v_ofk,
const Charge *const  chr,
K_Vectors p_kv,
ModulePW::PW_Basis_K wfc_basis,
const psi::Psi< std::complex< FPTYPE >, Device > *  psi_in 
)
Here is the call graph for this function:

◆ stress_nl()

template<typename FPTYPE , typename Device >
void Stress_Func< FPTYPE, Device >::stress_nl ( ModuleBase::matrix sigma,
const ModuleBase::matrix wg,
const ModuleBase::matrix ekb,
Structure_Factor p_sf,
K_Vectors p_kv,
ModuleSymmetry::Symmetry p_symm,
ModulePW::PW_Basis_K wfc_basis,
const psi::Psi< std::complex< FPTYPE >, Device > *  psi_in,
const pseudopot_cell_vnl nlpp_in,
const UnitCell ucell_in 
)

This routine computes the atomic force of non-local pseudopotential Stress^{NL}_{ij} = -1/\Omega \sum_{n,k}f_{nk}\sum_I \sum_{lm,l'm'}D_{l,l'}^{I} [ \sum_G \langle c_{nk}(\mathbf{G+K})|\beta_{lm}^I(\mathbf{G+K})\rangle * \sum_{G'}\langle \partial \beta_{lm}^I(\mathbf{G+K})/\partial \varepsilon_{ij} |c_{nk}(\mathbf{G+K})\rangle ] there would be three parts in the above equation: (1) sum over becp and dbecp with D_{l,l'}^{I} --— first line in the above equation (2) calculate becp = <psi | beta> --— second line in the above equation (3) calculate dbecp = <psi | dbeta> --— third line in the above equation.

Here is the call graph for this function:

◆ stress_onsite()

template<typename FPTYPE , typename Device >
void Stress_Func< FPTYPE, Device >::stress_onsite ( ModuleBase::matrix sigma,
const ModuleBase::matrix wg,
const ModulePW::PW_Basis_K wfc_basis,
const UnitCell ucell_in,
const psi::Psi< std::complex< FPTYPE >, Device > *  psi_in,
ModuleSymmetry::Symmetry p_symm 
)

This routine computes the stress contribution from the DFT+U and DeltaSpin calculations Stress^{NL}_{ij} = -1/\Omega \sum_{n,k}f_{nk}\sum_I \sum_{lm,l'm'}(V^U_{lmm'\sigma\sigma'} + f(\lambda,\sigma\sigma')) [ \sum_G \langle c_{nk}(\mathbf{G+K})|\alpha_{lm}^I(\mathbf{G+K})\rangle * \sum_{G'}\langle \partial \alpha_{lm}^I(\mathbf{G+K})/\partial \varepsilon_{ij} |c_{nk}(\mathbf{G+K})\rangle ] there would be three parts in the above equation: (1) sum over becp and dbecp with f(U+\lambda, \sigma\sigma', lmm')^{I} --— first line in the above equation (2) calculate becp = <psi | alpha> --— second line in the above equation (3) calculate dbecp = <psi | dalpha> --— third line in the above equation.

Here is the call graph for this function:

Member Data Documentation

◆ cpu_ctx

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
base_device::DEVICE_CPU* Stress_Func< FPTYPE, Device >::cpu_ctx = {}
protected

◆ ctx

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
Device* Stress_Func< FPTYPE, Device >::ctx = {}
protected

◆ device

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
base_device::AbacusDevice_t Stress_Func< FPTYPE, Device >::device = {}
protected

◆ nlpp

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
pseudopot_cell_vnl* Stress_Func< FPTYPE, Device >::nlpp = nullptr
protected

◆ ucell

template<typename FPTYPE , typename Device = base_device::DEVICE_CPU>
const UnitCell* Stress_Func< FPTYPE, Device >::ucell = nullptr
protected

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