forked from wannier-developers/wannier90
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathtransport.F90
4532 lines (3804 loc) · 176 KB
/
transport.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!-*- mode: F90 -*-!
!------------------------------------------------------------!
! This file is distributed as part of the Wannier90 code and !
! under the terms of the GNU General Public License. See the !
! file `LICENSE' in the root directory of the Wannier90 !
! distribution, or http://www.gnu.org/copyleft/gpl.txt !
! !
! The webpage of the Wannier90 code is www.wannier.org !
! !
! The Wannier90 code is hosted on GitHub: !
! !
! https://github.com/wannier-developers/wannier90 !
!------------------------------------------------------------!
!
!-----------------------------------------------------------------------!
! Based on !
! < dosqc_1.0 > !
! Density Of States and Quantum Conductance - Version 1.0 !
! Marco Buongiorno Nardelli, January 2000. !
! !
! Reference: !
! - M. Buongiorno Nardelli, "Electronic transport in extended systems: !
! application to carbon nanotubes", Phys. Rev. B, vol. 60(11), 7828 !
! (1999) !
!-----------------------------------------------------------------------!
!=======================================================================!
! Definition of parameters used in w90_transport !
!=======================================================================!
! !
! transport_mode = 'bulk' or 'lcr' !
! tran_win_min = minimum E !
! tran_win_max = maximum E !
! tran_energy_step = delta E !
! tran_num_bb = # of WFs in a principal layer of a perfectly !
! periodic bulk system !
! tran_num_ll = # of WFs in a principal layer of a left lead !
! tran_num_rr = # of WFs in a principal layer of a right lead !
! tran_num_cc = # of WFs in a disordered conducter cell !
! tran_num_lc = # of WFs in a disordered conducter cell that !
! are used to calculate interaction with a !
! left lead !
! tran_num_cr = # of WFs in a disordered conducter cell that !
! are used to calculate interaction with a !
! right lead !
! tran_num_bandc = width of band-diagonal hC matrix !
! tran_read_ht = .true. => read H matrix from h*.dat files !
! tran_write_ht = .true. => write H matrix from h*.dat files !
! tran_use_same_lead = .true. => in L-C-R construction, left and !
! right lead are the same kind !
! tran_num_cell_ll = # of unit cells in a PL of left lead !
! tran_num_cell_rr = # of unit cells in a PL of right lead !
! (equal to tran_num_cell_ll for now) !
! tran_group_threshold = distance defining the grouping of WFs !
!=======================================================================!
module w90_transport_mod
!! Module to handle ballistic transport.
!! Based on
!! < dosqc_1.0 >
!! Density Of States and Quantum Conductance - Version 1.0
!! Marco Buongiorno Nardelli, January 2000.
!! Reference:
!! - M. Buongiorno Nardelli, "Electronic transport in extended systems:
!! application to carbon nanotubes", Phys. Rev. B, vol. 60(11), 7828
!! (1999)
use w90_constants, only: dp
use w90_error
use w90_comms
implicit none
private
complex(kind=dp), parameter :: eta = (0.0_dp, 0.0005_dp) !! small complex number
integer, parameter :: nterx = 50 !! nterx = # of maximum iteration to calculate transfer matrix
public :: tran_main
contains
!================================================!
subroutine tran_main(atom_data, dis_manifold, fermi_energy_list, ham_logical, kpt_latt, &
output_file, real_space_ham, transport, print_output, wannier_data, &
ws_region, w90_calculation, ham_k, ham_r, u_matrix, u_matrix_opt, eigval, &
real_lattice, wannier_centres_translated, irvec, mp_grid, ndegen, &
shift_vec, nrpts, num_bands, num_kpts, num_wann, rpt_origin, &
bands_plot_mode, have_disentangled, lsitesymmetry, seedname, stdout, &
timer, error, comm)
!================================================!
!
!! Main transport subroutine
!
!================================================!
use w90_io, only: io_stopwatch_start, io_stopwatch_stop
use w90_error, only: w90_error_type, set_error_dealloc
use w90_hamiltonian, only: hamiltonian_get_hr, hamiltonian_write_hr, hamiltonian_setup
use w90_types, only: wannier_data_type, print_output_type, ws_region_type, &
atom_data_type, dis_manifold_type, timer_list_type
use w90_wannier90_types, only: w90_calculation_type, transport_type, output_file_type, &
real_space_ham_type, ham_logical_type
implicit none
! arguments
type(transport_type), intent(inout) :: transport
type(real_space_ham_type), intent(inout) :: real_space_ham
type(ws_region_type), intent(inout) :: ws_region
type(print_output_type), intent(in) :: print_output
type(w90_calculation_type), intent(in) :: w90_calculation
type(output_file_type), intent(in) :: output_file
type(wannier_data_type), intent(in) :: wannier_data
type(atom_data_type), intent(in) :: atom_data
type(dis_manifold_type), intent(in) :: dis_manifold
real(kind=dp), intent(in) :: kpt_latt(:, :)
real(kind=dp), allocatable, intent(in) :: fermi_energy_list(:)
type(ham_logical_type), intent(inout) :: ham_logical
type(timer_list_type), intent(inout) :: timer
type(w90_error_type), allocatable, intent(out) :: error
type(w90_comm_type), intent(in) :: comm
integer, intent(inout) :: rpt_origin
integer, intent(inout) :: nrpts
integer, intent(inout), allocatable :: ndegen(:)
integer, intent(inout), allocatable :: shift_vec(:, :)
integer, intent(inout), allocatable :: irvec(:, :)
integer, intent(in) :: num_wann
integer, intent(in) :: num_bands
integer, intent(in) :: num_kpts
integer, intent(in) :: stdout
integer, intent(in) :: mp_grid(3)
real(kind=dp), intent(inout), allocatable :: wannier_centres_translated(:, :)
real(kind=dp), intent(in) :: real_lattice(3, 3)
real(kind=dp), intent(in) :: eigval(:, :)
complex(kind=dp), intent(in) :: u_matrix(:, :, :)
complex(kind=dp), intent(in) :: u_matrix_opt(:, :, :)
complex(kind=dp), allocatable, intent(inout) :: ham_k(:, :, :)
complex(kind=dp), intent(inout), allocatable :: ham_r(:, :, :)
character(len=*), intent(in) :: bands_plot_mode
character(len=50), intent(in) :: seedname
logical, intent(in) :: have_disentangled
logical, intent(in) :: lsitesymmetry !YN:
! local variables
integer :: one_dim_vec
!! cartesian axis to which real_lattice(:,one_dim_vec) is parallel
integer :: nrpts_one_dim
integer :: num_pl
!! number of unit cell in a principal layer
integer :: coord(3)
!! coord : coord(1) defines the conduction direction according to
!1=x,2=y,3=z
!! coord(2),coord(3) define the other directions during sorting routines
integer, allocatable :: tran_sorted_idx(:)
!! index of sorted WF centres to unsorted
integer :: num_G
integer :: irvec_max ! size of hr_one_dim's last dimension (:, :, -irvec_max:irvec_max)
integer :: ierr
real(kind=dp), allocatable :: hr_one_dim(:, :, :)
real(kind=dp), allocatable :: signatures(:, :)
real(kind=dp), allocatable :: hB0(:, :)
real(kind=dp), allocatable :: hB1(:, :)
real(kind=dp), allocatable :: hC(:, :)
real(kind=dp), allocatable :: hCR(:, :)
real(kind=dp), allocatable :: hL0(:, :)
real(kind=dp), allocatable :: hL1(:, :)
real(kind=dp), allocatable :: hLC(:, :)
real(kind=dp), allocatable :: hR0(:, :)
real(kind=dp), allocatable :: hR1(:, :)
logical :: pl_warning
! fixme jj printout guards as elsewhere please (even if not yet parallel)
if (print_output%timing_level > 0) call io_stopwatch_start('tran: main', timer)
write (stdout, '(/1x,a)') '*---------------------------------------------------------------------------*'
write (stdout, '(1x,a)') '| TRANSPORT |'
write (stdout, '(1x,a)') '*---------------------------------------------------------------------------*'
write (stdout, *)
if (index(transport%mode, 'bulk') > 0) then
write (stdout, '(/1x,a/)') 'Calculation of Quantum Conductance and DoS: bulk mode'
if (.not. transport%read_ht) then
call hamiltonian_setup(ham_logical, print_output, ws_region, w90_calculation, ham_k, &
ham_r, real_lattice, wannier_centres_translated, irvec, mp_grid, &
ndegen, num_kpts, num_wann, nrpts, rpt_origin, bands_plot_mode, &
stdout, timer, error, transport%mode, comm)
if (allocated(error)) return
call hamiltonian_get_hr(atom_data, dis_manifold, ham_logical, real_space_ham, &
print_output, ham_k, ham_r, u_matrix, u_matrix_opt, eigval, &
kpt_latt, real_lattice, wannier_data%centres, &
wannier_centres_translated, irvec, shift_vec, nrpts, num_bands, &
num_kpts, num_wann, have_disentangled, stdout, timer, error, &
lsitesymmetry, comm)
if (allocated(error)) return
if (output_file%write_hr) then
! is this redundant after the plot call?
call hamiltonian_write_hr(ham_r, irvec, ndegen, nrpts, num_wann, &
print_output%timing_level, seedname, timer, error, comm)
endif
if (allocated(error)) return
call tran_reduce_hr(real_space_ham, ham_r, hr_one_dim, real_lattice, irvec, mp_grid, &
irvec_max, nrpts, nrpts_one_dim, num_wann, one_dim_vec, &
print_output%timing_level, stdout, timer, error, comm)
if (allocated(error)) return
call tran_cut_hr_one_dim(real_space_ham, transport, print_output, hr_one_dim, &
real_lattice, wannier_centres_translated, mp_grid, irvec_max, &
num_pl, num_wann, one_dim_vec, stdout, timer)
call tran_get_ht(fermi_energy_list, transport, hB0, hB1, hr_one_dim, irvec_max, num_pl, &
num_wann, print_output%timing_level, seedname, timer, error, comm)
if (allocated(error)) return
if (output_file%write_xyz) then
call tran_write_xyz(atom_data, transport, wannier_centres_translated, tran_sorted_idx, &
num_wann, seedname, stdout)
endif
end if
call tran_bulk(transport, hB0, hB1, print_output%timing_level, stdout, seedname, timer, error, comm)
if (allocated(error)) return
end if
if (index(transport%mode, 'lcr') > 0) then
write (stdout, '(/1x,a/)') 'Calculation of Quantum Conductance and DoS: lead-conductor-lead mode'
if (.not. transport%read_ht) then
call hamiltonian_setup(ham_logical, print_output, ws_region, w90_calculation, ham_k, &
ham_r, real_lattice, wannier_centres_translated, irvec, mp_grid, &
ndegen, num_kpts, num_wann, nrpts, rpt_origin, bands_plot_mode, &
stdout, timer, error, transport%mode, comm)
if (allocated(error)) return
call hamiltonian_get_hr(atom_data, dis_manifold, ham_logical, real_space_ham, &
print_output, ham_k, ham_r, u_matrix, u_matrix_opt, eigval, &
kpt_latt, real_lattice, wannier_data%centres, &
wannier_centres_translated, irvec, shift_vec, nrpts, num_bands, &
num_kpts, num_wann, have_disentangled, stdout, timer, error, &
lsitesymmetry, comm)
if (allocated(error)) return
if (output_file%write_hr) then
call hamiltonian_write_hr(ham_r, irvec, ndegen, nrpts, num_wann, &
print_output%timing_level, seedname, timer, error, comm)
if (allocated(error)) return
endif
call tran_reduce_hr(real_space_ham, ham_r, hr_one_dim, real_lattice, irvec, mp_grid, &
irvec_max, nrpts, nrpts_one_dim, num_wann, one_dim_vec, &
print_output%timing_level, stdout, timer, error, comm)
if (allocated(error)) return
call tran_cut_hr_one_dim(real_space_ham, transport, print_output, hr_one_dim, &
real_lattice, wannier_centres_translated, mp_grid, irvec_max, &
num_pl, num_wann, one_dim_vec, stdout, timer)
write (stdout, *) '------------------------- 2c2 Calculation Type: ------------------------------'
write (stdout, *) ' '
call tran_find_integral_signatures(signatures, num_G, print_output, real_lattice, &
u_matrix_opt, u_matrix, num_bands, num_wann, &
have_disentangled, wannier_centres_translated, stdout, &
seedname, timer, error, comm)
if (allocated(error)) return
call tran_lcr_2c2_sort(signatures, num_G, pl_warning, transport, atom_data, wannier_data, &
real_space_ham, print_output, real_lattice, num_wann, mp_grid, &
ham_r, irvec, nrpts, wannier_centres_translated, one_dim_vec, &
nrpts_one_dim, num_pl, coord, tran_sorted_idx, hr_one_dim, &
irvec_max, output_file%write_xyz, stdout, seedname, timer, error, comm)
if (allocated(error)) return
if (output_file%write_xyz) call tran_write_xyz(atom_data, transport, &
wannier_centres_translated, &
tran_sorted_idx, num_wann, seedname, stdout)
call tran_parity_enforce(signatures, print_output, transport, num_wann, tran_sorted_idx, &
hr_one_dim, irvec_max, stdout, timer)
call tran_lcr_2c2_build_ham(pl_warning, real_space_ham, fermi_energy_list, kpt_latt, &
num_wann, transport, print_output, real_lattice, mp_grid, &
ham_r, irvec, nrpts, wannier_centres_translated, one_dim_vec, &
nrpts_one_dim, num_pl, coord, tran_sorted_idx, hC, hCR, hL0, &
hL1, hLC, hR0, hR1, hr_one_dim, irvec_max, stdout, seedname, &
timer, error, comm)
if (allocated(error)) return
endif
call tran_lcr(transport, hC, hCR, hL0, hL1, hLC, hR0, hR1, print_output%timing_level, &
stdout, seedname, timer, error, comm)
if (allocated(error)) return
end if
if (print_output%timing_level > 0) call io_stopwatch_stop('tran: main', timer)
if (allocated(hR1)) then
deallocate (hR1, stat=ierr)
if (ierr /= 0) then
call set_error_dealloc(error, 'Error in deallocating hR1 in tran_main', comm)
return
endif
end if
if (allocated(hR0)) then
deallocate (hR0, stat=ierr)
if (ierr /= 0) then
call set_error_dealloc(error, 'Error in deallocating hR0 in tran_main', comm)
return
endif
end if
if (allocated(hL1)) then
deallocate (hL1, stat=ierr)
if (ierr /= 0) then
call set_error_dealloc(error, 'Error in deallocating hL1 in tran_main', comm)
return
endif
end if
if (allocated(hB1)) then
deallocate (hB1, stat=ierr)
if (ierr /= 0) then
call set_error_dealloc(error, 'Error in deallocating hB1 in tran_main', comm)
return
endif
end if
if (allocated(hB0)) then
deallocate (hB0, stat=ierr)
if (ierr /= 0) then
call set_error_dealloc(error, 'Error in deallocating hB0 in tran_main', comm)
return
endif
end if
if (allocated(hr_one_dim)) then
deallocate (hr_one_dim, stat=ierr)
if (ierr /= 0) then
call set_error_dealloc(error, 'Error in deallocating hr_one_dim in tran_main', comm)
return
endif
end if
end subroutine tran_main
!================================================!
subroutine tran_reduce_hr(real_space_ham, ham_r, hr_one_dim, real_lattice, irvec, mp_grid, &
irvec_max, nrpts, nrpts_one_dim, num_wann, one_dim_vec, timing_level, &
stdout, timer, error, comm)
!================================================!
!
! reduce ham_r from 3-d to 1-d
!
!================================================!
use w90_constants, only: dp, eps8
use w90_io, only: io_stopwatch_start, io_stopwatch_stop
use w90_wannier90_types, only: real_space_ham_type
use w90_error, only: w90_error_type, set_error_alloc, set_error_fatal
use w90_types, only: timer_list_type
implicit none
! passed vars
type(real_space_ham_type), intent(in) :: real_space_ham
type(timer_list_type), intent(inout) :: timer
type(w90_error_type), allocatable, intent(out) :: error
type(w90_comm_type), intent(in) :: comm
integer, intent(in) :: irvec(:, :)
integer, intent(inout) :: irvec_max ! limits of hr_one_dim final dim
integer, intent(in) :: mp_grid(3)
integer, intent(in) :: nrpts
integer, intent(in) :: num_wann
integer, intent(inout) :: nrpts_one_dim
integer, intent(inout) :: one_dim_vec
integer, intent(in) :: timing_level
integer, intent(in) :: stdout
real(kind=dp), allocatable, intent(inout) :: hr_one_dim(:, :, :)
real(kind=dp), intent(in) :: real_lattice(3, 3)
complex(kind=dp), intent(in) :: ham_r(:, :, :)
! local variables
integer :: ierr
integer :: irvec_tmp(3), two_dim_vec(2)
integer :: i, j
integer :: i1, i2, i3, n1, nrpts_tmp, loop_rpt
if (timing_level > 1) call io_stopwatch_start('tran: reduce_hr', timer)
! Find one_dim_vec which is parallel to one_dim_dir
! two_dim_vec - the other two lattice vectors
j = 0
do i = 1, 3
if (abs(abs(real_lattice(real_space_ham%one_dim_dir, i)) &
- sqrt(dot_product(real_lattice(:, i), real_lattice(:, i)))) .lt. eps8) then
one_dim_vec = i
j = j + 1
end if
end do
if (j .ne. 1) then
write (stdout, '(i3,a)') j, ' : 1-D LATTICE VECTOR NOT DEFINED'
call set_error_fatal(error, 'Error: 1-d lattice vector not defined in tran_reduce_hr', comm)
return
end if
j = 0
do i = 1, 3
if (i .ne. one_dim_vec) then
j = j + 1
two_dim_vec(j) = i
end if
end do
! starting H matrix should include all W-S supercell where
! the center of the cell spans the full space of the home cell
! adding one more buffer layer when mp_grid(one_dim_vec) is an odd number
!irvec_max = (mp_grid(one_dim_vec)+1)/2
irvec_tmp = maxval(irvec, DIM=2) + 1
irvec_max = irvec_tmp(one_dim_vec)
nrpts_one_dim = 2*irvec_max + 1
allocate (hr_one_dim(num_wann, num_wann, -irvec_max:irvec_max), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating hr_one_dim in tran_reduce_hr', comm)
return
endif
hr_one_dim = 0.0_dp
! check imaginary part
write (stdout, '(1x,a,F12.6)') 'Maximum imaginary part of the real-space Hamiltonian: ', &
maxval(abs(aimag(ham_r)))
! select a subset of ham_r, where irvec is 0 along the two other lattice vectors
nrpts_tmp = 0
loop_n1: do n1 = -irvec_max, irvec_max
do loop_rpt = 1, nrpts
i1 = mod(n1 - irvec(one_dim_vec, loop_rpt), mp_grid(one_dim_vec))
i2 = irvec(two_dim_vec(1), loop_rpt)
i3 = irvec(two_dim_vec(2), loop_rpt)
if (i1 .eq. 0 .and. i2 .eq. 0 .and. i3 .eq. 0) then
nrpts_tmp = nrpts_tmp + 1
hr_one_dim(:, :, n1) = real(ham_r(:, :, loop_rpt), dp)
cycle loop_n1
end if
end do
end do loop_n1
if (nrpts_tmp .ne. nrpts_one_dim) then
write (stdout, '(a)') 'FAILED TO EXTRACT 1-D HAMILTONIAN'
call set_error_fatal(error, 'Error: cannot extract 1d hamiltonian in tran_reduce_hr', comm)
return
end if
if (timing_level > 1) call io_stopwatch_stop('tran: reduce_hr', timer)
return
end subroutine tran_reduce_hr
!================================================!
subroutine tran_cut_hr_one_dim(real_space_ham, transport, print_output, hr_one_dim, &
real_lattice, wannier_centres_translated, mp_grid, irvec_max, &
num_pl, num_wann, one_dim_vec, stdout, timer)
!================================================!
use w90_constants, only: dp
use w90_io, only: io_stopwatch_start, io_stopwatch_stop
use w90_types, only: print_output_type, timer_list_type
use w90_wannier90_types, only: transport_type, real_space_ham_type
implicit none
! arguments
type(real_space_ham_type), intent(inout) :: real_space_ham
type(print_output_type), intent(in) :: print_output
type(transport_type), intent(inout) :: transport
type(timer_list_type), intent(inout) :: timer
integer, intent(in) :: mp_grid(3)
integer, intent(in) :: irvec_max
integer, intent(in) :: num_wann
integer, intent(in) :: one_dim_vec
integer, intent(inout) :: num_pl
integer, intent(in) :: stdout
real(kind=dp), intent(inout) :: hr_one_dim(:, :, -irvec_max:)
real(kind=dp), intent(in) :: real_lattice(3, 3)
real(kind=dp), intent(in) :: wannier_centres_translated(:, :)
! local variables
integer :: i, j, n1
real(kind=dp) :: hr_max
real(kind=dp) :: dist
real(kind=dp) :: dist_vec(3)
real(kind=dp) :: dist_ij_vec(3)
real(kind=dp) :: shift_vec(3, -irvec_max:irvec_max)
real(kind=dp) :: hr_tmp(num_wann, num_wann)
if (print_output%timing_level > 1) call io_stopwatch_start('tran: cut_hr_one_dim', timer)
!irvec_max = nrpts_one_dim/2 ! now passed as arg
! maximum possible dist_cutoff
dist = real(mp_grid(one_dim_vec), dp)*abs(real_lattice(real_space_ham%one_dim_dir, one_dim_vec)) &
/2.0_dp
if (real_space_ham%dist_cutoff .gt. dist) then
write (stdout, '(1x,a,1x,F10.5,1x,a)') 'dist_cutoff', real_space_ham%dist_cutoff, &
trim(print_output%length_unit), 'is too large'
real_space_ham%dist_cutoff = dist
! aam_2012-04-13
real_space_ham%dist_cutoff_hc = dist
write (stdout, '(4x,a,1x,F10.5,1x,a)') 'reset to', real_space_ham%dist_cutoff, &
trim(print_output%length_unit)
end if
do n1 = -irvec_max, irvec_max
shift_vec(:, n1) = real(n1, dp)*(real_lattice(:, one_dim_vec))
! write(stdout,'(a,3f10.6)') 'shift_vec', shift_vec(:,n1)
end do
! apply dist_cutoff first
if (index(real_space_ham%dist_cutoff_mode, 'one_dim') > 0) then
do i = 1, num_wann
do j = 1, num_wann
dist_ij_vec(real_space_ham%one_dim_dir) = wannier_centres_translated(real_space_ham%one_dim_dir, i) &
- wannier_centres_translated(real_space_ham%one_dim_dir, j)
do n1 = -irvec_max, irvec_max
dist_vec(real_space_ham%one_dim_dir) = dist_ij_vec(real_space_ham%one_dim_dir) &
+ shift_vec(real_space_ham%one_dim_dir, n1)
!MS: Add special case for lcr: We must not cut the elements that are within
! dist_cutoff under PBC's (single kpt assumed) in order to build
! hamiltonians correctly in tran_2c2_build_hams
if ((index(transport%mode, 'lcr') > 0) .and. &
!~ (transport%num_cell_ll .eq. 1) .and. &
(abs(dist_vec(real_space_ham%one_dim_dir)) .gt. real_space_ham%dist_cutoff)) then
! Move to right
dist_vec(real_space_ham%one_dim_dir) = dist_ij_vec(real_space_ham%one_dim_dir) &
+ real_lattice(real_space_ham%one_dim_dir, one_dim_vec)
! Move to left
if (abs(dist_vec(real_space_ham%one_dim_dir)) .gt. real_space_ham%dist_cutoff) &
dist_vec(real_space_ham%one_dim_dir) = dist_ij_vec(real_space_ham%one_dim_dir) &
- real_lattice(real_space_ham%one_dim_dir, one_dim_vec)
endif
!end MS
dist = abs(dist_vec(real_space_ham%one_dim_dir))
if (dist .gt. real_space_ham%dist_cutoff) hr_one_dim(j, i, n1) = 0.0_dp
end do
end do
end do
else
do i = 1, num_wann
do j = 1, num_wann
dist_ij_vec(:) = wannier_centres_translated(:, i) - wannier_centres_translated(:, j)
do n1 = -irvec_max, irvec_max
dist_vec(:) = dist_ij_vec(:) + shift_vec(:, n1)
dist = sqrt(dot_product(dist_vec, dist_vec))
! MS: Special case (as above) equivalent for alternate definition of cut off
if ((index(transport%mode, 'lcr') > 0) .and. &
!~ (transport%num_cell_ll .eq. 1) .and. &
(dist .gt. real_space_ham%dist_cutoff)) then
! Move to right
dist_vec(:) = dist_ij_vec(:) + real_lattice(:, one_dim_vec)
dist = sqrt(dot_product(dist_vec, dist_vec))
! Move to left
if (dist .gt. real_space_ham%dist_cutoff) then
dist_vec(:) = dist_ij_vec(:) - real_lattice(:, one_dim_vec)
dist = sqrt(dot_product(dist_vec, dist_vec))
endif
endif
! End MS
if (dist .gt. real_space_ham%dist_cutoff) hr_one_dim(j, i, n1) = 0.0_dp
end do
end do
end do
end if
! output maximum to check a decay of H as a function of lattice vector R
write (stdout, '(/1x,a78)') repeat('-', 78)
write (stdout, '(1x,4x,a)') &
'Maximum real part of the real-space Hamiltonian at each lattice point'
write (stdout, '(1x,8x,a62)') repeat('-', 62)
write (stdout, '(1x,11x,a,11x,a)') 'Lattice point R', 'Max |H_ij(R)|'
! calculate number of units inside a principal layer
num_pl = 0
do n1 = -irvec_max, irvec_max
hr_tmp(:, :) = abs(hr_one_dim(:, :, n1))
hr_max = maxval(hr_tmp)
if (hr_max .gt. real_space_ham%hr_cutoff) then
if (abs(n1) .gt. num_pl) num_pl = abs(n1)
else
hr_one_dim(:, :, n1) = 0.0_dp
end if
write (stdout, '(1x,9x,5x,I5,5x,12x,F12.6)') n1, hr_max
end do
write (stdout, '(1x,8x,a62)') repeat('-', 62)
if (index(transport%mode, 'lcr') > 0) then
write (stdout, '(/1x,a,I6)') 'Number of unit cells inside the principal layer:', &
transport%num_cell_ll
write (stdout, '(1x,a,I6)') 'Number of Wannier Functions inside the principal layer:', &
transport%num_ll
elseif (index(transport%mode, 'bulk') > 0) then
write (stdout, '(/1x,a,I6)') 'Number of unit cells inside the principal layer:', num_pl
write (stdout, '(1x,a,I6)') 'Number of Wannier Functions inside the principal layer:', &
num_pl*num_wann
endif
! apply hr_cutoff to each element inside the principal layer
do n1 = -num_pl, num_pl
do i = 1, num_wann
do j = 1, num_wann
if (abs(hr_one_dim(j, i, n1)) .lt. real_space_ham%hr_cutoff) hr_one_dim(j, i, n1) = 0.0_dp
end do
end do
end do
if (print_output%timing_level > 1) call io_stopwatch_stop('tran: cut_hr_one_dim', timer)
return
end subroutine tran_cut_hr_one_dim
!================================================!
subroutine tran_get_ht(fermi_energy_list, transport, hB0, hB1, hr_one_dim, irvec_max, num_pl, &
num_wann, timing_level, seedname, timer, error, comm)
!================================================!
!
!! Construct h00 and h01
!
!================================================!
use w90_constants, only: dp
use w90_io, only: io_stopwatch_start, io_stopwatch_stop, io_date
use w90_wannier90_types, only: transport_type
use w90_error, only: w90_error_type, set_error_alloc, set_error_fatal
use w90_types, only: timer_list_type
implicit none
! arguments
real(kind=dp), allocatable, intent(in) :: fermi_energy_list(:)
type(transport_type), intent(inout) :: transport
type(timer_list_type), intent(inout) :: timer
type(w90_error_type), allocatable, intent(out) :: error
type(w90_comm_type), intent(in) :: comm
integer, intent(in) :: num_pl
integer, intent(in) :: num_wann
integer, intent(in) :: irvec_max
integer, intent(in) :: timing_level
real(kind=dp), intent(in) :: hr_one_dim(:, :, -irvec_max:)
real(kind=dp), allocatable, intent(inout) :: hB0(:, :)
real(kind=dp), allocatable, intent(inout) :: hB1(:, :)
character(len=50), intent(in) :: seedname
! local variables
integer :: ierr, file_unit
integer :: i, j, n1, im, jm
integer :: fermi_n
character(len=9) :: cdate, ctime
if (timing_level > 1) call io_stopwatch_start('tran: get_ht', timer)
fermi_n = 0
if (allocated(fermi_energy_list)) fermi_n = size(fermi_energy_list)
if (fermi_n > 1) then
call set_error_fatal(error, "Error in tran_get_ht: nfermi>1. " &
//"Set the fermi level using the input parameter 'fermi_evel'", comm)
return
endif
transport%num_bb = num_pl*num_wann
allocate (hB0(transport%num_bb, transport%num_bb), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating hB0 in tran_get_ht', comm)
return
endif
allocate (hB1(transport%num_bb, transport%num_bb), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating hB1 in tran_get_ht', comm)
return
endif
hB0 = 0.0_dp
hB1 = 0.0_dp
! h00
do j = 0, num_pl - 1
do i = 0, num_pl - 1
n1 = i - j
im = i*num_wann
jm = j*num_wann
hB0(jm + 1:jm + num_wann, im + 1:im + num_wann) = hr_one_dim(:, :, n1)
end do
end do
! h01
do j = 1, num_pl
do i = 0, j - 1
n1 = i - (j - 1) + num_pl
im = i*num_wann
jm = (j - 1)*num_wann
hB1(jm + 1:jm + num_wann, im + 1:im + num_wann) = hr_one_dim(:, :, n1)
end do
end do
! shift by fermi_energy
do i = 1, transport%num_bb
hB0(i, i) = hB0(i, i) - fermi_energy_list(1)
end do
if (transport%write_ht) then
open (newunit=file_unit, file=trim(seedname)//'_htB.dat', status='unknown', form='formatted', &
action='write')
call io_date(cdate, ctime)
write (file_unit, *) 'written on '//cdate//' at '//ctime ! Date and time
write (file_unit, '(I6)') transport%num_bb
write (file_unit, '(6F12.6)') ((hB0(j, i), j=1, transport%num_bb), i=1, transport%num_bb)
write (file_unit, '(I6)') transport%num_bb
write (file_unit, '(6F12.6)') ((hB1(j, i), j=1, transport%num_bb), i=1, transport%num_bb)
close (file_unit)
end if
if (timing_level > 1) call io_stopwatch_stop('tran: get_ht', timer)
return
end subroutine tran_get_ht
!================================================!
subroutine tran_bulk(transport, hB0, hB1, timing_level, stdout, seedname, timer, error, comm)
!================================================!
use w90_constants, only: dp, cmplx_0, cmplx_1, cmplx_i, pi
use w90_io, only: io_stopwatch_start, io_stopwatch_stop, io_date
use w90_wannier90_types, only: transport_type
use w90_error, only: w90_error_type, set_error_alloc, set_error_dealloc
use w90_types, only: timer_list_type
implicit none
! arguments
integer, intent(in) :: stdout
integer, intent(in) :: timing_level
real(kind=dp), allocatable :: hB0(:, :)
real(kind=dp), allocatable :: hB1(:, :)
type(transport_type), intent(in) :: transport
type(timer_list_type), intent(inout) :: timer
type(w90_error_type), allocatable, intent(out) :: error
type(w90_comm_type), intent(in) :: comm
character(len=50), intent(in) :: seedname
! local variables
integer :: qc_unit, dos_unit
integer :: ierr
integer :: n_e, n, i
real(kind=dp) :: qc, dos
real(kind=dp) :: e_scan
complex(kind=dp) :: e_scan_cmp
complex(kind=dp), allocatable, dimension(:, :) :: tot, tott
complex(kind=dp), allocatable, dimension(:, :) :: g_B, gR, gL
complex(kind=dp), allocatable, dimension(:, :) :: sLr, sRr
complex(kind=dp), allocatable, dimension(:, :) :: s1, s2, c1
character(len=50) :: filename
character(len=9) :: cdate, ctime
if (timing_level > 1) call io_stopwatch_start('tran: bulk', timer)
allocate (tot(transport%num_bb, transport%num_bb), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating tot in tran_bulk', comm)
return
endif
allocate (tott(transport%num_bb, transport%num_bb), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating tott in tran_bulk', comm)
return
endif
allocate (g_B(transport%num_bb, transport%num_bb), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating g_B in tran_bulk', comm)
return
endif
allocate (gL(transport%num_bb, transport%num_bb), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating gL in tran_bulk', comm)
return
endif
allocate (gR(transport%num_bb, transport%num_bb), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating gR in tran_bulk', comm)
return
endif
allocate (sLr(transport%num_bb, transport%num_bb), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating sLr in tran_bulk', comm)
return
endif
allocate (sRr(transport%num_bb, transport%num_bb), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating sRr in tran_bulk', comm)
return
endif
allocate (s1(transport%num_bb, transport%num_bb), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating s1 in tran_bulk', comm)
return
endif
allocate (s2(transport%num_bb, transport%num_bb), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating s2 in tran_bulk', comm)
return
endif
allocate (c1(transport%num_bb, transport%num_bb), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating c1 in tran_bulk', comm)
return
endif
call io_date(cdate, ctime)
open (newunit=qc_unit, file=trim(seedname)//'_qc.dat', status='unknown', &
form='formatted', action='write')
write (qc_unit, *) '## written on '//cdate//' at '//ctime ! Date and time
open (newunit=dos_unit, file=trim(seedname)//'_dos.dat', status='unknown', &
form='formatted', action='write')
write (dos_unit, *) '## written on '//cdate//' at '//ctime ! Date and time
! set up the layer hamiltonians
if (transport%read_ht) then
allocate (hB0(transport%num_bb, transport%num_bb), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating hB0 in tran_bulk', comm)
return
endif
allocate (hB1(transport%num_bb, transport%num_bb), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating hB1 in tran_bulk', comm)
return
endif
filename = trim(seedname)//'_htB.dat'
call tran_read_htX(transport%num_bb, hB0, hB1, filename, stdout, error, comm)
if (allocated(error)) return
end if
! loop over the energies
n_e = floor((transport%win_max - transport%win_min)/transport%energy_step) + 1
write (stdout, '(/1x,a)', advance='no') 'Calculating quantum&
& conductance and density of states...'
do n = 1, n_e
e_scan = transport%win_min + real(n - 1, dp)*transport%energy_step
! if (mod(n,nint(0.1*n_e)).eq.0) write(stdout,'(a)',advance='no') '.'
! compute conductance according to Fisher and Lee
! retarded Green
e_scan_cmp = e_scan + eta
call tran_transfer(tot, tott, hB0, hB1, e_scan_cmp, transport%num_bb, stdout, error, comm)
if (allocated(error)) return
call tran_green(tot, tott, hB0, hB1, e_scan, g_B, 0, 1, transport%num_bb, stdout, error, comm)
if (allocated(error)) return
! compute S_Lr and S_Rr
c1(:, :) = cmplx(hB1(:, :), kind=dp)
! Self-energy (Sigma_L^r) : sLr = (hB1)^+ * tott
! Self-energy (Sigma_R^r) : sRr = (hB1) * tot
sLr = cmplx_0
sRr = cmplx_0
call ZGEMM('C', 'N', transport%num_bb, transport%num_bb, transport%num_bb, cmplx_1, c1, &
transport%num_bb, tott, transport%num_bb, cmplx_0, sLr, transport%num_bb)
call ZGEMM('N', 'N', transport%num_bb, transport%num_bb, transport%num_bb, cmplx_1, c1, &
transport%num_bb, tot, transport%num_bb, cmplx_0, sRr, transport%num_bb)
! Gamma_L = i(Sigma_L^r-Sigma_L^a)
gL = cmplx_i*(sLr - conjg(transpose(sLr)))
! Gamma_R = i(Sigma_R^r-Sigma_R^a)
gR = cmplx_i*(sRr - conjg(transpose(sRr)))
s1 = cmplx_0
s2 = cmplx_0
c1 = cmplx_0
! s1 = Gamma_L * g_B^r
call ZGEMM('N', 'N', transport%num_bb, transport%num_bb, transport%num_bb, cmplx_1, gL, &
transport%num_bb, g_B, transport%num_bb, cmplx_0, s1, transport%num_bb)
! s2 = Gamma_L * g_B^r * Gamma_R
call ZGEMM('N', 'N', transport%num_bb, transport%num_bb, transport%num_bb, cmplx_1, s1, &
transport%num_bb, gR, transport%num_bb, cmplx_0, s2, transport%num_bb)
! c1 = Gamma_L * g_B^r * Gamma_R * g_B^a
call ZGEMM('N', 'C', transport%num_bb, transport%num_bb, transport%num_bb, cmplx_1, s2, &
transport%num_bb, g_B, transport%num_bb, cmplx_0, c1, transport%num_bb)
qc = 0.0_dp
do i = 1, transport%num_bb
qc = qc + real(c1(i, i), dp)
end do
write (qc_unit, '(f15.9,f18.9)') e_scan, qc
dos = 0.0_dp
do i = 1, transport%num_bb
dos = dos - aimag(g_B(i, i))
end do
dos = dos/pi
write (dos_unit, '(f15.9,f18.9)') e_scan, dos
end do
write (stdout, '(a/)') ' done'
close (qc_unit)
close (dos_unit)
deallocate (c1, stat=ierr)
if (ierr /= 0) then
call set_error_dealloc(error, 'Error in deallocating c1 in tran_bulk', comm)
return
endif
deallocate (s2, stat=ierr)
if (ierr /= 0) then
call set_error_dealloc(error, 'Error in deallocating s2 in tran_bulk', comm)
return
endif
deallocate (s1, stat=ierr)
if (ierr /= 0) then
call set_error_dealloc(error, 'Error in deallocating s1 in tran_bulk', comm)
return
endif
deallocate (sRr, stat=ierr)
if (ierr /= 0) then
call set_error_dealloc(error, 'Error in deallocating sRr in tran_bulk', comm)
return
endif
deallocate (sLr, stat=ierr)
if (ierr /= 0) then
call set_error_dealloc(error, 'Error in deallocating sLr in tran_bulk', comm)
return
endif
deallocate (gR, stat=ierr)
if (ierr /= 0) then
call set_error_dealloc(error, 'Error in deallocating gR in tran_bulk', comm)
return
endif
deallocate (gL, stat=ierr)
if (ierr /= 0) then
call set_error_dealloc(error, 'Error in deallocating gL in tran_bulk', comm)
return
endif
deallocate (g_B, stat=ierr)
if (ierr /= 0) then
call set_error_dealloc(error, 'Error in deallocating g_B in tran_bulk', comm)
return
endif
deallocate (tott, stat=ierr)
if (ierr /= 0) then
call set_error_dealloc(error, 'Error in deallocating tott in tran_bulk', comm)
return
endif
deallocate (tot, stat=ierr)
if (ierr /= 0) then
call set_error_dealloc(error, 'Error in deallocating tot in tran_bulk', comm)
return
endif
if (timing_level > 1) call io_stopwatch_stop('tran: bulk', timer)
return
end subroutine tran_bulk