ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
Loading...
Searching...
No Matches
to_wannier90_pw.h
Go to the documentation of this file.
1#ifndef TOWannier90_PW_H
2#define TOWannier90_PW_H
3
4#include <iostream>
5#include <algorithm>
6#include <cmath>
7#include <cstdlib>
8#include <vector>
9
10#include "to_wannier90.h"
11
15#include "source_base/matrix.h"
16#include "source_base/matrix3.h"
17#include "source_cell/klist.h"
19#include "source_psi/psi.h"
20
22{
23 public:
25 const bool &out_wannier_mmn,
26 const bool &out_wannier_amn,
27 const bool &out_wannier_unk,
28 const bool &out_wannier_eig,
29 const bool &out_wannier_wvfn_formatted,
30 const std::string &nnkpfile,
31 const std::string &wannier_spin
32 );
34
35 void calculate(
36 const UnitCell& ucell,
37 const ModuleBase::matrix& ekb,
38 const ModulePW::PW_Basis_K* wfcpw,
39 const ModulePW::PW_Basis_Big* bigpw,
40 const K_Vectors& kv,
41 const psi::Psi<std::complex<double>>* psi
42 );
43
45 const UnitCell& ucell,
46 const ModuleBase::matrix& ekb,
47 const ModulePW::PW_Basis_K* wfcpw,
48 const ModulePW::PW_Basis_Big* bigpw,
49 const K_Vectors& kv,
51 )
52 {
53 throw std::logic_error("The wave function of toWannier90_PW is generally a std::complex<double> type.");
54 }
55
56 void cal_Amn(const psi::Psi<std::complex<double>>& psi_pw, const ModulePW::PW_Basis_K* wfcpw);
57 void cal_Mmn(const psi::Psi<std::complex<double>>& psi_pw, const ModulePW::PW_Basis_K* wfcpw);
58 void out_unk(
59 const psi::Psi<std::complex<double>>& psi_pw,
60 const ModulePW::PW_Basis_K* wfcpw,
61 const ModulePW::PW_Basis_Big* bigpw
62 );
63 void set_tpiba_omega(const double& tpiba, const double& omega);
64 protected:
65 // Radial section of trial orbitals
66 const int mesh_r = 333;
67 const double dx = 0.025;
68 const double x_min = -6.0;
69 double const *tpiba;
70 double const *omega;
71
72 void unkdotkb(
73 const psi::Psi<std::complex<double>>& psi_pw,
74 const ModulePW::PW_Basis_K* wfcpw,
75 const int& ik,
76 const int& ikb,
79 );
80
81 void gen_radial_function_in_q(std::vector<ModuleBase::matrix> &radial_in_q);
82
83 void integral(
84 const int meshr,
85 const double *psir,
86 const double *r,
87 const double *rab,
88 const int &l,
89 double *table
90 );
91
93 const psi::Psi<std::complex<double>>& psi_pw,
94 const int& ik,
95 const ModulePW::PW_Basis_K* wfcpw,
96 const std::vector<ModuleBase::matrix> &radial_in_q,
97 ModuleBase::ComplexMatrix& trial_orbitals_k
98 );
99
101 const int &orbital_L,
102 const int &orbital_m,
103 const ModuleBase::matrix &ylm,
105 const int &npw,
106 double *radial_in_q_single,
107 std::complex<double> *orbital_in_G_single
108 );
109
110 void unkdotW_A(
111 const psi::Psi<std::complex<double>>& psi_pw,
112 const ModulePW::PW_Basis_K* wfcpw,
113 const int& ik,
114 const std::vector<ModuleBase::matrix> &radial_in_q,
116 );
117
118};
119
120#endif
Definition klist.h:13
Definition complexmatrix.h:14
3 elements vector
Definition vector3.h:22
Definition matrix.h:19
Definition pw_basis_big.h:16
Special pw_basis class. It includes different k-points.
Definition pw_basis_k.h:57
Definition unitcell.h:16
Definition psi.h:37
Definition to_wannier90_pw.h:22
void unkdotW_A(const psi::Psi< std::complex< double > > &psi_pw, const ModulePW::PW_Basis_K *wfcpw, const int &ik, const std::vector< ModuleBase::matrix > &radial_in_q, ModuleBase::ComplexMatrix &Amn)
Definition to_wannier90_pw.cpp:1014
void produce_trial_in_pw(const psi::Psi< std::complex< double > > &psi_pw, const int &ik, const ModulePW::PW_Basis_K *wfcpw, const std::vector< ModuleBase::matrix > &radial_in_q, ModuleBase::ComplexMatrix &trial_orbitals_k)
Definition to_wannier90_pw.cpp:579
const double x_min
Definition to_wannier90_pw.h:68
void calculate(const UnitCell &ucell, const ModuleBase::matrix &ekb, const ModulePW::PW_Basis_K *wfcpw, const ModulePW::PW_Basis_Big *bigpw, const K_Vectors &kv, const psi::Psi< double > *psi)
Definition to_wannier90_pw.h:44
double const * omega
Definition to_wannier90_pw.h:70
void gen_radial_function_in_q(std::vector< ModuleBase::matrix > &radial_in_q)
Definition to_wannier90_pw.cpp:494
double const * tpiba
Definition to_wannier90_pw.h:69
void set_tpiba_omega(const double &tpiba, const double &omega)
Definition to_wannier90_pw.cpp:78
~toWannier90_PW()
Definition to_wannier90_pw.cpp:25
void unkdotkb(const psi::Psi< std::complex< double > > &psi_pw, const ModulePW::PW_Basis_K *wfcpw, const int &ik, const int &ikb, const ModuleBase::Vector3< double > G, ModuleBase::ComplexMatrix &Mmn)
Definition to_wannier90_pw.cpp:357
const int mesh_r
Definition to_wannier90_pw.h:66
void integral(const int meshr, const double *psir, const double *r, const double *rab, const int &l, double *table)
Definition to_wannier90_pw.cpp:972
void get_trial_orbitals_lm_k(const int &orbital_L, const int &orbital_m, const ModuleBase::matrix &ylm, const ModuleBase::Vector3< double > *gk, const int &npw, double *radial_in_q_single, std::complex< double > *orbital_in_G_single)
Definition to_wannier90_pw.cpp:946
const double dx
Definition to_wannier90_pw.h:67
Definition to_wannier90.h:23
void cal_Amn()
Definition to_wannier90.cpp:124
bool out_wannier_wvfn_formatted
Definition to_wannier90.h:76
std::string nnkpfile
Definition to_wannier90.h:78
std::string wannier_spin
Definition to_wannier90.h:80
void cal_Mmn()
Definition to_wannier90.cpp:128
void out_unk()
Definition to_wannier90.cpp:120
bool out_wannier_unk
Definition to_wannier90.h:74
bool out_wannier_mmn
Definition to_wannier90.h:72
void calculate()
Definition to_wannier90.cpp:69
bool out_wannier_amn
Definition to_wannier90.h:73
bool out_wannier_eig
Definition to_wannier90.h:75
Definition exx_lip.h:23