-
Notifications
You must be signed in to change notification settings - Fork 0
/
readHmnR.f90
executable file
·1199 lines (1024 loc) · 37 KB
/
readHmnR.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
subroutine readNormalHmnR()
!>> Read in the tight-binding model from wannier90_hr.dat
!> The format is defined by the wannier90 software
! Constructed by quansheng wu 4/2/2010
!
! Yifei Guan added the sparse hr file parsing June/2018
! License: GPL V3
use para
!> in: N of wann
!> out : nth atom
implicit none
character*4 :: c_temp
integer :: i, j, ir, ia, io
integer :: i1, i2, i3, i4, i5
integer :: n, m, ir0
integer :: add_electric_field
integer :: nwann, nwann_nsoc
real(dp) :: static_potential
real(dp) :: tot, rh, ih
real(dp) :: pos(Origin_cell%Num_atoms)
!> add a check for NumOccupied parameter
if (NumOccupied<=0 .or. NumOccupied>Num_wann) then
write(stdout, '(a, i6, a)')">>> ERROR: NumOccupied should be in [1, ", Num_wann, " ]"
write(stdout, '(a)')">>> Usually, it is the number of occupied Wannier bands."
stop
endif
if(cpuid.eq.0)write(stdout,*)' '
open(12, file=Hrfile, status='OLD')
if (index(Hrfile, 'HWR')==0) then
!> for Normal HmnR obtained from Wannier90 or sparse HmnR
!> skip a comment line
read(12, *)
!> number of Wannier orbitals in the hr file
nwann=0
read(12, *)nwann
if (nwann==0) then
stop "ERROR : num_wann is zero in hr file"
endif
nwann_nsoc=nwann
if (SOC>0) nwann_nsoc= nwann/2
if ((soc==0 .and. sum(Origin_cell%nprojs)/=nwann .and. .not.Add_Zeeman_Field) .or. &
(soc>0 .and. sum(Origin_cell%nprojs)/=nwann/2))then
print *, 'sum(Origin_cell%nprojs), num_wann, num_wann/2'
print *, sum(Origin_cell%nprojs), nwann, nwann/2
print *, "ERROR: Maybe the SOC tags in the SYSTEM is wrongly set"
stop "ERROR: the summation of all projectors times spin degeneracy is not equal to num_wann"
endif
!> number of lattice vectors taken into account
read(12, *)Nrpts
if (.not. allocated(HmnR)) allocate(HmnR(Num_wann,Num_wann,nrpts))
allocate(irvec(3,nrpts))
allocate(ndegen(nrpts))
irvec= 0
ndegen=1
!> The degeneracy of each R point, this is caused by the Fourier transform in the Wigner-Seitz cell
read(12, *)(ndegen(i), i=1, Nrpts)
ir=0
do ir=1,Nrpts
do n=1,nwann
do m=1,nwann
read(12,*,end=1001)i1, i2, i3, i4, i5, rh, ih
irvec(1,ir)=i1
irvec(2,ir)=i2
irvec(3,ir)=i3
HmnR(i4,i5,ir)=dcmplx(rh,ih)
end do
enddo
enddo
!> extract the fermi energy
do iR=1,Nrpts
if (Irvec(1,iR).eq.0.and.Irvec(2,iR).eq.0.and.Irvec(3,iR).eq.0)then
do i=1, Num_wann
HmnR(i,i,iR)=HmnR(i,i,iR)-E_fermi
enddo
endif
enddo
!> WannierTools codes use Hatree atomic units
if (index(Particle,'phonon')/=0) then
HmnR= HmnR*eV2Hartree*eV2Hartree ! from eV to Hartree
else
HmnR= HmnR*eV2Hartree ! from eV to Hartree
endif
else
!File *.HWR exist, We are using HmnR from WHM
! skip 8 lines
do i=1,8
read(12,*)
enddo
read(12,'(a11,f18.7)')c_temp,E_fermi
do iR=1,Nrpts
read(12,'(a3,3i5,a3,i4)')c_temp,irvec(1:3,iR),c_temp,ndegen(iR)
do i=1, Num_wann*Num_wann
read(12,*)n,m,rh,ih
HmnR(n,m,iR)=rh+ zi*ih ! in Hartree
enddo
if (sum(abs(irvec(:, ir)))==0) then
do i=1, Num_wann
HmnR(i,i,iR)= HmnR(i,i,iR)-E_fermi ! in Hartree
enddo
endif
enddo
if (cpuid==0) then
open(unit=105, file='wannier90_hr.dat')
write(105, *)'hr file transformed from HWR'
write(105, *)Num_wann
write(105, *)nrpts
write(105, '(15I5)')(ndegen(i), i=1, nrpts)
do ir=1, nrpts
do i=1, Num_wann
do j=1, Num_wann
write(105, '(5I5, 2f16.8)')irvec(:, ir), i, j, HmnR(i, j, ir)
enddo
enddo
enddo
close(105)
endif
endif ! HWR or not
1001 continue
close(12)
!call get_fermilevel
!check sum rule
tot= 0d0
do ir=1, Nrpts
tot= tot+ 1d0/ndegen(ir)
enddo
!> get the Cartesian coordinates for R points
allocate( crvec(3, nrpts))
!> get R coordinates
do iR=1, Nrpts
crvec(:, iR)= Origin_cell%Rua*irvec(1,iR) + Origin_cell%Rub*irvec(2,iR) + Origin_cell%Ruc*irvec(3,iR)
enddo
!> change from "up dn up dn" to "up up dn dn"
if (index( Package, 'QE')/=0.or.index( Package, 'quantumespresso')/=0 &
.or.index( Package, 'quantum-espresso')/=0.or.index(Package, 'VASP6')/=0.or.index( Package, 'pwscf')/=0) then
call reorder_wannierbasis
if (cpuid==0.and.export_newhr) then
!> write to new_hr.dat
outfileindex= outfileindex+ 1
open(unit=outfileindex, file='wannier90_hr_standard.dat')
write(outfileindex, '(a,1X,a,1X,a,1X, a, a)') &
'HmnR transformed from QE ', time_now, date_now, 'UTC', zone_now
write(outfileindex, *)Num_wann
write(outfileindex, *)nrpts
write(outfileindex, '(15I5)')ndegen
do ir=1, nrpts
do j=1, Num_wann
do i=1, Num_wann
write( outfileindex, '(5I5, 2f16.6)') &
irvec(:, ir), i, j, HmnR(i, j, ir)/eV2Hartree
end do
end do
end do
close(outfileindex)
endif
endif
!> Adding zeeman field
!> Bx=Bdirection(1)
!> By=Bdirection(2)
!> Bz=Bdirection(3)
!> Hz= Zeeman_energy_in_eV*(Bx*sx+By*sy+Bz*sz)/2d0
!> sx, sy, sz are Pauli matrices.
if (Add_Zeeman_Field) then
if (abs(Zeeman_energy_in_eV)<eps6)then
Zeeman_energy_in_eV= Bmagnitude*Effective_gfactor*Bohr_magneton
if (cpuid==0)then
write(stdout, '(1x, a, 3f16.6)')"Zeeman_energy_in_eV: ", Zeeman_energy_in_eV/eV2Hartree
endif
endif
!> After considering the Zeeman field, we already extended the spin space to spin-full.
SOC = 1
endif ! Add_Zeeman_Field
call get_stacking_direction_and_pos(add_electric_field, pos)
ir0=0
do ir=1, nrpts
if (irvec(1, ir)==0.and.irvec(2, ir)==0.and.irvec(3, ir)==0) ir0=ir
enddo
if (add_electric_field>0) then
io=0
do ia=1, Origin_cell%Num_atoms
!static_potential= pos(ia)*Origin_cell%cell_parameters(add_electric_field)*Electric_field_in_eVpA
!if (Inner_symmetrical_Electric_Field) then
! static_potential= abs(pos(ia)*Origin_cell%cell_parameters(add_electric_field))*Electric_field_in_eVpA
!endif
static_potential= abs(pos(ia)*Origin_cell%cell_parameters(add_electric_field))&
*Symmetrical_Electric_field_in_eVpA*eV2Hartree/Angstrom2atomic+ &
(pos(ia)*Origin_cell%cell_parameters(add_electric_field))*&
Electric_field_in_eVpA*eV2Hartree/Angstrom2atomic
do i=1, Origin_cell%nprojs(ia)
io=io+1
HmnR(io, io, ir0)= HmnR(io, io, ir0)+ static_potential
if (SOC>0) then
HmnR(io+Num_wann/2, io+Num_wann/2, ir0)= HmnR(io+Num_wann/2, io+Num_wann/2, ir0)+ static_potential
endif ! SOC
enddo ! nproj
enddo ! ia
endif ! add electric field or not
!> write out Hmn(R=0)
if (cpuid.eq.0 .and. Num_wann< 200)then
write(stdout, '(a)')" "
write(stdout, '(a)')" >> Hopping parameters in the home unit cell"
write(stdout, '(a)')" >> H00= Hmn(R=0) real part"
do i=1, Num_wann
write(stdout, '(50000f7.3)') real(HmnR(i, :, ir0))/eV2Hartree
enddo
write(stdout, '(a)')" "
write(stdout, '(a)')" >> H00= Hmn(R=0) imagary part"
do i=1, Num_wann
write(stdout, '(50000f7.3)') aimag(HmnR(i, :, ir0))/eV2Hartree
enddo
write(stdout, '(a)')" "
endif
call get_hmnr_cell(Cell_defined_by_surface)
return
end subroutine readNormalHmnR
subroutine read_valley_operator()
!>> Read in the valley operator from valley_operator.dat
! Constructed by quansheng wu 04 Nov. 2023
! License: GPL V3
use para
implicit none
character*4 :: c_temp
integer :: i, j, ir, ia, io
integer :: i1, i2, i3, i4, i5
integer :: n, m, ir0
integer :: add_electric_field
integer :: nwann, nwann_nsoc
logical :: exists
real(dp) :: static_potential
real(dp) :: tot, rh, ih
real(dp) :: pos(Origin_cell%Num_atoms)
if(cpuid.eq.0)write(stdout,*)' '
inquire (file ="valley_operator.dat", EXIST = exists)
if (exists)then
open(12, file="valley_operator.dat", status='OLD')
else
STOP ">> for valley projection , you have to prepare a file valley_operator.dat"
endif
!> skip a comment line
read(12, *)
!> number of Wannier orbitals in the hr file
nwann=0
read(12, *)nwann
if (nwann==0) then
stop "ERROR : num_wann is zero in hr file"
endif
nwann_nsoc=nwann
if (SOC>0) nwann_nsoc= nwann/2
!> number of lattice vectors taken into account
read(12, *)Nrpts_valley
!> The degeneracy of each R point, this is caused by the Fourier transform in the Wigner-Seitz cell
read(12, *)(j, i=1, Nrpts_valley)
allocate(irvec_valley(3, Nrpts_valley))
allocate(valley_operator_R(num_wann, num_wann, Nrpts_valley))
ir=0
do ir=1,Nrpts_valley
do n=1,nwann
do m=1,nwann
read(12,*,end=1001)i1, i2, i3, i4, i5, rh, ih
irvec_valley(1,ir)=i1
irvec_valley(2,ir)=i2
irvec_valley(3,ir)=i3
valley_operator_R(i4,i5,ir)=dcmplx(rh,ih)
end do
enddo
enddo
1001 continue
close(12)
if (cpuid.eq.0) write(stdout, *) ">>> finished reading of valley operator"
end subroutine read_valley_operator
subroutine get_hmnr_cell(cell)
!> Get new hmnr for a new cell with the same size as the previous one
use para
implicit none
type(cell_type) :: cell
!type(dense_tb_hr) :: cell_hr
integer :: ir, i, j, iter
real(dp) :: shift_vec_direct(3)
!> for newcell
real(dp) :: apos1d(3),apos2d(3)
!>count newcell nrpts
integer :: max_ir
integer :: nir1,nir2,nir3, ir_cell
integer :: nrpts_new, nrpts_max
real(dp) :: new_ia, new_ib, new_ic, max_val
!> all new irs
integer, allocatable :: rpts_array(:, :, :), rpts_map(:, :, :)
integer, allocatable :: allirs(:, :, :, :)
integer :: irn1(3),irn2(3)
max_ir=8
nrpts_max=(2*max_ir+1)**3
allocate( rpts_array(-max_ir:max_ir,-max_ir:max_ir,-max_ir:max_ir))
allocate( rpts_map(-max_ir:max_ir,-max_ir:max_ir,-max_ir:max_ir))
allocate(allirs(nrpts, num_wann, num_wann, 3))
call cart_direct_real(shift_to_topsurface_cart, shift_vec_direct, cell%lattice)
call date_and_time(DATE=date_now,ZONE=zone_now, TIME=time_now)
!> Get new Hrs
rpts_array=0
nrpts_new=0
rpts_map= 0
!> get number of R points for the new cell first
do ir=1, nrpts
do i=1, Num_wann
do j=1, Num_wann
apos1d= Origin_cell%wannier_centers_direct(:, i)- shift_vec_direct
apos2d= Origin_cell%wannier_centers_direct(:, j)- shift_vec_direct+irvec(:, ir)
call latticetransform(apos1d(1),apos1d(2),apos1d(3),new_ia,new_ib,new_ic)
irn1=floor([new_ia,new_ib,new_ic])
call latticetransform(apos2d(1),apos2d(2),apos2d(3),new_ia,new_ib,new_ic)
irn2=floor([new_ia,new_ib,new_ic])
nir1=irn2(1)-irn1(1)
nir2=irn2(2)-irn1(2)
nir3=irn2(3)-irn1(3)
if (abs(nir1)>max_ir .or. abs(nir2)>max_ir .or. abs(nir3)>max_ir) cycle
rpts_array(nir1,nir2,nir3)=1
enddo
enddo
enddo
!> find all irvec
nrpts_new= sum(rpts_array)
allocate(irvec_newcell(3, Nrpts_new))
iter= 0
do nir3=-max_ir,max_ir
do nir2=-max_ir,max_ir
do nir1=-max_ir,max_ir
if (rpts_array(nir1, nir2, nir3)==1) then
iter=iter+1
irvec_newcell(:, iter)=[nir1, nir2, nir3]
rpts_map(nir1, nir2, nir3)=iter
endif
enddo
enddo
enddo
!>> hamiltonian for the new cell
allocate(HmnR_newcell(Num_wann, Num_wann, Nrpts_new))
allocate(ndegen_newcell(Nrpts_new))
HmnR_newcell= 0d0
ndegen_newcell= 1
!> get number of R points for the new cell first
do ir=1, nrpts
do j=1, Num_wann
do i=1, Num_wann
!ia1=Origin_cell%spinorbital_to_atom_index(i)
!ia2=Origin_cell%spinorbital_to_atom_index(j)
!apos1d= Origin_cell%Atom_position_direct(:, ia1)- shift_vec_direct
!apos2d= Origin_cell%Atom_position_direct(:, ia2)- shift_vec_direct+irvec(:, ir)
apos1d= Origin_cell%wannier_centers_direct(:, i)- shift_vec_direct
apos2d= Origin_cell%wannier_centers_direct(:, j)- shift_vec_direct+irvec(:, ir)
call latticetransform(apos1d(1),apos1d(2),apos1d(3),new_ia,new_ib,new_ic)
irn1=floor([new_ia,new_ib,new_ic])
call latticetransform(apos2d(1),apos2d(2),apos2d(3),new_ia,new_ib,new_ic)
irn2=floor([new_ia,new_ib,new_ic])
nir1=irn2(1)-irn1(1)
nir2=irn2(2)-irn1(2)
nir3=irn2(3)-irn1(3)
if (abs(nir1)>max_ir .or. abs(nir2)>max_ir .or. abs(nir3)>max_ir) cycle
ir_cell= rpts_map(nir1, nir2, nir3)
HmnR_newcell(i,j,ir_cell)=HmnR(i,j,ir)/ndegen(ir)
enddo
enddo
enddo
!> do cut-off according to the hopping value
iter= 0
do ir=1, nrpts_new
max_val=maxval(abs(HmnR_newcell(:, :, ir)))/eV2Hartree
if (max_val<eps4) then
iter= iter+ 1
endif
enddo
if (cpuid==0.and.export_newhr) then
!> write to new_hr.dat
nrpts_new=sum(rpts_array)-iter
outfileindex= outfileindex+ 1
open(unit=outfileindex, file='wannier90_hr_newcell.dat')
write(outfileindex, '(a,1X,a,1X,a,1X, a, a)') &
'HmnR for new cell generated at ', time_now, date_now, 'UTC', zone_now
write(outfileindex, *)Num_wann
write(outfileindex, *)nrpts_new
write(outfileindex, '(15I5)')(1 , i=1, nrpts_new)
do ir=1, nrpts_new
max_val=maxval(abs(HmnR_newcell(:, :, ir)))/eV2Hartree
!if (max_val<eps4) cycle
do j=1, Num_wann
do i=1, Num_wann
write( outfileindex, '(5I5, 2f16.6)') &
irvec_newcell(:, ir), i, j, HmnR_newcell(i, j, ir)/eV2Hartree
end do
end do
end do
close(outfileindex)
endif
return
end subroutine get_hmnr_cell
!>-----------------------------------------------------------------------
! The sparse hamiltonian is stored as follows
! start with a comment line
! 10 ! splen: number of non-zero Hmn(R) lines
! 20 ! nwan: number of orbitals (including spin-degeneracy if SOC is included)
! 9 ! nrpts: number of R points
! 1 0 0 1 1 1.0 0.0 ! i1, i2, i3, n, m, rh, ih Hnm(i1, i2, i3)
! ... ! in total there are splen lines
!> Here the unit of energy is eV
!>-----------------------------------------------------------------------
subroutine readSparseHmnR
!> This subroutine not just read the sparse hr file, but also can read
!> the standard hr file defined in the Wannier90.
use para
implicit none
integer:: i,j,nwann,nwann_nsoc,i1,i2,i3,i4,i5,ir, ia, io
integer :: q,n,m, ir0
real(8) :: r1,r2, pos(Origin_cell%Num_atoms), static_potential
!> the direction which adding electric field which is also the stacking direction
integer :: add_electric_field
real(dp) :: Bx_in_au, By_in_au, Bz_in_au
complex(dp) :: h_value
open(12, file=Hrfile)
!> skip a comment line
read(12, *)
!> comparing with the standard hr file, we add another line to show howmany
!> lines that Hmn(R) is not zero.
if(Is_Sparse_Hr) then
read(12,*) splen_input !> number of non-zeros lines
end if
!> number of Wannier orbitals in the hr file
nwann=0
read(12, *)nwann
if (nwann==0) then
stop "ERROR : num_wann is zero in hr file"
endif
nwann_nsoc=nwann
if (SOC>0) nwann_nsoc= nwann/2
!> whether we need to add the electric field potential
add_electric_field=0
call get_stacking_direction_and_pos(add_electric_field, pos)
if (add_electric_field>0) then
!> with Electric_field
if (Add_Zeeman_Field) then
if (SOC==0) then
splen= splen_input*2+ 6*nwann_nsoc
else
splen= splen_input+ 6*nwann_nsoc
endif
else
if (SOC==0) then
splen= splen_input+ nwann_nsoc
else
splen= splen_input+ nwann_nsoc*2
endif
endif
else
!> without electric field
if (Add_Zeeman_Field) then
if (SOC==0) then
splen= splen_input*2+ 4*nwann_nsoc
else
splen= splen_input+ 4*nwann_nsoc
endif
else
!> no electrical field and Zeeman field
splen= splen_input
endif
end if
!> in order to include the Fermi level
!> for non-orthogonal_Basis, we can't include the fermi level in to H'
if (.not.Orthogonal_Basis) then
splen=splen
else
splen=splen+nwann
endif
allocate(hacoo(splen),hicoo(splen),hjcoo(splen),hirv(3, splen))
hacoo=(0d0, 0d0)
hicoo=0
hjcoo=0
hirv=0
!> will reread the above line
do j=1, splen_input
read(12,*,end=1001)i1, i2, i3, i4, i5, r1, r2
hirv (1, j)=i1
hirv (2, j)=i2
hirv (3, j)=i3
hicoo(j)=i4
hjcoo(j)=i5
if (trim(adjustl(Package))=="OPENMX") then
hacoo(j)=dcmplx(r1,r2)
else
hacoo(j)=dcmplx(r1,r2)*eV2Hartree
endif
enddo
j=j-1
1001 continue
!> Adding zeeman field
!> Bx=Bdirection(1)
!> By=Bdirection(2)
!> Bz=Bdirection(3)
!> Hz= Zeeman_energy_in_eV*(Bx*sx+By*sy+Bz*sz)/2d0
!> sx, sy, sz are Pauli matrices.
if (Add_Zeeman_Field) then
if (abs(Zeeman_energy_in_eV)<eps6)then
Zeeman_energy_in_eV= Bmagnitude*Effective_gfactor*Bohr_magneton
if (cpuid==0)then
write(stdout, '(1x, a, 3f16.6)')"Zeeman_energy_in_eV: ", Zeeman_energy_in_eV/eV2Hartree
endif
endif
!> for sparse hr format
!> extend from spinless to spinfull
if (SOC_in==0) then
do i=1, splen_input
j= i+ splen_input
hicoo(j)=hicoo(i)+nwann_nsoc
hjcoo(j)=hjcoo(i)+nwann_nsoc
hirv (1, j)=hirv (1, i)
hirv (2, j)=hirv (2, i)
hirv (3, j)=hirv (3, i)
hacoo(j)=hacoo(i)
enddo
endif
!> zeeman energy
if (SOC_in==0) then
j= splen_input*2
else
j= splen_input
endif
Bx_in_au= Bx
By_in_au= By
Bz_in_au= Bz
call add_zeeman_sparse_hr(Bx_in_au, By_in_au, Bz_in_au)
j=j+ nwann_nsoc*4
!> After adding zeeman term, SOC=1
SOC = 1
endif ! Add_Zeeman_Field
!> Adding electric field
if (add_electric_field>0) then
if (Add_Zeeman_Field) then
!> continue with the j used in the Add_Zeeman_Field
io=0
do ia=1, Origin_cell%Num_atoms
!static_potential= pos(ia)*Origin_cell%cell_parameters(add_electric_field)*Electric_field_in_eVpA
!if (Inner_symmetrical_Electric_Field) then
! static_potential= abs(pos(ia)*Origin_cell%cell_parameters(add_electric_field))*Electric_field_in_eVpA
!endif
static_potential= abs(pos(ia)*Origin_cell%cell_parameters(add_electric_field))&
*Symmetrical_Electric_field_in_eVpA*eV2Hartree/Angstrom2atomic+ &
(pos(ia)*Origin_cell%cell_parameters(add_electric_field))*&
Electric_field_in_eVpA*eV2Hartree/Angstrom2atomic
do i=1, Origin_cell%nprojs(ia)
io=io+1
!> spin up
j=j+1
hicoo(j)= io
hjcoo(j)= io
hirv (1, j)= 0
hirv (2, j)= 0
hirv (3, j)= 0
hacoo(j)= static_potential
!> spin down
j=j+1
hicoo(j)= io + nwann_nsoc
hjcoo(j)= io + nwann_nsoc
hirv (1, j)= 0
hirv (2, j)= 0
hirv (3, j)= 0
hacoo(j)= static_potential
enddo ! nproj
enddo ! ia
else
j=splen_input
io=0
do ia=1, Origin_cell%Num_atoms
!static_potential= pos(ia)*Origin_cell%cell_parameters(add_electric_field)*Electric_field_in_eVpA
!if (Inner_symmetrical_Electric_Field) then
! static_potential= abs(pos(ia)*Origin_cell%cell_parameters(add_electric_field))*Electric_field_in_eVpA
!endif
static_potential= abs(pos(ia)*Origin_cell%cell_parameters(add_electric_field))&
*Symmetrical_Electric_field_in_eVpA*eV2Hartree/Angstrom2atomic+ &
(pos(ia)*Origin_cell%cell_parameters(add_electric_field))*&
Electric_field_in_eVpA*eV2Hartree/Angstrom2atomic
do i=1, Origin_cell%nprojs(ia)
io=io+1
j=j+1
hicoo(j)= io
hjcoo(j)= io
hirv (1, j)= 0
hirv (2, j)= 0
hirv (3, j)= 0
hacoo(j)= static_potential
if (SOC>0) then ! with SOC
j=j+1
hicoo(j)= io + nwann_nsoc
hjcoo(j)= io + nwann_nsoc
hirv (1, j)= 0
hirv (2, j)= 0
hirv (3, j)= 0
hacoo(j)= static_potential
endif
enddo ! nproj
enddo ! ia
endif ! Add_Zeeman_Field or not
endif ! add_electric_field
!> the Fermi level can be added to the Hamiltonian directly only with in the orthogonal basis
if (Orthogonal_Basis) then
!> add Fermi level
do n=1,num_wann
j=j+1
hicoo(j)=n
hjcoo(j)=n
hirv (1, j)= 0
hirv (2, j)= 0
hirv (3, j)= 0
hacoo(j)= - E_FERMI*eV2Hartree
enddo
endif
!> transform it into dense format only when number of wannier orbitals is less than 100
outfileindex= outfileindex+ 1
if (cpuid.eq.0.and. nwann<500) then
! open (unit=outfileindex, file='hr.dat-dense')
! write(outfileindex, *) ' ! HmnR file from sparse hr file'
! write(outfileindex, '(I10, a)') nwann, ' ! Num_wann: number of wannier orbitals'
! write(outfileindex, '(I10, a)') Nrpts, ' ! NRPTS: number of R points'
! write(outfileindex, '(10I5)') ndegen(:)
! do iR=1, NRPTS
! i1= irvec(1, iR)
! i2= irvec(2, iR)
! i3= irvec(3, iR)
! do n=1, nwann
! do m=1, nwann
! h_value=0d0
! do i=1, splen
! if (hicoo(i)==n.and.hjcoo(i)==m.and.hirv(1, i)==i1.and.hirv(2, i)==i2.and.hirv(3, i)==i3)then
! h_value= h_value+ hacoo(i)
! endif
! enddo !i
! write(outfileindex, '(5I5, 2f12.6)')irvec(:, iR), n, m, h_value/eV2Hartree
! enddo !m
! enddo !n
! enddo !iR
! close(outfileindex)
endif
return
end subroutine readSparseHmnR
subroutine readsparse_overlap
!> This subroutine reads the overlap matrix S between the non-orthogonal basis
!> 2023 Dec. 4th
!> by QSWU ([email protected])
use para
implicit none
integer:: i,j,nwann,nwann_nsoc,i1,i2,i3,i4,i5,ir,n,m
real(dp) :: r1, r2
logical :: exists
if(cpuid.eq.0)write(stdout,*)' '
if(cpuid.eq.0)write(stdout,'(2a)')' >> Now we are reading the overlap files: ', trim(adjustl(Overlapfile))
inquire(file= Overlapfile, EXIST = exists)
if (exists)then
open(13, file= Overlapfile, status='OLD')
else
STOP ">> for non-orthogonal basis , you have to prepare a file like overlap.dat"
endif
!> skip a comment line
read(13, *)
!> comparing with the standard overlap file, we add another line to show howmany
!> lines that S is not zero.
! if (Is_Sparse_Hr) then
read(13,*) splen_overlap_input !> number of non-zeros lines
! end if
!> number of Wannier orbitals in the hr file
nwann=0
read(13, *)nwann
if (nwann==0.or.nwann.ne.Num_wann) then
print *, 'in readsparse_overlap'
print *, 'nwann, num_wann', nwann, num_wann
stop "ERROR : num_wann is zero or num_wann is not consistent with hr.dat"
endif
nwann_nsoc=nwann
if (SOC>0) nwann_nsoc= nwann/2
!> overlap matrix in sparse format
allocate(sacoo(splen_overlap_input))
allocate(sicoo(splen_overlap_input))
allocate(sjcoo(splen_overlap_input))
allocate(sirv(3, splen_overlap_input))
sacoo=(0d0, 0d0)
sicoo=0
sjcoo=0
sirv=0
!> will reread the above line
do i=1, splen_overlap_input
read(13,*,end=1003)i1, i2, i3, i4, i5, r1, r2
sirv(1, i)= i1
sirv(2, i)= i2
sirv(3, i)= i3
sicoo(i)=i4
sjcoo(i)=i5
sacoo(i)=dcmplx(r1,r2)
enddo
1003 continue
close(13)
if (cpuid.eq.0) write(stdout, '(a, i)')' >> Number of non-zeros splen_overlap_input', splen_overlap_input
if (cpuid.eq.0) write(stdout, '(a)')' >> Sparse overlap matrix reading finished '
return
end subroutine readsparse_overlap
subroutine readsparse_valley_operator
!> This subroutine not just read the sparse valley operator
!> 2023 Nov. 4
use para
implicit none
integer:: i,j,nwann,nwann_nsoc,i1,i2,i3,i4,i5,ir,n,m
real(dp) :: r1, r2
logical :: exists
!> the direction which adding electric field which is also the stacking direction
integer :: add_electric_field
real(dp) :: Bx_in_au, By_in_au, Bz_in_au
complex(dp) :: h_value
if(cpuid.eq.0)write(stdout,*)' '
inquire (file ="valley_operator.dat", EXIST = exists)
if (exists)then
open(13, file="valley_operator.dat", status='OLD')
else
STOP ">> for valley projection , you have to prepare a file valley_operator.dat"
endif
!> skip a comment line
read(13, *)
!> comparing with the standard hr file, we add another line to show howmany
!> lines that Hmn(R) is not zero.
if(Is_Sparse_Hr) then
read(13,*) splen_valley_input !> number of non-zeros lines
end if
!> number of Wannier orbitals in the hr file
nwann=0
read(13, *)nwann
if (nwann==0.or.nwann.ne.Num_wann) then
print *, 'nwann, num_wann', nwann, num_wann
stop "ERROR : num_wann is zero or num_wann is not consistent with hr.dat"
endif
nwann_nsoc=nwann
if (SOC>0) nwann_nsoc= nwann/2
!> number of lattice vectors taken into account
read(13, *)Nrpts_valley
allocate(irvec_valley(3, Nrpts_valley))
!> The degeneracy of each R point
read(13, *)(j, i=1, nrpts)
ir=0
allocate(valley_operator_acoo(splen_valley_input))
allocate(valley_operator_icoo(splen_valley_input))
allocate(valley_operator_jcoo(splen_valley_input))
allocate(valley_operator_irv(splen_valley_input))
valley_operator_acoo=(0d0, 0d0)
valley_operator_icoo=0
valley_operator_jcoo=0
valley_operator_irv=0
j=0
ir=1
irvec_valley=-9999
read(13,*,end=1002)i1, i2, i3, i4, i5, r1, r2
irvec_valley(1,ir)=i1
irvec_valley(2,ir)=i2
irvec_valley(3,ir)=i3
!> will reread the above line
backspace(13)
do i=1,Nrpts_valley
do n=1,nwann
do m=1,nwann
read(13,*,end=1002)i1, i2, i3, i4, i5, r1, r2
if (sum(abs(irvec_valley(:,ir)-[i1,i2,i3]))/=0) ir=ir+1
j=j+1
valley_operator_icoo(j)=i4
valley_operator_jcoo(j)=i5
valley_operator_irv (j)=ir
valley_operator_acoo(j)=dcmplx(r1,r2)
irvec_valley(1, ir)= i1
irvec_valley(2, ir)= i2
irvec_valley(3, ir)= i3
end do
enddo
enddo
1002 continue
close(13)
!> correct nrpts
Nrpts_valley=ir
if (cpuid.eq.0) write(stdout, '(a, i6)')' >> Nrpts_valley is ', Nrpts_valley
if (cpuid.eq.0) write(stdout, '(a)')' >> splen_valley_input', splen_valley_input, j
if (cpuid.eq.0) write(stdout, '(a)')' >> Valley operator reading finished '
return
end subroutine readsparse_valley_operator
!> return the maximum distance between two neighbour elements of array x
!> dx_max= mod(dx_max, 1)
!> x is in [0, 1)
subroutine get_max_nn_distance(n, x, dx_max)
use para, only: dp
implicit none
!> dimension of arr
integer, intent(in) :: n
!> an array
real(dp), intent(in) :: x(n)
real(dp) :: x1(n)
real(dp), intent(out) :: dx_max
integer :: i
real(dp) :: dx
dx_max=0d0
do i=1, n-1
dx= x1(i+1)-x1(i)
if (dx>dx_max) dx_max= dx
enddo
dx= 1d0+x(1)-x(n)
if (dx>dx_max) dx_max= dx
return
end subroutine get_max_nn_distance
!> output: add_electric_field is the direction which the E field is adding on.
!> if add_electric_field=0, we don't add E field.
!> pos are the direct coordinate of atoms along the direction where the Electric_field applied.
subroutine get_stacking_direction_and_pos(add_electric_field, pos)
use para
implicit none
integer, intent(inout) :: add_electric_field
real(dp), intent(out) :: pos(Origin_cell%Num_atoms)
real(dp), allocatable :: pos_selected_for_center(:, :)
integer :: i, ia, NumberofSelectedAtoms_center, iter, ia_g, ig
real(dp) :: dxmax(3), center
!> sum over all selected atoms
NumberofSelectedAtoms_center= sum(NumberofSelectedAtoms(:))
allocate(pos_selected_for_center(3, NumberofSelectedAtoms_center))
pos_selected_for_center= 0d0
iter= 0
do ig=1, NumberofSelectedAtoms_groups
do ia_g= 1, NumberofSelectedAtoms(ig)
iter=iter+1
ia= Selected_Atoms(ig)%iarray(ia_g)
pos_selected_for_center(:, iter)= Origin_cell%Atom_position_direct(:, ia)
enddo
enddo
if (cpuid.eq.0)then
write(stdout, '(a, 3f12.6)') ">>: Center position for zero electrical potential for adding electrical field", &
sum(pos_selected_for_center(:, :), dim=2)/NumberofSelectedAtoms_center
endif
!> Adding electric field
!> The electric field is only applied on the 2D system.
!> First we determine whether we have an vacuum and which direction
!> Here we only support adding electric field along the unit lattice vectors.
add_electric_field=0
pos=0d0
if (abs(Symmetrical_Electric_field_in_eVpA)>eps6.or.abs(Electric_field_in_eVpA)>eps6) then
do i=1, 3
pos= Origin_cell%Atom_position_direct(i, :)
!> shift all the values into [0, 1)
do ia=1, Origin_cell%Num_atoms
pos(ia)=mod(pos(ia), 1d0)
enddo
call sortheap(Origin_cell%Num_atoms, pos)
call get_max_nn_distance(Origin_cell%Num_atoms, pos, dxmax(i))
dxmax(i)= dxmax(i)*Origin_cell%cell_parameters(i)
if (dxmax(i)>9d0) then
add_electric_field=i
endif
enddo
endif ! add electric field or not
if (add_electric_field==0.and.cpuid==0) then
write(stdout, '(a)') " "
write(stdout, '(a)') " Note: This is not a 2D system or E=0, so we can't add electric field."
write(stdout, '(a)') " "
endif
!> shift all the positions together and centered at zero along the electric field
if (add_electric_field>0) then
pos=Origin_cell%Atom_position_direct(add_electric_field, :)
pos= mod(pos, 1d0)-0.5d0