ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
diago_bpcg.h
Go to the documentation of this file.
1#ifndef DIAGO_BPCG_H_
2#define DIAGO_BPCG_H_
3
11
12#include <ATen/core/tensor.h>
14#include <source_base/macros.h>
15
16namespace hsolver {
17
24template <typename T = std::complex<double>, typename Device = base_device::DEVICE_CPU>
26{
27 private:
28 // Note GetTypeReal<T>::type will
29 // return T if T is real type(float, double),
30 // otherwise return the real type of T(complex<float>, std::complex<double>)
31 using Real = typename GetTypeReal<T>::type;
32 // Column major psi in this class
33 public:
39 explicit DiagoBPCG(const Real* precondition);
40
44 ~DiagoBPCG();
45
57 void init_iter(const int nband, const int nband_l, const int nbasis, const int ndim);
58
59 using HPsiFunc = std::function<void(T*, T*, const int, const int)>;
60
71 void diag(const HPsiFunc& hpsi_func,
72 T* psi_in,
73 Real* eigenvalue_in,
74 const std::vector<double>& ethr_band);
75
76 private:
78 int n_band = 0;
80 int n_band_l = 0;
82 int n_basis = 0;
84 int n_dim = 0;
86 int nline = 4;
90
91
92 ct::DataType r_type = ct::DataType::DT_INVALID;
93 ct::DataType t_type = ct::DataType::DT_INVALID;
94 ct::DeviceType device_type = ct::DeviceType::UnKnown;
95
97
104
108 ct::Tensor psi = {}, hpsi = {};
109
111
114 ct::Tensor grad = {}, hgrad = {}, grad_old = {};
115
118
119 // These are for hsolver gemm_op use
121 Device * ctx = {};
122 // Pointer to objects of 1 and 0 for gemm
123 const T *one = nullptr, *zero = nullptr, *neg_one = nullptr;
124 const T one_ = static_cast<T>(1.0), zero_ = static_cast<T>(0.0), neg_one_ = static_cast<T>(-1.0);
125
142 void calc_prec();
143
162 const HPsiFunc& hpsi_func,
163 T *psi_in,
164 ct::Tensor& hpsi_out);
165
180 void diag_hsub(
181 const ct::Tensor& psi_in,
182 const ct::Tensor& hpsi_in,
183 ct::Tensor& hsub_out,
184 ct::Tensor& eigenvalue_out);
185
196 void rotate_wf(
197 const ct::Tensor& hsub_in,
198 ct::Tensor& psi_out,
199 ct::Tensor& workspace_in);
200
227 const ct::Tensor& prec_in,
228 ct::Tensor& err_out,
229 ct::Tensor& beta_out,
230 ct::Tensor& psi_in, ct::Tensor& hpsi_in,
231 ct::Tensor& grad_out, ct::Tensor& grad_old_out);
232
250 const HPsiFunc& hpsi_func,
251 T *psi_in,
252 ct::Tensor& psi_out, ct::Tensor& hpsi_out,
253 ct::Tensor& hsub_out, ct::Tensor& workspace_in,
254 ct::Tensor& eigenvalue_out);
255
273 ct::Tensor& psi_out,
274 ct::Tensor& hpsi_out,
275 ct::Tensor& hsub_out,
276 ct::Tensor& workspace_in,
277 ct::Tensor& eigenvalue_out);
278
289 void orth_projection(
290 const ct::Tensor& psi_in,
291 ct::Tensor& hsub_in,
292 ct::Tensor& grad_out);
293
307 void line_minimize(
308 ct::Tensor& grad_in,
309 ct::Tensor& hgrad_in,
310 ct::Tensor& psi_out,
311 ct::Tensor& hpsi_out);
312
321 void orth_cholesky(
322 ct::Tensor& workspace_in,
323 ct::Tensor& psi_out,
324 ct::Tensor& hpsi_out,
325 ct::Tensor& hsub_out);
326
334 bool test_error(const ct::Tensor& err_in, const std::vector<double>& ethr_band);
335
337 using setmem_var_op = ct::kernels::set_memory<Real, ct_Device>;
338 using resmem_var_op = ct::kernels::resize_memory<Real, ct_Device>;
339 using delmem_var_op = ct::kernels::delete_memory<Real, ct_Device>;
340 using syncmem_var_h2d_op = ct::kernels::synchronize_memory<Real, ct_Device, ct::DEVICE_CPU>;
341 using syncmem_var_d2h_op = ct::kernels::synchronize_memory<Real, ct::DEVICE_CPU, ct_Device>;
342
343 using setmem_complex_op = ct::kernels::set_memory<T, ct_Device>;
344 using delmem_complex_op = ct::kernels::delete_memory<T, ct_Device>;
345 using resmem_complex_op = ct::kernels::resize_memory<T, ct_Device>;
346 using syncmem_complex_op = ct::kernels::synchronize_memory<T, ct_Device, ct_Device>;
347
348 // note: these operators use template parameter base_device::Device_*
349 // defined in source_base/module_device/types.h
350 // different from ct_Device!
352
353};
354
355} // namespace hsolver
356#endif // DIAGO_BPCG_H_
this class is used to perform parallel matrix multiplication C = alpha * A^H * B + beta * C Here,...
Definition para_gemm.h:25
A multi-dimensional array of elements of a single data type.
Definition tensor.h:32
A class for diagonalization using the Blocked-PCG method.
Definition diago_bpcg.h:26
ct::Tensor hpsi
Definition diago_bpcg.h:108
int n_dim
valid dimension of psi
Definition diago_bpcg.h:84
ct::Tensor h_prec
Definition diago_bpcg.h:96
void init_iter(const int nband, const int nband_l, const int nbasis, const int ndim)
Initialize the class before diagonalization.
Definition diago_bpcg.cpp:37
void calc_hsub_with_block_exit(ct::Tensor &psi_out, ct::Tensor &hpsi_out, ct::Tensor &hsub_out, ct::Tensor &workspace_in, ct::Tensor &eigenvalue_out)
Apply the Hamiltonian operator to psi and obtain the hpsi matrix.
Definition diago_bpcg.cpp:244
ct::Tensor prec
Definition diago_bpcg.h:96
void calc_grad_with_block(const ct::Tensor &prec_in, ct::Tensor &err_out, ct::Tensor &beta_out, ct::Tensor &psi_in, ct::Tensor &hpsi_in, ct::Tensor &grad_out, ct::Tensor &grad_old_out)
Calculate the gradient for all bands used in CG method.
Definition diago_bpcg.cpp:137
int n_basis
the number of cols of the input psi
Definition diago_bpcg.h:82
void rotate_wf(const ct::Tensor &hsub_in, ct::Tensor &psi_out, ct::Tensor &workspace_in)
Inplace matrix multiplication to obtain the initial guessed wavefunction.
Definition diago_bpcg.cpp:180
void orth_projection(const ct::Tensor &psi_in, ct::Tensor &hsub_in, ct::Tensor &grad_out)
Orthogonalize column vectors in grad to column vectors in psi.
Definition diago_bpcg.cpp:165
const T zero_
Definition diago_bpcg.h:124
ct::kernels::set_memory< Real, ct_Device > setmem_var_op
Definition diago_bpcg.h:337
~DiagoBPCG()
Destructor for DiagoBPCG class.
Definition diago_bpcg.cpp:32
ct::kernels::synchronize_memory< Real, ct::DEVICE_CPU, ct_Device > syncmem_var_d2h_op
Definition diago_bpcg.h:341
void diag_hsub(const ct::Tensor &psi_in, const ct::Tensor &hpsi_in, ct::Tensor &hsub_out, ct::Tensor &eigenvalue_out)
Diagonalization of the subspace matrix.
Definition diago_bpcg.cpp:203
ct::kernels::synchronize_memory< T, ct_Device, ct_Device > syncmem_complex_op
Definition diago_bpcg.h:346
int nline
max iter steps for all-band cg loop
Definition diago_bpcg.h:86
void calc_hsub_with_block(const HPsiFunc &hpsi_func, T *psi_in, ct::Tensor &psi_out, ct::Tensor &hpsi_out, ct::Tensor &hsub_out, ct::Tensor &workspace_in, ct::Tensor &eigenvalue_out)
Apply the Hamiltonian operator to psi and obtain the hpsi matrix.
Definition diago_bpcg.cpp:219
const T * zero
Definition diago_bpcg.h:123
const T * one
Definition diago_bpcg.h:123
ct::Tensor grad
Definition diago_bpcg.h:114
int n_band
the number of bands of all processes
Definition diago_bpcg.h:78
ct::DataType t_type
Definition diago_bpcg.h:93
std::function< void(T *, T *, const int, const int)> HPsiFunc
Definition diago_bpcg.h:59
ct::kernels::resize_memory< Real, ct_Device > resmem_var_op
Definition diago_bpcg.h:338
ct::Tensor hgrad
Definition diago_bpcg.h:114
const T one_
Definition diago_bpcg.h:124
int n_band_l
the number of bands of current process
Definition diago_bpcg.h:80
Device * ctx
ctx is nothing but the devices used in gemm_op (Device * ctx = nullptr;),
Definition diago_bpcg.h:121
ct::kernels::set_memory< T, ct_Device > setmem_complex_op
Definition diago_bpcg.h:343
ModuleBase::PGemmCN< T, Device > pmmcn
parallel matrix multiplication
Definition diago_bpcg.h:88
typename ct::PsiToContainer< Device >::type ct_Device
Definition diago_bpcg.h:336
void diag(const HPsiFunc &hpsi_func, T *psi_in, Real *eigenvalue_in, const std::vector< double > &ethr_band)
Diagonalize the Hamiltonian using the BPCG method.
Definition diago_bpcg.cpp:262
ct::Tensor beta
The coefficient for mixing the current and previous step gradients, used in iterative methods.
Definition diago_bpcg.h:99
const T neg_one_
Definition diago_bpcg.h:124
void calc_hpsi_with_block(const HPsiFunc &hpsi_func, T *psi_in, ct::Tensor &hpsi_out)
Apply the H operator to psi and obtain the hpsi matrix.
Definition diago_bpcg.cpp:193
ct::kernels::synchronize_memory< Real, ct_Device, ct::DEVICE_CPU > syncmem_var_h2d_op
Definition diago_bpcg.h:340
void orth_cholesky(ct::Tensor &workspace_in, ct::Tensor &psi_out, ct::Tensor &hpsi_out, ct::Tensor &hsub_out)
Orthogonalize and normalize the column vectors in psi_out using Cholesky decomposition.
Definition diago_bpcg.cpp:114
void line_minimize(ct::Tensor &grad_in, ct::Tensor &hgrad_in, ct::Tensor &psi_out, ct::Tensor &hpsi_out)
Optimize psi as well as the hpsi.
Definition diago_bpcg.cpp:96
ct::Tensor err_st
Error state value, if it is smaller than the given threshold, then exit the iteration.
Definition diago_bpcg.h:101
ct::kernels::delete_memory< T, ct_Device > delmem_complex_op
Definition diago_bpcg.h:344
PLinearTransform< T, Device > plintrans
Definition diago_bpcg.h:89
ct::Tensor eigen
Calculated eigen.
Definition diago_bpcg.h:103
ct::Tensor grad_old
Definition diago_bpcg.h:114
ct::DataType r_type
Definition diago_bpcg.h:92
ct::Tensor hsub
Definition diago_bpcg.h:110
bool test_error(const ct::Tensor &err_in, const std::vector< double > &ethr_band)
Checks if the error satisfies the given threshold.
Definition diago_bpcg.cpp:70
const T * neg_one
Definition diago_bpcg.h:123
ct::kernels::resize_memory< T, ct_Device > resmem_complex_op
Definition diago_bpcg.h:345
void calc_prec()
Update the precondition array.
Definition diago_bpcg.cpp:159
typename GetTypeReal< T >::type Real
Definition diago_bpcg.h:31
ct::kernels::delete_memory< Real, ct_Device > delmem_var_op
Definition diago_bpcg.h:339
ct::Tensor work
work for some calculations within this class, including rotate_wf call
Definition diago_bpcg.h:117
ct::DeviceType device_type
Definition diago_bpcg.h:94
B = alpha * A * U + beta * B A and B are local matrice U can be a local matrix or a global matrix.
Definition para_linear_transform.h:20
#define T
Definition exp.cpp:237
DataType
Enumeration of data types for tensors. The DataType enum lists the supported data types for tensors....
Definition tensor_types.h:50
DeviceType
The type of memory used by an allocator.
Definition tensor_types.h:73
Definition diag_comm_info.h:9
Definition exx_lip.h:23
T type
Definition macros.h:8
Definition math_kernel_op.h:216
Definition tensor_types.h:113