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

A class that implements the block-Davidson algorithm for solving generalized eigenvalue problems. More...

#include <diago_david.h>

Collaboration diagram for hsolver::DiagoDavid< T, Device >:

Public Types

using HPsiFunc = std::function< void(T *, T *, const int, const int)>
 A function type representing the HX function.
 
using SPsiFunc = std::function< void(T *, T *, const int, const int)>
 A function type representing the SX function.
 

Public Member Functions

 DiagoDavid (const Real *precondition_in, const int nband_in, const int dim_in, const int david_ndim_in, const bool use_paw_in, const diag_comm_info &diag_comm_in)
 Constructor for the DiagoDavid class.
 
 ~DiagoDavid ()
 Destructor for the DiagoDavid class.
 
int diag (const HPsiFunc &hpsi_func, const SPsiFunc &spsi_func, const int ld_psi, T *psi_in, Real *eigenvalue_in, const std::vector< double > &ethr_band, const int david_maxiter, const int ntry_max=5, const int notconv_max=0)
 Performs iterative diagonalization using the David algorithm.
 

Private Types

using Real = typename GetTypeReal< T >::type
 
using resmem_complex_op = base_device::memory::resize_memory_op< T, Device >
 
using delmem_complex_op = base_device::memory::delete_memory_op< T, Device >
 
using setmem_complex_op = base_device::memory::set_memory_op< T, Device >
 
using resmem_var_op = base_device::memory::resize_memory_op< Real, Device >
 
using delmem_var_op = base_device::memory::delete_memory_op< Real, Device >
 
using setmem_var_op = base_device::memory::set_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_op = base_device::memory::synchronize_memory_op< T, Device, Device >
 
using castmem_complex_op = base_device::memory::cast_memory_op< std::complex< double >, T, Device, Device >
 
using syncmem_h2d_op = base_device::memory::synchronize_memory_op< T, Device, base_device::DEVICE_CPU >
 
using syncmem_d2h_op = base_device::memory::synchronize_memory_op< T, base_device::DEVICE_CPU, Device >
 

Private Member Functions

int diag_once (const HPsiFunc &hpsi_func, const SPsiFunc &spsi_func, const int dim, const int nband, const int ld_psi, T *psi_in, Real *eigenvalue_in, const std::vector< double > &ethr_band, const int david_maxiter)
 
void cal_grad (const HPsiFunc &hpsi_func, const SPsiFunc &spsi_func, const int &dim, const int &nbase, const int nbase_x, const int &notconv, T *hpsi, T *spsi, const T *vcc, const int *unconv, const Real *eigenvalue)
 
void cal_elem (const int &dim, int &nbase, const int nbase_x, const int &notconv, const T *hpsi, const T *spsi, T *hcc)
 
void refresh (const int &dim, const int &nband, int &nbase, const int nbase_x, const Real *eigenvalue, const T *psi_in, const int ld_psi, T *hpsi, T *spsi, T *hcc, T *vcc)
 
void SchmidtOrth (const int &dim, const int nband, const int m, const T *spsi, T *lagrange_m, const int mm_size, const int mv_size)
 
void planSchmidtOrth (const int nband, std::vector< int > &pre_matrix_mm_m, std::vector< int > &pre_matrix_mv_m)
 Plans the Schmidt orthogonalization for a given number of bands.
 
void diag_zhegvx (const int &nbase, const int &nband, const T *hcc, const int &nbase_x, Real *eigenvalue, T *vcc)
 
bool check_block_conv (const int &ntry, const int &notconv, const int &ntry_max, const int &notconv_max)
 Check the convergence of block eigenvectors in the Davidson iteration.
 

Private Attributes

bool use_paw = false
 
int test_david = 0
 
diag_comm_info diag_comm
 
const int nband
 number of required eigenpairs
 
const int dim
 dimension of the input matrix to be diagonalized
 
const int nbase_x
 maximum dimension of the reduced basis set
 
const int david_ndim = 4
 dimension of the subspace allowed in Davidson
 
int notconv = 0
 number of unconverged eigenvalues
 
const Realprecondition = nullptr
 precondition for diag, diagonal approximation of matrix A(i.e. Hamilt)
 
Reald_precondition = nullptr
 
Realeigenvalue = nullptr
 eigenvalue results
 
Tbasis = nullptr
 
Thpsi = nullptr
 pointer to basis set(dim, nbase_x), leading dimension = dim
 
Tspsi = nullptr
 the product of H and psi in the reduced basis set
 
Thcc = nullptr
 the Product of S and psi in the reduced basis set
 
Tvcc = nullptr
 Hamiltonian on the reduced basis.
 
Tlagrange_matrix = nullptr
 eigenvectors of hc
 
Device * ctx = {}
 device type of psi
 
base_device::DEVICE_CPU * cpu_ctx = {}
 
base_device::AbacusDevice_t device = {}
 
const Tone = nullptr
 
const Tzero = nullptr
 
const Tneg_one = nullptr
 
const T one_ = static_cast<T>(1.0)
 
const T zero_ = static_cast<T>(0.0)
 
const T neg_one_ = static_cast<T>(-1.0)
 

Detailed Description

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
class hsolver::DiagoDavid< T, Device >

A class that implements the block-Davidson algorithm for solving generalized eigenvalue problems.

The DiagoDavid class provides methods for performing iterative diagonalization using the Davidson algorithm. It supports both real and complex data types and can be executed on different devices (CPU or GPU).

Template Parameters
TThe data type of the matrices and arrays (e.g., float, double, std::complex<float>, std::complex<double>).
DeviceThe device type (e.g., base_device::DEVICE_CPU or DEVICE_GPU).

Member Typedef Documentation

◆ castmem_complex_op

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
using hsolver::DiagoDavid< T, Device >::castmem_complex_op = base_device::memory::cast_memory_op<std::complex<double>, T, Device, Device>
private

◆ delmem_complex_op

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
using hsolver::DiagoDavid< T, Device >::delmem_complex_op = base_device::memory::delete_memory_op<T, Device>
private

◆ delmem_var_op

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
using hsolver::DiagoDavid< T, Device >::delmem_var_op = base_device::memory::delete_memory_op<Real, Device>
private

◆ HPsiFunc

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
using hsolver::DiagoDavid< T, Device >::HPsiFunc = std::function<void(T*, T*, const int, const int)>

A function type representing the HX function.

This function type is used to define a matrix-blockvector operator H. For eigenvalue problem HX = λX or generalized eigenvalue problem HX = λSX, this function computes the product of the Hamiltonian matrix H and a blockvector X.

Called as follows: hpsi(X, HX, ld, nvec) where X and HX are (ld, nvec)-shaped blockvectors. Result HX = H * X is stored in HX.

Parameters
[out]XHead address of input blockvector of type T*.
[in]HXHead address of output blockvector of type T*.
[in]ldLeading dimension of blockvector.
[in]nvecNumber of vectors in a block.
Warning
X and HX are the exact address to read input X and store output H*X,
both of size ld * nvec.

◆ Real

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
using hsolver::DiagoDavid< T, Device >::Real = typename GetTypeReal<T>::type
private

◆ resmem_complex_op

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
using hsolver::DiagoDavid< T, Device >::resmem_complex_op = base_device::memory::resize_memory_op<T, Device>
private

◆ resmem_var_op

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
using hsolver::DiagoDavid< T, Device >::resmem_var_op = base_device::memory::resize_memory_op<Real, Device>
private

◆ setmem_complex_op

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
using hsolver::DiagoDavid< T, Device >::setmem_complex_op = base_device::memory::set_memory_op<T, Device>
private

◆ setmem_var_op

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
using hsolver::DiagoDavid< T, Device >::setmem_var_op = base_device::memory::set_memory_op<Real, Device>
private

◆ SPsiFunc

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
using hsolver::DiagoDavid< T, Device >::SPsiFunc = std::function<void(T*, T*, const int, const int)>

A function type representing the SX function.

nrow is leading dimension of spsi, npw is leading dimension of psi, nbands is number of vecs

This function type is used to define a matrix-blockvector operator S. For generalized eigenvalue problem HX = λSX, this function computes the product of the overlap matrix S and a blockvector X.

Parameters
[in]XPointer to the input blockvector.
[out]SXPointer to the output blockvector.
[in]ld_psiLeading dimension of psi and spsi. Dimension of X&SX: ld * nvec.
[in]nvecNumber of vectors.

◆ syncmem_complex_op

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
using hsolver::DiagoDavid< T, Device >::syncmem_complex_op = base_device::memory::synchronize_memory_op<T, Device, Device>
private

◆ syncmem_d2h_op

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
using hsolver::DiagoDavid< T, Device >::syncmem_d2h_op = base_device::memory::synchronize_memory_op<T, base_device::DEVICE_CPU, Device>
private

◆ syncmem_h2d_op

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
using hsolver::DiagoDavid< T, Device >::syncmem_h2d_op = base_device::memory::synchronize_memory_op<T, Device, base_device::DEVICE_CPU>
private

◆ syncmem_var_d2h_op

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
using hsolver::DiagoDavid< T, Device >::syncmem_var_d2h_op = base_device::memory::synchronize_memory_op<Real, base_device::DEVICE_CPU, Device>
private

◆ syncmem_var_h2d_op

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
using hsolver::DiagoDavid< T, Device >::syncmem_var_h2d_op = base_device::memory::synchronize_memory_op<Real, Device, base_device::DEVICE_CPU>
private

Constructor & Destructor Documentation

◆ DiagoDavid()

template<typename T , typename Device >
hsolver::DiagoDavid< T, Device >::DiagoDavid ( const Real precondition_in,
const int  nband_in,
const int  dim_in,
const int  david_ndim_in,
const bool  use_paw_in,
const diag_comm_info diag_comm_in 
)

Constructor for the DiagoDavid class.

Parameters
[in]precondition_inPointer to the preconditioning matrix.
[in]nband_inNumber of eigenpairs required(i.e. bands).
[in]dim_inDimension of the matrix.
[in]david_ndim_inDimension of the reduced basis set of Davidson. david_ndim_in * nband_in is the maximum allowed size of the reduced basis set before restart of Davidson.
[in]use_paw_inFlag indicating whether to use PAW.
[in]diag_comm_inCommunication information for diagonalization.
Template Parameters
TThe data type of the matrices and arrays.
DeviceThe device type (base_device::DEVICE_CPU or DEVICE_GPU).
Note
Auxiliary memory is allocated in the constructor and deallocated in the destructor.

initialize variables k_first = 0 means that nks is more like a dimension of "basis" to be contracted in "HPsi".In LR-TDDFT the formula writes : $$\sum_{ jb\mathbf{k}'}A^I_{ia\mathbf{k}, jb\mathbf{k}' }X ^ I_{ jb\mathbf{k}'}$$ In the code :

  • "H" means A
  • "Psi" means X
  • "band" means the superscript I : the number of excited states to be solved
  • k : k-points, the same meaning as the ground state
  • "basis" : number of occupied ks-orbitals(subscripts i,j) * number of unoccupied ks-orbitals(subscripts a,b), corresponding to "bands" of the ground state

◆ ~DiagoDavid()

template<typename T , typename Device >
hsolver::DiagoDavid< T, Device >::~DiagoDavid ( )

Destructor for the DiagoDavid class.

This destructor releases the dynamically allocated memory used by the class members. It deletes the basis, hpsi, spsi, hcc, vcc, lagrange_matrix, and eigenvalue arrays.

Member Function Documentation

◆ cal_elem()

template<typename T , typename Device >
void DiagoDavid::cal_elem ( const int &  dim,
int &  nbase,
const int  nbase_x,
const int &  notconv,
const T hpsi,
const T spsi,
T hcc 
)
private

Calculates the elements of the diagonalization matrix for the DiagoDavid class.

Parameters
dimThe dimension of the problem.
nbaseThe current dimension of the reduced basis.
nbase_xThe maximum dimension of the reduced basis set.
notconvThe number of newly added basis vectors.
hpsiThe output array for the Hamiltonian H times blockvector psi.
spsiThe output array for the overlap matrix S times blockvector psi.
hccPointer to the array where the calculated Hamiltonian matrix elements will be stored.
Here is the call graph for this function:

◆ cal_grad()

template<typename T , typename Device >
void DiagoDavid::cal_grad ( const HPsiFunc hpsi_func,
const SPsiFunc spsi_func,
const int &  dim,
const int &  nbase,
const int  nbase_x,
const int &  notconv,
T hpsi,
T spsi,
const T vcc,
const int *  unconv,
const Real eigenvalue 
)
private

Calculates the preconditioned gradient of the eigenvectors in Davidson method.

Parameters
hpsi_funcThe function to calculate the matrix-blockvector product H * psi.
spsi_funcThe function to calculate the matrix-blockvector product overlap S * psi.
dimThe dimension of the blockvector.
nbaseThe current dimension of the reduced basis.
nbase_xThe maximum dimension of the reduced basis set.
notconvThe number of unconverged eigenpairs.
hpsiThe output array for the Hamiltonian H times blockvector psi.
spsiThe output array for the overlap matrix S times blockvector psi.
vccThe input array for the eigenvector coefficients.
unconvThe array of indices for the unconverged eigenpairs.
eigenvalueThe array of eigenvalues.
Here is the call graph for this function:

◆ check_block_conv()

template<typename T , typename Device >
bool DiagoDavid::check_block_conv ( const int &  ntry,
const int &  notconv,
const int &  ntry_max,
const int &  notconv_max 
)
inlineprivate

Check the convergence of block eigenvectors in the Davidson iteration.

This function determines whether the block eigenvectors have reached convergence during the iterative diagonalization process. Convergence is judged based on the number of eigenvectors that have not converged and the maximum allowed number of such eigenvectors.

Template Parameters
TThe data type for the eigenvalues and eigenvectors (e.g., float, double).
DeviceThe device type (e.g., base_device::DEVICE_CPU).
Parameters
ntryThe current number of tries for diagonalization.
notconvThe current number of eigenvectors that have not converged.
ntry_maxThe maximum allowed number of tries for diagonalization.
notconv_maxThe maximum allowed number of eigenvectors that can fail to converge.
Returns
true if the eigenvectors are considered converged or the maximum number of tries has been reached, false otherwise.
Note
Exits the diagonalization loop if either the convergence criteria are met or the maximum number of tries is exceeded.

◆ diag()

template<typename T , typename Device >
int hsolver::DiagoDavid< T, Device >::diag ( const HPsiFunc hpsi_func,
const SPsiFunc spsi_func,
const int  ld_psi,
T psi_in,
Real eigenvalue_in,
const std::vector< double > &  ethr_band,
const int  david_maxiter,
const int  ntry_max = 5,
const int  notconv_max = 0 
)

Performs iterative diagonalization using the David algorithm.

Warning
Please see docs of HPsiFunc for more information about the hpsi mat-vec interface.
Template Parameters
TThe type of the elements in the matrix.
DeviceThe device type (CPU or GPU).
Parameters
hpsi_funcThe function object that computes the matrix-blockvector product H * psi.
spsi_funcThe function object that computes the matrix-blockvector product overlap S * psi.
ld_psiThe leading dimension of the psi_in array.
psi_inThe input wavefunction.
eigenvalue_inThe array to store the eigenvalues.
david_diag_thrThe convergence threshold for the diagonalization.
david_maxiterThe maximum number of iterations for the diagonalization.
ntry_maxThe maximum number of attempts for the diagonalization restart.
notconv_maxThe maximum number of bands unconverged allowed.
Returns
The total number of iterations performed during the diagonalization.
Note
ntry_max is an empirical parameter that should be specified in external routine, default 5 notconv_max is determined by the accuracy required for the calculation, default 0

record the times of trying iterative diagonalization

Here is the caller graph for this function:

◆ diag_once()

template<typename T , typename Device >
int DiagoDavid::diag_once ( const HPsiFunc hpsi_func,
const SPsiFunc spsi_func,
const int  dim,
const int  nband,
const int  ld_psi,
T psi_in,
Real eigenvalue_in,
const std::vector< double > &  ethr_band,
const int  david_maxiter 
)
private
Here is the call graph for this function:

◆ diag_zhegvx()

template<typename T , typename Device >
void DiagoDavid::diag_zhegvx ( const int &  nbase,
const int &  nband,
const T hcc,
const int &  nbase_x,
Real eigenvalue,
T vcc 
)
private
Here is the call graph for this function:

◆ planSchmidtOrth()

template<typename T , typename Device >
void DiagoDavid::planSchmidtOrth ( const int  nband,
std::vector< int > &  pre_matrix_mm_m,
std::vector< int > &  pre_matrix_mv_m 
)
private

Plans the Schmidt orthogonalization for a given number of bands.

Template Parameters
TThe type of the elements in the vectors.
DeviceThe device on which the computation will be performed.
Parameters
nbandThe number of bands.
pre_matrix_mm_mThe vector to store the matrix sizes.
pre_matrix_mv_mThe vector to store the number of matrix-vector multiplications.

◆ refresh()

template<typename T , typename Device >
void DiagoDavid::refresh ( const int &  dim,
const int &  nband,
int &  nbase,
const int  nbase_x,
const Real eigenvalue,
const T psi_in,
const int  ld_psi,
T hpsi,
T spsi,
T hcc,
T vcc 
)
private

Refreshes the diagonalization solver by updating the basis and the reduced Hamiltonian.

Parameters
dimThe dimension of the problem.
nbandThe number of bands.
nbaseThe number of basis states.
nbase_xThe maximum dimension of the reduced basis set.
eigenvalue_inPointer to the array of eigenvalues.
psi_inPointer to the array of wavefunctions.
ld_psiThe leading dimension of the wavefunction array.
hpsiPointer to the output array for the updated basis set.
spsiPointer to the output array for the updated basis set (nband-th column).
hccPointer to the output array for the updated reduced Hamiltonian.
vccPointer to the output array for the updated eigenvector matrix.
Here is the call graph for this function:

◆ SchmidtOrth()

template<typename T , typename Device >
void DiagoDavid::SchmidtOrth ( const int &  dim,
const int  nband,
const int  m,
const T spsi,
T lagrange_m,
const int  mm_size,
const int  mv_size 
)
private

SchmidtOrth function performs orthogonalization of the starting eigenfunction to those already calculated. It takes the dimension of the basis, number of bands, index of the current band, starting eigenfunction psi_m, lagrange_m array, mm_size, and mv_size as input parameters.

Parameters
dimThe dimension of the basis.
nbandThe number of bands.
mThe index of the current band.
spsiPointer to the starting eigenfunction psi_m.
lagrange_mPointer to the lagrange_m array.
mm_sizeThe size of the square matrix for future lagranges.
mv_sizeThe size of the lagrange_m array.
Here is the call graph for this function:

Member Data Documentation

◆ basis

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
T* hsolver::DiagoDavid< T, Device >::basis = nullptr
private

◆ cpu_ctx

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
base_device::DEVICE_CPU* hsolver::DiagoDavid< T, Device >::cpu_ctx = {}
private

◆ ctx

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
Device* hsolver::DiagoDavid< T, Device >::ctx = {}
private

device type of psi

◆ d_precondition

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
Real* hsolver::DiagoDavid< T, Device >::d_precondition = nullptr
private

◆ david_ndim

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
const int hsolver::DiagoDavid< T, Device >::david_ndim = 4
private

dimension of the subspace allowed in Davidson

◆ device

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
base_device::AbacusDevice_t hsolver::DiagoDavid< T, Device >::device = {}
private

◆ diag_comm

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
diag_comm_info hsolver::DiagoDavid< T, Device >::diag_comm
private

◆ dim

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
const int hsolver::DiagoDavid< T, Device >::dim
private

dimension of the input matrix to be diagonalized

◆ eigenvalue

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
Real* hsolver::DiagoDavid< T, Device >::eigenvalue = nullptr
private

eigenvalue results

◆ hcc

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
T* hsolver::DiagoDavid< T, Device >::hcc = nullptr
private

the Product of S and psi in the reduced basis set

◆ hpsi

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
T* hsolver::DiagoDavid< T, Device >::hpsi = nullptr
private

pointer to basis set(dim, nbase_x), leading dimension = dim

◆ lagrange_matrix

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
T* hsolver::DiagoDavid< T, Device >::lagrange_matrix = nullptr
private

eigenvectors of hc

◆ nband

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
const int hsolver::DiagoDavid< T, Device >::nband
private

number of required eigenpairs

◆ nbase_x

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
const int hsolver::DiagoDavid< T, Device >::nbase_x
private

maximum dimension of the reduced basis set

◆ neg_one

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
const T * hsolver::DiagoDavid< T, Device >::neg_one = nullptr
private

◆ neg_one_

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
const T hsolver::DiagoDavid< T, Device >::neg_one_ = static_cast<T>(-1.0)
private

◆ notconv

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
int hsolver::DiagoDavid< T, Device >::notconv = 0
private

number of unconverged eigenvalues

◆ one

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
const T* hsolver::DiagoDavid< T, Device >::one = nullptr
private

◆ one_

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
const T hsolver::DiagoDavid< T, Device >::one_ = static_cast<T>(1.0)
private

◆ precondition

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
const Real* hsolver::DiagoDavid< T, Device >::precondition = nullptr
private

precondition for diag, diagonal approximation of matrix A(i.e. Hamilt)

◆ spsi

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
T* hsolver::DiagoDavid< T, Device >::spsi = nullptr
private

the product of H and psi in the reduced basis set

◆ test_david

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
int hsolver::DiagoDavid< T, Device >::test_david = 0
private

◆ use_paw

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
bool hsolver::DiagoDavid< T, Device >::use_paw = false
private

◆ vcc

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
T* hsolver::DiagoDavid< T, Device >::vcc = nullptr
private

Hamiltonian on the reduced basis.

◆ zero

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
const T * hsolver::DiagoDavid< T, Device >::zero = nullptr
private

◆ zero_

template<typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
const T hsolver::DiagoDavid< T, Device >::zero_ = static_cast<T>(0.0)
private

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