Skip to main content

Rustb/
conductivity.rs

1//! Calculation of linear and nonlinear conductivity tensors using Kubo formalism.
2//!
3//! This module implements various conductivity calculations including:
4//! - Anomalous Hall conductivity
5//! - Spin Hall conductivity  
6//! - Nonlinear Hall conductivity
7//! - Berry curvature and orbital magnetization
8//!
9//! The implementations are based on the Kubo formula and semiclassical wave-packet
10//! dynamics, providing both intrinsic and extrinsic contributions to transport.
11
12//! Derivation of nonlinear Hall effect using Niu-Qian equations
13//! This section derives formulas for linear and nonlinear Hall conductivities using the Niu-Qian formalism.
14//! Starting from the current density formula:
15//!$$\bm J=-e\int_\tx{BZ}\dd\bm k\sum_n f_n\bm v_n$$
16//! Here $n$ labels bands, $f_n$ is the Fermi-Dirac distribution. The velocity operator according to Niu-Qian is:
17//! $$\bm v=\f{1}{\hbar}\f{\p\ve_n}{\p\bm k}-\f{e}{\hbar}\bm E\times\bm\Og_n$$
18//! The $n$-th order Hall conductivity is defined as:
19//! $$\sg_{\ap_1,\ap_2,\cdots,\ap_n;d}=\f{1}{n!}\left\.\f{\p^n J_d}{\p E_{\ap_1}\cdots\p E_{\ap_n}}\right\vert_{\bm E=0}$$
20//! To obtain its expression, we define series expansions:
21//! $$\lt\\\{\\begin{aligned}
22//! f_n=f_n^{(0)}+f_n^{(1)}+f_n^{(2)}\cdots\\\\
23//! \bm v_n=\bm v_n^{(0)}+\bm v_n^{(1)}+\bm v_n^{(2)}\cdots\\\\
24//! \\end{aligned}\rt\.$$
25//! This gives:
26//! $$ \\begin{aligned}\bm J^{(0)}&=-e\int_\tx{BZ}\dd\bm k\sum_n f_n^{(0)}\bm v_n^{(0)}\\\\
27//! \bm J^{(1)}&=-e\int_\tx{BZ}\dd\bm k\sum_n f_n^{(1)}\bm v_n^{(0)}+f_n^{(0)}\bm v_n^{(1)}\\\\
28//! \bm J^{(2)}&=-e\int_\tx{BZ}\dd\bm k\sum_n f_n^{(2)}\bm v_n^{(0)}+f_n^{(1)}\bm v_n^{(1)}+f_n^{(0)}\bm v_n^{(2)}\\\\
29//! \\end{aligned}$$
30//!
31//! Now consider the corrections to $f$. Using the Boltzmann equation:
32//! $$\p_t f-\f{e}{\hbar}\bm E\cdot\nb_{\bm k} f=-\f{f-f_0}{\tau}$$
33//! Setting $f=\sum_{s=1}e^{is\og t} f_n^{(s)}$, we have:
34//! $$\\begin{aligned} is\og\sum_{s=1}f_n^{(s)}-\f{e}{\hbar}\bm E\cdot\nb_{\bm k}\sum_{s=0} f_n^{(s)}=-\f{1}{\tau}\sum_{s=1} f_n^{(s)}\\\\
35//! \Rightarrow (is\og+\f{1}{\tau})\sum_{s=1} f_n^{(i)}-\f{e}{\hbar}\bm E\cdot\nb_{\bm k}\sum_{i=0} f_n^{(i)}=0\\\\
36//! \\end{aligned}$$
37//! Finally, we obtain the higher-order Fermi distribution:
38//! $$f_n^{(l)}=\f{e}{\hbar} \f{\bm E\nb_{\bm k} f_n^{(l-1)}}{i l \og+1/\tau}=\lt(\f{e/\hbar}{i\og+1/\tau}\rt)\bm E^l\nb^l_{\bm k} f_n^{(0)}$$
39//! Taking the zero-frequency limit: $$\lim_{\og\to 0} f_n^{(l)}\approx \lt(\f{e\tau}{\hbar}\rt)^l \bm E^l\nb^l_{\bm k} f_n^{(0)}$$
40//!
41//! For the Fermi velocity $\bm v_n=\f{1}{\hbar}\pdv{\ve_n}{\bm k}+\f{e}{\hbar}\bm E\times\bm \Og_n$,
42//! we define order-by-order expansions:
43//! $$\\begin{aligned}
44//! \bm v_n^{(0)}&=\f{1}{\hbar}\pdv{\ve_n^{(0)}}{\bm k}\\\\
45//! \bm v_n^{(1)}&=\f{1}{\hbar}\pdv{\ve_n^{(1)}}{\bm k}+\f{e}{\hbar}\bm E\times\bm \Og_n^{(0)}\\\\
46//! \bm v_n^{(2)}&=\f{1}{\hbar}\pdv{\ve_n^{(2)}}{\bm k}+\f{e}{\hbar}\bm E\times\bm \Og_n^{(1)}\\\\
47//! \\end{aligned}$$
48//! Next, starting from the Hamiltonian under electric field:
49//! $$H_{\bm k}=\sum_{mn}\lt(\ve_n^{(0)}\dt_{nm}-e\bm E\cdot\bra{\psi_n}\bm r\ket{\psi_n}\rt)\ket{\psi_n}\bra{\psi_m}$$
50//!我们将其拆成两部分, 对角部分和非对角部分
51//!$$\\begin{aligned}
52//!H_{\bm k}^{(0)}&=\sum_{n}\lt(\ve_{n\bm k}^{(0)}-e\bm E\cdot\bm A_n\rt)\dyad{\psi_n}\\\\
53//!H_{\bm k}^{(1)}&=\sum_{n=\not m}\lt(-e\bm E\cdot\bm A_{mn}\rt)\ket{\psi_m}\bra{\psi_n}\\\\
54//!\\end{aligned}$$
55//!这里 $\bm A_{mn}=\bra{\psi_m}\bm r\ket{\psi_n}=i\bra{\psi_m}\p_{\bm k}\ket{\psi_n}$
56//!
57//!显然, 我们知道公式
58//!$$e^{\hat S}\hat{\mathcal{O}}e^{-\hat S}=\mathcal{O}+\lt[\hat S,\hat{\mcl{O}}\rt]+\f{1}{2}\lt[\hat S,\lt[\hat S,\hat{\mcl{O}}\rt]\rt]+\f{1}{6}\lt[\hat S,\lt[\hat S,\lt[\hat S,\hat{\mcl{O}}\rt]\rt]\rt]\cdots$$
59//! 为了方便计算, 我们可以选择一个 $\hat S$, 让 $H_{\bm k}^{(1)}+\lt[\hat S,\hat H_{\bm k}^{(0)}\rt]=0$, 我们有
60//!$$\\begin{aligned}
61//!H^\prime_{\bm k}&=e^{\hat S}H_{\bm k} e^{-\hat S}=H_{\bm k}^{(0)}+\lt(H_{\bm k}^{(1)}+\lt[\hat S,\hat H_{\bm k}^{(0)}\rt]\rt)+\lt(\lt[\hat S,\hat H_{\bm k}^{(1)}\rt]+\f{1}{2}\lt[\hat S,\lt[\hat S,\hat H_{\bm k}^{(0)}\rt]\rt]\rt)\cdots\\\\
62//!&=H_{\bm k}^{(0)}+\f{1}{2}\lt[S,H_{\bm k}^{(1)}\rt]+\f{1}{3}\lt[S,\lt[S,H_{\bm k}^{(1)}\rt]\rt]\cdots
63//!\\end{aligned}$$
64//!为了满足条件, 我们选择 $$S_{nn}=0,\ S_{nm}=\f{-e\bm E\cdot \bm A_{nm}}{\ve_{nm}-e\bm E\cdot \bm A_{nm}}$$
65//!
66//!因为我们有
67//!$$\\begin{aligned} \lt[S,H_{\bm k}^{(0)}\rt]&=SH_{\bm k}^{(0)}-H_{\bm k}^{(0)}S=\sum_{j=\not m} S_{mj}H_{\bm k,jn}^{(0)}-\sum_{j=\not n }H_{\bm k,mj}^{(0)}S_{jn}\\\\
68//!&=\sum_{j=\not m}\f{-e\bm E\cdot \bm A_{mj}\lt(\ve_j^{(0)}-e\bm E\cdot\bm A_j\rt)\dt_{jn}}{\ve_{mj}-e\bm E\cdot\lt(\bm A_m-\bm A_j\rt)}-\sum_{j=\not n}\f{-e\lt(\ve_j^{(0)}-e\bm E\cdot\bm A_j\rt)\lt(\bm E\cdot \bm A_{jn}\rt)\dt_{mj}}{\ve_{jn}-e\bm E\cdot\lt(\bm A_j-\bm A_n\rt)}\\\\
69//!&=\f{e\lt(\bm E\cdot\bm A_{mn}\rt)\lt[\ve_{mn}- e\bm E\cdot\lt(\bm A_m-\bm A_n\rt)\rt]}{\ve_{mn}-e\bm E\cdot(\bm A_m-\bm A_n)}=-H_{\bm k}^{(1)}
70//!\\end{aligned}$$
71//!这样我们就验证了我们的结果, 我们将 $\hat S$ 进行化简和展开有
72//!$$S_{nm}\approx \f{-e\bm E\cdot\bm A_{nm}}{\ve_n^{(0)}-\ve_m^{(0)}}-\f{ e^2\lt(\bm E\cdot\bm A_{nm}\rt)\lt(\bm E\cdot\lt(\bm A_n-\bm A_m\rt)\rt)}{\lt(\ve_n^{(0)}-\ve_m^{(0)}\rt)^2}$$
73//!这样我们就能得到能带的各阶扰动
74//!$$\\begin{aligned}
75//!\ve_n^{(1)}&=-e\bm E\cdot\bm A_n\\\\
76//!\ve_n^{(2)}&=\f{e^2}{2}E_a E_b \sum_{m=\not n}\f{A_{nm}^a A_{mn}^b+A_{mn}^a A_{nm}^b}{\ve_n-\ve_m}=e^2 G_n^{ab}E_a E_b\\\\
77//!\ve_n^{(3)}&=-e^3E_a E_b E_c \lt( \sum_{m=\not n}\sum_{l=\not m,n}\f{A_{nl}^a A_{lm}^b A_{mn}^c}{(\ve_n-\ve_m)(\ve_n-\ve_l)}\rt)+e^3 E_a E_b E_c\lt( \sum_{m=\not n}\sum_{l=\not m,n}\f{A_{nm}^a A_{mn}^b (A_n^c-A_m^c)}{(\ve_n-\ve_m)^2}\rt)\\\\
78//!\\end{aligned}$$
79//!这里 $$G_n^{ab}=\sum_{m=\not n}\f{A_{nm}^a A_{mn}^b+A_{mn}^a A_{nm}^b}{\ve_n-\ve_m}=\sum_{m=\not n} 2\tx{Re}\f{v_{nm}^a v_{mn}^b}{(\ve_n-\ve_m)^3}$$
80//!到这里, 我们将能带的扰动得到了. 但是有一个问题, 就是 $\bm A$ 是一个规范变换, 所以并不是唯一的,
81//!同时, 对于带内的贡献 $\bm A_{n}=i\bra{\psi_{n\bm k}}\p_{\bm k}\ket{\psi_{n\bm k}}$, 没有好的求解方法,因为 $\bm A=-e\bra{\psi_n}\bm r\ket{\psi_n}$ 破坏了平移对称性. 但是我们总是可以选择一个规范.
82//!在这里我们选择 $-e\bm E\cdot\bm A_n=0$, 这个规范令 $\ve_n^{(1)}=0$. 这种规范在物理的意义上, 可以理解为贝利联络和电场的方向垂直. 对于贝利曲率的高阶项, 利用 $\bm A\to\bm A^\prime=A+\lt[\hat S,\bm A\rt]+\f{1}{2}\lt[\hat S,\lt[\hat S,\bm A\rt]\rt]\cdots$, 我们有
83//!$$\\begin{aligned}
84//!\lt(A_n^b\rt)^{(1)}&=-e\bm E_a G_n^{ab}\\\\
85//!\lt(A_n^c\rt)^{(2)}&=e^2 E_a E_b \lt( \sum_{m=\not n}\sum_{l=\not m,n}\f{A_{nl}^a A_{lm}^b A_{mn}^c}{(\ve_n-\ve_m)(\ve_n-\ve_l)}\rt)+e^2 E_a E_b\lt( \sum_{m=\not n}\f{A_{nm}^a A_{mn}^b (A_n^c-A_m^c)}{(\ve_n-\ve_m)^2}\rt)\\\\
86//!&=e^2 E_a E_b\lt(S_n^{abc}-F_n^{abc}\rt)
87//!\\end{aligned}$$
88//!
89//!这样利用贝利曲率公式 $\Og_n^{ab}=\p_a A_n^b -\p_b A_n^a$
90//!我们有 $$\\begin{aligned}
91//!\lt(\Og_n^{ab}\rt)^{(1)}&=-e E_c\lt(\p_a G_n^{bc}-\p_b G_n^{ac}\rt)\\\\
92//!\lt(\Og_n^{ab}\rt)^{(2)}&=e^2 E_{\ap}E_{\bt}\lt(\p_a S^{\ap\bt b}-\p_b S^{\ap\bt a}-\p_a F^{\ap\bt b}+\p_b F^{\ap\bt a}\rt)
93//!\\end{aligned}$$
94//!最终我们带入到电导率公式中, 有
95//!$$\begin{aligned}
96//!\sigma_{ab}=&-\f{e^2}{\hbar}\int_\tx{BZ} \f{\dd\bm k}{(2\pi)^3}\sum_n f_n\Og_n^{ab}+\f{e^2\tau}{\hbar^2}\sum_n \int_\tx{BZ}\f{\dd\bm k}{(2\pi)^3}\f{\p^2\ve_n}{\p k_a\p k_b}\\\\
97//!\sigma_{abc}=&-\f{e^3\tau^2}{\hbar^3}\sum_n\int_\tx{BZ}\f{\dd\bm k}{(2\pi)^3}\f{\p^3\ve_n}{\p k_a \p k_b \p k_c}
98//!+\f{e^3\tau}{\hbar^2}\sum_n \int_\tx{BZ}\f{\dd\bm k}{(2\pi)^3} \f{1}{2} f_n \lt(\p_a\Og_n^{bc}+\p_b\Og_n^{ac}\rt)\\\\
99//!&-\f{e^3}{\hbar}\sum_n\int_\tx{BZ}\f{\dd\bm k}{(2\pi)^3} f_n\lt(2\p_c G_n^{ab}-\f{1}{2}\lt(\p_a G_n^{bc}+\p_b G_n^{ac}\rt)\rt)
100//!\end{aligned}$$
101//!
102//! ## Berry connection 的化简
103//!
104//!为了实际的计算, 我们需要将 Berry connection 的形式修改一下, 我们首先按照微分的定理有 $$\p_{\bm k}\lt(H_{\bm k}\ket{\psi_{n\bm k}}\rt)=\lt(\p_{\bm k}H_{\bm k}+H_{\bm k}\p_{\bm k}\rt)\ket{\psi_{n\bm k}}$$
105//!然后我们又因为 $H_{\bm k}\ket{\psi_{n\bm k}}=\ve_{n\bm k}\ket{\psi_{n\bm k}}$, 所以 $$\p_{\bm k}\lt(H_{\bm k}\ket{\psi_{n\bm k}}\rt)=\p_{\bm k}\ve_{n\bm k}\ket{\psi_{n\bm k}}+\ve_{n\bm k}\p_{\bm k}\ket{\psi_{n\bm k}}$$
106//!所以我们有 $$\\begin{aligned}
107//!\p_{\bm k}H_{\bm k}\ket{\psi_{n\bm k}}+H_{\bm k}\p_{\bm k}\ket{\psi_{n\bm k}}=\p_{\bm k}\ve_{n\bm k}\ket{\psi_{n\bm k}}+\ve_{n\bm k}\p_{\bm k}\ket{\psi_{n\bm k}}
108//!\\end{aligned}$$
109//!显然我们将上式等号两边的左侧插入一个完备算符 $\sum_m \dyad{\psi_{m\bm k}}$ 有 $$\sum_m\lt[\bra{\psi_{m\bm k}}\p_{\bm k}H_{\bm k}\ket{\psi_{n\bm k}}+\lt(\ve_{m\bm k}-\ve_{n\bm k}\rt)\bra{\psi_{m\bm k}}\p_{\bm k}\ket{\psi_{n\bm k}}\rt]\ket{\psi_{m\bm k}}=\bra{\psi_{n\bm k}}\p_{\bm k}\ket{\psi_{n\bm k}}\ket{\psi_{n\bm k}} $$
110//!根据上面的式子, 我们很容易得到当 $m=\not n$ 时 $$\bra{\psi_{m\bm k}}\p_{\bm k}\ket{\psi_{n\bm k}}=\f{\bra{\psi_{m\bm k}}\p_{\bm k}\ket{\psi_{n\bm k}}}{\ve_{n\bm k}-\ve_{m\bm k}}$$
111//!也就是说, 我们能够最终得到 $$\bm A_{mn}=i\f{\bra{\psi_{m\bm k}}\p_{\bm k}\ket{\psi_{n\bm k}}}{\ve_{n\bm k}-\ve_{m\bm k}}$$
112
113use crate::error::{Result, TbError};
114use crate::kpoints::{gen_kmesh, gen_krange};
115use crate::math::*;
116use crate::phy_const::mu_B;
117use crate::solve_ham::solve;
118use crate::velocity::*;
119use crate::{Gauge, Model};
120use ndarray::linalg::kron;
121use ndarray::prelude::*;
122use ndarray::*;
123use ndarray_linalg::conjugate;
124use ndarray_linalg::*;
125use num_complex::Complex;
126use rayon::prelude::*;
127use std::f64::consts::PI;
128use std::ops::AddAssign;
129use std::ops::MulAssign;
130
131/**
132这个函数是用来做自适应积分算法的
133
134对于任意维度的积分 n, 我们的将区域刨分成 n+1面体的小块, 然后用线性插值来近似这个n+1的积分结果
135
136设被积函数为 $f(x_1,x_2,...,x_n)$, 存在 $n+1$ 个点 $(y_{01},y_{02},\cdots y_{0n})\cdots(y_{n1},y_{n2}\cdots y_{nn})$, 对应的值为 $z_0,z_1,...,z_n$
137
138这样我们就能得到这一块积分的近似值为 $$ \f{1}{(n+1)!}\times\sum_{i=0}^n z_i *\dd V.$$ 其中$\dd V$ 是正 $n+1$ 面体的体积.
139
140在这里, 对于一维体系, 线性插值积分等价于梯形积分. 在两个相邻的数据点 ($x_1$,$f_1$) 和 ($x_2,$f_2$), 其积分结果为$\Delta=\f{f_1+f_2}{2}*(x_2-x_2)$.
141
142对于二维系统, 用三角形进行近似, 对于任意的小三角形得到的积分结果为 $\Delta=S\sum_{i=1}^3 f_i/3!$
143
144对于三维系统, 线性插值的结果为 $\Delta=S\sum_{i=1}^4 f_i/4!$
145*/
146
147#[inline(always)]
148pub fn adapted_integrate_quick(
149    f0: &dyn Fn(&Array1<f64>) -> f64,
150    k_range: &Array2<f64>,
151    re_err: f64,
152    ab_err: f64,
153) -> f64 {
154    let dim = k_range.len_of(Axis(0));
155    match dim {
156        1 => {
157            //对于一维情况, 我们就是用梯形算法的 (a+b)*h/2, 这里假设的是函数的插值为线性插值.
158            let mut use_range = vec![(k_range.clone(), re_err, ab_err)];
159            let mut result = 0.0;
160            while let Some((k_range, re_err, ab_err)) = use_range.pop() {
161                let kvec_l: Array1<f64> = arr1(&[k_range[[0, 0]]]);
162                let kvec_r: Array1<f64> = arr1(&[k_range[[0, 1]]]);
163                let kvec_m: Array1<f64> = arr1(&[(k_range[[0, 1]] + k_range[[0, 0]]) / 2.0]);
164                let dk: f64 = k_range[[0, 1]] - k_range[[0, 0]];
165                let y_l: f64 = f0(&kvec_l);
166                let y_r: f64 = f0(&kvec_r);
167                let y_m: f64 = f0(&kvec_m);
168                let all: f64 = (y_l + y_r) * dk / 2.0;
169                let all_1 = (y_l + y_m) * dk / 4.0;
170                let all_2 = (y_r + y_m) * dk / 4.0;
171                let err = all_1 + all_2 - all;
172                let abs_err = if ab_err > all * re_err {
173                    ab_err
174                } else {
175                    all * re_err
176                };
177                if err < abs_err {
178                    result += all_1 + all_2;
179                } else {
180                    let k_range_l = arr2(&[[kvec_l[[0]], kvec_m[[0]]]]);
181                    let k_range_r = arr2(&[[kvec_m[[0]], kvec_r[[0]]]]);
182                    use_range.push((k_range_l, re_err, ab_err / 2.0));
183                    use_range.push((k_range_r, re_err, ab_err / 2.0));
184                }
185            }
186            return result;
187        }
188        2 => {
189            //对于二维, 我们依旧假设线性插值, 这样我们考虑的就是二维平面上的三角形上的任意一点的值是到其余三个点的距离的加权系数的平均值, 我们将四边形变成两个三角形来考虑.
190            let area_1: Array2<f64> = arr2(&[
191                [k_range.row(0)[0], k_range.row(1)[0]],
192                [k_range.row(0)[1], k_range.row(1)[0]],
193                [k_range.row(0)[0], k_range.row(1)[1]],
194            ]); //第一个三角形
195            let area_2: Array2<f64> = arr2(&[
196                [k_range.row(0)[1], k_range.row(1)[1]],
197                [k_range.row(0)[1], k_range.row(1)[0]],
198                [k_range.row(0)[0], k_range.row(1)[1]],
199            ]); //第二个三角形
200            #[inline(always)]
201            fn adapt_integrate_triangle(
202                f0: &dyn Fn(&Array1<f64>) -> f64,
203                kvec: &Array2<f64>,
204                re_err: f64,
205                ab_err: f64,
206                s1: f64,
207                s2: f64,
208                s3: f64,
209                S: f64,
210            ) -> f64 {
211                //这个函数是用来进行自适应算法的
212                let mut result = 0.0;
213                let mut use_kvec = vec![(kvec.clone(), re_err, ab_err, s1, s2, s3, S)];
214                while let Some((kvec, re_err, ab_err, s1, s2, s3, S)) = use_kvec.pop() {
215                    let kvec_m = kvec.mean_axis(Axis(0)).unwrap();
216                    let sm: f64 = f0(&kvec_m.to_owned());
217
218                    let mut new_kvec = kvec.to_owned();
219                    new_kvec.push_row(kvec_m.view());
220                    let kvec_1 = new_kvec.select(Axis(0), &[0, 1, 3]);
221                    let kvec_2 = new_kvec.select(Axis(0), &[0, 3, 2]);
222                    let kvec_3 = new_kvec.select(Axis(0), &[3, 1, 2]);
223                    let all: f64 = (s1 + s2 + s3) * S / 6.0;
224                    let all_new: f64 = all / 3.0 * 2.0 + sm * S / 6.0;
225                    let abs_err: f64 = if ab_err > all * re_err {
226                        ab_err
227                    } else {
228                        all * re_err
229                    };
230                    if (all_new - all).abs() > abs_err && S > 1e-8 {
231                        let S1 = S / 3.0;
232                        use_kvec.push((kvec_1, re_err, ab_err / 3.0, s1, s2, sm, S1));
233                        use_kvec.push((kvec_2, re_err, ab_err / 3.0, s1, sm, s3, S1));
234                        use_kvec.push((kvec_3, re_err, ab_err / 3.0, sm, s2, s3, S1));
235                    } else {
236                        result += all_new;
237                    }
238                }
239                result
240            }
241            let S = (k_range[[0, 1]] - k_range[[0, 0]]) * (k_range[[1, 1]] - k_range[[1, 0]]);
242            let s1 = f0(&arr1(&[k_range.row(0)[0], k_range.row(1)[0]]));
243            let s2 = f0(&arr1(&[k_range.row(0)[1], k_range.row(1)[0]]));
244            let s3 = f0(&arr1(&[k_range.row(0)[0], k_range.row(1)[1]]));
245            let s4 = f0(&arr1(&[k_range.row(0)[1], k_range.row(1)[1]]));
246            let all_1 = adapt_integrate_triangle(f0, &area_1, re_err, ab_err / 2.0, s1, s2, s3, S);
247            let all_2 = adapt_integrate_triangle(f0, &area_2, re_err, ab_err / 2.0, s4, s2, s3, S);
248            return all_1 + all_2;
249        }
250        3 => {
251            //对于三位情况, 需要用到四面体, 所以需要先将6面体变成6个四面体
252            #[inline(always)]
253            fn adapt_integrate_tetrahedron(
254                f0: &dyn Fn(&Array1<f64>) -> f64,
255                kvec: &Array2<f64>,
256                re_err: f64,
257                ab_err: f64,
258                s1: f64,
259                s2: f64,
260                s3: f64,
261                s4: f64,
262                S: f64,
263            ) -> f64 {
264                //这个函数是用来进行自适应算法的
265                let mut result = 0.0;
266                let mut use_kvec = vec![(kvec.clone(), re_err, ab_err, s1, s2, s3, s4, S)];
267                while let Some((kvec, re_err, ab_err, s1, s2, s3, s4, S)) = use_kvec.pop() {
268                    let kvec_m = kvec.mean_axis(Axis(0)).unwrap();
269                    let sm = f0(&kvec_m.to_owned());
270                    let mut new_kvec = kvec.to_owned();
271                    new_kvec.push_row(kvec_m.view());
272                    let kvec_1 = new_kvec.select(Axis(0), &[0, 1, 2, 4]);
273                    let kvec_2 = new_kvec.select(Axis(0), &[0, 1, 4, 3]);
274                    let kvec_3 = new_kvec.select(Axis(0), &[0, 4, 2, 3]);
275                    let kvec_4 = new_kvec.select(Axis(0), &[4, 1, 2, 3]);
276
277                    let all = (s1 + s2 + s3 + s4) * S / 24.0;
278                    let all_new = all * 0.75 + sm * S / 24.0;
279                    let S1 = S * 0.25;
280                    let abs_err = if ab_err > all * re_err {
281                        ab_err
282                    } else {
283                        all * re_err
284                    };
285                    if (all_new - all).abs() > abs_err && S > 1e-9 {
286                        use_kvec.push((kvec_1, re_err, ab_err * 0.25, s1, s2, s3, sm, S1));
287                        use_kvec.push((kvec_2, re_err, ab_err * 0.25, s1, s2, sm, s4, S1));
288                        use_kvec.push((kvec_3, re_err, ab_err * 0.25, s1, sm, s3, s4, S1));
289                        use_kvec.push((kvec_4, re_err, ab_err * 0.25, sm, s2, s3, s4, S1));
290                    } else {
291                        result += all_new;
292                    }
293                }
294                result
295            }
296
297            let k_points: Array2<f64> = arr2(&[
298                [k_range.row(0)[0], k_range.row(1)[0], k_range.row(2)[0]],
299                [k_range.row(0)[1], k_range.row(1)[0], k_range.row(2)[0]],
300                [k_range.row(0)[0], k_range.row(1)[1], k_range.row(2)[0]],
301                [k_range.row(0)[1], k_range.row(1)[1], k_range.row(2)[0]],
302                [k_range.row(0)[0], k_range.row(1)[0], k_range.row(2)[1]],
303                [k_range.row(0)[1], k_range.row(1)[0], k_range.row(2)[1]],
304                [k_range.row(0)[0], k_range.row(1)[1], k_range.row(2)[1]],
305                [k_range.row(0)[1], k_range.row(1)[1], k_range.row(2)[1]],
306            ]); //六面体的顶点
307
308            let area_1 = k_points.select(Axis(0), &[0, 1, 2, 5]);
309            let area_2 = k_points.select(Axis(0), &[0, 2, 4, 5]);
310            let area_3 = k_points.select(Axis(0), &[6, 2, 4, 5]);
311            let area_4 = k_points.select(Axis(0), &[1, 2, 3, 5]);
312            let area_5 = k_points.select(Axis(0), &[7, 2, 3, 5]);
313            let area_6 = k_points.select(Axis(0), &[7, 2, 6, 5]);
314            let s0 = f0(&k_points.row(0).to_owned());
315            let s1 = f0(&k_points.row(1).to_owned());
316            let s2 = f0(&k_points.row(2).to_owned());
317            let s3 = f0(&k_points.row(3).to_owned());
318            let s4 = f0(&k_points.row(4).to_owned());
319            let s5 = f0(&k_points.row(5).to_owned());
320            let s6 = f0(&k_points.row(6).to_owned());
321            let s7 = f0(&k_points.row(7).to_owned());
322            let V = (k_range[[0, 1]] - k_range[[0, 0]])
323                * (k_range[[1, 1]] - k_range[[1, 0]])
324                * (k_range[[2, 1]] - k_range[[2, 0]]);
325            let all_1 =
326                adapt_integrate_tetrahedron(f0, &area_1, re_err, ab_err / 6.0, s0, s1, s2, s5, V);
327            let all_2 =
328                adapt_integrate_tetrahedron(f0, &area_2, re_err, ab_err / 6.0, s0, s2, s4, s5, V);
329            let all_3 =
330                adapt_integrate_tetrahedron(f0, &area_3, re_err, ab_err / 6.0, s6, s2, s4, s5, V);
331            let all_4 =
332                adapt_integrate_tetrahedron(f0, &area_4, re_err, ab_err / 6.0, s1, s2, s3, s5, V);
333            let all_5 =
334                adapt_integrate_tetrahedron(f0, &area_5, re_err, ab_err / 6.0, s7, s2, s3, s5, V);
335            let all_6 =
336                adapt_integrate_tetrahedron(f0, &area_5, re_err, ab_err / 6.0, s7, s2, s6, s5, V);
337            return all_1 + all_2 + all_3 + all_4 + all_5 + all_6;
338        }
339        _ => {
340            panic!(
341                "wrong, the row_dim if k_range must be 1,2 or 3, but you's give {}",
342                dim
343            );
344        }
345    }
346}
347
348pub trait BerryCurvature: Velocity {
349    //! 这个trait是用来提供Berry Curvature 的
350    //!
351    /// 产生第n条能带的贝利曲率
352    fn berry_curvature_n_onek<S: Data<Elem = f64>>(
353        &self,
354        k_vec: &ArrayBase<S, Ix1>,
355        dir_1: &Array1<f64>,
356        dir_2: &Array1<f64>,
357        spin: usize,
358        eta: f64,
359    ) -> (Array1<f64>, Array1<f64>);
360
361    ///给定一个 k 点, 指定 dir_1=$\alpha$, dir_2=$\beta$, T 代表温度, og= $\og$,
362    ///mu=$\mu$ 为费米能级, spin=0,1,2,3 为$\sg_0,\sg_x,\sg_y,\sg_z$,
363    ///当体系不存在自旋的时候无论如何输入spin都默认 spin=0
364    ///eta=$\eta$ 是一个小量
365    /// 这个函数返回的是
366    /// $$ \sum_n f_n\Omega_{n,\ap\bt}^\gm(\bm k)=\sum_n \f{1}{e^{(\ve_{n\bm k}-\mu)/T/k_B}+1} \sum_{m=\not n}\f{J_{\ap,nm}^\gm v_{\bt,mn}}{(\ve_{n\bm k}-\ve_{m\bm k})^2-(\og+i\eta)^2}$$
367    /// 其中 $J_\ap^\gm=\\{s_\gm,v_\ap\\}$
368    fn berry_curvature_onek<S: Data<Elem = f64>>(
369        &self,
370        k_vec: &ArrayBase<S, Ix1>,
371        dir_1: &Array1<f64>,
372        dir_2: &Array1<f64>,
373        mu: f64,
374        T: f64,
375        spin: usize,
376        eta: f64,
377    ) -> f64;
378
379    ///这个是用来并行计算大量k点的贝利曲率
380    ///这个可以用来画能带上的贝利曲率, 或者画一个贝利曲率的热图
381    #[allow(non_snake_case)]
382    fn berry_curvature<S: Data<Elem = f64>>(
383        &self,
384        k_vec: &ArrayBase<S, Ix2>,
385        dir_1: &Array1<f64>,
386        dir_2: &Array1<f64>,
387        mu: f64,
388        T: f64,
389        spin: usize,
390        eta: f64,
391    ) -> Array1<f64>;
392}
393
394impl BerryCurvature for Model {
395    #[allow(non_snake_case)]
396    #[inline(always)]
397    fn berry_curvature_n_onek<S: Data<Elem = f64>>(
398        &self,
399        k_vec: &ArrayBase<S, Ix1>,
400        dir_1: &Array1<f64>,
401        dir_2: &Array1<f64>,
402        spin: usize,
403        eta: f64,
404    ) -> (Array1<f64>, Array1<f64>) {
405        let li: Complex<f64> = 1.0 * Complex::i();
406        //let (band, evec) = self.solve_onek(&k_vec);
407        let (mut v, hamk): (Array3<Complex<f64>>, Array2<Complex<f64>>) =
408            self.gen_v(&k_vec, Gauge::Atom); //这是速度算符
409        let (band, evec) = if let Ok((eigvals, eigvecs)) = hamk.eigh(UPLO::Lower) {
410            (eigvals, eigvecs)
411        } else {
412            todo!()
413        };
414        let mut J = v.view();
415        let (J, v) = if self.spin {
416            let J = match spin {
417                0 => {
418                    let J = J
419                        .outer_iter()
420                        .zip(dir_1.iter())
421                        .fold(Array2::zeros((self.nsta(), self.nsta())), |acc, (x, d)| {
422                            acc + &x * (*d + 0.0 * li)
423                        });
424                    J
425                }
426                1 => {
427                    let pauli = arr2(&[
428                        [0.0 + 0.0 * li, 1.0 + 0.0 * li],
429                        [1.0 + 0.0 * li, 0.0 + 0.0 * li],
430                    ]) / 2.0;
431                    let mut X: Array2<Complex<f64>> = Array2::eye(self.nsta());
432                    X = kron(&pauli, &Array2::eye(self.norb()));
433                    let J = J
434                        .outer_iter()
435                        .zip(dir_1.iter())
436                        .fold(Array2::zeros((self.nsta(), self.nsta())), |acc, (x, d)| {
437                            acc + &anti_comm(&X, &x) * (*d * 0.5 + 0.0 * li)
438                        });
439                    J
440                }
441                2 => {
442                    let pauli = arr2(&[
443                        [0.0 + 0.0 * li, 0.0 - 1.0 * li],
444                        [0.0 + 1.0 * li, 0.0 + 0.0 * li],
445                    ]) / 2.0;
446                    let mut X: Array2<Complex<f64>> = Array2::eye(self.nsta());
447                    X = kron(&pauli, &Array2::eye(self.norb()));
448                    let J = J
449                        .outer_iter()
450                        .zip(dir_1.iter())
451                        .fold(Array2::zeros((self.nsta(), self.nsta())), |acc, (x, d)| {
452                            acc + &anti_comm(&X, &x) * (*d * 0.5 + 0.0 * li)
453                        });
454                    J
455                }
456                3 => {
457                    let pauli = arr2(&[
458                        [1.0 + 0.0 * li, 0.0 + 0.0 * li],
459                        [0.0 + 0.0 * li, -1.0 + 0.0 * li],
460                    ]) / 2.0;
461                    let mut X: Array2<Complex<f64>> = Array2::eye(self.nsta());
462                    X = kron(&pauli, &Array2::eye(self.norb()));
463                    let J = J
464                        .outer_iter()
465                        .zip(dir_1.iter())
466                        .fold(Array2::zeros((self.nsta(), self.nsta())), |acc, (x, d)| {
467                            acc + &anti_comm(&X, &x) * (*d * 0.5 + 0.0 * li)
468                        });
469                    J
470                }
471                _ => panic!("Wrong, spin should be 0, 1, 2, 3, but you input {}", spin),
472            };
473            let v = v
474                .outer_iter()
475                .zip(dir_2.iter())
476                .fold(Array2::zeros((self.nsta(), self.nsta())), |acc, (x, d)| {
477                    acc + &x * (*d + 0.0 * li)
478                });
479            (J, v)
480        } else {
481            if spin != 0 {
482                println!("Warning, the model haven't got spin, so the spin input will be ignord");
483            }
484
485            let J = J
486                .outer_iter()
487                .zip(dir_1.iter())
488                .fold(Array2::zeros((self.nsta(), self.nsta())), |acc, (x, d)| {
489                    acc + &x * (*d + 0.0 * li)
490                });
491            let v = v
492                .outer_iter()
493                .zip(dir_2.iter())
494                .fold(Array2::zeros((self.nsta(), self.nsta())), |acc, (x, d)| {
495                    acc + &x * (*d + 0.0 * li)
496                });
497            (J, v)
498        };
499
500        let evec_conj = evec.t();
501        let evec = evec.mapv(|x| x.conj());
502        let A1 = J.dot(&evec);
503        let A1 = &evec_conj.dot(&A1);
504        let A2 = v.dot(&evec);
505        let A2 = evec_conj.dot(&A2);
506        let A2 = A2.reversed_axes();
507        let AA = A1 * A2;
508        let Complex { re, im } = AA.view().split_complex();
509        let im = im.mapv(|x| -2.0 * x);
510        assert_eq!(
511            band.len(),
512            self.nsta(),
513            "this is strange for band's length is not equal to self.nsta()"
514        );
515        let mut UU = Array2::<f64>::zeros((self.nsta(), self.nsta()));
516        for i in 0..self.nsta() {
517            for j in 0..self.nsta() {
518                let a = band[[i]] - band[[j]];
519                //这里用η进行展宽
520                UU[[i, j]] = 1.0 / (a.powi(2) + eta.powi(2));
521                /*
522                if a.abs() < 1e-8 {
523                    UU[[i, j]] = 0.0;
524                } else {
525                    UU[[i, j]] = 1.0 / (a.powi(2)+eta.powi(2));
526                }
527                */
528            }
529        }
530        let omega_n = im
531            .outer_iter()
532            .zip(UU.outer_iter())
533            .map(|(a, b)| a.dot(&b))
534            .collect();
535        let omega_n = Array1::from_vec(omega_n);
536        (omega_n, band)
537    }
538
539    #[allow(non_snake_case)]
540    fn berry_curvature_onek<S: Data<Elem = f64>>(
541        &self,
542        k_vec: &ArrayBase<S, Ix1>,
543        dir_1: &Array1<f64>,
544        dir_2: &Array1<f64>,
545        mu: f64,
546        T: f64,
547        spin: usize,
548        eta: f64,
549    ) -> f64 {
550        //!给定一个 k 点, 指定 dir_1=$\alpha$, dir_2=$\beta$, T 代表温度, og= $\og$,
551        //!mu=$\mu$ 为费米能级, spin=0,1,2,3 为$\sg_0,\sg_x,\sg_y,\sg_z$,
552        //!当体系不存在自旋的时候无论如何输入spin都默认 spin=0
553        //!eta=$\eta$ 是一个小量
554        //! 这个函数返回的是
555        //! $$ \sum_n f_n\Omega_{n,\ap\bt}^\gm(\bm k)=\sum_n \f{1}{e^{(\ve_{n\bm k}-\mu)/T/k_B}+1} \sum_{m=\not n}\f{J_{\ap,nm}^\gm v_{\bt,mn}}{(\ve_{n\bm k}-\ve_{m\bm k})^2-(\og+i\eta)^2}$$
556        //! 其中 $J_\ap^\gm=\\{s_\gm,v_\ap\\}$
557        let (omega_n, band) = self.berry_curvature_n_onek(&k_vec, &dir_1, &dir_2, spin, eta);
558        let mut omega: f64 = 0.0;
559        let fermi_dirac = if T == 0.0 {
560            band.mapv(|x| if x > mu { 0.0 } else { 1.0 })
561        } else {
562            let beta = 1.0 / T / 8.617e-5;
563            band.mapv(|x| ((beta * (x - mu)).exp() + 1.0).recip())
564        };
565        let omega = omega_n.dot(&fermi_dirac);
566        omega
567    }
568    #[allow(non_snake_case)]
569    fn berry_curvature<S: Data<Elem = f64>>(
570        &self,
571        k_vec: &ArrayBase<S, Ix2>,
572        dir_1: &Array1<f64>,
573        dir_2: &Array1<f64>,
574        mu: f64,
575        T: f64,
576        spin: usize,
577        eta: f64,
578    ) -> Array1<f64> {
579        //!这个是用来并行计算大量k点的贝利曲率
580        //!这个可以用来画能带上的贝利曲率, 或者画一个贝利曲率的热图
581        if dir_1.len() != self.dim_r() || dir_2.len() != self.dim_r() {
582            panic!(
583                "Wrong, the dir_1 or dir_2 you input has wrong length, it must equal to dim_r={}, but you input {},{}",
584                self.dim_r(),
585                dir_1.len(),
586                dir_2.len()
587            )
588        }
589        let nk = k_vec.len_of(Axis(0));
590        let omega: Vec<f64> = k_vec
591            .axis_iter(Axis(0))
592            .into_par_iter()
593            .map(|x| {
594                let omega_one =
595                    self.berry_curvature_onek(&x.to_owned(), &dir_1, &dir_2, mu, T, spin, eta);
596                omega_one
597            })
598            .collect();
599        let omega = arr1(&omega);
600        omega
601    }
602}
603
604#[allow(non_snake_case)]
605impl Model {
606    //! 这个模块是用来提供电导率张量的, 包括自旋霍尔电导率和霍尔电导率, 以及非线性霍尔电导率.
607    //!
608    //!
609    //!
610    ///这个是计算某个费米能级, 某个温度下的Hall conductivity 的, 输出的单位为 $e^2/\hbar/\AA$.
611    #[allow(non_snake_case)]
612    pub fn Hall_conductivity(
613        &self,
614        k_mesh: &Array1<usize>,
615        dir_1: &Array1<f64>,
616        dir_2: &Array1<f64>,
617        mu: f64,
618        T: f64,
619        spin: usize,
620        eta: f64,
621    ) -> Result<f64> {
622        //!这个是用来计算霍尔电导的.
623        //!这里采用的是均匀撒点的方法, 利用 berry_curvature, 我们有
624        //!$$\sg_{\ap\bt}^\gm=\f{1}{N(2\pi)^r V}\sum_{\bm k} \Og_{\ap\bt}(\bm k),$$ 其中 $N$ 是 k 点数目,
625        let kvec: Array2<f64> = gen_kmesh(&k_mesh)?;
626        let nk: usize = kvec.len_of(Axis(0));
627        let omega = self.berry_curvature(&kvec, &dir_1, &dir_2, mu, T, spin, eta);
628        //目前求积分的方法上, 还是直接求和最有用, 其他的典型积分方法, 如gauss 法等,
629        //都因为存在间断点而效率不高.
630        //对于非零温的, 使用梯形法应该效果能好一些.
631        let conductivity: f64 = omega.sum() / (nk as f64) / self.lat.det().unwrap();
632        Ok(conductivity)
633    }
634    #[allow(non_snake_case)]
635    ///这个是采用自适应积分算法来计算霍尔电导的, 一般来说, 我们建议 re_err 设置为 1, 而 ab_err 设置为 0.01
636    pub fn Hall_conductivity_adapted(
637        &self,
638        k_mesh: &Array1<usize>,
639        dir_1: &Array1<f64>,
640        dir_2: &Array1<f64>,
641        mu: f64,
642        T: f64,
643        spin: usize,
644        eta: f64,
645        re_err: f64,
646        ab_err: f64,
647    ) -> Result<f64> {
648        let mut k_range = gen_krange(k_mesh)?; //将要计算的区域分成小块
649        let n_range = k_range.len_of(Axis(0));
650        let ab_err = ab_err / (n_range as f64);
651        let use_fn =
652            |k0: &Array1<f64>| self.berry_curvature_onek(k0, &dir_1, &dir_2, mu, T, spin, eta);
653        let inte = |k_range| adapted_integrate_quick(&use_fn, &k_range, re_err, ab_err);
654        let omega: Vec<f64> = k_range
655            .axis_iter(Axis(0))
656            .into_par_iter()
657            .map(|x| inte(x.to_owned()))
658            .collect();
659        let omega: Array1<f64> = arr1(&omega);
660        let conductivity: f64 = omega.sum() / self.lat.det().unwrap();
661        Ok(conductivity)
662    }
663    ///用来计算多个 $\mu$ 值的, 这个函数是先求出 $\Omega_n$, 然后再分别用不同的费米能级来求和, 这样速度更快, 因为避免了重复求解 $\Omega_n$, 但是相对来说更耗内存, 而且不能做到自适应积分算法.
664    pub fn Hall_conductivity_mu(
665        &self,
666        k_mesh: &Array1<usize>,
667        dir_1: &Array1<f64>,
668        dir_2: &Array1<f64>,
669        mu: &Array1<f64>,
670        T: f64,
671        spin: usize,
672        eta: f64,
673    ) -> Result<Array1<f64>> {
674        let kvec: Array2<f64> = gen_kmesh(&k_mesh)?;
675        let nk: usize = kvec.len_of(Axis(0));
676        let (omega_n, band): (Vec<_>, Vec<_>) = kvec
677            .axis_iter(Axis(0))
678            .into_par_iter()
679            .map(|x| {
680                let (omega_n, band) =
681                    self.berry_curvature_n_onek(&x.to_owned(), &dir_1, &dir_2, spin, eta);
682                (omega_n, band)
683            })
684            .collect();
685        let omega_n = Array2::<f64>::from_shape_vec(
686            (nk, self.nsta()),
687            omega_n.into_iter().flatten().collect(),
688        )
689        .unwrap();
690        let band =
691            Array2::<f64>::from_shape_vec((nk, self.nsta()), band.into_iter().flatten().collect())
692                .unwrap();
693        let n_mu: usize = mu.len();
694        let conductivity = if T == 0.0 {
695            let conductivity_new: Vec<f64> = mu
696                .into_par_iter()
697                .map(|x| {
698                    let mut omega = Array1::<f64>::zeros(nk);
699                    for k in 0..nk {
700                        for i in 0..self.nsta() {
701                            omega[[k]] += if band[[k, i]] > *x {
702                                0.0
703                            } else {
704                                omega_n[[k, i]]
705                            };
706                        }
707                    }
708                    omega.sum() / self.lat.det().unwrap() / (nk as f64)
709                })
710                .collect();
711            Array1::<f64>::from_vec(conductivity_new)
712        } else {
713            let beta = 1.0 / (T * 8.617e-5);
714            let conductivity_new: Vec<f64> = mu
715                .into_par_iter()
716                .map(|x| {
717                    let fermi_dirac = band.mapv(|x0| 1.0 / ((beta * (x0 - x)).exp() + 1.0));
718                    let omega: Vec<f64> = omega_n
719                        .axis_iter(Axis(0))
720                        .zip(fermi_dirac.axis_iter(Axis(0)))
721                        .map(|(a, b)| (&a * &b).sum())
722                        .collect();
723                    let omega: Array1<f64> = arr1(&omega);
724                    omega.sum() / self.lat.det().unwrap() / (nk as f64)
725                })
726                .collect();
727            Array1::<f64>::from_vec(conductivity_new)
728        };
729        Ok(conductivity)
730    }
731
732    pub fn berry_curvature_dipole_n_onek(
733        &self,
734        k_vec: &Array1<f64>,
735        dir_1: &Array1<f64>,
736        dir_2: &Array1<f64>,
737        dir_3: &Array1<f64>,
738        og: f64,
739        spin: usize,
740        eta: f64,
741    ) -> (Array1<f64>, Array1<f64>) {
742        //! 这个是用来计算 $$\pdv{\ve_{n\bm k}}{k_\gm}\Og_{n,\ap\bt}$$
743        //!
744        //!这里需要注意的一点是, 一般来说对于 $\p_\ap\ve_{\bm k}$, 需要用差分法来求解, 我这里提供了一个算法.
745        //!$$ \ve_{\bm k}=U^\dag H_{\bm k} U\Rightarrow \pdv{\ve_{\bm k}}{\bm k}=U^\dag\pdv{H_{\bm k}}{\bm k}U+\pdv{U^\dag}{\bm k} H_{\bm k}U+U^\dag H_{\bm k}\pdv{U}{\bm k}$$
746        //!因为 $U^\dag U=1\Rightarrow \p_{\bm k}U^\dag U=-U^\dag\p_{\bm k}U$, $\p_{\bm k}H_{\bm k}=v_{\bm k}$我们有
747        //!$$\pdv{\ve_{\bm k}}{\bm k}=v_{\bm k}+\lt[\ve_{\bm k},U^\dag\p_{\bm k}U\rt]$$
748        //!而这里面唯一比较难求的项是 $D_{\bm k}=U^\dag\p_{\bm k}U$. 按照 vanderbilt 2008 年的论文中的公式, 用微扰论有
749        //!$$D_{mn,\bm k}=\left\\{\\begin{aligned}\f{v_{mn,\bm k}}{\ve_n-\ve_m} \quad &\text{if}\\ m\\ =\not n\\\ 0 \quad \quad &\text{if}\\ m\\ = n\\end{aligned}\right\.$$
750        //!我们观察到第二项对对角部分没有贡献, 所以我们可以直接设置为
751        //!$$\pdv{\ve_{\bm k}}{\bm k}=\text{diag}\lt(v_{\bm k}\rt)$$
752        //我们首先求解 omega_n 和 U^\dag j
753
754        let li: Complex<f64> = 1.0 * Complex::i();
755        //let (band, evec) = self.solve_onek(&k_vec);
756        let (mut v, hamk): (Array3<Complex<f64>>, Array2<Complex<f64>>) =
757            self.gen_v(&k_vec, Gauge::Atom); //这是速度算符
758        let mut J: Array3<Complex<f64>> = v.clone();
759        let mut v0 = Array2::<Complex<f64>>::zeros((self.nsta(), self.nsta())); //这个是速度项, 对应的dir_3 的速度
760        for r in 0..self.dim_r() {
761            v0 = v0 + v.slice(s![r, .., ..]).to_owned() * dir_3[[r]];
762        }
763        if self.spin {
764            let mut X: Array2<Complex<f64>> = Array2::eye(self.nsta());
765            let pauli: Array2<Complex<f64>> = match spin {
766                0 => arr2(&[
767                    [1.0 + 0.0 * li, 0.0 + 0.0 * li],
768                    [0.0 + 0.0 * li, 1.0 + 0.0 * li],
769                ]),
770                1 => {
771                    arr2(&[
772                        [0.0 + 0.0 * li, 1.0 + 0.0 * li],
773                        [1.0 + 0.0 * li, 0.0 + 0.0 * li],
774                    ]) / 2.0
775                }
776                2 => {
777                    arr2(&[
778                        [0.0 + 0.0 * li, 0.0 - 1.0 * li],
779                        [0.0 + 1.0 * li, 0.0 + 0.0 * li],
780                    ]) / 2.0
781                }
782                3 => {
783                    arr2(&[
784                        [1.0 + 0.0 * li, 0.0 + 0.0 * li],
785                        [0.0 + 0.0 * li, -1.0 + 0.0 * li],
786                    ]) / 2.0
787                }
788                _ => panic!("Wrong, spin should be 0, 1, 2, 3, but you input {}", spin),
789            };
790            X = kron(&pauli, &Array2::eye(self.norb()));
791            for i in 0..self.dim_r() {
792                let j = J.slice(s![i, .., ..]).to_owned();
793                let j = anti_comm(&X, &j) / 2.0; //这里做反对易
794                J.slice_mut(s![i, .., ..]).assign(&(j * dir_1[[i]]));
795                v.slice_mut(s![i, .., ..])
796                    .mul_assign(Complex::new(dir_2[[i]], 0.0));
797            }
798        } else {
799            if spin != 0 {
800                println!("Warning, the model haven't got spin, so the spin input will be ignord");
801            }
802            for i in 0..self.dim_r() {
803                J.slice_mut(s![i, .., ..])
804                    .mul_assign(Complex::new(dir_1[[i]], 0.0));
805                v.slice_mut(s![i, .., ..])
806                    .mul_assign(Complex::new(dir_2[[i]], 0.0));
807            }
808        };
809
810        let J: Array2<Complex<f64>> = J.sum_axis(Axis(0));
811        let v: Array2<Complex<f64>> = v.sum_axis(Axis(0));
812
813        let (band, evec) = if let Ok((eigvals, eigvecs)) = hamk.eigh(UPLO::Lower) {
814            (eigvals, eigvecs)
815        } else {
816            todo!()
817        };
818        let evec_conj = evec.t();
819        let evec = evec.mapv(|x| x.conj());
820
821        let v0 = v0.dot(&evec.t());
822        let v0 = &evec_conj.dot(&v0);
823        let partial_ve = v0.diag().map(|x| x.re);
824        let A1 = J.dot(&evec.t());
825        let A1 = &evec_conj.dot(&A1);
826        let A2 = v.dot(&evec.t());
827        let A2 = &evec_conj.dot(&A2);
828        let mut U0 = Array2::<Complex<f64>>::zeros((self.nsta(), self.nsta()));
829        for i in 0..self.nsta() {
830            for j in 0..self.nsta() {
831                if i != j {
832                    U0[[i, j]] = 1.0 / ((band[[i]] - band[[j]]).powi(2) - (og + li * eta).powi(2));
833                } else {
834                    U0[[i, j]] = Complex::new(0.0, 0.0);
835                }
836            }
837        }
838        //let omega_n:Array1::<f64>=(-Complex::new(2.0,0.0)*(A1*U0).dot(&A2)).diag().map(|x| x.im).to_owned();
839        let mut omega_n = Array1::<f64>::zeros(self.nsta());
840        let A1 = A1 * U0;
841        for i in 0..self.nsta() {
842            omega_n[[i]] = -2.0 * A1.slice(s![i, ..]).dot(&A2.slice(s![.., i])).im;
843        }
844
845        //let (omega_n,band)=self.berry_curvature_n_onek(&k_vec,&dir_1,&dir_2,og,spin,eta);
846        let omega_n: Array1<f64> = omega_n * partial_ve;
847        (omega_n, band) //最后得到的 D
848    }
849    pub fn berry_curvature_dipole_n(
850        &self,
851        k_vec: &Array2<f64>,
852        dir_1: &Array1<f64>,
853        dir_2: &Array1<f64>,
854        dir_3: &Array1<f64>,
855        og: f64,
856        spin: usize,
857        eta: f64,
858    ) -> (Array2<f64>, Array2<f64>) {
859        //这个是在 onek的基础上进行并行计算得到一系列k点的berry curvature dipole
860        //!This function performs parallel computation based on the onek function to obtain a series of Berry curvature dipoles at different k-points.
861        //!这个方法用的是对费米分布的修正, 因为高阶的dipole 修正导致的非线性霍尔电导为 $$\sg_{\ap\bt\gm}=\tau\int\dd\bm k\sum_n\p_\gm\ve_{n\bm k}\Og_{n,\ap\bt}\lt\.\pdv{f_{\bm k}}{\ve}\rt\rvert_{E=\ve_{n\bm k}}.$$ 所以我们这里输出的是
862        //!$$\p_\gm\ve_{n\bm k}\Og_{n,\ap\bt}.$$
863        if dir_1.len() != self.dim_r() || dir_2.len() != self.dim_r() || dir_3.len() != self.dim_r()
864        {
865            panic!(
866                "Wrong, the dir_1 or dir_2 you input has wrong length, it must equal to dim_r={}, but you input {},{}",
867                self.dim_r(),
868                dir_1.len(),
869                dir_2.len()
870            )
871        }
872        let nk = k_vec.len_of(Axis(0));
873        let (omega, band): (Vec<_>, Vec<_>) = k_vec
874            .axis_iter(Axis(0))
875            .into_par_iter()
876            .map(|x| {
877                let (omega_one, band) = self.berry_curvature_dipole_n_onek(
878                    &x.to_owned(),
879                    &dir_1,
880                    &dir_2,
881                    &dir_3,
882                    og,
883                    spin,
884                    eta,
885                );
886                (omega_one, band)
887            })
888            .collect();
889        let omega =
890            Array2::<f64>::from_shape_vec((nk, self.nsta()), omega.into_iter().flatten().collect())
891                .unwrap();
892        let band =
893            Array2::<f64>::from_shape_vec((nk, self.nsta()), band.into_iter().flatten().collect())
894                .unwrap();
895        (omega, band)
896    }
897    pub fn Nonlinear_Hall_conductivity_Extrinsic(
898        &self,
899        k_mesh: &Array1<usize>,
900        dir_1: &Array1<f64>,
901        dir_2: &Array1<f64>,
902        dir_3: &Array1<f64>,
903        mu: &Array1<f64>,
904        T: f64,
905        og: f64,
906        spin: usize,
907        eta: f64,
908    ) -> Result<Array1<f64>> {
909        //这个是用 berry curvature dipole 对整个布里渊去做积分得到非线性霍尔电导, 是extrinsic 的
910        //!This function calculates the extrinsic nonlinear Hall conductivity by integrating the Berry curvature dipole over the entire Brillouin zone. The Berry curvature dipole is first computed at a series of k-points using parallel computation based on the onek function.
911
912        //! 我们基于 berry_curvature_n_dipole 来并行得到所有 k 点的 $\p_\gm\ve_{n\bm k}\Og_{n,\ap\bt}$,
913        //! 但是我们最后的公式为
914        //! $$\\mathcal D_{\ap\bt\gm}=\int \dd\bm k \sum_n\lt(-\pdv{f_{n}}{\ve}\rt)\p_\gm\ve_{n\bm k}\Og_{n,\ap\bt}$$
915        //! 然而,
916        //! $$-\pdv{f_{n}}{\ve}=\beta\f{e^{beta(\ve_n-\mu)}}{(e^{beta(\ve_n-\mu)}+1)^2}=\beta f_n(1-f_n)$$
917        //! 对于 T=0 的情况, 我们将采用四面体积分来替代, 这个需要很高的k点密度, 不建议使用
918        //! 对于 T!=0 的情况, 我们会采用类似 Dos 的方法来计算
919
920        if dir_1.len() != self.dim_r() || dir_2.len() != self.dim_r() || dir_3.len() != self.dim_r()
921        {
922            panic!(
923                "Wrong, the dir_1 or dir_2 you input has wrong length, it must equal to dim_r={}, but you input {},{}",
924                self.dim_r(),
925                dir_1.len(),
926                dir_2.len()
927            )
928        }
929        let kvec: Array2<f64> = gen_kmesh(&k_mesh)?;
930        let nk: usize = kvec.len_of(Axis(0));
931        //为了节省内存, 本来是可以直接算完求和, 但是为了方便, 我是先存下来再算, 让程序结构更合理
932        let (omega, band) =
933            self.berry_curvature_dipole_n(&kvec, &dir_1, &dir_2, &dir_3, og, spin, eta);
934        let omega = omega.into_raw_vec();
935        let band = band.into_raw_vec();
936        let n_e = mu.len();
937        let mut conductivity = Array1::<f64>::zeros(n_e);
938        if T != 0.0 {
939            let beta = 1.0 / T / (8.617e-5);
940            let use_iter = band.iter().zip(omega.iter()).par_bridge();
941            conductivity = use_iter
942                .fold(
943                    || Array1::<f64>::zeros(n_e),
944                    |acc, (energy, omega0)| {
945                        let f = 1.0 / (beta * (mu - *energy)).mapv(|x| x.exp() + 1.0);
946                        acc + &f * (1.0 - &f) * beta * *omega0
947                    },
948                )
949                .reduce(|| Array1::<f64>::zeros(n_e), |acc, x| acc + x);
950            conductivity = conductivity.clone() / (nk as f64) / self.lat.det().unwrap();
951        } else {
952            //采用四面体积分法, 或者对于二维体系, 采用三角形积分法
953            //积分的思路是, 通过将一个六面体变成5个四面体, 然后用线性插值的方法, 得到费米面,
954            //以及费米面上的数, 最后, 通过积分算出来结果
955            panic!("When T=0, the algorithm have not been writed, please wait for next version");
956        }
957        Ok(conductivity)
958    }
959
960    pub fn berry_connection_dipole_onek(
961        &self,
962        k_vec: &Array1<f64>,
963        dir_1: &Array1<f64>,
964        dir_2: &Array1<f64>,
965        dir_3: &Array1<f64>,
966        spin: usize,
967    ) -> (Array1<f64>, Array1<f64>, Option<Array1<f64>>) {
968        //!这个是根据 Nonlinear_Hall_conductivity_intrinsic 的注释, 当不存在自旋的时候提供
969        //!$$v_\ap G_{\bt\gm}-v_\bt G_{\ap\gm}$$
970        //!其中 $$ G_{ij}=-2\text{Re}\sum_{m=\not n}\f{v_{i,nm}v_{j,mn}}{\lt(\ve_n-\ve_m\rt)^3} $$
971        //!如果存在自旋, 即spin不等于0, 则还存在 $\p_{h_i} G_{jk}$ 项, 具体请看下面的非线性霍尔部分
972        //!我们这里暂时不考虑磁场, 只考虑电场
973        let (mut v, hamk): (Array3<Complex<f64>>, Array2<Complex<f64>>) =
974            self.gen_v(&k_vec, Gauge::Atom); //这是速度算符
975        let mut J = v.clone();
976        //let (band, evec) = self.solve_onek(&k_vec); //能带和本征值
977        let (band, evec) = if let Ok((eigvals, eigvecs)) = hamk.eigh(UPLO::Lower) {
978            (eigvals, eigvecs)
979        } else {
980            todo!()
981        };
982        let evec_conj = evec.t();
983        let evec = evec.mapv(|x| x.conj());
984        for i in 0..self.dim_r() {
985            let v_s = v.slice(s![i, .., ..]).to_owned();
986            let v_s = evec_conj
987                .clone()
988                .dot(&(v_s.dot(&evec.clone().reversed_axes()))); //变换到本征态基函数
989            v.slice_mut(s![i, .., ..]).assign(&v_s); //将 v 变换到以本征态为基底
990        }
991        //现在速度算符已经是以本征态为基函数
992        let mut v_1 = Array2::<Complex<f64>>::zeros((self.nsta(), self.nsta())); //三个方向的速度算符
993        let mut v_2 = Array2::<Complex<f64>>::zeros((self.nsta(), self.nsta()));
994        let mut v_3 = Array2::<Complex<f64>>::zeros((self.nsta(), self.nsta()));
995        for i in 0..self.dim_r() {
996            v_1 = v_1.clone() + v.slice(s![i, .., ..]).to_owned() * Complex::new(dir_1[[i]], 0.0);
997            v_2 = v_2.clone() + v.slice(s![i, .., ..]).to_owned() * Complex::new(dir_2[[i]], 0.0);
998            v_3 = v_3.clone() + v.slice(s![i, .., ..]).to_owned() * Complex::new(dir_3[[i]], 0.0);
999        }
1000        //三个方向的速度算符都得到了
1001        let mut U0 = Array2::<f64>::zeros((self.nsta(), self.nsta()));
1002        for i in 0..self.nsta() {
1003            for j in 0..self.nsta() {
1004                if (band[[i]] - band[[j]]).abs() < 1e-5 {
1005                    U0[[i, j]] = 0.0;
1006                } else {
1007                    U0[[i, j]] = 1.0 / (band[[i]] - band[[j]]);
1008                }
1009            }
1010        }
1011        //这样U0[[i,j]]=1/(E_i-E_j), 这样就可以省略判断, 减少计算量
1012
1013        //开始计算能带的导数, 详细的公式请看 berry_curvature_dipole_onek 的公式
1014        //其实就是速度算符的对角项
1015        //开始计算速度的偏导项, 这里偏导来自实空间
1016        let partial_ve_1 = v_1.diag().map(|x| x.re);
1017        let partial_ve_2 = v_2.diag().map(|x| x.re);
1018        let partial_ve_3 = v_3.diag().map(|x| x.re);
1019
1020        //开始最后的计算
1021        if self.spin {
1022            //如果考虑自旋, 我们就计算 \partial_h G_{ij}
1023            let mut S: Array2<Complex<f64>> = Array2::eye(self.nsta());
1024            let li = Complex::<f64>::new(0.0, 1.0);
1025            let pauli: Array2<Complex<f64>> = match spin {
1026                0 => Array2::<Complex<f64>>::eye(2),
1027                1 => {
1028                    arr2(&[
1029                        [0.0 + 0.0 * li, 1.0 + 0.0 * li],
1030                        [1.0 + 0.0 * li, 0.0 + 0.0 * li],
1031                    ]) / 2.0
1032                }
1033                2 => {
1034                    arr2(&[
1035                        [0.0 + 0.0 * li, 0.0 - 1.0 * li],
1036                        [0.0 + 1.0 * li, 0.0 + 0.0 * li],
1037                    ]) / 2.0
1038                }
1039                3 => {
1040                    arr2(&[
1041                        [1.0 + 0.0 * li, 0.0 + 0.0 * li],
1042                        [0.0 + 0.0 * li, -1.0 + 0.0 * li],
1043                    ]) / 2.0
1044                }
1045                _ => panic!("Wrong, spin should be 0, 1, 2, 3, but you input {}", spin),
1046            };
1047            let X = kron(&pauli, &Array2::eye(self.norb()));
1048            let mut S = Array3::<Complex<f64>>::zeros((self.dim_r(), self.nsta(), self.nsta()));
1049            for i in 0..self.dim_r() {
1050                let v0 = J.slice(s![i, .., ..]).to_owned();
1051                let v0 = anti_comm(&X, &v0) / 2.0;
1052                let v0 = evec_conj
1053                    .clone()
1054                    .dot(&(v0.dot(&evec.clone().reversed_axes()))); //变换到本征态基函数
1055                S.slice_mut(s![i, .., ..]).assign(&v0);
1056            }
1057            let mut s_1 = Array2::<Complex<f64>>::zeros((self.nsta(), self.nsta())); //三个方向的速度算符
1058            let mut s_2 = Array2::<Complex<f64>>::zeros((self.nsta(), self.nsta()));
1059            let mut s_3 = Array2::<Complex<f64>>::zeros((self.nsta(), self.nsta()));
1060            for i in 0..self.dim_r() {
1061                s_1 =
1062                    s_1.clone() + S.slice(s![i, .., ..]).to_owned() * Complex::new(dir_1[[i]], 0.0);
1063                s_2 =
1064                    s_2.clone() + S.slice(s![i, .., ..]).to_owned() * Complex::new(dir_2[[i]], 0.0);
1065                s_3 =
1066                    s_3.clone() + S.slice(s![i, .., ..]).to_owned() * Complex::new(dir_3[[i]], 0.0);
1067            }
1068            let G_23: Array1<f64> = {
1069                //用来计算  beta gamma 的 G
1070                let A = &v_2 * (U0.map(|x| Complex::<f64>::new(x.powi(3), 0.0)));
1071                let mut G = Array1::<f64>::zeros(self.nsta());
1072                for i in 0..self.nsta() {
1073                    G[[i]] = A.slice(s![i, ..]).dot(&v_3.slice(s![.., i])).re * 2.0
1074                }
1075                G
1076            };
1077            let G_13_h: Array1<f64> = {
1078                //用来计算 alpha gamma 的 G
1079                let A = &s_1 * (U0.map(|x| Complex::<f64>::new(x.powi(3), 0.0)));
1080                let mut G = Array1::<f64>::zeros(self.nsta());
1081                for i in 0..self.nsta() {
1082                    G[[i]] = A.slice(s![i, ..]).dot(&v_3.slice(s![.., i])).re * 2.0
1083                }
1084                G
1085            };
1086            //开始计算partial_s
1087            let partial_s_1 = s_1.clone().diag().map(|x| x.re);
1088            let partial_s_2 = s_2.clone().diag().map(|x| x.re);
1089            let partial_s_3 = s_3.clone().diag().map(|x| x.re);
1090            let mut partial_s = Array2::<f64>::zeros((self.dim_r(), self.nsta()));
1091            for r in 0..self.dim_r() {
1092                let s0 = S.slice(s![r, .., ..]).to_owned();
1093                partial_s
1094                    .slice_mut(s![r, ..])
1095                    .assign(&s0.diag().map(|x| x.re)); //\p_i s算符的对角项
1096            }
1097            //开始计算partial G
1098            let partial_G: Array1<f64> = {
1099                let mut A = Array1::<Complex<f64>>::zeros(self.nsta()); //第一项
1100                for i in 0..self.nsta() {
1101                    for j in 0..self.nsta() {
1102                        A[[i]] += 3.0
1103                            * (partial_s_1[[i]] - partial_s_1[[j]])
1104                            * v_2[[i, j]]
1105                            * v_3[[j, i]]
1106                            * U0[[i, j]].powi(4);
1107                    }
1108                }
1109                let mut B = Array1::<Complex<f64>>::zeros(self.nsta()); //第二项
1110                for n in 0..self.nsta() {
1111                    for n1 in 0..self.nsta() {
1112                        for n2 in 0..self.nsta() {
1113                            B[[n]] += s_1[[n, n2]]
1114                                * (v_2[[n2, n1]] * v_3[[n1, n]] + v_3[[n2, n1]] * v_2[[n1, n]])
1115                                * U0[[n, n1]].powi(3)
1116                                * U0[[n, n2]];
1117                        }
1118                    }
1119                }
1120                let mut C = Array1::<Complex<f64>>::zeros(self.nsta()); //第三项
1121                for n in 0..self.nsta() {
1122                    for n1 in 0..self.nsta() {
1123                        for n2 in 0..self.nsta() {
1124                            C[[n]] += s_1[[n1, n2]]
1125                                * (v_2[[n2, n]] * v_3[[n, n1]] + v_3[[n2, n]] * v_2[[n, n1]])
1126                                * U0[[n, n1]].powi(3)
1127                                * U0[[n1, n2]];
1128                        }
1129                    }
1130                }
1131                2.0 * (A - B - C).map(|x| x.re)
1132            };
1133            //计算结束
1134            //开始最后的输出
1135            return (
1136                (partial_s_1 * G_23 - partial_ve_2 * G_13_h),
1137                band,
1138                Some(partial_G),
1139            );
1140        } else {
1141            //开始计算 G_{ij}
1142            //G_{ij}=2Re\sum_{m\neq n} v_{i,nm}v_{j,mn}/(E_n-E_m)^3
1143            let G_23: Array1<f64> = {
1144                //用来计算  beta gamma 的 G
1145                let A = &v_2 * (U0.map(|x| Complex::<f64>::new(x.powi(3), 0.0)));
1146                let mut G = Array1::<f64>::zeros(self.nsta());
1147                for i in 0..self.nsta() {
1148                    G[[i]] = A.slice(s![i, ..]).dot(&v_3.slice(s![.., i])).re * 2.0
1149                }
1150                G
1151            };
1152            let G_13: Array1<f64> = {
1153                //用来计算 alpha gamma 的 G
1154                let A = &v_1 * (U0.map(|x| Complex::<f64>::new(x.powi(3), 0.0)));
1155                let mut G = Array1::<f64>::zeros(self.nsta());
1156                for i in 0..self.nsta() {
1157                    G[[i]] = A.slice(s![i, ..]).dot(&v_3.slice(s![.., i])).re * 2.0
1158                }
1159                G
1160            };
1161            return (partial_ve_1 * G_23 - partial_ve_2 * G_13, band, None);
1162        }
1163    }
1164    pub fn berry_connection_dipole(
1165        &self,
1166        k_vec: &Array2<f64>,
1167        dir_1: &Array1<f64>,
1168        dir_2: &Array1<f64>,
1169        dir_3: &Array1<f64>,
1170        spin: usize,
1171    ) -> (Array2<f64>, Array2<f64>, Option<Array2<f64>>) {
1172        //! 这个是基于 onek 的, 进行关于 k 点并行求解
1173        if dir_1.len() != self.dim_r() || dir_2.len() != self.dim_r() || dir_3.len() != self.dim_r()
1174        {
1175            panic!(
1176                "Wrong, the dir_1 or dir_2 you input has wrong length, it must equal to dim_r={}, but you input {},{}",
1177                self.dim_r(),
1178                dir_1.len(),
1179                dir_2.len()
1180            )
1181        }
1182        let nk = k_vec.len_of(Axis(0));
1183
1184        if self.spin {
1185            let ((omega, band), partial_G): ((Vec<_>, Vec<_>), Vec<_>) = k_vec
1186                .axis_iter(Axis(0))
1187                .into_par_iter()
1188                .map(|x| {
1189                    let (omega_one, band, partial_G) = self.berry_connection_dipole_onek(
1190                        &x.to_owned(),
1191                        &dir_1,
1192                        &dir_2,
1193                        &dir_3,
1194                        spin,
1195                    );
1196                    let partial_G = partial_G.unwrap();
1197                    ((omega_one, band), partial_G)
1198                })
1199                .collect();
1200
1201            let omega = Array2::<f64>::from_shape_vec(
1202                (nk, self.nsta()),
1203                omega.into_iter().flatten().collect(),
1204            )
1205            .unwrap();
1206            let band = Array2::<f64>::from_shape_vec(
1207                (nk, self.nsta()),
1208                band.into_iter().flatten().collect(),
1209            )
1210            .unwrap();
1211            let partial_G = Array2::<f64>::from_shape_vec(
1212                (nk, self.nsta()),
1213                partial_G.into_iter().flatten().collect(),
1214            )
1215            .unwrap();
1216
1217            return (omega, band, Some(partial_G));
1218        } else {
1219            let (omega, band): (Vec<_>, Vec<_>) = k_vec
1220                .axis_iter(Axis(0))
1221                .into_par_iter()
1222                .map(|x| {
1223                    let (omega_one, band, partial_G) = self.berry_connection_dipole_onek(
1224                        &x.to_owned(),
1225                        &dir_1,
1226                        &dir_2,
1227                        &dir_3,
1228                        spin,
1229                    );
1230                    (omega_one, band)
1231                })
1232                .collect();
1233            let omega = Array2::<f64>::from_shape_vec(
1234                (nk, self.nsta()),
1235                omega.into_iter().flatten().collect(),
1236            )
1237            .unwrap();
1238            let band = Array2::<f64>::from_shape_vec(
1239                (nk, self.nsta()),
1240                band.into_iter().flatten().collect(),
1241            )
1242            .unwrap();
1243            return (omega, band, None);
1244        }
1245    }
1246    pub fn Nonlinear_Hall_conductivity_Intrinsic(
1247        &self,
1248        k_mesh: &Array1<usize>,
1249        dir_1: &Array1<f64>,
1250        dir_2: &Array1<f64>,
1251        dir_3: &Array1<f64>,
1252        mu: &Array1<f64>,
1253        T: f64,
1254        spin: usize,
1255    ) -> Result<Array1<f64>> {
1256        //! The Intrinsic Nonlinear Hall Conductivity arises from the correction of the Berry connection by the electric and magnetic fields [PRL 112, 166601 (2014)]. The formula employed is:
1257        //!$$\tilde\bm\Og_{\bm k}=\nb_{\bm k}\times\lt(\bm A_{\bm k}+\bm A_{\bm k}^\prime\rt)$$
1258        //!and the $\bm A_{i,\bm k}^\prime=F_{ij}B_j+G_{ij}E_j$, where
1259        //!$$
1260        //!\\begin{aligned}
1261        //!F_{ij}&=\text{Im}\sum_{m=\not n}\f{v_{i,nm}\og_{j,mn}}{\lt(\ve_{n}-\ve_m\rt)^2}\\\\
1262        //!G_{ij}&=2\text{Re}\sum_{m=\not n}\f{v_{i,nm}v_{j,mn}}{\lt(\ve_n-\ve_m\rt)^3}\\\\
1263        //!\og_{\ap,mn}&=-i\ep_{\ap\bt\gm}\sum_{l=\not n}\f{\lt(v_{\bt,ml}+\p_\bt \ve_{\bm k}\dt_{ml}\rt)v_{\gm,ln}}{\ve_l-\ve_n}
1264        //!\\end{aligned}
1265        //!$$
1266        //!最后我们有
1267        //!$$
1268        //!\bm j^\prime=\bm E\times\int\f{\dd\bm k}{(2\pi)^3}\lt[\p_{\bm k}\ve_{\bm k}\times\bm A^\prime+\bm\Og\lt(\bm B\cdot\bm m\rt)\rt]\pdv{f_{\bm k}}{\ve}
1269        //!$$
1270        //!对其对电场和磁场进行偏导, 有
1271        //!$$
1272        //!\\begin{aligned}
1273        //!\f{\p^2 j_{\ap}^\prime}{\p E_\bt\p E_\gm}&=\int\f{\dd\bm k}{(2\pi)^3}\lt(\p_\ap\ve_{\bm k} G_{\bt\gm}-\p_\bt\ve_{\bm k} G_{\ap\gm}\rt)\pdv{f_{\bm k}}{\ve}\\\\
1274        //!\f{\p^2 j_{\ap}^\prime}{\p E_\bt\p B_\gm}&=\int\f{\dd\bm k}{(2\pi)^3}\lt(\p_\ap\ve_{\bm k} F_{\bt\gm}-\p_\bt\ve_{\bm k} F_{\ap\gm}+\ep_{\ap\bt\ell}\Og_{\ell} m_\gm\rt)\pdv{f_{\bm k}}{\ve}
1275        //!\\end{aligned}
1276        //!$$
1277        //!由于存在 $\pdv{f_{\bm k}}{\ve}$, 不建议将温度 T=0
1278        //!
1279        //!可以考虑当 T=0 时候, 利用高斯公式, 将费米面内的部分进行积分, 得到精确解. 但是我现在还没办法很好的求解费米面, 所以暂时不考虑这个算法.而且对于二维体系, 公式还不一样, 还得分步讨论, 后面有时间再考虑这个程序.
1280        //!
1281        //!对于自旋霍尔效应, 按照文章 [PRL 112, 166601 (2014)], 非线性自旋霍尔电导为
1282        //!$$\sg_{\ap\bt\gm}^i=-\int\dd\bm k \lt[\f{1}{2}f_{\bm k}\pdv{G_{\bt\gm}}{h_\ap}+\pdv{f_{\bm k}}{\ve}\lt(\p_{\ap}s_{\bm k}^i G_{\bt\gm}-\p_\bt\ve_{\bm k}G_{\ap\gm}^h\rt)\rt]$$
1283        //!其中
1284        //!$$\f{\p G_{\bt\gm,n}}{\p h_\ap}=2\text{Re}\sum_{n^\pr =\not n}\f{3\lt(s^i_{\ap,n}-s^i_{\ap,n_1}\rt)v_{\bt,nn_1} v_{\gm,n^\pr n}}{\lt(\ve_n-\ve_{n^\pr}\rt)^4}-2\text{Re}\sum_{n_1=\not n}\sum_{n_2=\not n}\lt[\f{s^i_{\ap,nn_2} v_{\bt,n_2n_1} v_{\gm,n_1 n}}{\lt(\ve_n-\ve_{n_1}\rt)^3(\ve_n-\ve_{n_2})}+(\bt \leftrightarrow \gm)\rt]-2\text{Re}\sum_{n_1=\not n}\sum_{n_2=\not n_1}\lt[\f{s^i_{\ap,n_1n_2} v_{\bt,n_2n} v_{\gm,n n_1}}{\lt(\ve_n-\ve_{n_1}\rt)^3(\ve_{n_1}-\ve_{n_2})}+(\bt \leftrightarrow \gm)\rt]$$
1285        //!以及
1286        //!$$
1287        //!\lt\\\{\\begin{aligned}
1288        //!G_{\ap\bt}&=2\text{Re}\sum_{m=\not n}\f{v_{\ap,nm}v_{\bt,mn}}{\lt(\ve_n-\ve_m\rt)^3}\\\\
1289        //!G_{\ap\bt}^h&=2\text{Re}\sum_{m=\not n}\f{s^i_{\ap,nm}v_{\bt,mn}}{\lt(\ve_n-\ve_m\rt)^3}\\\\
1290        //!\\end{aligned}\rt\.
1291        //!$$
1292        //!
1293        //!这里 $s^i_{\ap,mn}$ 的具体形式, 原文中没有明确给出, 但是我根据霍尔效应的类比, 我猜是
1294        //!$\\\{\hat s^i,v_\ap\\\}$
1295
1296        let kvec: Array2<f64> = gen_kmesh(&k_mesh)?;
1297        let nk: usize = kvec.len_of(Axis(0));
1298        let (omega, band, mut partial_G): (Array2<f64>, Array2<f64>, Option<Array2<f64>>) =
1299            self.berry_connection_dipole(&kvec, &dir_1, &dir_2, &dir_3, spin);
1300        let omega = omega.into_raw_vec();
1301        let omega = Array1::from(omega);
1302        let band0 = band.clone();
1303        let band = band.into_raw_vec();
1304        let band = Array1::from(band);
1305        let n_e = mu.len();
1306        let mut conductivity = Array1::<f64>::zeros(n_e);
1307        if T != 0.0 {
1308            let beta = 1.0 / T / 8.617e-5;
1309            let use_iter = band.iter().zip(omega.iter()).par_bridge();
1310            conductivity = use_iter
1311                .fold(
1312                    || Array1::<f64>::zeros(n_e),
1313                    |acc, (energy, omega0)| {
1314                        let f = 1.0 / ((beta * (mu - *energy)).mapv(|x| x.exp() + 1.0));
1315                        acc + &f * (1.0 - &f) * beta * *omega0
1316                    },
1317                )
1318                .reduce(|| Array1::<f64>::zeros(n_e), |acc, x| acc + x);
1319            if self.spin {
1320                let partial_G = partial_G.unwrap();
1321                let conductivity_new: Vec<f64> = mu
1322                    .into_par_iter()
1323                    .map(|x| {
1324                        let f = band0.map(|x0| 1.0 / ((beta * (x - x0)).exp() + 1.0));
1325                        let mut omega = Array1::<f64>::zeros(nk);
1326                        for i in 0..nk {
1327                            omega[[i]] = (partial_G.row(i).to_owned() * f.row(i).to_owned()).sum();
1328                        }
1329                        omega.sum() / 2.0
1330                    })
1331                    .collect();
1332                let conductivity_new = Array1::<f64>::from_vec(conductivity_new);
1333                conductivity = conductivity.clone() + conductivity_new;
1334            }
1335            conductivity = conductivity.clone() / (nk as f64) / self.lat.det().unwrap();
1336        } else {
1337            //采用四面体积分法, 或者对于二维体系, 采用三角形积分法
1338            //积分的思路是, 通过将一个六面体变成5个四面体, 然后用线性插值的方法, 得到费米面,
1339            //以及费米面上的数, 最后, 通过积分算出来结果
1340            panic!("the code can not support for T=0");
1341        }
1342        Ok(conductivity)
1343    }
1344}