Skip to main content

Rustb/
lib.rs

1#![allow(warnings)]
2pub mod SKmodel;
3pub mod atom_struct;
4pub mod conductivity;
5pub mod cut;
6pub mod error;
7pub mod generics;
8pub mod geometry;
9pub mod io;
10pub mod kpath;
11pub mod kpoints;
12pub mod magnetic_field;
13pub mod math;
14pub mod model;
15pub mod model_build;
16pub mod model_physics;
17pub mod model_utils;
18pub mod ndarray_lapack;
19pub mod optical_conductivity;
20pub mod orbital_angular;
21pub mod output;
22pub mod phy_const;
23pub mod solve_ham;
24pub mod surfgreen;
25pub mod unfold;
26pub mod velocity;
27pub mod wannier90;
28pub use crate::SKmodel::{SkAtom, SkParams, SlaterKosterModel, ToTbModel};
29pub use crate::atom_struct::{Atom, OrbProj};
30pub use crate::conductivity::*;
31pub use crate::cut::CutModel;
32pub use crate::error::{Result, TbError};
33use crate::generics::usefloat;
34pub use crate::geometry::*;
35pub use crate::io::*;
36pub use crate::kpath::*;
37pub use crate::kpoints::{gen_kmesh, gen_krange};
38pub use crate::magnetic_field::*;
39pub use crate::math::*;
40pub use crate::model::*;
41pub use crate::model_physics::*;
42pub use crate::optical_conductivity::*;
43pub use crate::output::*;
44pub use crate::solve_ham::solve;
45pub use crate::surfgreen::surf_Green;
46pub use crate::unfold::Unfold;
47pub use crate::velocity::*;
48pub use crate::wannier90::*;
49
50/// # Rustb - Tight-Binding Model Library
51///
52/// A comprehensive Rust library for tight-binding model calculations with support for:
53/// - Band structure calculations
54/// - Surface state computations using Green's functions
55/// - Linear and nonlinear transport properties
56/// - Berry phase and curvature calculations
57/// - Slater-Koster parameterized models
58///
59/// ## Key Features
60///
61/// - **Band Structure**: Solve eigenvalue problems $H(\mathbf{k}) \psi_n(\mathbf{k}) = E_n(\mathbf{k}) \psi_n(\mathbf{k})$
62/// - **Surface States**: Compute surface Green's functions $G^s(\omega, \mathbf{k}_\parallel)$
63/// - **Transport Properties**:
64///   - Anomalous Hall conductivity $\sigma_{xy} = \frac{e^2}{\hbar} \int \frac{d^2k}{(2\pi)^2} \Omega_z(\mathbf{k})$
65///   - Spin Hall conductivity
66///   - Nonlinear conductivity tensors
67/// - **Topological Invariants**: Chern numbers, Wilson loops, and Wannier centers
68/// - **Slater-Koster Models**: Parameterized tight-binding models with two-center integrals
69/// - **File I/O**: Utilities for writing calculation results to text files
70///
71/// ## Mathematical Foundation
72///
73/// The library implements the tight-binding Hamiltonian:
74/// $$
75/// H = \sum_{i,j} t_{ij} c_i^\dagger c_j + \sum_i \epsilon_i c_i^\dagger c_i
76/// $$
77/// where $t_{ij}$ are hopping parameters and $\epsilon_i$ are on-site energies.
78///
79/// For transport calculations, we compute the Berry curvature:
80/// $$
81/// \Omega_n(\mathbf{k}) = -2\,\text{Im}\sum_{m\neq n} \frac{\bra{n}\partial_{k_x} H\ket{m}\bra{m}\partial_{k_y} H\ket{n}}{(E_n - E_m)^2}
82/// $$
83///
84/// ## File I/O Utilities
85///
86/// The library provides file output utilities through the `io` module:
87/// - `write_txt`: Write 2D arrays to text files with formatted output
88/// - `write_txt_1`: Write 1D arrays to text files with formatted output
89///
90/// These functions automatically handle number formatting and spacing for scientific data.
91
92///An example
93///
94///```
95///use gnuplot::{Color,Figure, AxesCommon, AutoOption::Fix,HOT};
96///use gnuplot::Major;
97///use ndarray::*;
98///use ndarray::prelude::*;
99///use num_complex::Complex;
100///use Rustb::*;
101///
102///fn graphene(){
103///    let li:Complex<f64>=1.0*Complex::i();
104///    let t1=1.0+0.0*li;
105///    let t2=0.1+0.0*li;
106///    let t3=0.0+0.0*li;
107///    let delta=0.5;
108///    let dim_r:usize=2;
109///    let norb:usize=2;
110///    let lat=arr2(&[[3.0_f64.sqrt(),-1.0],[3.0_f64.sqrt(),1.0]]);
111///    let orb=arr2(&[[0.0,0.0],[1.0/3.0,1.0/3.0]]);
112///    let mut model=Model::tb_model(dim_r,lat,orb,false,None);
113///    model.set_onsite(&arr1(&[delta,-delta]),SpinDirection::None);
114///    model.add_hop(t1,0,1,&array![0,0],SpinDirection::None);
115///    model.add_hop(t1,0,1,&array![-1,0],SpinDirection::None);
116///    model.add_hop(t1,0,1,&array![0,-1],SpinDirection::None);
117///    model.add_hop(t2,0,0,&array![1,0],SpinDirection::None);
118///    model.add_hop(t2,1,1,&array![1,0],SpinDirection::None);
119///    model.add_hop(t2,0,0,&array![0,1],SpinDirection::None);
120///    model.add_hop(t2,1,1,&array![0,1],SpinDirection::None);
121///    model.add_hop(t2,0,0,&array![1,-1],SpinDirection::None);
122///    model.add_hop(t2,1,1,&array![1,-1],SpinDirection::None);
123///    model.add_hop(t3,0,1,&array![1,-1],SpinDirection::None);
124///    model.add_hop(t3,0,1,&array![-1,1],SpinDirection::None);
125///    model.add_hop(t3,0,1,&array![-1,-1],SpinDirection::None);
126///    let nk:usize=1001;
127///    let path=[[0.0,0.0],[2.0/3.0,1.0/3.0],[0.5,0.5],[0.0,0.0]];
128///    let path=arr2(&path);
129///    let (k_vec,k_dist,k_node)=model.k_path(&path,nk);
130///    let (eval,evec)=model.solve_all_parallel(&k_vec);
131///    let label=vec!["G","K","M","G"];
132///    let (k_vec,k_dist,k_node)=model.k_path(&path,nk); //generate the k vector
133///    let eval=model.solve_band_all_parallel(&k_vec);  //calculate the bands
134///    let mut fg = Figure::new();
135///    let x:Vec<f64>=k_dist.to_vec();
136///    let axes=fg.axes2d();
137///    for i in 0..model.nsta{
138///        let y:Vec<f64>=eval.slice(s![..,i]).to_owned().to_vec();
139///        axes.lines(&x, &y, &[Color("black")]);
140///    }
141///    let axes=axes.set_x_range(Fix(0.0), Fix(k_node[[k_node.len()-1]]));
142///    let label=label.clone();
143///    let mut show_ticks=Vec::new();
144///    for i in 0..k_node.len(){
145///        let A=k_node[[i]];
146///        let B=label[i];
147///        show_ticks.push(Major(A,Fix(B)));
148///    }
149///    axes.set_x_ticks_custom(show_ticks.into_iter(),&[],&[]);
150///    let k_node=k_node.to_vec();
151///    let mut jpg_name=String::new();
152///    jpg_name.push_str("band.jpg");
153///    fg.set_terminal("jpeg", &jpg_name);
154///    fg.show();
155///
156///    //start to draw the band structure
157///    //Starting to calculate the edge state, first is the zigzag state
158///    let nk:usize=501;
159///    let U=arr2(&[[1.0,1.0],[-1.0,1.0]]);
160///    let super_model=model.make_supercell(&U);
161///    let zig_model=super_model.cut_piece(100,0);
162///    let path=[[0.0,0.0],[0.0,0.5],[0.0,1.0]];
163///    let path=arr2(&path);
164///    let label=vec!["G","M","G"];
165///    zig_model.show_band(&path,&label,nk,"graphene_zig");
166///    //Starting to calculate the DOS of graphene
167///    let nk:usize=101;
168///    let kmesh=arr1(&[nk,nk]);
169///    let E_min=-3.0;
170///    let E_max=3.0;
171///    let E_n=1000;
172///    let (E0,dos)=model.dos(&kmesh,E_min,E_max,E_n,1e-2);
173///    //start to show DOS
174///    let mut fg = Figure::new();
175///    let x:Vec<f64>=E0.to_vec();
176///    let axes=fg.axes2d();
177///    let y:Vec<f64>=dos.to_vec();
178///    axes.lines(&x, &y, &[Color("black")]);
179///    let mut show_ticks=Vec::<String>::new();
180///    let mut pdf_name=String::new();
181///    pdf_name.push_str("dos.jpg");
182///    fg.set_terminal("pdfcairo", &pdf_name);
183///    fg.show();
184///}
185///```
186///
187///
188
189#[cfg(test)]
190mod tests {
191    use super::*;
192    use gnuplot::{
193        AutoOption, AxesCommon, Color, Figure, Fix, Font, LineStyle, Major, PointSymbol, Rotate,
194        Solid, TextOffset,
195    };
196    use ndarray::concatenate;
197    use ndarray::linalg::kron;
198    use ndarray::prelude::*;
199    use ndarray::*;
200    use ndarray_linalg::conjugate;
201    use ndarray_linalg::*;
202    use ndarray_linalg::{Eigh, UPLO};
203    use num_complex::Complex;
204    use rayon::prelude::*;
205    use std::f64::consts::PI;
206    use std::fs::File;
207    use std::io::Write;
208    use std::time::{Duration, Instant};
209    use std::fs::create_dir_all;
210
211    fn write_txt(data: Array2<f64>, output: &str) -> std::io::Result<()> {
212        let mut file = File::create(output).expect("Unable to BAND.dat");
213        let n = data.len_of(Axis(0));
214        let s = data.len_of(Axis(1));
215        let mut s0 = String::new();
216        for i in 0..n {
217            for j in 0..s {
218                if data[[i, j]] >= 0.0 {
219                    s0.push_str("     ");
220                } else {
221                    s0.push_str("    ");
222                }
223                let aa = format!("{:.6}", data[[i, j]]);
224                s0.push_str(&aa);
225            }
226            s0.push_str("\n");
227        }
228        writeln!(file, "{}", s0)?;
229        Ok(())
230    }
231
232    fn write_txt_1(data: Array1<f64>, output: &str) -> std::io::Result<()> {
233        use std::fs::File;
234        use std::io::Write;
235        let mut file = File::create(output).expect("Unable to BAND.dat");
236        let n = data.len_of(Axis(0));
237        let mut s0 = String::new();
238        for i in 0..n {
239            if data[[i]] >= 0.0 {
240                s0.push_str(" ");
241            }
242            let aa = format!("{:.6}\n", data[[i]]);
243            s0.push_str(&aa);
244        }
245        writeln!(file, "{}", s0)?;
246        Ok(())
247    }
248    #[test]
249    fn test_gen_v() {
250        //判断两个Array1<f64> 是否足够接近
251        fn are_arrays_close(a: &Array1<f64>, b: &Array1<f64>, tolerance: f64) -> bool {
252            a.iter()
253                .zip(b.iter())
254                .all(|(&x, &y)| (x - y).abs() < tolerance)
255        }
256
257        //判断两个Array2<Compelx<f64>> 是否足够接近
258        fn are_complex_arrays_close(
259            a: &Array2<Complex<f64>>,
260            b: &Array2<Complex<f64>>,
261            tolerance: f64,
262        ) -> bool {
263            a.iter()
264                .zip(b.iter())
265                .all(|(&x, &y)| (x.re - y.re).abs() < tolerance && (x.im - y.im).abs() < tolerance)
266        }
267        let li: Complex<f64> = 1.0 * Complex::i();
268        let t = 1.0;
269        let delta = 0.0;
270        let dim_r: usize = 2;
271        let norb: usize = 2;
272        let lat = arr2(&[[1.0, 0.0], [0.5, 3.0_f64.sqrt() / 2.0]]);
273        let orb = arr2(&[[1.0 / 3.0, 1.0 / 3.0], [2.0 / 3.0, 2.0 / 3.0]]);
274        let mut model = Model::tb_model(dim_r, lat, orb, false, None).unwrap();
275        model.set_onsite(&arr1(&[-delta, delta]), SpinDirection::None);
276        let R0: Array2<isize> = arr2(&[[0, 0], [-1, 0], [0, -1]]);
277        for (i, R) in R0.axis_iter(Axis(0)).enumerate() {
278            let R = R.to_owned();
279            model.set_hop(t, 0, 1, &R, SpinDirection::None);
280        }
281        assert_eq!(model.solve_band_onek(&array![0.0, 0.0]), array![-3.0, 3.0]);
282        let result = model.solve_band_onek(&array![1.0 / 3.0, 2.0 / 3.0]);
283        assert!(
284            are_arrays_close(&result, &array![0.0, 0.0], 1e-5),
285            "wrong!, the solve_band_onek get wrong result! please check it!"
286        );
287        let (result, _) = model.gen_v(&array![1.0 / 3.0, 1.0 / 3.0], Gauge::Atom);
288        let resulty = array![
289            [0.0 * li, -0.4698463103929542 - 0.17101007166283436 * li],
290            [-0.4698463103929542 + 0.17101007166283436 * li, 0.0 * li]
291        ];
292        let resultx = array![
293            [0.0 * li, -0.8137976813493737 - 0.2961981327260237 * li],
294            [-0.8137976813493737 + 0.2961981327260237 * li, 0.0 * li]
295        ];
296        println!("result={}", result);
297        assert!(
298            are_complex_arrays_close(&result.slice(s![0, .., ..]).to_owned(), &resultx, 1e-8),
299            "Wrong! the gen_v is get wrong results! please check it!"
300        );
301        assert!(
302            are_complex_arrays_close(&result.slice(s![1, .., ..]).to_owned(), &resulty, 1e-8),
303            "Wrong! the gen_v is get wrong results! please check it!"
304        );
305
306        let (result, _) = model.gen_v(&array![1.0 / 3.0, 1.0 / 3.0], Gauge::Lattice);
307        let resultx = array![
308            [
309                0.0 * li,
310                -3.0 * 3.0_f64.sqrt() / 4.0 * t + 3.0 / 4.0 * t * li
311            ],
312            [
313                -3.0 * 3.0_f64.sqrt() / 4.0 * t - 3.0 / 4.0 * t * li,
314                0.0 * li
315            ]
316        ];
317        println!("result={}", &result - &resultx);
318        assert!(
319            are_complex_arrays_close(&result.slice(s![0, .., ..]).to_owned(), &resultx, 1e-8),
320            "Wrong! the gen_v is get wrong results! please check it!"
321        );
322
323        let kvec = array![1.0 / 3.0, 1.0 / 3.0];
324        let (band, evec) = model.solve_onek(&kvec);
325        let ham = model.gen_ham(&kvec, Gauge::Atom);
326        let evec_conj = evec.map(|x| x.conj());
327        let evec = evec.t();
328        let ham = ham.dot(&evec);
329        let ham = evec_conj.dot(&ham);
330        let new_band = ham.diag().map(|x| x.re);
331        assert!(
332            are_arrays_close(&new_band, &band, 1e-5),
333            "wrong!, the solve_onek get wrong result! please check it!"
334        );
335    }
336    #[test]
337    fn conductivity_test() {
338        //这个是用 Haldan 模型来测试
339        let li: Complex<f64> = 1.0 * Complex::i();
340        let t = -1.0 + 0.0 * li;
341        let t2 = -1.0 + 0.0 * li;
342        let delta = 0.7;
343        let dim_r: usize = 2;
344        let norb: usize = 2;
345        let lat = arr2(&[[1.0, 0.0], [0.5, 3.0_f64.sqrt() / 2.0]]);
346        let orb = arr2(&[[1.0 / 3.0, 1.0 / 3.0], [2.0 / 3.0, 2.0 / 3.0]]);
347        let mut model = Model::tb_model(dim_r, lat, orb, false, None).unwrap();
348        model.set_onsite(&arr1(&[-delta, delta]), SpinDirection::None);
349        let R0: Array2<isize> = arr2(&[[0, 0], [-1, 0], [0, -1]]);
350        for (i, R) in R0.axis_iter(Axis(0)).enumerate() {
351            let R = R.to_owned();
352            model.add_hop(t, 0, 1, &R, SpinDirection::None);
353        }
354        let R0: Array2<isize> = arr2(&[[1, 0], [-1, 1], [0, -1]]);
355        for (i, R) in R0.axis_iter(Axis(0)).enumerate() {
356            let R = R.to_owned();
357            model.add_hop(t2 * li, 0, 0, &R, SpinDirection::None);
358        }
359        let R0: Array2<isize> = arr2(&[[-1, 0], [1, -1], [0, 1]]);
360        for (i, R) in R0.axis_iter(Axis(0)).enumerate() {
361            let R = R.to_owned();
362            model.add_hop(t2 * li, 1, 1, &R, SpinDirection::None);
363        }
364        let k_vec = array![1.0 / 3.0, 2.0 / 3.0];
365        let dir_1 = array![1.0, 0.0];
366        let dir_2 = array![0.0, 1.0];
367        let mu = 0.0;
368        let T = 0.0;
369        let og = 0.0;
370        let spin = 0;
371        let eta = 1e-3;
372        let result1 =
373            model.berry_curvature_onek(&k_vec, &dir_1, &dir_2, mu, T, spin, eta) * (2.0 * PI);
374
375        let mut k_list = Array2::zeros((9, 2));
376        let dk = 0.0001;
377        k_list.row_mut(0).assign(&(&k_vec + dk * &dir_1));
378        k_list
379            .row_mut(1)
380            .assign(&(&k_vec + dk * &dir_1 + dk * &dir_2));
381        k_list.row_mut(2).assign(&(&k_vec + dk * &dir_2));
382        k_list
383            .row_mut(3)
384            .assign(&(&k_vec - dk * &dir_1 + dk * &dir_2));
385        k_list.row_mut(4).assign(&(&k_vec - dk * &dir_1));
386        k_list
387            .row_mut(5)
388            .assign(&(&k_vec - dk * &dir_1 - dk * &dir_2));
389        k_list.row_mut(6).assign(&(&k_vec - dk * &dir_2));
390        k_list
391            .row_mut(7)
392            .assign(&(&k_vec + dk * &dir_1 - dk * &dir_2));
393        k_list.row_mut(8).assign(&(&k_vec + dk * &dir_1));
394        let result2 = model.berry_loop(&k_list, &vec![0]);
395        let result2 = result2[[0]] / (dk.powi(2)) / 4.0 / (2.0 * PI) * 3_f64.sqrt() / 2.0;
396        println!("result2={},result1={}", result2, result1);
397        assert!(
398            (result2 - result1).abs() < 1e-4,
399            "Wrong!, the berry_curvature or berry_flux mut be false"
400        );
401    }
402    #[test]
403    fn gen_v_speed_test() {
404        println!("开始测试各个函数的运行速度, 用次近邻的石墨烯模型");
405        let li: Complex<f64> = 1.0 * Complex::i();
406        let t = 2.0 + 0.0 * li;
407        let t2 = -1.0 + 0.0 * li;
408        let delta = 0.7;
409        let dim_r: usize = 2;
410        let norb: usize = 2;
411        let lat = arr2(&[[1.0, 0.0], [0.5, 3.0_f64.sqrt() / 2.0]]);
412        let orb = arr2(&[[1.0 / 3.0, 1.0 / 3.0], [2.0 / 3.0, 2.0 / 3.0]]);
413        let mut model = Model::tb_model(dim_r, lat, orb, false, None).unwrap();
414        model.set_onsite(&arr1(&[-delta, delta]), SpinDirection::None);
415        let R0: Array2<isize> = arr2(&[[0, 0], [-1, 0], [0, -1]]);
416        for (i, R) in R0.axis_iter(Axis(0)).enumerate() {
417            let R = R.to_owned();
418            model.add_hop(t, 0, 1, &R, SpinDirection::None);
419        }
420        let R0: Array2<isize> = arr2(&[[1, 0], [-1, 1], [0, -1]]);
421        for (i, R) in R0.axis_iter(Axis(0)).enumerate() {
422            let R = R.to_owned();
423            model.add_hop(t2 * li, 0, 0, &R, SpinDirection::None);
424        }
425        let R0: Array2<isize> = arr2(&[[-1, 0], [1, -1], [0, 1]]);
426        for (i, R) in R0.axis_iter(Axis(0)).enumerate() {
427            let R = R.to_owned();
428            model.add_hop(t2 * li, 1, 1, &R, SpinDirection::None);
429        }
430        println!("{:?}", model.atom_list());
431        let U = array![[3.0, 0.0], [0.0, 3.0]];
432        let model = model.make_supercell(&U).unwrap();
433
434        let nk = 101;
435        let k_mesh = array![nk, nk];
436        let kvec = gen_kmesh(&k_mesh).unwrap();
437
438        {
439            println!("开始计算 gen_v 的耗时速度, 为了平均, 我们单线程求解gen_v");
440            let start = Instant::now(); // 开始计时
441            let A: Vec<_> = kvec
442                .outer_iter()
443                .into_par_iter()
444                .map(|x| {
445                    let (a, _) = model.gen_v(&x.to_owned(), Gauge::Atom);
446                    a
447                })
448                .collect();
449            let end = Instant::now(); // 结束计时
450            let duration = end.duration_since(start); // 计算执行时间
451            println!(
452                "run gen_v {} times took {} seconds",
453                kvec.nrows(),
454                duration.as_secs_f64()
455            ); // 输出执行时间
456        }
457    }
458    #[test]
459    fn Haldan_model() {
460        let li: Complex<f64> = 1.0 * Complex::i();
461        let t = -1.0 + 0.0 * li;
462        let t2 = -1.0 + 0.0 * li;
463        let delta = 0.7;
464        let dim_r: usize = 2;
465        let norb: usize = 2;
466        let lat = arr2(&[[1.0, 0.0], [0.5, 3.0_f64.sqrt() / 2.0]]);
467        let orb = arr2(&[[1.0 / 3.0, 1.0 / 3.0], [2.0 / 3.0, 2.0 / 3.0]]);
468        let mut model = Model::tb_model(dim_r, lat, orb, false, None).unwrap();
469        model.set_onsite(&arr1(&[-delta, delta]), SpinDirection::None);
470        let R0: Array2<isize> = arr2(&[[0, 0], [-1, 0], [0, -1]]);
471        for (i, R) in R0.axis_iter(Axis(0)).enumerate() {
472            let R = R.to_owned();
473            model.add_hop(t, 0, 1, &R, SpinDirection::None);
474        }
475        let R0: Array2<isize> = arr2(&[[1, 0], [-1, 1], [0, -1]]);
476        for (i, R) in R0.axis_iter(Axis(0)).enumerate() {
477            let R = R.to_owned();
478            model.add_hop(t2 * li, 0, 0, &R, SpinDirection::None);
479        }
480        let R0: Array2<isize> = arr2(&[[-1, 0], [1, -1], [0, 1]]);
481        for (i, R) in R0.axis_iter(Axis(0)).enumerate() {
482            let R = R.to_owned();
483            model.add_hop(t2 * li, 1, 1, &R, SpinDirection::None);
484        }
485        let nk: usize = 101;
486        let path = [
487            [0.0, 0.0],
488            [2.0 / 3.0, 1.0 / 3.0],
489            [0.5, 0.5],
490            [1.0 / 3.0, 2.0 / 3.0],
491            [0.0, 0.0],
492        ];
493        let path = arr2(&path);
494        let (k_vec, k_dist, k_node) = model.k_path(&path, nk).unwrap();
495        let (eval, evec) = model.solve_all_parallel(&k_vec);
496        let label = vec!["G", "K", "M", "K'", "G"];
497        model.show_band(&path, &label, nk, "tests/Haldan").unwrap();
498        /////开始计算体系的霍尔电导率//////
499        let nk: usize = 31;
500        let T: f64 = 0.0;
501        let eta: f64 = 0.001;
502        let og: f64 = 0.0;
503        let mu: f64 = 0.0;
504        let dir_1 = arr1(&[1.0, 0.0]);
505        let dir_2 = arr1(&[0.0, 1.0]);
506        let dir_3 = arr1(&[0.0, 1.0]);
507        let spin: usize = 0;
508        let kmesh = arr1(&[nk, nk]);
509
510        let start = Instant::now(); // 开始计时
511        let conductivity = model
512            .Hall_conductivity(&kmesh, &dir_1, &dir_2, mu, T, spin, eta)
513            .unwrap();
514        let end = Instant::now(); // 结束计时
515        let duration = end.duration_since(start); // 计算执行时间
516        println!("quantom_Hall_effect={}", conductivity * (2.0 * PI));
517        assert!(
518            (conductivity * (2.0 * PI) - 1.0).abs() < 1e-3,
519            "Wrong!, the Hall conductivity is wrong!"
520        );
521        println!("function_a took {} seconds", duration.as_secs_f64()); // 输出执行时间
522
523        let mu = Array1::linspace(-2.0, 2.0, 101);
524        let start = Instant::now(); // 开始计时
525        let conductivity_mu = model
526            .Hall_conductivity_mu(&kmesh, &dir_1, &dir_2, &mu, T, spin, eta)
527            .unwrap();
528        let end = Instant::now(); // 结束计时
529        let duration = end.duration_since(start); // 计算执行时间
530        println!("quantom_Hall_effect={}", conductivity_mu[[50]] * (2.0 * PI));
531        assert!(
532            (conductivity_mu[[50]] - conductivity).abs() < 1e-3,
533            "Wrong!, the Hall conductivity is wrong!, Hall_mu's result is {}, but Hall conductivity is {}",
534            conductivity_mu[[50]],
535            conductivity
536        );
537        println!("function_a took {} seconds", duration.as_secs_f64()); // 输出执行时间
538        let conductivity = model
539            .Hall_conductivity(&kmesh, &dir_1, &dir_2, -2.0, T, spin, eta)
540            .unwrap();
541        assert!(
542            (conductivity_mu[[0]] - conductivity).abs() < 1e-3,
543            "Wrong!, the Hall conductivity is wrong!, Hall_mu's result is {}, but Hall conductivity is {}",
544            conductivity_mu[[0]],
545            conductivity
546        );
547        let conductivity = model
548            .Hall_conductivity(&kmesh, &dir_1, &dir_2, 2.0, T, spin, eta)
549            .unwrap();
550        assert!(
551            (conductivity_mu[[100]] - conductivity).abs() < 1e-3,
552            "Wrong!, the Hall conductivity is wrong!, Hall_mu's result is {}, but Hall conductivity is {}",
553            conductivity_mu[[100]],
554            conductivity
555        );
556        //开始绘图
557        let mut fg = Figure::new();
558        let x: Vec<f64> = mu.to_vec();
559        let axes = fg.axes2d();
560        let y: Vec<f64> = (conductivity_mu * 2.0 * PI).to_vec();
561        axes.lines(&x, &y, &[Color("black")]);
562        let mut show_ticks = Vec::<String>::new();
563        let mut pdf_name = String::new();
564        pdf_name.push_str("tests/Haldan");
565        pdf_name.push_str("/hall_mu.pdf");
566        fg.set_terminal("pdfcairo", &pdf_name);
567        fg.show();
568
569        let mu = 0.0;
570        let nk: usize = 31;
571        let kmesh = arr1(&[nk, nk]);
572        let start = Instant::now(); // 开始计时
573        let conductivity = model
574            .Hall_conductivity_adapted(&kmesh, &dir_1, &dir_2, mu, T, spin, eta, 0.01, 0.0001)
575            .unwrap();
576        let end = Instant::now(); // 结束计时
577        let duration = end.duration_since(start); // 计算执行时间
578        println!("霍尔电导率{}", conductivity * (2.0 * PI));
579        assert!(
580            (conductivity * (2.0 * PI) - 1.0).abs() < 1e-3,
581            "Wrong!, the Hall conductivity is wrong!"
582        );
583        println!("function_a took {} seconds", duration.as_secs_f64()); // 输出执行时间
584        //画一下3000k的时候的费米导数分布
585        let T = 100.0;
586        let nk: usize = 101;
587        let kmesh = arr1(&[nk, nk]);
588        println!("{}", kmesh);
589        let E_min = -3.0;
590        let E_max = 3.0;
591        let E_n = 1000;
592        let mu = Array1::linspace(E_min, E_max, E_n);
593        let beta: f64 = 1.0 / T / (8.617e-5);
594        let f: Array1<f64> = 1.0 / ((beta * &mu).mapv(f64::exp) + 1.0);
595        let par_f = beta * &f * (1.0 - &f);
596        let mut fg = Figure::new();
597        let x: Vec<f64> = mu.to_vec();
598        let axes = fg.axes2d();
599        let y: Vec<f64> = par_f.to_vec();
600        axes.lines(&x, &y, &[Color("black")]);
601        let mut show_ticks = Vec::<String>::new();
602        let mut pdf_name = String::new();
603        pdf_name.push_str("tests/Haldan");
604        pdf_name.push_str("/par_f.pdf");
605        fg.set_terminal("pdfcairo", &pdf_name);
606        fg.show();
607
608        //画一下omega_n 随能量的分布
609        let kvec: Array2<f64> = gen_kmesh(&kmesh).unwrap();
610        let nk: usize = kvec.len_of(Axis(0));
611        let (omega, band) =
612            model.berry_curvature_dipole_n(&kvec, &dir_1, &dir_2, &dir_3, og, spin, eta);
613        let omega = omega.into_raw_vec();
614        let omega = Array1::from(omega);
615        let band = band.into_raw_vec();
616        let band = Array1::from(band);
617        let mut fg = Figure::new();
618        let x: Vec<f64> = band.to_vec();
619        let axes = fg.axes2d();
620        let y: Vec<f64> = omega.to_vec();
621        axes.points(
622            x.iter(),
623            y.iter(),
624            &[Color("black"), PointSymbol((".").chars().next().unwrap())],
625        );
626        let mut show_ticks = Vec::<String>::new();
627        let mut pdf_name = String::new();
628        pdf_name.push_str("tests/Haldan");
629        pdf_name.push_str("/omega_energy.pdf");
630        fg.set_terminal("pdfcairo", &pdf_name);
631        fg.show();
632
633        //画一下表面态
634        let nk = 101;
635        let green = surf_Green::from_Model(&model, 0, 1e-3, None).unwrap();
636        let E_min = -3.0;
637        let E_max = 3.0;
638        let E_n = 101;
639        let path = [[0.0], [0.5], [1.0]];
640        let path = arr2(&path);
641        let label = vec!["G", "M", "G"];
642        green.show_surf_state("tests/Haldan/surf", &path, &label, nk, E_min, E_max, E_n, 0);
643
644        //-----算一下wilson loop 的结果-----------------------
645        let dir_1 = arr1(&[1.0, 0.0]);
646        let dir_2 = arr1(&[0.0, 1.0]);
647        let occ = vec![0];
648        let wcc = model.wannier_centre(&occ, &array![0.0, 0.0], &dir_1, &dir_2, 101, 101);
649        let nocc = occ.len();
650
651        let mut fg = Figure::new();
652        let x: Vec<f64> = Array1::<f64>::linspace(0.0, 1.0, 101).to_vec();
653        let axes = fg.axes2d();
654        for j in -1..2 {
655            for i in 0..nocc {
656                let a = wcc.row(i).to_owned() + (j as f64) * 2.0 * PI;
657                let y: Vec<f64> = a.to_vec();
658                axes.points(
659                    &x,
660                    &y,
661                    &[
662                        Color("black"),
663                        gnuplot::PointSymbol('O'),
664                        gnuplot::PointSize(0.2),
665                    ],
666                );
667            }
668        }
669        let axes = axes.set_x_range(Fix(0.0), Fix(1.0));
670        let axes = axes.set_y_range(Fix(0.0), Fix(2.0 * PI));
671        let show_ticks = vec![
672            Major(0.0, Fix("0")),
673            Major(0.5, Fix("π")),
674            Major(1.0, Fix("2π")),
675        ];
676        axes.set_x_ticks_custom(
677            show_ticks.into_iter(),
678            &[],
679            &[Font("Times New Roman", 32.0)],
680        );
681        let show_ticks = vec![
682            Major(0.0, Fix("0")),
683            Major(PI, Fix("π")),
684            Major(2.0 * PI, Fix("2π")),
685        ];
686        axes.set_y_ticks_custom(
687            show_ticks.into_iter(),
688            &[],
689            &[Font("Times New Roman", 32.0)],
690        );
691        axes.set_x_label(
692            "k_x",
693            &[Font("Times New Roman", 32.0), TextOffset(0.0, -0.5)],
694        );
695        axes.set_y_label(
696            "WCC",
697            &[
698                Font("Times New Roman", 32.0),
699                Rotate(90.0),
700                TextOffset(-1.0, 0.0),
701            ],
702        );
703        let mut pdf_name = String::new();
704        pdf_name.push_str("tests/Haldan/wcc.pdf");
705        fg.set_terminal("pdfcairo", &pdf_name);
706        fg.show();
707        //-----------用 berry_flux 算一下
708        let C = model
709            .berry_flux(
710                &occ,
711                &array![0.0, 0.0],
712                &array![1.0, 0.0],
713                &array![0.0, 1.0],
714                101,
715                101,
716            )
717            .sum()
718            / PI
719            / 2.0;
720        println!("The Chern number of Haldan model is {}", C);
721    }
722    #[test]
723    fn graphene() {
724        let li: Complex<f64> = 1.0 * Complex::i();
725        let t1 = 1.0 + 0.0 * li;
726        let t2 = 0.0 + 0.0 * li;
727        let t3 = 0.0 + 0.0 * li;
728        let delta = 0.0;
729        let dim_r: usize = 2;
730        let norb: usize = 2;
731        let lat = arr2(&[[3.0_f64.sqrt(), -1.0], [3.0_f64.sqrt(), 1.0]]);
732        let orb = arr2(&[[0.0, 0.0], [1.0 / 3.0, 1.0 / 3.0]]);
733        let mut model = Model::tb_model(dim_r, lat, orb, false, None).unwrap();
734        model.set_onsite(&arr1(&[delta, -delta]), SpinDirection::None);
735        model.add_hop(t1, 0, 1, &array![0, 0], SpinDirection::None);
736        model.add_hop(t1, 0, 1, &array![-1, 0], SpinDirection::None);
737        model.add_hop(t1, 0, 1, &array![0, -1], SpinDirection::None);
738        model.add_hop(t2, 0, 0, &array![1, 0], SpinDirection::None);
739        model.add_hop(t2, 1, 1, &array![1, 0], SpinDirection::None);
740        model.add_hop(t2, 0, 0, &array![0, 1], SpinDirection::None);
741        model.add_hop(t2, 1, 1, &array![0, 1], SpinDirection::None);
742        model.add_hop(t2, 0, 0, &array![1, -1], SpinDirection::None);
743        model.add_hop(t2, 1, 1, &array![1, -1], SpinDirection::None);
744        model.add_hop(t3, 0, 1, &array![1, -1], SpinDirection::None);
745        model.add_hop(t3, 0, 1, &array![-1, 1], SpinDirection::None);
746        model.add_hop(t3, 0, 1, &array![-1, -1], SpinDirection::None);
747        let nk: usize = 101;
748        let path = [[0.0, 0.0], [2.0 / 3.0, 1.0 / 3.0], [0.5, 0.5], [0.0, 0.0]];
749        let path = arr2(&path);
750        let (k_vec, k_dist, k_node) = model.k_path(&path, nk).unwrap();
751        let (eval, evec) = model.solve_all_parallel(&k_vec);
752        let label = vec!["G", "K", "M", "G"];
753        model
754            .show_band(&path, &label, nk, "tests/graphene")
755            .unwrap();
756
757        // 开始计算两个本征态
758        let k1 = array![1.0 / 3.0 - 0.002, 2.0 / 3.0];
759        let k2 = array![1.0 / 3.0 + 0.001, 2.0 / 3.0];
760        let (eval1, evec1) = model.solve_onek(&k1);
761        let (eval2, evec2) = model.solve_onek(&k2);
762        let evec1 = evec1.reversed_axes();
763        let evec2 = evec2.mapv(|x| x.conj());
764        println!("{},{}", eval1, eval2);
765        println!("{}", evec2.dot(&evec1).mapv(|x| x.norm().round()));
766
767        /////开始计算体系的霍尔电导率//////
768        let nk: usize = 11;
769        let T: f64 = 0.0;
770        let eta: f64 = 0.001;
771        let og: f64 = 0.0;
772        let mu: f64 = 0.0;
773        let dir_1 = arr1(&[1.0, 0.0]);
774        let dir_2 = arr1(&[0.0, 1.0]);
775        let spin: usize = 0;
776        let kmesh = arr1(&[nk, nk]);
777        let (eval, evec) = model.solve_onek(&arr1(&[0.3, 0.5]));
778        let conductivity = model.Hall_conductivity(&kmesh, &dir_1, &dir_2, mu, T, spin, eta);
779        //println!("{}",conductivity/(2.0*PI));
780        //开始计算边缘态, 首先是zigsag态
781        let nk: usize = 501;
782        let U = arr2(&[[1.0, 1.0], [-1.0, 1.0]]);
783        let super_model = model.make_supercell(&U).unwrap();
784        let zig_model = super_model.cut_piece(100, 0).unwrap();
785        let path = [[0.0, 0.0], [0.0, 0.5], [0.0, 1.0]];
786        //let path=[[0.0,0.0],[0.5,0.0],[1.0,0.0]];
787        //let path=[[0.0,0.0],[0.5,0.0],[0.5,0.5],[0.0,0.5],[0.0,0.0]];
788        let path = arr2(&path);
789        let (k_vec, k_dist, k_node) = super_model.k_path(&path, nk).unwrap();
790        let (eval, evec) = super_model.solve_all_parallel(&k_vec);
791        //let label=vec!["G","X","M","Y","G"];
792        let label = vec!["G", "M", "G"];
793        zig_model.show_band(&path, &label, nk, "tests/graphene_zig");
794
795        //开始计算石墨烯的态密度
796        let nk: usize = 51;
797        let kmesh = arr1(&[nk, nk]);
798        let E_min = -3.0;
799        let E_max = 3.0;
800        let E_n = 1000;
801        let (E0, dos) = model.dos(&kmesh, E_min, E_max, E_n, 1e-2).unwrap();
802        //开始绘制dos
803        let mut fg = Figure::new();
804        let x: Vec<f64> = E0.to_vec();
805        let axes = fg.axes2d();
806        let y: Vec<f64> = dos.to_vec();
807        axes.lines(&x, &y, &[Color("black")]);
808        let mut show_ticks = Vec::<String>::new();
809        let mut pdf_name = String::new();
810        pdf_name.push_str("tests/graphene");
811        pdf_name.push_str("/dos.pdf");
812        fg.set_terminal("pdfcairo", &pdf_name);
813        fg.show();
814
815        //开始计算非线性霍尔电导
816        let dir_1 = arr1(&[1.0, 0.0]);
817        let dir_2 = arr1(&[0.0, 1.0]);
818        let dir_3 = arr1(&[1.0, 0.0]);
819        let og = 0.0;
820        let mu = Array1::linspace(E_min, E_max, E_n);
821        let T = 300.0;
822        let sigma: Array1<f64> = model
823            .Nonlinear_Hall_conductivity_Extrinsic(
824                &kmesh, &dir_1, &dir_2, &dir_3, &mu, T, og, 0, 1e-5,
825            )
826            .unwrap();
827
828        //开始绘制非线性电导
829        let mut fg = Figure::new();
830        let x: Vec<f64> = mu.to_vec();
831        let axes = fg.axes2d();
832        let y: Vec<f64> = sigma.to_vec();
833        axes.lines(&x, &y, &[Color("black")]);
834        let mut show_ticks = Vec::<String>::new();
835        let mut pdf_name = String::new();
836        pdf_name.push_str("tests/graphene");
837        pdf_name.push_str("/nonlinear_ex.pdf");
838        fg.set_terminal("pdfcairo", &pdf_name);
839        fg.show();
840    }
841
842    #[test]
843    fn kane_mele() {
844        let li: Complex<f64> = 1.0 * Complex::i();
845        let t = -1.0;
846        let delta = 0.0;
847        let alter = 0.0 + 0.0 * li;
848        let soc = 0.06 * t;
849        let rashba = 0.0 * t;
850        let dim_r: usize = 2;
851        let norb: usize = 2;
852        let lat = arr2(&[[1.0, 0.0], [0.5, 3.0_f64.sqrt() / 2.0]]);
853        let orb = arr2(&[[1.0 / 3.0, 1.0 / 3.0], [2.0 / 3.0, 2.0 / 3.0]]);
854        let mut model = Model::tb_model(dim_r, lat, orb, true, None).unwrap();
855        model.set_onsite(&arr1(&[delta, -delta]), SpinDirection::None);
856        let R0: Array2<isize> = arr2(&[[0, 0], [-1, 0], [0, -1]]);
857        for (i, R) in R0.axis_iter(Axis(0)).enumerate() {
858            let R = R.to_owned();
859            model.set_hop(t, 0, 1, &R, SpinDirection::None);
860        }
861        let R0: Array2<isize> = arr2(&[[1, 0], [-1, 1], [0, -1]]);
862        for (i, R) in R0.axis_iter(Axis(0)).enumerate() {
863            let R = R.to_owned();
864            model.set_hop(soc * li, 0, 0, &R, SpinDirection::z);
865        }
866        let R0: Array2<isize> = arr2(&[[-1, 0], [1, -1], [0, 1]]);
867        for (i, R) in R0.axis_iter(Axis(0)).enumerate() {
868            let R = R.to_owned();
869            model.set_hop(soc * li, 1, 1, &R, SpinDirection::z);
870        }
871        //加入rashba项
872        let R0: Array2<isize> = arr2(&[[1, 0], [-1, 1], [0, -1]]);
873        for (i, R) in R0.axis_iter(Axis(0)).enumerate() {
874            let R = R.to_owned();
875            let r0 = R.map(|x| *x as f64).dot(&model.lat);
876            model.add_hop(rashba * li * r0[[1]], 0, 0, &R, SpinDirection::x);
877            model.add_hop(rashba * li * r0[[0]], 0, 0, &R, SpinDirection::y);
878        }
879
880        let R0: Array2<isize> = arr2(&[[-1, 0], [1, -1], [0, 1]]);
881        for (i, R) in R0.axis_iter(Axis(0)).enumerate() {
882            let R = R.to_owned();
883            let r0 = R.map(|x| *x as f64).dot(&model.lat);
884            model.add_hop(-rashba * li * r0[[1]], 1, 1, &R, SpinDirection::x);
885            model.add_hop(-rashba * li * r0[[0]], 1, 1, &R, SpinDirection::y);
886        }
887        let nk: usize = 101;
888        let path = [
889            [0.0, 0.0],
890            [2.0 / 3.0, 1.0 / 3.0],
891            [0.5, 0.5],
892            [1.0 / 3.0, 2.0 / 3.0],
893            [0.0, 0.0],
894        ];
895        let path = arr2(&path);
896        let (k_vec, k_dist, k_node) = model.k_path(&path, nk).unwrap();
897        let (eval, evec) = model.solve_all_parallel(&k_vec);
898        let label = vec!["G", "K", "M", "K'", "G"];
899        model.show_band(&path, &label, nk, "tests/kane");
900        //开始计算超胞
901
902        let super_model = model.cut_piece(50, 0).unwrap();
903        let path = [[0.0, 0.0], [0.0, 0.5], [0.0, 1.0]];
904        let path = arr2(&path);
905        let label = vec!["G", "M", "G"];
906        super_model
907            .show_band(&path, &label, nk, "tests/kane_super")
908            .unwrap();
909        //开始计算表面态
910        let nk = 101;
911        let green = surf_Green::from_Model(&model, 0, 1e-3, None).unwrap();
912        let E_min = -1.0;
913        let E_max = 1.0;
914        let E_n = 101;
915        let path = [[0.0], [0.5], [1.0]];
916        let path = arr2(&path);
917        let label = vec!["G", "M", "G"];
918        green.show_surf_state("tests/kane", &path, &label, nk, E_min, E_max, E_n, 0);
919
920        //-----算一下wilson loop 结果-----------------------
921        let n = 51;
922        let dir_1 = arr1(&[1.0, 0.0]);
923        let dir_2 = arr1(&[0.0, 1.0]);
924        let occ = vec![0, 1];
925        let wcc = model.wannier_centre(&occ, &array![0.0, 0.0], &dir_1, &dir_2, n, n);
926        let nocc = occ.len();
927        let mut fg = Figure::new();
928        let x: Vec<f64> = Array1::<f64>::linspace(0.0, 1.0, n).to_vec();
929        let axes = fg.axes2d();
930        for j in -1..2 {
931            for i in 0..nocc {
932                let a = wcc.row(i).to_owned() + (j as f64) * 2.0 * PI;
933                let y: Vec<f64> = a.to_vec();
934                axes.points(&x, &y, &[Color("black"), gnuplot::PointSymbol('O')]);
935            }
936        }
937        let axes = axes.set_x_range(Fix(0.0), Fix(1.0));
938        let axes = axes.set_y_range(Fix(0.0), Fix(2.0 * PI));
939        let show_ticks = vec![
940            Major(0.0, Fix("0")),
941            Major(0.5, Fix("π")),
942            Major(1.0, Fix("2π")),
943        ];
944        axes.set_x_ticks_custom(
945            show_ticks.into_iter(),
946            &[],
947            &[Font("Times New Roman", 32.0)],
948        );
949        let show_ticks = vec![
950            Major(0.0, Fix("0")),
951            Major(PI, Fix("π")),
952            Major(2.0 * PI, Fix("2π")),
953        ];
954        axes.set_y_ticks_custom(
955            show_ticks.into_iter(),
956            &[],
957            &[Font("Times New Roman", 32.0)],
958        );
959        axes.set_x_label(
960            "k_x",
961            &[Font("Times New Roman", 32.0), TextOffset(0.0, -0.5)],
962        );
963        axes.set_y_label(
964            "WCC",
965            &[
966                Font("Times New Roman", 32.0),
967                Rotate(90.0),
968                TextOffset(-1.0, 0.0),
969            ],
970        );
971        let mut pdf_name = String::new();
972        pdf_name.push_str("tests/kane/wcc.pdf");
973        fg.set_terminal("pdfcairo", &pdf_name);
974        fg.show();
975
976        /////开始计算体系的霍尔电导率//////
977        let nk: usize = 31;
978        let T: f64 = 0.0;
979        let eta: f64 = 0.001;
980        let og: f64 = 0.0;
981        let mu: f64 = 0.0;
982        //let dir_1=arr1(&[3.0_f64.sqrt()/2.0,-0.5]);
983        let dir_1 = arr1(&[1.0, 0.0]);
984        let dir_2 = arr1(&[0.0, 1.0]);
985        let spin: usize = 3;
986        let kmesh = arr1(&[nk, nk]);
987        let start = Instant::now(); // 开始计时
988        let conductivity = model
989            .Hall_conductivity(&kmesh, &dir_1, &dir_2, mu, T, spin, eta)
990            .unwrap();
991        let end = Instant::now(); // 结束计时
992        let duration = end.duration_since(start); // 计算执行时间
993        println!("{}", conductivity * (2.0 * PI));
994        println!("function_a took {} seconds", duration.as_secs_f64()); // 输出执行时间
995        let nk: usize = 21;
996        let kmesh = arr1(&[nk, nk]);
997        let start = Instant::now(); // 开始计时
998        let conductivity = model
999            .Hall_conductivity_adapted(&kmesh, &dir_1, &dir_2, mu, T, spin, eta, 0.01, 0.01)
1000            .unwrap();
1001        let end = Instant::now(); // 结束计时
1002        let duration = end.duration_since(start); // 计算执行时间
1003        println!("{}", conductivity * (2.0 * PI));
1004        println!("function_a took {} seconds", duration.as_secs_f64()); // 输出执行时间
1005
1006        let (E0, dos) = model.dos(&kmesh, E_min, E_max, E_n, 1e-2).unwrap();
1007        //开始绘制dos
1008        let mut fg = Figure::new();
1009        let x: Vec<f64> = E0.to_vec();
1010        let axes = fg.axes2d();
1011        let y: Vec<f64> = dos.to_vec();
1012        axes.lines(&x, &y, &[Color("black")]);
1013        let mut show_ticks = Vec::<String>::new();
1014        let mut pdf_name = String::new();
1015        pdf_name.push_str("tests/kane");
1016        pdf_name.push_str("/dos.pdf");
1017        fg.set_terminal("pdfcairo", &pdf_name);
1018        fg.show();
1019        //绘制非线性霍尔电导的平面图
1020
1021        //画一下贝利曲率的分布
1022        let nk: usize = 31;
1023        let kmesh = arr1(&[nk, nk]);
1024        let kvec = gen_kmesh(&kmesh).unwrap();
1025        //let kvec=kvec-0.5;
1026        let kvec = kvec * 2.0;
1027        let kvec = model.lat.dot(&(kvec.reversed_axes()));
1028        let kvec = kvec.reversed_axes();
1029        let berry_curv = model.berry_curvature(&kvec, &dir_1, &dir_2, T, 0.0, 1, 1e-3);
1030        let data = berry_curv.into_shape((nk, nk)).unwrap();
1031        draw_heatmap(
1032            &(-data).map(|x| (x + 1.0).log(10.0)),
1033            "./tests/kane/berry_curvature_distribution.pdf",
1034        );
1035
1036        //开始考虑磁场, 加入磁性
1037        let B = 0.1 + 0.0 * li;
1038        let tha = 0.0 / 180.0 * PI;
1039
1040        model.add_hop(B * tha.cos(), 0, 0, &array![0, 0], SpinDirection::x);
1041        model.add_hop(B * tha.cos(), 1, 1, &array![0, 0], SpinDirection::x);
1042        model.add_hop(B * tha.sin(), 0, 0, &array![0, 0], SpinDirection::y);
1043        model.add_hop(B * tha.sin(), 1, 1, &array![0, 0], SpinDirection::y);
1044        //考虑添加onsite 项破坏空间反演和mirror
1045
1046        let green = surf_Green::from_Model(&model, 0, 1e-3, None).unwrap();
1047        let E_min = -1.0;
1048        let E_max = 1.0;
1049        let E_n = nk;
1050        let path = [[0.0], [0.5], [1.0]];
1051        let path = arr2(&path);
1052        let label = vec!["G", "M", "G"];
1053        green.show_surf_state(
1054            "tests/kane/magnetic",
1055            &path,
1056            &label,
1057            nk,
1058            E_min,
1059            E_max,
1060            E_n,
1061            0,
1062        );
1063
1064        //-----算一下wilson loop 结果-----------------------
1065        let n = 51;
1066        let dir_1 = arr1(&[1.0, 0.0]);
1067        let dir_2 = arr1(&[0.0, 1.0]);
1068        let occ = vec![0, 1];
1069        let wcc = model.wannier_centre(&occ, &array![0.0, 0.0], &dir_1, &dir_2, n, n);
1070        let nocc = occ.len();
1071        let mut fg = Figure::new();
1072        let x: Vec<f64> = Array1::<f64>::linspace(0.0, 1.0, n).to_vec();
1073        let axes = fg.axes2d();
1074        for j in -1..2 {
1075            for i in 0..nocc {
1076                let a = wcc.row(i).to_owned() + (j as f64) * 2.0 * PI;
1077                let y: Vec<f64> = a.to_vec();
1078                axes.points(&x, &y, &[Color("black"), gnuplot::PointSymbol('O')]);
1079            }
1080        }
1081        let axes = axes.set_x_range(Fix(0.0), Fix(1.0));
1082        let axes = axes.set_y_range(Fix(0.0), Fix(2.0 * PI));
1083        let show_ticks = vec![
1084            Major(0.0, Fix("0")),
1085            Major(0.5, Fix("π")),
1086            Major(1.0, Fix("2π")),
1087        ];
1088        axes.set_x_ticks_custom(
1089            show_ticks.into_iter(),
1090            &[],
1091            &[Font("Times New Roman", 32.0)],
1092        );
1093        let show_ticks = vec![
1094            Major(0.0, Fix("0")),
1095            Major(PI, Fix("π")),
1096            Major(2.0 * PI, Fix("2π")),
1097        ];
1098        axes.set_y_ticks_custom(
1099            show_ticks.into_iter(),
1100            &[],
1101            &[Font("Times New Roman", 32.0)],
1102        );
1103        axes.set_x_label(
1104            "k_x",
1105            &[Font("Times New Roman", 32.0), TextOffset(0.0, -0.5)],
1106        );
1107        axes.set_y_label(
1108            "WCC",
1109            &[
1110                Font("Times New Roman", 32.0),
1111                Rotate(90.0),
1112                TextOffset(-1.0, 0.0),
1113            ],
1114        );
1115        let mut pdf_name = String::new();
1116        pdf_name.push_str("tests/kane/magnetic/wcc.pdf");
1117        fg.set_terminal("pdfcairo", &pdf_name);
1118        fg.show();
1119
1120        //开始计算角态
1121        let model = model
1122            .make_supercell(&array![[0.0, -1.0], [1.0, 0.0]])
1123            .unwrap();
1124        let num = 19;
1125        /*
1126        let model_1=model.cut_piece(num,0).unwrap();
1127        let new_model=model_1.cut_piece(num,1);
1128        */
1129        let new_model = model.cut_dot(num, 6, None).unwrap();
1130        let mut s = 0;
1131        let start = Instant::now();
1132        let (band, evec) = new_model.solve_range_onek(&arr1(&[0.0, 0.0]), (-0.3, 0.3), 1e-5);
1133        let end = Instant::now(); // 结束计时
1134        let duration = end.duration_since(start); // 计算执行时间
1135        println!("solve_band_all took {} seconds", duration.as_secs_f64()); // 输出执行时间
1136        let nresults = band.len();
1137        let show_evec = evec.to_owned().map(|x| x.norm_sqr());
1138        let mut size = Array2::<f64>::zeros((new_model.nsta(), new_model.natom()));
1139        let norb = new_model.norb();
1140        for i in 0..nresults {
1141            let mut s = 0;
1142            for j in 0..new_model.natom() {
1143                for k in 0..new_model.atoms[j].norb() {
1144                    size[[i, j]] += show_evec[[i, s]] + show_evec[[i, s + new_model.norb()]];
1145                    s += 1;
1146                }
1147            }
1148        }
1149
1150        let show_str = new_model.atom_position().dot(&model.lat);
1151        let show_str = show_str.slice(s![.., 0..2]).to_owned();
1152        let show_size = size.row(new_model.norb()).to_owned();
1153        create_dir_all("tests/kane/magnetic").expect("can't creat the file");
1154        write_txt_1(band, "tests/kane/magnetic/band.txt");
1155        write_txt(size, "tests/kane/magnetic/evec.txt");
1156        write_txt(show_str, "tests/kane/magnetic/structure.txt");
1157        //开始绘制角态
1158    }
1159
1160    #[test]
1161    fn Enonlinear() {
1162        //!arxiv 1706.07702
1163        //!用来测试非线性外在的霍尔电导
1164        let li: Complex<f64> = 1.0 * Complex::i();
1165        let delta = 0.;
1166        let t1 = 1.0 + 0.0 * li;
1167        let t2 = 0.2 * t1;
1168        let t3 = 0.2 * t1;
1169        let dim_r: usize = 3;
1170        let norb: usize = 2;
1171        let lat = arr2(&[
1172            [1.0, 0.0, 0.0],
1173            [0.5, 3.0_f64.sqrt() / 2.0, 0.0],
1174            [0.0, 0.0, 1.0],
1175        ]);
1176        let orb = arr2(&[[1.0 / 3.0, 1.0 / 3.0, 0.0], [2.0 / 3.0, 2.0 / 3.0, 0.0]]);
1177        let mut model = Model::tb_model(dim_r, lat, orb, false, None).unwrap();
1178        model.set_onsite(&arr1(&[delta, -delta]), SpinDirection::None);
1179        let R0: Array2<isize> = arr2(&[[0, 0, 0], [-1, 0, 0], [0, -1, 0]]);
1180        for (i, R) in R0.axis_iter(Axis(0)).enumerate() {
1181            let R = R.to_owned();
1182            model.set_hop(t1, 0, 1, &R, SpinDirection::None);
1183        }
1184        let R0: Array2<isize> = arr2(&[[1, 0, 1], [-1, 1, 1], [0, -1, 1]]);
1185        for (i, R) in R0.axis_iter(Axis(0)).enumerate() {
1186            let R = R.to_owned();
1187            model.set_hop(t2, 0, 0, &R, SpinDirection::None);
1188        }
1189        let R0: Array2<isize> = arr2(&[[1, 0, -1], [-1, 1, -1], [0, -1, -1]]);
1190        for (i, R) in R0.axis_iter(Axis(0)).enumerate() {
1191            let R = R.to_owned();
1192            model.set_hop(t2, 1, 1, &R, SpinDirection::None);
1193        }
1194        let R = arr1(&[0, 0, 1]);
1195        model.set_hop(t3, 0, 0, &R, SpinDirection::None);
1196        model.set_hop(t3, 1, 1, &R, SpinDirection::None);
1197        let path = array![
1198            [0.0, 0.0, 0.0],
1199            [1.0 / 3.0, 2.0 / 3.0, 0.0],
1200            [0.5, 0.5, 0.0],
1201            [0.0, 0.0, 0.0],
1202            [1.0 / 3.0, 2.0 / 3.0, 0.0],
1203            [1.0 / 3.0, 2.0 / 3.0, 0.5],
1204            [0.0, 0.0, 0.0],
1205            [0.0, 0.0, 0.5],
1206            [1.0 / 3.0, 2.0 / 3.0, 0.5],
1207            [0.5, 0.5, 0.5],
1208            [0.0, 0.0, 0.5]
1209        ];
1210        let label = vec!["G", "K", "M", "G", "K", "H", "G", "A", "H", "L", "A"];
1211        let nk = 101;
1212        model.show_band(&path, &label, nk, "tests/Enonlinear");
1213
1214        //开始计算非线性霍尔电导
1215        let dir_1 = arr1(&[1.0, 0.0, 0.0]);
1216        let dir_2 = arr1(&[0.0, 1.0, 0.0]);
1217        let dir_3 = arr1(&[0.0, 0.0, 1.0]);
1218        let nk: usize = 21;
1219        let kmesh = arr1(&[nk, nk, nk]);
1220        let E_min = -3.0;
1221        let E_max = 3.0;
1222        let E_n = 1000;
1223        let og = 0.0;
1224        let mu = Array1::linspace(E_min, E_max, E_n);
1225        let T = 30.0;
1226        let sigma: Array1<f64> = model
1227            .Nonlinear_Hall_conductivity_Extrinsic(
1228                &kmesh, &dir_1, &dir_2, &dir_3, &mu, T, og, 0, 1e-5,
1229            )
1230            .unwrap();
1231
1232        //开始绘制非线性电导
1233        let mut fg = Figure::new();
1234        let x: Vec<f64> = mu.to_vec();
1235        let axes = fg.axes2d();
1236        let y: Vec<f64> = sigma.to_vec();
1237        axes.lines(&x, &y, &[Color("black")]);
1238        axes.set_y_range(Fix(-10.0), Fix(10.0));
1239        axes.set_x_range(Fix(E_min), Fix(E_max));
1240        let mut show_ticks = Vec::<String>::new();
1241        let mut pdf_name = String::new();
1242        pdf_name.push_str("tests/Enonlinear");
1243        pdf_name.push_str("/nonlinear_ex.pdf");
1244        fg.set_terminal("pdfcairo", &pdf_name);
1245        fg.show();
1246
1247        let sigma: Array1<f64> = model
1248            .Nonlinear_Hall_conductivity_Intrinsic(&kmesh, &dir_1, &dir_2, &dir_3, &mu, T, 3)
1249            .unwrap();
1250        //开始绘制非线性电导
1251        let mut fg = Figure::new();
1252        let x: Vec<f64> = mu.to_vec();
1253        let axes = fg.axes2d();
1254        let y: Vec<f64> = sigma.to_vec();
1255        axes.lines(&x, &y, &[Color("black")]);
1256        axes.set_y_range(Fix(-10.0), Fix(10.0));
1257        axes.set_x_range(Fix(E_min), Fix(E_max));
1258        let mut show_ticks = Vec::<String>::new();
1259        let mut pdf_name = String::new();
1260        pdf_name.push_str("tests/Enonlinear");
1261        pdf_name.push_str("/nonlinear_in.pdf");
1262        fg.set_terminal("pdfcairo", &pdf_name);
1263        fg.show();
1264
1265        let (E0, dos) = model.dos(&kmesh, E_min, E_max, E_n, 1e-2).unwrap();
1266        //开始绘制dos
1267        let mut fg = Figure::new();
1268        let x: Vec<f64> = E0.to_vec();
1269        let axes = fg.axes2d();
1270        let y: Vec<f64> = dos.to_vec();
1271        axes.lines(&x, &y, &[Color("black")]);
1272        let mut show_ticks = Vec::<String>::new();
1273        let mut pdf_name = String::new();
1274        pdf_name.push_str("tests/Enonlinear");
1275        pdf_name.push_str("/dos.pdf");
1276        fg.set_terminal("pdfcairo", &pdf_name);
1277        fg.show();
1278    }
1279    #[test]
1280    fn kagome() {
1281        let li: Complex<f64> = 1.0 * Complex::i();
1282        let t1 = 1.0 + 0.0 * li;
1283        let t2 = 0.1 + 0.0 * li;
1284        let dim_r: usize = 2;
1285        let norb: usize = 2;
1286        let lat = arr2(&[[3.0_f64.sqrt(), -1.0], [3.0_f64.sqrt(), 1.0]]);
1287        let orb = arr2(&[[0.0, 0.0], [1.0 / 3.0, 0.0], [0.0, 1.0 / 3.0]]);
1288        let mut model = Model::tb_model(dim_r, lat, orb, false, None).unwrap();
1289        //最近邻hopping
1290        model.add_hop(t1, 0, 1, &array![0, 0], SpinDirection::None);
1291        model.add_hop(t1, 2, 0, &array![0, 0], SpinDirection::None);
1292        model.add_hop(t1, 1, 2, &array![0, 0], SpinDirection::None);
1293        model.add_hop(t1, 0, 2, &array![0, -1], SpinDirection::None);
1294        model.add_hop(t1, 0, 1, &array![-1, 0], SpinDirection::None);
1295        model.add_hop(t1, 2, 1, &array![-1, 1], SpinDirection::None);
1296        let nk: usize = 101;
1297        let path = [[0.0, 0.0], [2.0 / 3.0, 1.0 / 3.0], [0.5, 0.], [0.0, 0.0]];
1298        let path = arr2(&path);
1299        let label = vec!["G", "K", "M", "G"];
1300        model.show_band(&path, &label, nk, "tests/kagome/");
1301        //start to draw the band structure
1302        //Starting to calculate the edge state, first is the zigzag state
1303        let nk: usize = 101;
1304        let U = arr2(&[[1.0, 1.0], [-1.0, 1.0]]);
1305        let super_model = model.make_supercell(&U).unwrap();
1306        let zig_model = super_model.cut_piece(30, 0).unwrap();
1307        let path = [[0.0, 0.0], [0.0, 0.5], [0.0, 1.0]];
1308        let path = arr2(&path);
1309        let (k_vec, k_dist, k_node) = super_model.k_path(&path, nk).unwrap();
1310        let (eval, evec) = super_model.solve_all_parallel(&k_vec);
1311        let label = vec!["G", "M", "G"];
1312        zig_model.show_band(&path, &label, nk, "tests/kagome_zig/");
1313
1314        let green = surf_Green::from_Model(&super_model, 0, 1e-3, None).unwrap();
1315        let E_min = -2.0;
1316        let E_max = 4.0;
1317        let E_n = nk;
1318        let path = [[0.0], [0.5], [1.0]];
1319        let path = arr2(&path);
1320        let label = vec!["G", "M", "G"];
1321        green.show_surf_state("tests/kagome_zig", &path, &label, nk, E_min, E_max, E_n, 0);
1322
1323        //Starting to calculate the DOS of kagome
1324        let nk: usize = 51;
1325        let kmesh = arr1(&[nk, nk]);
1326        let E_min = -3.0;
1327        let E_max = 3.0;
1328        let E_n = 1000;
1329        let (E0, dos) = model.dos(&kmesh, E_min, E_max, E_n, 1e-2).unwrap();
1330        //start to show DOS
1331        let mut fg = Figure::new();
1332        let x: Vec<f64> = E0.to_vec();
1333        let axes = fg.axes2d();
1334        let y: Vec<f64> = dos.to_vec();
1335        axes.lines(&x, &y, &[Color("black")]);
1336        let mut show_ticks = Vec::<String>::new();
1337        let mut pdf_name = String::new();
1338        pdf_name.push_str("tests/kagome/");
1339        pdf_name.push_str("dos.pdf");
1340        fg.set_terminal("pdfcairo", &pdf_name);
1341        fg.show();
1342    }
1343
1344    #[test]
1345    fn SSH() {
1346        let li: Complex<f64> = 1.0 * Complex::i();
1347        let t1 = 1.0 + 0.0 * li;
1348        let t2 = 0.5 + 0.0 * li;
1349        let Delta = 0.0;
1350        let dim_r: usize = 1;
1351        let norb: usize = 2;
1352        let lat = arr2(&[[1.0]]);
1353        let orb = arr2(&[[0.3], [0.5]]);
1354        let mut model = Model::tb_model(dim_r, lat, orb, false, None).unwrap();
1355        model.add_hop(t1, 0, 1, &array![0], SpinDirection::None);
1356        model.add_hop(t2, 0, 1, &array![-1], SpinDirection::None);
1357        model.add_onsite(&array![Delta, -Delta], 0);
1358
1359        let nk: usize = 101;
1360        let path = [[0.0], [0.5], [1.0]];
1361        let path = arr2(&path);
1362        let label = vec!["G", "M", "G"];
1363        model.show_band(&path, &label, nk, "tests/SSH/");
1364        let mut super_model = model.cut_piece(5, 0).unwrap();
1365
1366        let (band, evec) = super_model.solve_onek(&array![0.0]);
1367        println!("{}", band);
1368    }
1369    #[test]
1370    fn BBH_model() {
1371        let li: Complex<f64> = 1.0 * Complex::i();
1372        let t1 = 0.1 + 0.0 * li;
1373        let t2 = 1.0 + 0.0 * li;
1374        let i0 = -1.0;
1375        let dim_r: usize = 2;
1376        let norb: usize = 2;
1377        let lat = arr2(&[[1.0, 0.0], [0.0, 1.0]]);
1378        let orb = arr2(&[[0.0, 0.0], [0.5, 0.0], [0.5, 0.5], [0.0, 0.5]]);
1379        let mut model = Model::tb_model(dim_r, lat, orb, false, None).unwrap();
1380        model.add_hop(t1, 0, 1, &array![0, 0], SpinDirection::None);
1381        model.add_hop(t1, 1, 2, &array![0, 0], SpinDirection::None);
1382        model.add_hop(t1, 2, 3, &array![0, 0], SpinDirection::None);
1383        model.add_hop(i0 * t1, 3, 0, &array![0, 0], SpinDirection::None);
1384        model.add_hop(t2, 0, 1, &array![-1, 0], SpinDirection::None);
1385        model.add_hop(i0 * t2, 0, 3, &array![0, -1], SpinDirection::None);
1386        model.add_hop(t2, 2, 3, &array![1, 0], SpinDirection::None);
1387        model.add_hop(t2, 2, 1, &array![0, 1], SpinDirection::None);
1388        let nk: usize = 101;
1389        let path = [[0.0, 0.0], [0.5, 0.0], [0.5, 0.5], [0.0, 0.0]];
1390        let path = arr2(&path);
1391        let label = vec!["G", "X", "M", "G"];
1392        model.show_band(&path, &label, nk, "tests/BBH/").unwrap();
1393        model.output_hr("tests/BBH/", "wannier90");
1394
1395        //算一下wilson loop
1396        let n = 51;
1397        let dir_1 = arr1(&[1.0, 0.0]);
1398        let dir_2 = arr1(&[0.0, 1.0]);
1399        let occ = vec![0, 1];
1400        let wcc = model.wannier_centre(&occ, &array![0.0, 0.0], &dir_1, &dir_2, n, n);
1401        let nocc = occ.len();
1402        let mut fg = Figure::new();
1403        let x: Vec<f64> = Array1::<f64>::linspace(0.0, 1.0, n).to_vec();
1404        let axes = fg.axes2d();
1405        for j in -1..2 {
1406            for i in 0..nocc {
1407                let a = wcc.row(i).to_owned() + (j as f64) * 2.0 * PI;
1408                let y: Vec<f64> = a.to_vec();
1409                axes.points(&x, &y, &[Color("black"), gnuplot::PointSymbol('O')]);
1410            }
1411        }
1412        let axes = axes.set_x_range(Fix(0.0), Fix(1.0));
1413        let axes = axes.set_y_range(Fix(0.0), Fix(2.0 * PI));
1414        let show_ticks = vec![
1415            Major(0.0, Fix("0")),
1416            Major(0.5, Fix("π")),
1417            Major(1.0, Fix("2π")),
1418        ];
1419        axes.set_x_ticks_custom(show_ticks.into_iter(), &[], &[]);
1420        let show_ticks = vec![
1421            Major(0.0, Fix("0")),
1422            Major(PI, Fix("π")),
1423            Major(2.0 * PI, Fix("2π")),
1424        ];
1425        axes.set_y_ticks_custom(show_ticks.into_iter(), &[], &[]);
1426        let mut pdf_name = String::new();
1427        pdf_name.push_str("tests/BBH/wcc.pdf");
1428        fg.set_terminal("pdfcairo", &pdf_name);
1429        fg.show();
1430        //算一下边界态
1431        let green = surf_Green::from_Model(&model, 0, 1e-3, None).unwrap();
1432        let E_min = -2.0;
1433        let E_max = 2.0;
1434        let E_n = nk;
1435        let path = [[0.0], [0.5], [1.0]];
1436        let path = arr2(&path);
1437        let label = vec!["G", "X", "G"];
1438        green.show_surf_state("tests/BBH", &path, &label, nk, E_min, E_max, E_n, 0);
1439
1440        //算一下corner state
1441        let num = 10;
1442        let model_1 = model.cut_piece(num, 0).unwrap();
1443        let new_model = model_1.cut_piece(2 * num, 1).unwrap();
1444        let mut s = 0;
1445        let start = Instant::now();
1446        let (band, evec) = new_model.solve_onek(&arr1(&[0.0, 0.0]));
1447        println!(
1448            "band shape is {:?}, evec shape is {:?}",
1449            band.shape(),
1450            evec.shape()
1451        );
1452        let end = Instant::now(); // 结束计时
1453        let duration = end.duration_since(start); // 计算执行时间
1454        println!("solve_band_all took {} seconds", duration.as_secs_f64()); // 输出执行时间
1455        let nresults = band.len();
1456        let show_evec = evec.to_owned().map(|x| x.norm_sqr());
1457        let norb = new_model.norb();
1458        let size = show_evec;
1459        let show_str = new_model.atom_position().dot(&model.lat);
1460        create_dir_all("tests/BBH/corner").expect("can't creat the file");
1461        write_txt_1(band, "tests/BBH/corner/band.txt");
1462        write_txt(size, "tests/BBH/corner/evec.txt");
1463        write_txt(show_str, "tests/BBH/corner/structure.txt");
1464    }
1465
1466    #[test]
1467    fn graphene_magnetic_field() {
1468        use crate::{MagneticField, Model, SpinDirection};
1469        use ndarray::{Axis, arr1, arr2};
1470        use num_complex::Complex;
1471        // 如果你在其他地方定义了画图函数,请确保 use 进来,例如:
1472        // use crate::draw_heatmap;
1473
1474        // 1. 设置模型基本参数
1475        let t = Complex::new(-1.0, 0.0);
1476        let delta = 0.0;
1477        let dim_r: usize = 2;
1478        let norb: usize = 2;
1479
1480        // 石墨烯晶格:a1 = (1, 0), a2 = (1/2, √3/2)
1481        let lat = arr2(&[[1.0, 0.0], [0.5, 3.0_f64.sqrt() / 2.0]]);
1482        // 轨道的相对分数坐标 (Fractional Coordinates)
1483        let orb = arr2(&[[1.0 / 3.0, 1.0 / 3.0], [2.0 / 3.0, 2.0 / 3.0]]);
1484
1485        let mut model = Model::tb_model(dim_r, lat, orb, false, None).unwrap();
1486        model.set_onsite(&arr1(&[-delta, delta]), SpinDirection::None);
1487
1488        // 添加最近邻跃迁
1489        let r0: ndarray::Array2<isize> = arr2(&[[0, 0], [-1, 0], [0, -1]]);
1490        for r in r0.axis_iter(Axis(0)) {
1491            model.add_hop(t, 0, 1, &r.to_owned(), SpinDirection::None);
1492        }
1493
1494        // 2. 施加磁场
1495        // 重要:二维中,面外磁场垂直于 xy 平面,对应的索引必须为 z 轴 (即 mag_dir = 2)
1496        // 扩胞 [3, 3] 表示 a1 和 a2 两个方向各扩胞 3 倍
1497        // 磁通 8 表示整体超胞含有 8 个磁通量子 (每个原胞 8/9 个磁通)
1498        let magnetic_model = model.add_magnetic_field(2, [9, 9], 40).unwrap();
1499
1500        // 3. 高对称路径 (注意:六角晶格的 K 点在分数倒格矢坐标下是 1/3, 1/3)
1501        let path = arr2(&[[0.0, 0.0], [0.5, 0.0], [1.0 / 3.0, 1.0 / 3.0], [0.0, 0.0]]);
1502        let label = vec!["Γ", "M", "K", "Γ"];
1503        let nk = 1001;
1504
1505        // 4. 绘制折叠态下的超胞能带 (Hofstadter 蝴蝶状能带切片)
1506        magnetic_model
1507            .show_band(&path, &label, nk, "tests/graphene_magnetic")
1508            .unwrap();
1509
1510        // 5. 展开能带 (Unfold) 回到原胞 Brillouin 区
1511        let u_matrix = arr2(&[[9.0, 0.0], [0.0, 9.0]]);
1512
1513        // 生成对应的谱函数 (Spectral Weight)
1514        let a_spectral = magnetic_model
1515            .unfold(&u_matrix, &path, nk, -3.0, 3.0, nk, 1e-3, 1e-5)
1516            .unwrap();
1517
1518        // 将谱函数输出为热力图
1519        // (假定 draw_heatmap 接收二维热力矩阵以及保存路径)
1520        draw_heatmap(
1521            &a_spectral.reversed_axes(),
1522            "./tests/graphene_magnetic/unfold_band.pdf",
1523        );
1524    }
1525
1526    #[test]
1527    fn test_hofstadter_butterfly_gnuplot() {
1528        use gnuplot::AutoOption::Fix;
1529        use gnuplot::{AxesCommon, Color, Figure, PointSize, PointSymbol};
1530
1531        // 1. 初始化一个标准的 2D 正方晶格
1532        let t = Complex::new(-1.0, 0.0);
1533        let dim_r = 2;
1534        // 正方晶格基矢 a1=(1,0), a2=(0,1)
1535        let lat = arr2(&[[1.0, 0.0], [0.0, 1.0]]);
1536        // 单轨道坐标位于原点
1537        let orb = arr2(&[[0.0, 0.0]]);
1538
1539        let mut model = Model::tb_model(dim_r, lat, orb, false, None).unwrap();
1540        model.set_onsite(&arr1(&[0.0]), SpinDirection::None);
1541
1542        // 加上最近邻跃迁 (上下左右4个方向)
1543        let r0: ndarray::Array2<isize> = arr2(&[[1, 0], [-1, 0], [0, 1], [0, -1]]);
1544        for r in r0.axis_iter(Axis(0)) {
1545            model.add_hop(t, 0, 0, &r.to_owned(), SpinDirection::None);
1546        }
1547
1548        // 2. 准备扫描磁通并记录坐标
1549        let q = 81; // 分母 q 设定为 49 形成的蝴蝶分形已经足够好看
1550
1551        // 我们用两个动态数组分别记录散点图的 X 轴 (磁通) 和 Y 轴 (能量)
1552        let mut x_data = Vec::new();
1553        let mut y_data = Vec::new();
1554
1555        println!("开始计算 Hofstadter 蝴蝶能谱,进度: ");
1556
1557        for p in 0..=q {
1558            // 在 y 轴方向扩胞 q 倍,即 [1, q],形成一个一维超胞
1559            let mag_model = model.add_magnetic_field(2, [1, q], p as isize).unwrap();
1560            let norb = mag_model.norb();
1561
1562            // 提取 Gamma 点 k = (0,0) 的哈密顿量矩阵
1563            let mut h_k0 = Array2::<Complex<f64>>::zeros((norb, norb));
1564            for iR in 0..mag_model.hamR.nrows() {
1565                for i in 0..norb {
1566                    for j in 0..norb {
1567                        h_k0[[i, j]] += mag_model.ham[[iR, i, j]];
1568                    }
1569                }
1570            }
1571
1572            // --- 求解本征值 ---
1573            let evals = h_k0.eigvalsh(UPLO::Upper).expect("矩阵对角化失败");
1574
1575            // 收集坐标点
1576            let flux_ratio = (p as f64) / (q as f64);
1577            for &e in evals.iter() {
1578                x_data.push(flux_ratio);
1579                y_data.push(e);
1580            }
1581
1582            if p % 10 == 0 {
1583                println!("已完成 {}/{}", p, q);
1584            }
1585        }
1586
1587        println!(
1588            "计算完成,共有 {} 个能级点,正在使用 gnuplot 绘图...",
1589            x_data.len()
1590        );
1591
1592        // 3. 使用 Gnuplot 直接输出图像
1593        // 确保目录存在
1594        create_dir_all("tests").expect("无法创建 tests 文件夹");
1595
1596        let mut fg = Figure::new();
1597        let axes = fg.axes2d();
1598
1599        axes.set_title("Hofstadter's Butterfly", &[]);
1600        axes.set_x_label("Magnetic Flux (\\Phi / \\Phi_0)", &[]);
1601        axes.set_y_label("Energy (E/t)", &[]);
1602
1603        // 固定 x 和 y 的坐标范围
1604        let axes = axes.set_x_range(Fix(0.0), Fix(1.0));
1605        let axes = axes.set_y_range(Fix(-10.0), Fix(10.0));
1606
1607        // 使用 .points 绘制散点图
1608        // PointSymbol('.') 表示画极小的像素点,最适合用来展示密集的分形结构
1609        // Color("navy") 使用深蓝色
1610        axes.points(
1611            &x_data,
1612            &y_data,
1613            &[Color("navy"), PointSymbol('.'), PointSize(0.6)],
1614        );
1615
1616        // 将图像渲染为 PDF
1617        fg.set_terminal("pdfcairo", "tests/hofstadter_butterfly.pdf");
1618        fg.show().expect("Gnuplot 画图失败");
1619
1620        println!("完美!图像已保存至 tests/hofstadter_butterfly.pdf");
1621    }
1622}