ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
rdmft_tools.h
Go to the documentation of this file.
1//==========================================================
2// Author: Jingang Han
3// DATE : 2024-03-11
4//==========================================================
5#ifndef RDMFT_TOOLS_H
6#define RDMFT_TOOLS_H
7
8#include "source_cell/klist.h"
9#include "source_psi/psi.h"
10#include "source_base/matrix.h"
20
24
25
26#ifdef __EXX
28// there are some operator reload to print data in different formats
29#endif
30
31#include <iostream>
32#include <type_traits>
33#include <complex>
34#include <vector>
35#include <iomanip>
36
37
38
39namespace rdmft
40{
41
43// for the dft-xc-functional part of xc-functional, just use the default is right! Or don't use the function
44double occNum_func(double eta, int symbol = 0, const std::string XC_func_rdmft = "hf", const double alpha_power = 1.0);
45
46
47template <typename TK>
49{
50 TK* pwfc = &wfc(0, 0, 0);
51 for(int i=0; i<wfc.size(); ++i) { pwfc[i] = std::conj( pwfc[i] );
52}
53}
54
55
56template <>
58
59
60// wfc and H_wfc need to be k_firest and provide wfc(ik, 0, 0) and H_wfc(ik, 0, 0)
62template <typename TK>
63void HkPsi(const Parallel_Orbitals* ParaV, const TK& HK, const TK& wfc, TK& H_wfc)
64{
65
66 const int one_int = 1;
67 //const double one_double = 1.0, zero_double = 0.0;
68 const std::complex<double> one_complex = {1.0, 0.0};
69 const std::complex<double> zero_complex = {0.0, 0.0};
70 const char N_char = 'N';
71 const char C_char = 'C'; // Using 'C' is consistent with the formula
72
73#ifdef __MPI
74 const int nbasis = ParaV->desc[2];
75 const int nbands = ParaV->desc_wfc[3];
76
77 //because wfc(bands, basis'), H(basis, basis'), we do wfc*H^T(in the perspective of cpp, not in fortran). And get H_wfc(bands, basis) is correct.
78 ScalapackConnector::gemm( C_char, N_char, nbasis, nbands, nbasis, one_complex, &HK, one_int, one_int, ParaV->desc,
79 &wfc, one_int, one_int, ParaV->desc_wfc, zero_complex, &H_wfc, one_int, one_int, ParaV->desc_wfc );
80#endif
81}
82
83
84template <>
85void HkPsi<double>(const Parallel_Orbitals* ParaV, const double& HK, const double& wfc, double& H_wfc);
86
87
89template <typename TK>
90void cal_bra_op_ket(const Parallel_Orbitals* ParaV, const Parallel_2D& para_Eij_in, const TK& wfc, const TK& H_wfc, std::vector<TK>& Dmn)
91{
92 const int one_int = 1;
93 const std::complex<double> one_complex = {1.0, 0.0};
94 const std::complex<double> zero_complex = {0.0, 0.0};
95 const char N_char = 'N';
96 const char C_char = 'C';
97
98 const int nrow_bands = para_Eij_in.get_row_size();
99 const int ncol_bands = para_Eij_in.get_col_size();
100
101#ifdef __MPI
102 const int nbasis = ParaV->desc[2];
103 const int nbands = ParaV->desc_wfc[3];
104
105 ScalapackConnector::gemm( C_char, N_char, nbands, nbands, nbasis, one_complex, &wfc, one_int, one_int, ParaV->desc_wfc,
106 &H_wfc, one_int, one_int, ParaV->desc_wfc, zero_complex, &Dmn[0], one_int, one_int, para_Eij_in.desc );
107#endif
108}
109
110
111template <>
112void cal_bra_op_ket<double>(const Parallel_Orbitals* ParaV, const Parallel_2D& para_Eij_in,
113 const double& wfc, const double& H_wfc, std::vector<double>& Dmn);
114
115
117template <typename TK>
118void _diagonal_in_serial(const Parallel_2D& para_Eij_in, const std::vector<TK>& Dmn, double* wfcHwfc)
119{
120 const int nrow_bands = para_Eij_in.get_row_size();
121 const int ncol_bands = para_Eij_in.get_col_size();
122
123 for(int i=0; i<nrow_bands; ++i)
124 {
125 int i_global = para_Eij_in.local2global_row(i);
126 for(int j=0; j<ncol_bands; ++j)
127 {
128 int j_global = para_Eij_in.local2global_col(j);
129 if(i_global==j_global)
130 {
131 // because the Dmn obtained from pzgemm_() is stored column-major
132 wfcHwfc[j_global] = std::real( Dmn[i+j*nrow_bands] );
133 }
134 }
135 }
136}
137
138
140template <typename TK>
141void occNum_MulPsi(const Parallel_Orbitals* ParaV, const ModuleBase::matrix& occ_number, psi::Psi<TK>& wfc, int symbol = 0,
142 const std::string XC_func_rdmft = "hf", const double alpha = 1.0)
143{
144 const int nk_local = wfc.get_nk();
145 const int nbands_local = wfc.get_nbands();
146 const int nbasis_local = wfc.get_nbasis();
147
148 // const int nbasis = ParaV->desc[2]; // need to be deleted
149 // const int nbands = ParaV->desc_wfc[3];
150
151 for (int ik = 0; ik < nk_local; ++ik)
152 {
153 for (int ib_local = 0; ib_local < nbands_local; ++ib_local) // ib_local < nbands_local , some problem, ParaV->ncol_bands
154 {
155 const double occNum_local = occNum_func( occ_number(ik, ParaV->local2global_col(ib_local)), symbol, XC_func_rdmft, alpha);
156 TK* wfc_pointer = &(wfc(ik, ib_local, 0));
157 BlasConnector::scal(nbasis_local, occNum_local, wfc_pointer, 1);
158 }
159 }
160}
161
162
164template <typename TK>
165void add_psi(const Parallel_Orbitals* ParaV,
166 const K_Vectors* kv,
167 const ModuleBase::matrix& occ_number,
168 psi::Psi<TK>& psi_TV,
169 psi::Psi<TK>& psi_hartree,
170 psi::Psi<TK>& psi_dft_XC,
171 psi::Psi<TK>& psi_exx_XC,
172 psi::Psi<TK>& occNum_Hpsi,
173 const std::string XC_func_rdmft = "hf",
174 const double alpha = 1.0)
175{
176 const int nk = psi_TV.get_nk();
177 const int nbn_local = psi_TV.get_nbands();
178 const int nbs_local = psi_TV.get_nbasis();
179 occNum_MulPsi(ParaV, occ_number, psi_TV);
180 occNum_MulPsi(ParaV, occ_number, psi_hartree);
181 occNum_MulPsi(ParaV, occ_number, psi_dft_XC);
182 occNum_MulPsi(ParaV, occ_number, psi_exx_XC, 2, XC_func_rdmft, alpha);
183
184 // const int nbasis = ParaV->desc[2];
185 // const int nbands = ParaV->desc_wfc[3];
186
187 for(int ik=0; ik<nk; ++ik)
188 {
189 for(int inbn=0; inbn<nbn_local; ++inbn)
190 {
191 TK* p_occNum_Hpsi = &( occNum_Hpsi(ik, inbn, 0) );
192 for(int inbs=0; inbs<nbs_local; ++inbs)
193 {
194 p_occNum_Hpsi[inbs] = psi_TV(ik, inbn, inbs) + psi_hartree(ik, inbn, inbs) + psi_dft_XC(ik, inbn, inbs) + psi_exx_XC(ik, inbn, inbs);
195 }
196
197 // test, consider the wk into psi or dE/d(wfc)
198 BlasConnector::scal(nbs_local, kv->wk[ik], p_occNum_Hpsi, 1);
199 }
200 }
201
202}
203
209void occNum_Mul_wfcHwfc(const ModuleBase::matrix& occ_number,
210 const ModuleBase::matrix& wfcHwfc,
211 ModuleBase::matrix& occNum_wfcHwfc,
212 int symbol = 0,
213 const std::string XC_func_rdmft = "hf",
214 const double alpha = 1.0);
215
216
221void add_occNum(const K_Vectors& kv,
222 const ModuleBase::matrix& occ_number,
223 const ModuleBase::matrix& wfcHwfc_TV_in,
224 const ModuleBase::matrix& wfcHwfc_hartree_in,
225 const ModuleBase::matrix& wfcHwfc_dft_XC_in,
226 const ModuleBase::matrix& wfcHwfc_exx_XC_in,
227 ModuleBase::matrix& occNum_wfcHwfc,
228 const std::string XC_func_rdmft = "hf",
229 const double alpha = 1.0);
230
231
233void add_wfcHwfc(const ModuleBase::matrix& wg,
234 const ModuleBase::matrix& wk_fun_occNum,
235 const ModuleBase::matrix& wfcHwfc_TV_in,
236 const ModuleBase::matrix& wfcHwfc_hartree_in,
237 const ModuleBase::matrix& wfcHwfc_XC_in,
238 ModuleBase::matrix& occNum_wfcHwfc,
239 const std::string XC_func_rdmft,
240 const double alpha);
241
242
244double getEnergy(const ModuleBase::matrix& occNum_wfcHwfc);
245
246
247
248
249
250
252template <typename TK, typename TR>
253class Veff_rdmft : public hamilt::OperatorLCAO<TK, TR>
254{
255 public:
260 const std::vector<ModuleBase::Vector3<double>>& kvec_d_in,
261 elecstate::Potential* pot_in,
263 const UnitCell* ucell_in,
264 const std::vector<double>& orb_cutoff,
265 const Grid_Driver* GridD_in,
266 const int& nspin,
267 const Charge* charge_in,
268 const ModulePW::PW_Basis* rho_basis_in,
269 const ModuleBase::matrix* vloc_in,
270 const ModuleBase::ComplexMatrix* sf_in,
271 const std::string potential_in,
272 double* etxc_in = nullptr,
273 double* vtxc_in = nullptr)
274 : orb_cutoff_(orb_cutoff), pot(pot_in), ucell(ucell_in),
275 gd(GridD_in), hamilt::OperatorLCAO<TK, TR>(hsk_in, kvec_d_in, hR_in), charge_(charge_in),
276 rho_basis_(rho_basis_in), vloc_(vloc_in), sf_(sf_in), potential_(potential_in), etxc(etxc_in), vtxc(vtxc_in)
277 {
279
280 this->initialize_HR(ucell_in, GridD_in);
281 }
282
284
291 virtual void contributeHR() override;
292
293 const UnitCell* ucell = nullptr;
294
295 const Grid_Driver* gd = nullptr;
296
297 private:
298
299 std::vector<double> orb_cutoff_;
300
301 // Charge calculating method in LCAO base and contained grid base calculation: DM_R, DM, pvpR_reduced
302
304
305 int nspin = 1;
307
313 void initialize_HR(const UnitCell* ucell_in, const Grid_Driver* GridD_in);
314
315 // added by jghan
316
317 const Charge* charge_ = nullptr;
318
319 std::string potential_;
320
322
324
326
327 double* etxc = nullptr;
328
329 double* vtxc = nullptr;
330
331};
332
333}
334
335#endif
const std::complex< double > i
Definition cal_pLpR.cpp:46
static void scal(const int n, const float alpha, float *X, const int incX, base_device::AbacusDevice_t device_type=base_device::AbacusDevice_t::CpuDevice)
Definition blas_connector_vector.cpp:80
Definition charge.h:17
Definition sltk_grid_driver.h:40
Definition klist.h:12
std::vector< double > wk
Definition klist.h:18
Definition complexmatrix.h:13
3 elements vector
Definition vector3.h:24
Definition matrix.h:18
A class which can convert a function of "r" to the corresponding linear superposition of plane waves ...
Definition pw_basis.h:56
This class packs the basic information of 2D-block-cyclic parallel distribution of an arbitrary matri...
Definition parallel_2d.h:12
int local2global_col(const int ilc) const
get the global index of a local index (col)
Definition parallel_2d.h:63
int get_row_size() const
number of local rows
Definition parallel_2d.h:21
int local2global_row(const int ilr) const
get the global index of a local index (row)
Definition parallel_2d.h:57
int desc[9]
ScaLAPACK descriptor.
Definition parallel_2d.h:103
int get_col_size() const
number of local columns
Definition parallel_2d.h:27
Definition parallel_orbitals.h:9
int desc_wfc[9]
Definition parallel_orbitals.h:37
static void gemm(const char transa, const char transb, const int M, const int N, const int K, const double alpha, const double *A, const int IA, const int JA, const int *DESCA, const double *B, const int IB, const int JB, const int *DESCB, const double beta, double *C, const int IC, const int JC, const int *DESCC)
Definition scalapack_connector.h:261
Definition unitcell.h:15
Definition potential_new.h:48
Definition hcontainer.h:144
Definition hs_matrix_k.hpp:11
Definition operator_lcao.h:12
OperatorLCAO(HS_Matrix_K< TK > *hsk_in, const std::vector< ModuleBase::Vector3< double > > &kvec_d_in, HContainer< TR > *hR_in)
H(R) matrix, R is the Bravis lattice vector.
Definition operator_lcao.h:15
enum calculation_type cal_type
Definition operator.h:105
Definition psi.h:37
const int & get_nbands() const
Definition psi.cpp:341
const int & get_nk() const
Definition psi.cpp:335
size_t size() const
Definition psi.cpp:353
const int & get_nbasis() const
Definition psi.cpp:347
this part of the code is copying from class Veff and do some modifications.
Definition rdmft_tools.h:254
double * etxc
Definition rdmft_tools.h:327
double * vtxc
Definition rdmft_tools.h:329
std::string potential_
Definition rdmft_tools.h:319
const ModuleBase::ComplexMatrix * sf_
Definition rdmft_tools.h:325
const Grid_Driver * gd
Definition rdmft_tools.h:295
elecstate::Potential * pot
Definition rdmft_tools.h:303
Veff_rdmft(hamilt::HS_Matrix_K< TK > *hsk_in, const std::vector< ModuleBase::Vector3< double > > &kvec_d_in, elecstate::Potential *pot_in, hamilt::HContainer< TR > *hR_in, const UnitCell *ucell_in, const std::vector< double > &orb_cutoff, const Grid_Driver *GridD_in, const int &nspin, const Charge *charge_in, const ModulePW::PW_Basis *rho_basis_in, const ModuleBase::matrix *vloc_in, const ModuleBase::ComplexMatrix *sf_in, const std::string potential_in, double *etxc_in=nullptr, double *vtxc_in=nullptr)
Construct a new Veff object for multi-kpoint calculation.
Definition rdmft_tools.h:259
int current_spin
Definition rdmft_tools.h:306
virtual void contributeHR() override
contributeHR() is used to calculate the HR matrix <phi_{\mu, 0}|V_{eff}|phi_{\nu, R}> the contributio...
const ModulePW::PW_Basis * rho_basis_
Definition rdmft_tools.h:321
const UnitCell * ucell
Definition rdmft_tools.h:293
const Charge * charge_
Definition rdmft_tools.h:317
int nspin
Definition rdmft_tools.h:305
void initialize_HR(const UnitCell *ucell_in, const Grid_Driver *GridD_in)
initialize HR, search the nearest neighbor atoms HContainer is used to store the electronic kinetic m...
Definition rdmft_tools.cpp:191
const ModuleBase::matrix * vloc_
Definition rdmft_tools.h:323
std::vector< double > orb_cutoff_
Definition rdmft_tools.h:299
Definition hamilt.h:13
Reduced Density Matrix Functional Theory (RDMFT)
Definition rdmft.cpp:22
void cal_bra_op_ket< double >(const Parallel_Orbitals *ParaV, const Parallel_2D &para_Eij_in, const double &wfc, const double &H_wfc, std::vector< double > &Dmn)
Definition rdmft_tools.cpp:54
void HkPsi< double >(const Parallel_Orbitals *ParaV, const double &HK, const double &wfc, double &H_wfc)
Definition rdmft_tools.cpp:30
void occNum_Mul_wfcHwfc(const ModuleBase::matrix &occ_number, const ModuleBase::matrix &wfcHwfc, ModuleBase::matrix &occNum_wfcHwfc, int symbol, const std::string XC_func_rdmft, const double alpha)
occNum_wfcHwfc = occNum*wfcHwfc + occNum_wfcHwfc
Definition rdmft_tools.cpp:81
void add_occNum(const K_Vectors &kv, const ModuleBase::matrix &occ_number, const ModuleBase::matrix &wfcHwfc_TV_in, const ModuleBase::matrix &wfcHwfc_hartree_in, const ModuleBase::matrix &wfcHwfc_dft_XC_in, const ModuleBase::matrix &wfcHwfc_exx_XC_in, ModuleBase::matrix &occNum_wfcHwfc, const std::string XC_func_rdmft, const double alpha)
Default symbol = 0 for the gradient of Etotal with respect to occupancy, symbol = 1 for the relevant ...
Definition rdmft_tools.cpp:99
void conj_psi< double >(psi::Psi< double > &wfc)
Definition rdmft_tools.cpp:27
void cal_bra_op_ket(const Parallel_Orbitals *ParaV, const Parallel_2D &para_Eij_in, const TK &wfc, const TK &H_wfc, std::vector< TK > &Dmn)
implement matrix multiplication of sum_mu conj(wfc(ik, m ,mu)) * op_wfc(ik, n, mu)
Definition rdmft_tools.h:90
double getEnergy(const ModuleBase::matrix &occNum_wfcHwfc)
give certain occNum_wfcHwfc, get the corresponding energy
Definition rdmft_tools.cpp:145
void HkPsi(const Parallel_Orbitals *ParaV, const TK &HK, const TK &wfc, TK &H_wfc)
implement matrix multiplication of Hk^dagger and psi
Definition rdmft_tools.h:63
void add_psi(const Parallel_Orbitals *ParaV, const K_Vectors *kv, const ModuleBase::matrix &occ_number, psi::Psi< TK > &psi_TV, psi::Psi< TK > &psi_hartree, psi::Psi< TK > &psi_dft_XC, psi::Psi< TK > &psi_exx_XC, psi::Psi< TK > &occNum_Hpsi, const std::string XC_func_rdmft="hf", const double alpha=1.0)
add psi with eta and g(eta)
Definition rdmft_tools.h:165
void _diagonal_in_serial(const Parallel_2D &para_Eij_in, const std::vector< TK > &Dmn, double *wfcHwfc)
for Dmn that conforms to the 2d-block rule, get its diagonal elements
Definition rdmft_tools.h:118
double occNum_func(const double eta, const int symbol, const std::string XC_func_rdmft, double alpha)
now support XC_func_rdmft = "hf", "muller", "power", "pbe", "pbe0". "wp22" and "cwp22" is realizing.
Definition rdmft_tools.cpp:162
void occNum_MulPsi(const Parallel_Orbitals *ParaV, const ModuleBase::matrix &occ_number, psi::Psi< TK > &wfc, int symbol=0, const std::string XC_func_rdmft="hf", const double alpha=1.0)
realize occNum_wfc = occNum * wfc. Calling this function and we can get wfc = occNum*wfc.
Definition rdmft_tools.h:141
void add_wfcHwfc(const ModuleBase::matrix &wg, const ModuleBase::matrix &wk_fun_occNum, const ModuleBase::matrix &wfcHwfc_TV_in, const ModuleBase::matrix &wfcHwfc_hartree_in, const ModuleBase::matrix &wfcHwfc_XC_in, ModuleBase::matrix &occNum_wfcHwfc, const std::string XC_func_rdmft, const double alpha)
do wk*g(occNum)*wfcHwfc and add for TV, hartree, XC. This function just use once, so it can be replac...
Definition rdmft_tools.cpp:128
void conj_psi(psi::Psi< TK > &wfc)
Definition rdmft_tools.h:48