forked from wannier-developers/wannier90
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathplot.F90
3145 lines (2795 loc) · 129 KB
/
plot.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
!
!------------------------------------------------------------!
! 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 !
!------------------------------------------------------------!
! !
! w90_plot: plot various data !
! !
!------------------------------------------------------------!
module w90_plot_mod
!! This module handles various plots
use w90_comms, only: comms_reduce, w90_comm_type, mpisize, mpirank, comms_array_split
implicit none
private
public :: plot_main
contains
!================================================!
subroutine plot_main(atom_data, band_plot, dis_manifold, fermi_energy_list, fermi_surface_plot, &
ham_logical, kmesh_info, kpt_latt, output_file, wvfn_read, real_space_ham, &
kpoint_path, print_output, wannier_data, wannier_plot, ws_region, &
w90_calculation, ham_k, ham_r, m_matrix, u_matrix, u_matrix_opt, eigval, &
real_lattice, wannier_centres_translated, bohr, irvec, mp_grid, ndegen, &
shift_vec, nrpts, num_bands, num_kpts, num_wann, rpt_origin, &
transport_mode, have_disentangled, lsitesymmetry, w90_system, seedname, &
stdout, timer, dist_k, error, comm)
!================================================!
!
!! Main plotting routine
!
!================================================!
use w90_constants, only: eps6, dp
use w90_hamiltonian, only: hamiltonian_get_hr, hamiltonian_write_hr, hamiltonian_setup, &
hamiltonian_write_tb
use w90_io, only: io_stopwatch_start, io_stopwatch_stop
use w90_types, only: kmesh_info_type, wannier_data_type, atom_data_type, dis_manifold_type, &
kpoint_path_type, print_output_type, ws_region_type, ws_distance_type, timer_list_type, &
w90_system_type
use w90_utility, only: utility_recip_lattice_base
use w90_wannier90_types, only: w90_calculation_type, wvfn_read_type, output_file_type, &
fermi_surface_plot_type, band_plot_type, wannier_plot_type, real_space_ham_type, &
ham_logical_type
use w90_ws_distance, only: ws_translate_dist, ws_write_vec
use w90_error, only: w90_error_type
implicit none
! arguments
type(real_space_ham_type), intent(inout) :: real_space_ham
type(w90_error_type), allocatable, intent(out) :: error
type(timer_list_type), intent(inout) :: timer
type(ham_logical_type), intent(inout) :: ham_logical
complex(kind=dp), intent(inout), allocatable :: ham_r(:, :, :)
complex(kind=dp), intent(inout), allocatable :: ham_k(:, :, :)
real(kind=dp), intent(inout), allocatable :: wannier_centres_translated(:, :)
integer, intent(inout), allocatable :: irvec(:, :)
integer, intent(inout), allocatable :: ndegen(:)
integer, intent(inout), allocatable :: shift_vec(:, :)
integer, intent(inout) :: nrpts
integer, intent(inout) :: rpt_origin
type(atom_data_type), intent(in) :: atom_data
type(band_plot_type), intent(in) :: band_plot
type(dis_manifold_type), intent(in) :: dis_manifold
type(fermi_surface_plot_type), intent(in) :: fermi_surface_plot
type(kmesh_info_type), intent(in) :: kmesh_info
type(kpoint_path_type), intent(in) :: kpoint_path
type(output_file_type), intent(in) :: output_file
type(print_output_type), intent(in) :: print_output
type(w90_calculation_type), intent(in) :: w90_calculation
type(w90_comm_type), intent(in) :: comm
type(w90_system_type), intent(in) :: w90_system
type(wannier_data_type), intent(in) :: wannier_data
type(wannier_plot_type), intent(in) :: wannier_plot
type(ws_region_type), intent(in) :: ws_region
type(wvfn_read_type), intent(in) :: wvfn_read
complex(kind=dp), intent(in) :: m_matrix(:, :, :, :)
complex(kind=dp), intent(in) :: u_matrix_opt(:, :, :)
complex(kind=dp), intent(in) :: u_matrix(:, :, :)
real(kind=dp), intent(in), allocatable :: fermi_energy_list(:)
real(kind=dp), intent(in) :: bohr
real(kind=dp), intent(in) :: eigval(:, :)
real(kind=dp), intent(in) :: kpt_latt(:, :)
real(kind=dp), intent(in) :: real_lattice(3, 3)
integer, intent(in) :: mp_grid(3)
integer, intent(in) :: num_bands
integer, intent(in) :: num_kpts
integer, intent(in) :: num_wann
integer, intent(in) :: stdout
integer, intent(in) :: dist_k(:)
character(len=20), intent(in) :: transport_mode
character(len=50), intent(in) :: seedname
logical, intent(in) :: have_disentangled
logical, intent(in) :: lsitesymmetry
! local variables
type(ws_distance_type) :: ws_distance
real(kind=dp) :: recip_lattice(3, 3), volume
integer :: nkp, bands_num_spec_points, my_node_id, i
logical :: have_gamma
logical :: on_root = .false.
my_node_id = mpirank(comm)
if (my_node_id == 0) on_root = .true.
call utility_recip_lattice_base(real_lattice, recip_lattice, volume)
! setup RS cluster for calculations that need it
if (output_file%write_hr .or. &
output_file%write_hr_diag .or. &
output_file%write_r2mn .or. &
output_file%write_rmn .or. &
output_file%write_tb .or. &
w90_calculation%bands_plot .or. &
w90_calculation%fermi_surface_plot) 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, band_plot%mode, stdout, &
timer, error, transport_mode, comm)
if (allocated(error)) return
endif
! setup RS Hamilton eqn for calculations that need it
if (output_file%write_hr .or. &
output_file%write_hr_diag .or. &
output_file%write_tb .or. &
w90_calculation%bands_plot .or. &
w90_calculation%fermi_surface_plot) then
! Check if the kmesh includes the gamma point
have_gamma = .false.
do nkp = 1, num_kpts
if (all(abs(kpt_latt(:, nkp)) < eps6)) have_gamma = .true.
end do
if (.not. have_gamma) then
write (stdout, '(1x,a)') '!!!! Kpoint grid does not include Gamma. '// &
' Interpolation may be incorrect. !!!!'
endif
! Transform Hamiltonian to WF basis
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
endif
if (on_root) then
if (print_output%timing_level > 0) call io_stopwatch_start('plot: main', timer)
if (print_output%iprint > 0) then
! Print the header only if there is something to plot
if (output_file%write_hr .or. &
output_file%write_r2mn .or. &
output_file%write_rmn .or. &
output_file%write_tb .or. &
output_file%write_u_matrices .or. &
w90_calculation%bands_plot .or. &
w90_calculation%fermi_surface_plot .or. &
w90_calculation%wannier_plot) then
write (stdout, '(1x,a)') '*---------------------------------------------------------------------------*'
write (stdout, '(1x,a)') '| PLOTTING |'
write (stdout, '(1x,a)') '*---------------------------------------------------------------------------*'
write (stdout, *)
endif
endif
if (w90_calculation%fermi_surface_plot) then
! plot_fermi_surface can be trivially parallelised--fixme
call plot_fermi_surface(fermi_energy_list, recip_lattice, fermi_surface_plot, num_wann, &
ham_r, irvec, ndegen, nrpts, print_output%timing_level, stdout, &
seedname, timer, error, comm)
if (allocated(error)) return
endif
if (output_file%write_hr .or. output_file%write_tb) then
call ws_translate_dist(ws_distance, ws_region, num_wann, &
wannier_data%centres, real_lattice, mp_grid, nrpts, irvec, &
error, comm, force_recompute=.false.)
if (allocated(error)) return
call ws_write_vec(ws_distance, nrpts, irvec, num_wann, ws_region%use_ws_distance, &
seedname, error, comm)
if (allocated(error)) return
end if
! calculate and write projection of WFs on original bands in outer window
! only meaningful in disentanglement case (and otherwise lwindow is not available)
if (output_file%write_proj .and. num_bands > num_wann) then
call plot_calc_projection(num_bands, num_wann, num_kpts, u_matrix_opt, eigval, &
dis_manifold%lwindow, print_output%timing_level, &
print_output%iprint, stdout, timer)
if (allocated(error)) return
endif
if (output_file%write_bvec) then
call plot_bvec(kmesh_info, num_kpts, seedname, error, comm)
if (allocated(error)) return
endif
if (output_file%write_hr) then
! this is a trivial matrix write; no need to parallelize
call hamiltonian_write_hr(ham_r, irvec, ndegen, nrpts, num_wann, &
print_output%timing_level, seedname, timer, error, comm)
if (allocated(error)) return
endif
if (output_file%write_hr_diag) then
if (print_output%iprint > 0) then
write (stdout, *)
write (stdout, '(1x,a)') 'On-site Hamiltonian matrix elements'
write (stdout, '(3x,a)') ' n <0n|H|0n> (eV)'
write (stdout, '(3x,a)') '-------------------------'
do i = 1, num_wann
write (stdout, '(3x,i3,5x,f12.6)') i, real(ham_r(i, i, rpt_origin), kind=dp)
enddo
write (stdout, *)
endif
endif
if (output_file%write_u_matrices) then
call plot_u_matrices(u_matrix_opt, u_matrix, kpt_latt, dis_manifold, have_disentangled, &
num_wann, num_kpts, num_bands, seedname)
endif
if (output_file%write_vdw_data) then
! aam: write data required for vdW utility
call plot_write_vdw_data(num_wann, wannier_data, real_lattice, u_matrix, u_matrix_opt, &
have_disentangled, w90_system, error, comm, stdout, seedname)
if (allocated(error)) return
endif
if (output_file%write_xyz) then
call plot_write_xyz(real_space_ham%translate_home_cell, num_wann, wannier_data%centres, &
real_lattice, atom_data, print_output, error, comm, stdout, seedname)
if (allocated(error)) return
endif
end if !on_root
if (w90_calculation%bands_plot) then
bands_num_spec_points = 0
if (allocated(kpoint_path%labels)) bands_num_spec_points = size(kpoint_path%labels)
call plot_interpolate_bands(mp_grid, real_lattice, band_plot, kpoint_path, real_space_ham, &
ws_region, print_output, recip_lattice, num_wann, wannier_data, &
ham_r, irvec, ndegen, nrpts, wannier_centres_translated, &
ws_distance, bands_num_spec_points, stdout, seedname, timer, &
error, comm)
if (allocated(error)) return
endif
if (output_file%svd_omega) then
call plot_svd_omega_i(num_wann, num_kpts, kmesh_info, m_matrix, print_output, timer, dist_k, &
error, comm, stdout)
if (allocated(error)) return
endif
if (output_file%write_rmn) then
! parallel write_rmn
call plot_write_rmn(kmesh_info, m_matrix, kpt_latt, irvec, nrpts, num_kpts, num_wann, &
seedname, dist_k, error, comm)
if (allocated(error)) return
endif
if (output_file%write_r2mn) then
! write matrix elements <m|r^2|n> to file
call plot_write_r2mn(num_kpts, num_wann, kmesh_info, m_matrix, seedname, dist_k, error, comm)
if (allocated(error)) return
endif
if (output_file%write_tb) then
call hamiltonian_write_tb(kmesh_info, ham_r, m_matrix, kpt_latt, real_lattice, irvec, &
ndegen, nrpts, num_kpts, num_wann, print_output%timing_level, &
seedname, timer, dist_k, error, comm)
if (allocated(error)) return
endif
if (w90_calculation%wannier_plot) then
call plot_wannier(wannier_plot, wvfn_read, wannier_data, print_output, u_matrix_opt, &
dis_manifold, real_lattice, atom_data, kpt_latt, u_matrix, num_kpts, &
num_bands, num_wann, have_disentangled, w90_system%spinors, bohr, stdout, &
seedname, timer, dist_k, error, comm)
if (allocated(error)) return
endif
if (on_root .and. print_output%timing_level > 0) call io_stopwatch_stop('plot: main', timer)
end subroutine plot_main
!-----------------------------------------------------------------
!----------------- Private Routines ------------------------------
!-----------------------------------------------------------------
!================================================!
subroutine plot_interpolate_bands(mp_grid, real_lattice, band_plot, kpoint_path, real_space_ham, &
ws_region, print_output, recip_lattice, num_wann, &
wannier_data, ham_r, irvec, ndegen, nrpts, &
wannier_centres_translated, ws_distance, &
bands_num_spec_points, stdout, seedname, timer, error, comm)
!================================================!
! !
!! Plots the interpolated band structure
! !
!================================================!
use w90_constants, only: dp, cmplx_0, twopi
use w90_io, only: io_time, io_stopwatch_start, io_stopwatch_stop
use w90_ws_distance, only: ws_translate_dist
use w90_utility, only: utility_metric
use w90_types, only: wannier_data_type, kpoint_path_type, print_output_type, ws_region_type, &
ws_distance_type, timer_list_type
use w90_wannier90_types, only: band_plot_type, real_space_ham_type
use w90_error, only: w90_error_type, set_error_alloc, set_error_dealloc, set_error_fatal, &
set_error_warn
implicit none
! arguments
type(band_plot_type), intent(in) :: band_plot
type(kpoint_path_type), intent(in) :: kpoint_path
type(print_output_type), intent(in) :: print_output
type(real_space_ham_type), intent(in) :: real_space_ham
type(wannier_data_type), intent(in) :: wannier_data
type(ws_distance_type), intent(inout) :: ws_distance
type(ws_region_type), intent(in) :: ws_region
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) :: nrpts
integer, intent(in) :: ndegen(:)
integer, intent(in) :: irvec(:, :)
integer, intent(in) :: mp_grid(3)
integer, intent(in) :: num_wann
integer, intent(in) :: bands_num_spec_points
integer, intent(in) :: stdout
real(kind=dp), intent(in) :: real_lattice(3, 3)
real(kind=dp), intent(in) :: recip_lattice(3, 3)
real(kind=dp), intent(in) :: wannier_centres_translated(:, :)
complex(kind=dp), intent(in) :: ham_r(:, :, :)
character(len=50), intent(in) :: seedname
! local variables
integer, allocatable :: irvec_cut(:, :)
integer :: irvec_max(3)
integer :: nrpts_cut
integer, allocatable :: iwork(:), ifail(:)
integer :: info, i, j
integer :: irpt, nfound, loop_kpt, counter
integer :: loop_spts, total_pts, loop_i, nkp, ideg
integer :: num_paths, num_spts, ierr
integer :: bndunit, gnuunit, loop_w, loop_p
integer, allocatable :: kpath_pts(:)
integer, allocatable :: idx_special_points(:)
integer, allocatable :: label_idx_special_points(:)
real(kind=dp) :: rdotk, vec(3), emin, emax, time0
real(kind=dp), allocatable :: kpath_len(:)
real(kind=dp), allocatable :: rwork(:)
real(kind=dp), allocatable :: xval(:)
real(kind=dp), allocatable :: eig_int(:, :), plot_kpoint(:, :)
real(kind=dp), allocatable :: bands_proj(:, :)
real(kind=dp), allocatable :: xval_special_points(:)
real(kind=dp) :: recip_metric(3, 3)
complex(kind=dp) :: fac
complex(kind=dp), allocatable :: ham_r_cut(:, :, :)
complex(kind=dp), allocatable :: ham_pack(:)
complex(kind=dp), allocatable :: ham_kprm(:, :)
complex(kind=dp), allocatable :: U_int(:, :)
complex(kind=dp), allocatable :: cwork(:)
logical, allocatable :: kpath_print_first_point(:)
character(len=20), allocatable :: glabel(:)
character(len=20), allocatable :: xlabel(:)
character(len=20), allocatable :: ctemp(:)
! mpi variables
integer :: my_node_id, num_nodes
logical :: on_root
integer, allocatable :: counts(:)
integer, allocatable :: displs(:)
num_nodes = mpisize(comm)
my_node_id = mpirank(comm)
on_root = .false.
if (my_node_id == 0) on_root = .true.
allocate (counts(0:num_nodes - 1))
allocate (displs(0:num_nodes - 1))
!
if (on_root) then
if (print_output%timing_level > 1) then
call io_stopwatch_start('plot: interpolate_bands', timer)
endif
!
time0 = io_time()
write (stdout, *)
write (stdout, '(1x,a)') 'Calculating interpolated band-structure'
write (stdout, *)
endif ! on_root
call utility_metric(recip_lattice, recip_metric)
!
allocate (ham_pack((num_wann*(num_wann + 1))/2), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating ham_pack in plot_interpolate_bands', comm)
return
endif
allocate (ham_kprm(num_wann, num_wann), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating ham_kprm in plot_interpolate_bands', comm)
return
endif
allocate (U_int(num_wann, num_wann), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating U_int in plot_interpolate_bands', comm)
return
endif
allocate (cwork(2*num_wann), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating cwork in plot_interpolate_bands', comm)
return
endif
allocate (rwork(7*num_wann), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating rwork in plot_interpolate_bands', comm)
return
endif
allocate (iwork(5*num_wann), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating iwork in plot_interpolate_bands', comm)
return
endif
allocate (ifail(num_wann), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating ifail in plot_interpolate_bands', comm)
return
endif
if (.not. kpoint_path%bands_kpt_explicit) then
allocate (kpath_len(bands_num_spec_points/2), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating kpath_len in plot_interpolate_bands', comm)
return
endif
allocate (kpath_pts(bands_num_spec_points/2), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating kpath_pts in plot_interpolate_bands', comm)
return
endif
allocate (kpath_print_first_point(bands_num_spec_points/2), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating kpath_print_first_point in plot_interpolate_bands', comm)
return
endif
allocate (idx_special_points(bands_num_spec_points), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating idx_special_points in plot_interpolate_bands', comm)
return
endif
allocate (xval_special_points(bands_num_spec_points), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating xval_special_points in plot_interpolate_bands', comm)
return
endif
idx_special_points = -1
xval_special_points = -1._dp
end if
!
! Work out how many points in the total path and the positions of the special points
!
if (kpoint_path%bands_kpt_explicit) then
total_pts = size(kpoint_path%bands_kpt_frac, 2)
! Count the total number of special points
num_spts = 0
do i = 1, total_pts
do j = 1, bands_num_spec_points
if (sum((kpoint_path%bands_kpt_frac(:, i) - kpoint_path%points(:, j))**2) <= 1.e-6) then
num_spts = num_spts + 1
exit
end if
end do
end do
else
num_paths = bands_num_spec_points/2
kpath_print_first_point = .false.
! Loop over paths, set to False print_first_point if the starting point
! is the same as the ending point of the previous path.
! I skip the first path for which I always want to print the first point.
kpath_print_first_point(1) = .true.
do i = 2, num_paths
! If either the coordinates are different or the label is different, compute again the point
! (it will end up at the same x coordinate)
if ((SUM((kpoint_path%points(:, (i - 1)*2) - &
kpoint_path%points(:, (i - 1)*2 + 1))**2) > 1.e-6) .or. &
(TRIM(kpoint_path%labels((i - 1)*2)) .ne. &
TRIM(kpoint_path%labels((i - 1)*2 + 1)))) then
kpath_print_first_point(i) = .true.
end if
enddo
! Count the total number of special points
num_spts = num_paths
do i = 1, num_paths
if (kpath_print_first_point(i)) num_spts = num_spts + 1
end do
do loop_spts = 1, num_paths
vec = kpoint_path%points(:, 2*loop_spts) - kpoint_path%points(:, 2*loop_spts - 1)
kpath_len(loop_spts) = sqrt(dot_product(vec, (matmul(recip_metric, vec))))
if (loop_spts == 1) then
kpath_pts(loop_spts) = kpoint_path%num_points_first_segment
else
kpath_pts(loop_spts) = nint(real(kpoint_path%num_points_first_segment, dp) &
*kpath_len(loop_spts)/kpath_len(1))
! At least 1 point
!if (kpath_pts(loop_spts) .eq. 0) kpath_pts(loop_spts) = 1
end if
end do
total_pts = sum(kpath_pts)
do i = 1, num_paths
if (kpath_print_first_point(i)) total_pts = total_pts + 1
end do
end if
allocate (plot_kpoint(3, total_pts), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating plot_kpoint in plot_interpolate_bands', comm)
return
endif
allocate (xval(total_pts), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating xval in plot_interpolate_bands', comm)
return
endif
allocate (eig_int(num_wann, total_pts), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating eig_int in plot_interpolate_bands', comm)
return
endif
allocate (bands_proj(num_wann, total_pts), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating bands_proj in plot_interpolate_bands', comm)
return
endif
allocate (glabel(num_spts), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating glabel in plot_interpolate_bands', comm)
return
endif
allocate (xlabel(num_spts), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating xlabel in plot_interpolate_bands', comm)
return
endif
allocate (ctemp(bands_num_spec_points), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating ctemp in plot_interpolate_bands', comm)
return
endif
eig_int = 0.0_dp; bands_proj = 0.0_dp
!
! Find the position of each kpoint in the path
!
if (kpoint_path%bands_kpt_explicit) then
allocate (idx_special_points(num_spts), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating idx_special_points in plot_interpolate_bands', comm)
return
endif
allocate (xval_special_points(num_spts), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating xval_special_points in plot_interpolate_bands', comm)
return
endif
idx_special_points = -1
xval_special_points = -1._dp
allocate (label_idx_special_points(num_spts), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating label_idx_special_points in plot_interpolate_bands', comm)
return
endif
plot_kpoint(:, :) = kpoint_path%bands_kpt_frac(:, :)
xval = 0.0_dp
counter = 0
do i = 1, total_pts
if (i > 1) then
vec = plot_kpoint(:, i) - plot_kpoint(:, i - 1)
xval(i) = xval(i - 1) + sqrt(dot_product(vec, (matmul(recip_metric, vec))))
end if
do j = 1, bands_num_spec_points
if (sum((kpoint_path%bands_kpt_frac(:, i) - kpoint_path%points(:, j))**2) <= 1.e-6) then
counter = counter + 1
idx_special_points(counter) = i
label_idx_special_points(counter) = j
xval_special_points(counter) = xval(i)
if (counter > 1) then
if (idx_special_points(counter) == idx_special_points(counter - 1) + 1) then
! If the two points are consecutive, it means that the x coordinate should be the same
xval(i) = xval(i - 1)
xval_special_points(counter) = xval(i)
end if
end if
exit
end if
end do
end do
else
counter = 0
do loop_spts = 1, num_paths
if (kpath_print_first_point(loop_spts)) then
counter = counter + 1
if (counter == 1) then
xval(counter) = 0.0_dp
else
! If we are printing the first point in a path,
! It means that the coordinate did not change (otherwise
! we would not be printing it). Therefore I do not move
! on the x axis, there was a jump in the path here.
xval(counter) = xval(counter - 1)
endif
plot_kpoint(:, counter) = kpoint_path%points(:, 2*loop_spts - 1)
idx_special_points(2*loop_spts - 1) = counter
xval_special_points(2*loop_spts - 1) = xval(counter)
end if
! This is looping on all points but the first (1 is the first point
! after the first in the path)
do loop_i = 1, kpath_pts(loop_spts)
counter = counter + 1
! Set xval, the x position on the path of the current path
if (counter == 1) then
! This case should never happen but I keep it in for "safety"
xval(counter) = 0.0_dp
else
xval(counter) = xval(counter - 1) + kpath_len(loop_spts)/real(kpath_pts(loop_spts), dp)
endif
plot_kpoint(:, counter) = kpoint_path%points(:, 2*loop_spts - 1) + &
(kpoint_path%points(:, 2*loop_spts) &
- kpoint_path%points(:, 2*loop_spts - 1))* &
(real(loop_i, dp)/real(kpath_pts(loop_spts), dp))
end do
idx_special_points(2*loop_spts) = counter
xval_special_points(2*loop_spts) = xval(counter)
end do
!xval(total_pts)=sum(kpath_len)
plot_kpoint(:, total_pts) = kpoint_path%points(:, bands_num_spec_points)
end if
!
! Write out the kpoints in the path
!
if (on_root) then
open (newunit=bndunit, file=trim(seedname)//'_band.kpt', form='formatted')
write (bndunit, *) total_pts
do loop_spts = 1, total_pts
write (bndunit, '(3f12.6,3x,a)') (plot_kpoint(loop_i, loop_spts), loop_i=1, 3), "1.0"
end do
close (bndunit)
!
! Write out information on high-symmetry points in the path
!
open (newunit=bndunit, file=trim(seedname)//'_band.labelinfo.dat', form='formatted')
if (kpoint_path%bands_kpt_explicit) then
do loop_spts = 1, num_spts
write (bndunit, '(a,3x,I10,3x,4f18.10)') &
kpoint_path%labels(label_idx_special_points(loop_spts)), &
idx_special_points(loop_spts), &
xval_special_points(loop_spts), &
(plot_kpoint(loop_i, idx_special_points(loop_spts)), loop_i=1, 3)
end do
else
do loop_spts = 1, bands_num_spec_points
if ((MOD(loop_spts, 2) .eq. 1) .and. &
(kpath_print_first_point((loop_spts + 1)/2) .eqv. .false.)) cycle
write (bndunit, '(a,3x,I10,3x,4f18.10)') &
kpoint_path%labels(loop_spts), &
idx_special_points(loop_spts), &
xval_special_points(loop_spts), &
(plot_kpoint(loop_i, idx_special_points(loop_spts)), loop_i=1, 3)
end do
end if
close (bndunit)
endif ! on_root
!
! Cut H matrix in real-space
!
if (index(band_plot%mode, 'cut') .ne. 0) then
call plot_cut_hr(real_space_ham, real_lattice, mp_grid, num_wann, &
wannier_centres_translated, stdout, error)
if (allocated(error)) return
endif
!
! Interpolate the Hamiltonian at each kpoint
!
if (ws_region%use_ws_distance) then
if (index(band_plot%mode, 's-k') .ne. 0) then
call ws_translate_dist(ws_distance, ws_region, num_wann, &
wannier_data%centres, real_lattice, mp_grid, nrpts, &
irvec, error, comm, force_recompute=.true.)
if (allocated(error)) return
elseif (index(band_plot%mode, 'cut') .ne. 0) then
call ws_translate_dist(ws_distance, ws_region, num_wann, &
wannier_data%centres, real_lattice, mp_grid, nrpts_cut, &
irvec_cut, error, comm, force_recompute=.true.)
if (allocated(error)) return
else
call set_error_warn(error, 'Error in plot_interpolate bands: value of bands_plot_mode not recognised', comm)
return
endif
endif
if (on_root .and. print_output%timing_level > 2) then
call io_stopwatch_start('plot: interpolate_bands: loop_kpoints', timer)
endif
! [lp] the s-k and cut codes are very similar when use_ws_distance is used, a complete
! merge after this point is not impossible
call comms_array_split(total_pts, counts, displs, comm)
! Don't worry about serial run, it is OK to call comms_array_split!
do loop_kpt = displs(my_node_id) + 1, displs(my_node_id) + counts(my_node_id)
! do loop_kpt = 1, total_pts
ham_kprm = cmplx_0
!
if (index(band_plot%mode, 's-k') .ne. 0) then
do irpt = 1, nrpts
! [lp] Shift the WF to have the minimum distance IJ, see also ws_distance.F90
if (ws_region%use_ws_distance) then
do j = 1, num_wann
do i = 1, num_wann
do ideg = 1, ws_distance%ndeg(i, j, irpt)
rdotk = twopi*dot_product(plot_kpoint(:, loop_kpt), &
real(ws_distance%irdist(:, ideg, i, j, irpt), dp))
fac = cmplx(cos(rdotk), sin(rdotk), dp) &
/real(ndegen(irpt)*ws_distance%ndeg(i, j, irpt), dp)
ham_kprm(i, j) = ham_kprm(i, j) + fac*ham_r(i, j, irpt)
enddo
enddo
enddo
else
! [lp] Original code, without IJ-dependent shift:
rdotk = twopi*dot_product(plot_kpoint(:, loop_kpt), irvec(:, irpt))
fac = cmplx(cos(rdotk), sin(rdotk), dp)/real(ndegen(irpt), dp)
ham_kprm = ham_kprm + fac*ham_r(:, :, irpt)
endif
end do
! end of s-k mode
elseif (index(band_plot%mode, 'cut') .ne. 0) then
do irpt = 1, nrpts_cut
! [lp] Shift the WF to have the minimum distance IJ, see also ws_distance.F90
if (ws_region%use_ws_distance) then
do j = 1, num_wann
do i = 1, num_wann
do ideg = 1, ws_distance%ndeg(i, j, irpt)
rdotk = twopi*dot_product(plot_kpoint(:, loop_kpt), &
real(ws_distance%irdist(:, ideg, i, j, irpt), dp))
fac = cmplx(cos(rdotk), sin(rdotk), dp)/real(ws_distance%ndeg(i, j, irpt), dp)
ham_kprm(i, j) = ham_kprm(i, j) + fac*ham_r_cut(i, j, irpt)
enddo
enddo
enddo
! [lp] Original code, without IJ-dependent shift:
else
rdotk = twopi*dot_product(plot_kpoint(:, loop_kpt), irvec_cut(:, irpt))
!~[aam] check divide by ndegen?
fac = cmplx(cos(rdotk), sin(rdotk), dp)
ham_kprm = ham_kprm + fac*ham_r_cut(:, :, irpt)
endif ! end of use_ws_distance
end do
endif ! end of "cut" mode
!
! Diagonalise H_k (->basis of eigenstates)
do j = 1, num_wann
do i = 1, j
ham_pack(i + ((j - 1)*j)/2) = ham_kprm(i, j)
enddo
enddo
call ZHPEVX('V', 'A', 'U', num_wann, ham_pack, 0.0_dp, 0.0_dp, 0, 0, -1.0_dp, &
nfound, eig_int(1, loop_kpt), U_int, num_wann, cwork, rwork, iwork, ifail, info)
if (info < 0) then
write (stdout, '(a,i3,a)') 'THE ', -info, ' ARGUMENT OF ZHPEVX HAD AN ILLEGAL VALUE'
call set_error_fatal(error, 'Error in plot_interpolate_bands', comm)
return
endif
if (info > 0) then
write (stdout, '(i3,a)') info, ' EIGENVECTORS FAILED TO CONVERGE'
call set_error_warn(error, 'Error in plot_interpolate_bands', comm)
return
endif
! Compute projection onto WF if requested
if (allocated(band_plot%project)) then
do loop_w = 1, num_wann
do loop_p = 1, num_wann
if (any(band_plot%project == loop_p)) then
bands_proj(loop_w, loop_kpt) = bands_proj(loop_w, loop_kpt) + &
abs(U_int(loop_p, loop_w))**2
end if
end do
end do
end if
!
end do
call comms_reduce(eig_int(1, 1), num_wann*total_pts, 'SUM', error, comm)
call comms_reduce(bands_proj(1, 1), num_wann*total_pts, 'SUM', error, comm)
if (on_root .and. print_output%timing_level > 2) then
call io_stopwatch_stop('plot: interpolate_bands: loop_kpoints', timer)
endif
!
! Interpolation Finished!
! Now we write plotting files
!
if (on_root) then
emin = minval(eig_int) - 1.0_dp
emax = maxval(eig_int) + 1.0_dp
if (index(band_plot%format, 'gnu') > 0) then
call plot_interpolate_gnuplot(band_plot, kpoint_path, bands_num_spec_points, num_wann)
endif
if (index(band_plot%format, 'xmgr') > 0) then
call plot_interpolate_xmgrace(kpoint_path, bands_num_spec_points, num_wann, error)
endif
write (stdout, '(1x,a,f11.3,a)') &
'Time to calculate interpolated band structure ', io_time() - time0, ' (sec)'
write (stdout, *)
endif
if (allocated(ham_r_cut)) then
deallocate (ham_r_cut, stat=ierr)
if (ierr /= 0) then
call set_error_dealloc(error, 'Error in deallocating ham_r_cut in plot_interpolate_bands', comm)
return
endif
endif
if (allocated(irvec_cut)) then
deallocate (irvec_cut, stat=ierr)
if (ierr /= 0) then
call set_error_dealloc(error, 'Error in deallocating irvec_cut in plot_interpolate_bands', comm)
return
endif
endif
if (print_output%timing_level > 1) call io_stopwatch_stop('plot: interpolate_bands', timer)
if (allocated(idx_special_points)) then
deallocate (idx_special_points, stat=ierr)
if (ierr /= 0) then
call set_error_dealloc(error, 'Error in deallocating idx_special_points in &
&plot_interpolate_bands', comm)
return
endif
endif
if (allocated(xval_special_points)) then
deallocate (xval_special_points, stat=ierr)
if (ierr /= 0) then
call set_error_dealloc(error, 'Error in deallocating xval_special_points in &
&plot_interpolate_bands', comm)
return
endif
endif
if (allocated(label_idx_special_points)) then
deallocate (label_idx_special_points, stat=ierr)
if (ierr /= 0) then
call set_error_dealloc(error, 'Error in deallocating label_idx_special_points in &
&plot_interpolate_bands', comm)
return
endif
endif
contains
!================================================!
subroutine plot_cut_hr(real_space_ham, real_lattice, mp_grid, num_wann, &
wannier_centres_translated, stdout, error)
!================================================!
!
!! In real-space picture, ham_r(j,i,k) is an interaction between
!! j_th WF at 0 and i_th WF at the lattice point translated
!! by matmul(real_lattice(:,:),irvec(:,k))
!! We truncate Hamiltonian matrix when
!! 1) | r_i(0) - r_j (R) | > dist_cutoff
!! 2) | ham_r(i,j,k) | < hr_cutoff
!! while the condition 1) is essential to get a meaningful band structure,
!! ( dist_cutoff must be smaller than the shortest distance from
!! the center of W-S supercell to the points at the cell boundaries )
!! the condition 2) is optional.
!!
!! limitation: when bands_plot_dim .ne. 3
!! one_dim_vec must be parallel to one of the cartesian axis
!! and perpendicular to the other two primitive lattice vectors
!================================================!
use w90_constants, only: dp, cmplx_0, eps8
use w90_wannier90_types, only: band_plot_type, real_space_ham_type
use w90_error, only: w90_error_type, set_error_alloc, set_error_warn
implicit none
! arguments
type(real_space_ham_type), intent(in) :: real_space_ham
type(w90_error_type), allocatable, intent(out) :: error
real(kind=dp), intent(in) :: real_lattice(3, 3)
real(kind=dp), intent(in) :: wannier_centres_translated(:, :)
integer, intent(in) :: mp_grid(3)
integer, intent(in) :: num_wann
integer, intent(in) :: stdout
! local variables
integer :: nrpts_tmp
integer :: one_dim_vec, two_dim_vec(2)
integer :: i, j, n1, n2, n3, i1, i2, i3
real(kind=dp), allocatable :: ham_r_tmp(:, :)
real(kind=dp), allocatable :: shift_vec(:, :)
real(kind=dp) :: dist_ij_vec(3)
real(kind=dp) :: dist_vec(3)
real(kind=dp) :: dist
allocate (ham_r_tmp(num_wann, num_wann), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating ham_r_tmp in plot_cut_hr', comm)
return
endif
irvec_max = maxval(irvec, DIM=2) + 1
if (real_space_ham%system_dim .ne. 3) then
! Find one_dim_vec which is parallel to one_dim_dir
! two_dim_vec - the other two lattice vectors
! Along the confined directions, take only irvec=0
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
call set_error_warn(error, 'Error: 1-d lattice vector not defined in plot_cut_hr', comm)
return
endif
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
if (real_space_ham%system_dim .eq. 1) then
irvec_max(two_dim_vec(1)) = 0
irvec_max(two_dim_vec(2)) = 0
end if
if (real_space_ham%system_dim .eq. 2) irvec_max(one_dim_vec) = 0
end if
nrpts_cut = (2*irvec_max(1) + 1)*(2*irvec_max(2) + 1)*(2*irvec_max(3) + 1)
allocate (ham_r_cut(num_wann, num_wann, nrpts_cut), stat=ierr)
if (ierr /= 0) then
call set_error_alloc(error, 'Error in allocating ham_r_cut in plot_cut_hr', comm)
return
endif