ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
berryphase.h
Go to the documentation of this file.
1#ifndef MODULE_IO_BERRYPHASE_H
2#define MODULE_IO_BERRYPHASE_H
3#include "unk_overlap_pw.h"
4#ifdef __LCAO
5#include "unk_overlap_lcao.h"
6#endif
8#include "source_cell/klist.h"
9#include "source_psi/psi.h"
11
13{
14
15 public:
16 berryphase(); // for pw-line
17#ifdef __LCAO
18 berryphase(const Parallel_Orbitals* paraV_in); // for lcao-line
19#endif
21
22 // mohan add 2021-02-16
23 static bool berry_phase_flag;
25#ifdef __LCAO
26 unkOverlap_lcao lcao_method;
27 const Parallel_Orbitals* paraV;
28#endif
29
31 std::vector<std::vector<int>> k_index;
32 int nppstr=0;
33 int direction=0;
35 int GDIR;
36
38#ifdef __LCAO
39 void lcao_init(const UnitCell& ucell,
40 const Grid_Driver& gd,
41 const K_Vectors& kv,
42 const LCAO_Orbitals& orb);
43#endif
44 void set_kpoints(const K_Vectors& kv, const int direction);
45
46 double stringPhase(const UnitCell& ucell,
47 int index_str,
48 int nbands,
49 const int npwx,
50 const psi::Psi<std::complex<double>>* psi_in,
51 const ModulePW::PW_Basis* rhopw,
52 const ModulePW::PW_Basis_K* wfcpw,
53 const K_Vectors& kv);
54
55 void Berry_Phase(const UnitCell& ucell,
56 int nbands,
57 double& pdl_elec_tot,
58 int& mod_elec_tot,
59 const int npwx,
60 const psi::Psi<std::complex<double>>* psi_in,
61 const ModulePW::PW_Basis* rhopw,
62 const ModulePW::PW_Basis_K* wfcpw,
63 const K_Vectors& kv);
64
66 const int npwx,
67 const psi::Psi<double>* psi_in,
68 const ModulePW::PW_Basis* rhopw,
69 const ModulePW::PW_Basis_K* wfcpw,
70 const K_Vectors& kv)
71 {
72 throw std::logic_error("berry phase supports only multi-k");
73 };
74 void Macroscopic_polarization(const UnitCell& ucell,
75 const int npwx,
76 const psi::Psi<std::complex<double>>* psi_in,
77 const ModulePW::PW_Basis* rhopw,
78 const ModulePW::PW_Basis_K* wfcpw,
79 const K_Vectors& kv);
80
81 std::string outFormat(const double polarization, const double modulus, const ModuleBase::Vector3<double> project);
82};
83
84#endif
Definition sltk_grid_driver.h:43
Definition klist.h:13
Definition ORB_read.h:19
3 elements vector
Definition vector3.h:22
Special pw_basis class. It includes different k-points.
Definition pw_basis_k.h:57
A class which can convert a function of "r" to the corresponding linear superposition of plane waves ...
Definition pw_basis.h:56
Definition parallel_orbitals.h:9
Definition unitcell.h:17
Definition berryphase.h:13
void set_kpoints(const K_Vectors &kv, const int direction)
Definition berryphase.cpp:58
int GDIR
Definition berryphase.h:35
static bool berry_phase_flag
Definition berryphase.h:23
berryphase()
Definition berryphase.cpp:11
void Macroscopic_polarization(const UnitCell &ucell, const int npwx, const psi::Psi< double > *psi_in, const ModulePW::PW_Basis *rhopw, const ModulePW::PW_Basis_K *wfcpw, const K_Vectors &kv)
Definition berryphase.h:65
void get_occupation_bands()
Definition berryphase.cpp:27
~berryphase()
Definition berryphase.cpp:23
int occ_nbands
Definition berryphase.h:34
void Berry_Phase(const UnitCell &ucell, int nbands, double &pdl_elec_tot, int &mod_elec_tot, const int npwx, const psi::Psi< std::complex< double > > *psi_in, const ModulePW::PW_Basis *rhopw, const ModulePW::PW_Basis_K *wfcpw, const K_Vectors &kv)
Definition berryphase.cpp:346
std::vector< std::vector< int > > k_index
Definition berryphase.h:31
double stringPhase(const UnitCell &ucell, int index_str, int nbands, const int npwx, const psi::Psi< std::complex< double > > *psi_in, const ModulePW::PW_Basis *rhopw, const ModulePW::PW_Basis_K *wfcpw, const K_Vectors &kv)
Definition berryphase.cpp:220
unkOverlap_pw pw_method
Definition berryphase.h:24
int total_string
Definition berryphase.h:30
int nppstr
Definition berryphase.h:32
int direction
Definition berryphase.h:33
std::string outFormat(const double polarization, const double modulus, const ModuleBase::Vector3< double > project)
Definition berryphase.cpp:680
Definition psi.h:37
Definition unk_overlap_lcao.h:22
Definition unk_overlap_pw.h:18