1use 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#[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 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 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 ]); 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 ]); #[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 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 #[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 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 ]); 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 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 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 #[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 (mut v, hamk): (Array3<Complex<f64>>, Array2<Complex<f64>>) =
408 self.gen_v(&k_vec, Gauge::Atom); 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 UU[[i, j]] = 1.0 / (a.powi(2) + eta.powi(2));
521 }
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 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 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 #[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 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 let conductivity: f64 = omega.sum() / (nk as f64) / self.lat.det().unwrap();
632 Ok(conductivity)
633 }
634 #[allow(non_snake_case)]
635 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)?; 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 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 let li: Complex<f64> = 1.0 * Complex::i();
755 let (mut v, hamk): (Array3<Complex<f64>>, Array2<Complex<f64>>) =
757 self.gen_v(&k_vec, Gauge::Atom); let mut J: Array3<Complex<f64>> = v.clone();
759 let mut v0 = Array2::<Complex<f64>>::zeros((self.nsta(), self.nsta())); 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; 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 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: Array1<f64> = omega_n * partial_ve;
847 (omega_n, band) }
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 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 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 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 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 let (mut v, hamk): (Array3<Complex<f64>>, Array2<Complex<f64>>) =
974 self.gen_v(&k_vec, Gauge::Atom); let mut J = v.clone();
976 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()))); v.slice_mut(s![i, .., ..]).assign(&v_s); }
991 let mut v_1 = Array2::<Complex<f64>>::zeros((self.nsta(), self.nsta())); 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 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 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 if self.spin {
1022 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()))); S.slice_mut(s![i, .., ..]).assign(&v0);
1056 }
1057 let mut s_1 = Array2::<Complex<f64>>::zeros((self.nsta(), self.nsta())); 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 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 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 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)); }
1097 let partial_G: Array1<f64> = {
1099 let mut A = Array1::<Complex<f64>>::zeros(self.nsta()); 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()); 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()); 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 return (
1136 (partial_s_1 * G_23 - partial_ve_2 * G_13_h),
1137 band,
1138 Some(partial_G),
1139 );
1140 } else {
1141 let G_23: Array1<f64> = {
1144 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 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 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 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 panic!("the code can not support for T=0");
1341 }
1342 Ok(conductivity)
1343 }
1344}