-
Notifications
You must be signed in to change notification settings - Fork 145
/
Copy pathpostw90_common.F90
1508 lines (1332 loc) · 60.3 KB
/
postw90_common.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 !
!------------------------------------------------------------!
module w90_postw90_common
!==============================================================================
!! This contains the common variables and procedures needed to set up a Wannier
!! interpolatation calculation for any physical property
!==============================================================================
! Should we remove this 'use w90_comms' and invoke in individual routines
! when needed?
!
use w90_comms
use w90_constants, only: dp
implicit none
private
public :: pw90common_wanint_setup, pw90common_wanint_get_kpoint_file, pw90common_wanint_param_dist
public :: pw90common_wanint_data_dist, pw90common_get_occ
public :: pw90common_fourier_R_to_k, pw90common_fourier_R_to_k_new, pw90common_fourier_R_to_k_vec
public :: nrpts, rpt_origin, v_matrix, ndegen, irvec, crvec
public :: num_int_kpts_on_node, int_kpts, weight
public :: pw90common_kmesh_spacing
public :: pw90common_fourier_R_to_k_new_second_d, pw90common_fourier_R_to_k_new_second_d_TB_conv, &
pw90common_fourier_R_to_k_vec_dadb, pw90common_fourier_R_to_k_vec_dadb_TB_conv
! AAM PROBABLY REMOVE THIS
! This 'save' statement could probably be ommited, since this module
! is USEd by the main program 'wannier_parint'
!
save
! AAM REMOVE THIS
! Default accessibility is PUBLIC
!
! private :: wigner_seitz
!
! private :: kmesh_spacing_singleinteger, kmesh_spacing_mesh
interface pw90common_kmesh_spacing
module procedure kmesh_spacing_singleinteger
module procedure kmesh_spacing_mesh
end interface pw90common_kmesh_spacing
! Parameters describing the direct lattice points R on a
! Wigner-Seitz supercell
!
integer, allocatable :: irvec(:, :)
real(kind=dp), allocatable :: crvec(:, :)
integer, allocatable :: ndegen(:)
integer :: nrpts
integer :: rpt_origin
integer :: max_int_kpts_on_node, num_int_kpts
integer, allocatable :: num_int_kpts_on_node(:)
real(kind=dp), allocatable :: int_kpts(:, :), weight(:)
complex(kind=dp), allocatable :: v_matrix(:, :, :)
contains
!===========================================================!
! PUBLIC PROCEDURES !
!===========================================================!
! Public procedures have names starting with wanint_
subroutine pw90common_wanint_setup
!! Setup data ready for interpolation
use w90_constants, only: dp, cmplx_0
use w90_io, only: io_error, io_file_unit, stdout, seedname
use w90_utility, only: utility_cart_to_frac
use w90_parameters, only: real_lattice, effective_model, num_wann
integer :: ierr, i, j, k, ikpt, ir, file_unit, num_wann_loc
! Find nrpts, the number of points in the Wigner-Seitz cell
!
if (effective_model) then
if (on_root) then
! nrpts is read from file, together with num_wann
file_unit = io_file_unit()
open (file_unit, file=trim(seedname)//'_HH_R.dat', form='formatted', &
status='old', err=101)
read (file_unit, *) !header
read (file_unit, *) num_wann_loc
if (num_wann_loc /= num_wann) &
call io_error('Inconsistent values of num_wann in ' &
//trim(seedname)//'_HH_R.dat and '//trim(seedname)//'.win')
read (file_unit, *) nrpts
close (file_unit)
endif
call comms_bcast(nrpts, 1)
else
call wigner_seitz(count_pts=.true.)
endif
! Now can allocate several arrays
!
allocate (irvec(3, nrpts), stat=ierr)
if (ierr /= 0) call io_error('Error in allocating irvec in pw90common_wanint_setup')
irvec = 0
allocate (crvec(3, nrpts), stat=ierr)
if (ierr /= 0) call io_error('Error in allocating crvec in pw90common_wanint_setup')
crvec = 0.0_dp
allocate (ndegen(nrpts), stat=ierr)
if (ierr /= 0) call io_error('Error in allocating ndegen in pw90common_wanint_setup')
ndegen = 0
!
! Also rpt_origin, so that when effective_model=.true it is not
! passed to get_HH_R without being initialized.
rpt_origin = 0
! If effective_model, this is done in get_HH_R
if (.not. effective_model) then
!
! Set up the lattice vectors on the Wigner-Seitz supercell
! where the Wannier functions live
!
call wigner_seitz(count_pts=.false.)
!
! Convert from reduced to Cartesian coordinates
!
do ir = 1, nrpts
! Note that 'real_lattice' stores the lattice vectors as *rows*
crvec(:, ir) = matmul(transpose(real_lattice), irvec(:, ir))
end do
endif
return
101 call io_error('Error in pw90common_wanint_setup: problem opening file '// &
trim(seedname)//'_HH_R.dat')
end subroutine pw90common_wanint_setup
!===========================================================!
subroutine pw90common_wanint_get_kpoint_file
!===========================================================!
! !
!! read kpoints from kpoint.dat and distribute
! !
!===========================================================!
use w90_constants, only: dp
use w90_io, only: io_error, io_file_unit, &
io_date, io_time, io_stopwatch
integer :: k_unit
integer :: loop_nodes, loop_kpt, i, ierr
real(kind=dp) :: sum
k_unit = io_file_unit()
if (on_root) then
open (unit=k_unit, file='kpoint.dat', status='old', form='formatted', err=106)
read (k_unit, *) num_int_kpts
end if
call comms_bcast(num_int_kpts, 1)
allocate (num_int_kpts_on_node(0:num_nodes - 1))
num_int_kpts_on_node(:) = num_int_kpts/num_nodes
max_int_kpts_on_node = num_int_kpts - (num_nodes - 1)*(num_int_kpts/num_nodes)
num_int_kpts_on_node(0) = max_int_kpts_on_node
! if(my_node_id < num_int_kpts- num_int_kpts_on_node*num_nodes) num_int_kpts_on_node= num_int_kpts_on_node+1
allocate (int_kpts(3, max_int_kpts_on_node), stat=ierr)
if (ierr /= 0) call io_error('Error allocating max_int_kpts_on_node in param_read_um')
int_kpts = 0.0_dp
allocate (weight(max_int_kpts_on_node), stat=ierr)
if (ierr /= 0) call io_error('Error allocating weight in param_read_um')
weight = 0.0_dp
sum = 0.0_dp
if (on_root) then
do loop_nodes = 1, num_nodes - 1
do loop_kpt = 1, num_int_kpts_on_node(loop_nodes)
read (k_unit, *) (int_kpts(i, loop_kpt), i=1, 3), weight(loop_kpt)
sum = sum + weight(loop_kpt)
end do
call comms_send(int_kpts(1, 1), 3*num_int_kpts_on_node(loop_nodes), loop_nodes)
call comms_send(weight(1), num_int_kpts_on_node(loop_nodes), loop_nodes)
end do
do loop_kpt = 1, num_int_kpts_on_node(0)
read (k_unit, *) (int_kpts(i, loop_kpt), i=1, 3), weight(loop_kpt)
sum = sum + weight(loop_kpt)
end do
! print*,'rsum',sum
end if
if (.not. on_root) then
call comms_recv(int_kpts(1, 1), 3*num_int_kpts_on_node(my_node_id), root_id)
call comms_recv(weight(1), num_int_kpts_on_node(my_node_id), root_id)
end if
return
106 call io_error('Error: Problem opening file kpoint.dat in pw90common_wanint_get_kpoint_file')
end subroutine pw90common_wanint_get_kpoint_file
!===========================================================!
subroutine pw90common_wanint_param_dist
!===========================================================!
! !
!! distribute the parameters across processors
!! NOTE: we only send the ones postw90 uses, not all in w90
! !
!===========================================================!
use w90_constants, only: dp, cmplx_0, cmplx_i, twopi
use w90_io, only: io_error, io_file_unit, io_date, io_time, &
io_stopwatch
use w90_parameters
integer :: ierr
call comms_bcast(effective_model, 1)
call comms_bcast(eig_found, 1)
if (.not. effective_model) then
call comms_bcast(mp_grid(1), 3)
call comms_bcast(num_kpts, 1)
call comms_bcast(num_bands, 1)
endif
call comms_bcast(num_wann, 1)
call comms_bcast(timing_level, 1)
call comms_bcast(iprint, 1)
call comms_bcast(ws_distance_tol, 1)
call comms_bcast(ws_search_size(1), 3)
! call comms_bcast(num_atoms,1) ! Ivo: not used in postw90, right?
! call comms_bcast(num_species,1) ! Ivo: not used in postw90, right?
call comms_bcast(real_lattice(1, 1), 9)
call comms_bcast(recip_lattice(1, 1), 9)
call comms_bcast(real_metric(1, 1), 9)
call comms_bcast(recip_metric(1, 1), 9)
call comms_bcast(cell_volume, 1)
call comms_bcast(dos_energy_step, 1)
call comms_bcast(dos_adpt_smr, 1)
call comms_bcast(dos_smr_index, 1)
call comms_bcast(dos_kmesh_spacing, 1)
call comms_bcast(dos_kmesh(1), 3)
call comms_bcast(dos_adpt_smr_max, 1)
call comms_bcast(dos_smr_fixed_en_width, 1)
call comms_bcast(dos_adpt_smr_fac, 1)
call comms_bcast(num_dos_project, 1)
call comms_bcast(berry, 1)
call comms_bcast(berry_task, len(berry_task))
call comms_bcast(berry_kmesh_spacing, 1)
call comms_bcast(berry_kmesh(1), 3)
call comms_bcast(berry_curv_adpt_kmesh, 1)
call comms_bcast(berry_curv_adpt_kmesh_thresh, 1)
call comms_bcast(berry_curv_unit, len(berry_curv_unit))
! Tsirkin
call comms_bcast(gyrotropic, 1)
call comms_bcast(gyrotropic_task, len(gyrotropic_task))
call comms_bcast(gyrotropic_kmesh_spacing, 1)
call comms_bcast(gyrotropic_kmesh(1), 3)
call comms_bcast(gyrotropic_smr_fixed_en_width, 1)
call comms_bcast(gyrotropic_smr_index, 1)
call comms_bcast(gyrotropic_eigval_max, 1)
call comms_bcast(gyrotropic_nfreq, 1)
call comms_bcast(gyrotropic_degen_thresh, 1)
call comms_bcast(gyrotropic_num_bands, 1)
call comms_bcast(gyrotropic_box(1, 1), 9)
call comms_bcast(gyrotropic_box_corner(1), 3)
call comms_bcast(gyrotropic_smr_max_arg, 1)
call comms_bcast(gyrotropic_smr_fixed_en_width, 1)
call comms_bcast(gyrotropic_smr_index, 1)
call comms_bcast(spinors, 1)
call comms_bcast(shc_freq_scan, 1)
call comms_bcast(shc_alpha, 1)
call comms_bcast(shc_beta, 1)
call comms_bcast(shc_gamma, 1)
call comms_bcast(shc_bandshift, 1)
call comms_bcast(shc_bandshift_firstband, 1)
call comms_bcast(shc_bandshift_energyshift, 1)
call comms_bcast(kubo_adpt_smr, 1)
call comms_bcast(kubo_adpt_smr_fac, 1)
call comms_bcast(kubo_adpt_smr_max, 1)
call comms_bcast(kubo_smr_fixed_en_width, 1)
call comms_bcast(kubo_smr_index, 1)
call comms_bcast(kubo_eigval_max, 1)
call comms_bcast(kubo_nfreq, 1)
call comms_bcast(nfermi, 1)
call comms_bcast(dos_energy_min, 1)
call comms_bcast(dos_energy_max, 1)
call comms_bcast(spin_kmesh_spacing, 1)
call comms_bcast(spin_kmesh(1), 3)
call comms_bcast(wanint_kpoint_file, 1)
call comms_bcast(dis_win_min, 1)
call comms_bcast(dis_win_max, 1)
call comms_bcast(sc_eta, 1)
call comms_bcast(sc_w_thr, 1)
call comms_bcast(sc_phase_conv, 1)
! ----------------------------------------------
!
! New input variables in development
!
call comms_bcast(devel_flag, len(devel_flag))
call comms_bcast(spin_moment, 1)
call comms_bcast(spin_axis_polar, 1)
call comms_bcast(spin_axis_azimuth, 1)
call comms_bcast(spin_decomp, 1)
call comms_bcast(use_degen_pert, 1)
call comms_bcast(degen_thr, 1)
call comms_bcast(num_valence_bands, 1)
call comms_bcast(dos, 1)
call comms_bcast(dos_task, len(dos_task))
call comms_bcast(kpath, 1)
call comms_bcast(kpath_task, len(kpath_task))
call comms_bcast(kpath_bands_colour, len(kpath_bands_colour))
call comms_bcast(kslice, 1)
call comms_bcast(kslice_task, len(kslice_task))
call comms_bcast(kslice_corner(1), 3)
call comms_bcast(kslice_b1(1), 3)
call comms_bcast(kslice_b2(1), 3)
call comms_bcast(kslice_2dkmesh(1), 2)
call comms_bcast(kslice_fermi_lines_colour, len(kslice_fermi_lines_colour))
call comms_bcast(transl_inv, 1)
call comms_bcast(num_elec_per_state, 1)
call comms_bcast(scissors_shift, 1)
!
! Do these have to be broadcasted? (Plots done on root node only)
!
! call comms_bcast(bands_num_points,1)
! call comms_bcast(bands_num_spec_points,1)
! if(allocated(bands_spec_points)) &
! call comms_bcast(bands_spec_points(1,1),3*bands_num_spec_points)
! if(allocated(bands_label)) &
! call comms_bcast(bands_label(:),len(bands_label(1))*bands_num_spec_points)
! ----------------------------------------------
call comms_bcast(geninterp, 1)
call comms_bcast(geninterp_alsofirstder, 1)
call comms_bcast(geninterp_single_file, 1)
! [gp-begin, Apr 12, 2012]
! BoltzWann variables
call comms_bcast(boltzwann, 1)
call comms_bcast(boltz_calc_also_dos, 1)
call comms_bcast(boltz_2d_dir_num, 1)
call comms_bcast(boltz_dos_energy_step, 1)
call comms_bcast(boltz_dos_energy_min, 1)
call comms_bcast(boltz_dos_energy_max, 1)
call comms_bcast(boltz_dos_adpt_smr, 1)
call comms_bcast(boltz_dos_smr_fixed_en_width, 1)
call comms_bcast(boltz_dos_adpt_smr_fac, 1)
call comms_bcast(boltz_dos_adpt_smr_max, 1)
call comms_bcast(boltz_mu_min, 1)
call comms_bcast(boltz_mu_max, 1)
call comms_bcast(boltz_mu_step, 1)
call comms_bcast(boltz_temp_min, 1)
call comms_bcast(boltz_temp_max, 1)
call comms_bcast(boltz_temp_step, 1)
call comms_bcast(boltz_kmesh_spacing, 1)
call comms_bcast(boltz_kmesh(1), 3)
call comms_bcast(boltz_tdf_energy_step, 1)
call comms_bcast(boltz_relax_time, 1)
call comms_bcast(boltz_TDF_smr_fixed_en_width, 1)
call comms_bcast(boltz_TDF_smr_index, 1)
call comms_bcast(boltz_dos_smr_index, 1)
call comms_bcast(boltz_bandshift, 1)
call comms_bcast(boltz_bandshift_firstband, 1)
call comms_bcast(boltz_bandshift_energyshift, 1)
! [gp-end]
call comms_bcast(use_ws_distance, 1)
! These variables are different from the ones above in that they are
! allocatable, and in param_read they were allocated on the root node only
!
if (.not. on_root) then
allocate (fermi_energy_list(nfermi), stat=ierr)
if (ierr /= 0) call io_error( &
'Error allocating fermi_energy_read in postw90_param_dist')
allocate (kubo_freq_list(kubo_nfreq), stat=ierr)
if (ierr /= 0) call io_error( &
'Error allocating kubo_freq_list in postw90_param_dist')
allocate (gyrotropic_band_list(gyrotropic_num_bands), stat=ierr)
if (ierr /= 0) call io_error( &
'Error allocating gyrotropic_band_list in postw90_param_dist')
allocate (gyrotropic_freq_list(gyrotropic_nfreq), stat=ierr)
if (ierr /= 0) call io_error( &
'Error allocating gyrotropic_freq_list in postw90_param_dist')
allocate (dos_project(num_dos_project), stat=ierr)
if (ierr /= 0) &
call io_error('Error allocating dos_project in postw90_param_dist')
if (.not. effective_model) then
if (eig_found) then
allocate (eigval(num_bands, num_kpts), stat=ierr)
if (ierr /= 0) &
call io_error('Error allocating eigval in postw90_param_dist')
end if
allocate (kpt_latt(3, num_kpts), stat=ierr)
if (ierr /= 0) &
call io_error('Error allocating kpt_latt in postw90_param_dist')
endif
end if
if (nfermi > 0) call comms_bcast(fermi_energy_list(1), nfermi)
call comms_bcast(gyrotropic_freq_list(1), gyrotropic_nfreq)
call comms_bcast(gyrotropic_band_list(1), gyrotropic_num_bands)
call comms_bcast(kubo_freq_list(1), kubo_nfreq)
call comms_bcast(dos_project(1), num_dos_project)
if (.not. effective_model) then
if (eig_found) then
call comms_bcast(eigval(1, 1), num_bands*num_kpts)
end if
call comms_bcast(kpt_latt(1, 1), 3*num_kpts)
endif
! kmesh: only nntot,wb, and bk are needed to evaluate the WF matrix
! elements of the position operator in reciprocal space. For the
! extra matrix elements entering the orbital magnetization, also
! need nnlist. In principle could only broadcast those four variables
if (.not. effective_model) then
call comms_bcast(nnh, 1)
call comms_bcast(nntot, 1)
if (.not. on_root) then
allocate (nnlist(num_kpts, nntot), stat=ierr)
if (ierr /= 0) &
call io_error('Error in allocating nnlist in pw90common_wanint_param_dist')
allocate (neigh(num_kpts, nntot/2), stat=ierr)
if (ierr /= 0) &
call io_error('Error in allocating neigh in pw90common_wanint_param_dist')
allocate (nncell(3, num_kpts, nntot), stat=ierr)
if (ierr /= 0) &
call io_error('Error in allocating nncell in pw90common_wanint_param_dist')
allocate (wb(nntot), stat=ierr)
if (ierr /= 0) &
call io_error('Error in allocating wb in pw90common_wanint_param_dist')
allocate (bka(3, nntot/2), stat=ierr)
if (ierr /= 0) &
call io_error('Error in allocating bka in pw90common_wanint_param_dist')
allocate (bk(3, nntot, num_kpts), stat=ierr)
if (ierr /= 0) &
call io_error('Error in allocating bk in pw90common_wanint_param_dist')
end if
call comms_bcast(nnlist(1, 1), num_kpts*nntot)
call comms_bcast(neigh(1, 1), num_kpts*nntot/2)
call comms_bcast(nncell(1, 1, 1), 3*num_kpts*nntot)
call comms_bcast(wb(1), nntot)
call comms_bcast(bka(1, 1), 3*nntot/2)
call comms_bcast(bk(1, 1, 1), 3*nntot*num_kpts)
endif
end subroutine pw90common_wanint_param_dist
!===========================================================!
subroutine pw90common_wanint_data_dist
!===========================================================!
! !
!! Distribute the um and chk files
! !
!===========================================================!
use w90_constants, only: dp, cmplx_0, cmplx_i, twopi
use w90_io, only: io_error, io_file_unit, &
io_date, io_time, io_stopwatch
use w90_parameters, only: num_wann, num_kpts, num_bands, have_disentangled, &
u_matrix_opt, u_matrix, m_matrix, &
ndimwin, lwindow, nntot, wannier_centres, &
num_valence_bands, scissors_shift
implicit none
integer :: ierr, loop_kpt, m, i, j
if (.not. on_root) then
! wannier_centres is allocated in param_read, so only on root node
! It is then read in param_read_chpkt
! Therefore, now we need to allocate it on all nodes, and then broadcast it
allocate (wannier_centres(3, num_wann), stat=ierr)
if (ierr /= 0) call io_error('Error allocating wannier_centres in pw90common_wanint_data_dist')
end if
call comms_bcast(wannier_centres(1, 1), 3*num_wann)
! -------------------
! Ivo: added 8april11
! -------------------
!
! Calculate the matrix that describes the combined effect of
! disentanglement and maximal localization. This is the combination
! that is most often needed for interpolation purposes
!
! Allocate on all nodes
allocate (v_matrix(num_bands, num_wann, num_kpts), stat=ierr)
if (ierr /= 0) &
call io_error('Error allocating v_matrix in pw90common_wanint_data_dist')
! u_matrix and u_matrix_opt are stored on root only
if (on_root) then
if (.not. have_disentangled) then
v_matrix = u_matrix
else
v_matrix = cmplx_0
do loop_kpt = 1, num_kpts
do j = 1, num_wann
do m = 1, ndimwin(loop_kpt)
do i = 1, num_wann
v_matrix(m, j, loop_kpt) = v_matrix(m, j, loop_kpt) &
+ u_matrix_opt(m, i, loop_kpt)*u_matrix(i, j, loop_kpt)
enddo
enddo
enddo
enddo
endif
if (allocated(u_matrix_opt)) deallocate (u_matrix_opt)
if (.not. (num_valence_bands > 0 .and. abs(scissors_shift) > 1.0e-7_dp)) then
if (allocated(u_matrix)) deallocate (u_matrix)
endif
endif
call comms_bcast(v_matrix(1, 1, 1), num_bands*num_wann*num_kpts)
if (num_valence_bands > 0 .and. abs(scissors_shift) > 1.0e-7_dp) then
if (.not. on_root .and. .not. allocated(u_matrix)) then
allocate (u_matrix(num_wann, num_wann, num_kpts), stat=ierr)
if (ierr /= 0) &
call io_error('Error allocating u_matrix in pw90common_wanint_data_dist')
endif
call comms_bcast(u_matrix(1, 1, 1), num_wann*num_wann*num_kpts)
endif
! if (.not.on_root .and. .not.allocated(m_matrix)) then
! allocate(m_matrix(num_wann,num_wann,nntot,num_kpts),stat=ierr)
! if (ierr/=0)&
! call io_error('Error allocating m_matrix in pw90common_wanint_data_dist')
! endif
! call comms_bcast(m_matrix(1,1,1,1),num_wann*num_wann*nntot*num_kpts)
call comms_bcast(have_disentangled, 1)
if (have_disentangled) then
if (.not. on_root) then
! Do we really need these 'if not allocated'? Didn't use them for
! eigval and kpt_latt above!
! if (.not.allocated(u_matrix_opt)) then
! allocate(u_matrix_opt(num_bands,num_wann,num_kpts),stat=ierr)
! if (ierr/=0)&
! call io_error('Error allocating u_matrix_opt in pw90common_wanint_data_dist')
! endif
if (.not. allocated(lwindow)) then
allocate (lwindow(num_bands, num_kpts), stat=ierr)
if (ierr /= 0) &
call io_error('Error allocating lwindow in pw90common_wanint_data_dist')
endif
if (.not. allocated(ndimwin)) then
allocate (ndimwin(num_kpts), stat=ierr)
if (ierr /= 0) &
call io_error('Error allocating ndimwin in pw90common_wanint_data_dist')
endif
end if
! call comms_bcast(u_matrix_opt(1,1,1),num_bands*num_wann*num_kpts)
call comms_bcast(lwindow(1, 1), num_bands*num_kpts)
call comms_bcast(ndimwin(1), num_kpts)
end if
end subroutine pw90common_wanint_data_dist
!=======================================================================
subroutine pw90common_get_occ(eig, occ, ef)
!! Compute the electronic occupancy
use w90_constants, only: dp !,eps7
use w90_parameters, only: num_wann !,smear_temp
! use w90_constants, only : elem_charge_SI,k_B_SI
! Arguments
!
real(kind=dp), intent(in) :: eig(num_wann)
!! Eigenvalues
real(kind=dp), intent(in) :: ef
!! Fermi level
real(kind=dp), intent(out) :: occ(num_wann)
!! Occupancy of states
! Misc/Dummy
!
integer :: i
! real(kind=dp) :: kt
! State occupancies
!
! if(smear_temp < eps7) then
!
! Use a step function occupancy (T=0)
!
occ(:) = 0.0_dp
do i = 1, num_wann
if (eig(i) < ef) occ(i) = 1.0_dp
end do
! else
!
! Use a Fermi-Dirac occupancy (T=smear_temp, in Kelvin)
!
! k_B.T in electron-volts
!
! kt=k_B_SI*smear_temp/elem_charge_SI
! do i=1,num_wann
! occ(i)=1.0_dp/(exp((eig(i)-ef)/kt)+1.0_dp)
! end do
! end if
end subroutine pw90common_get_occ
!=======================================================================
function kmesh_spacing_singleinteger(num_points)
!! Set up the value of the interpolation mesh spacing, needed for
!! adaptive smearing [see Eqs. (34-35) YWVS07]. Choose it as the largest of
!! the three Delta_k's for each of the primitive translations b1, b2, and b3
use w90_parameters, only: recip_lattice
integer, intent(in) :: num_points
real(kind=dp) :: kmesh_spacing_singleinteger
integer :: i
real(kind=dp) :: Delta_k_i(3)
! NOTE: The vectors b_i are stored as *rows* in recip_lattice (argh!).
! Hence I believe Jonathan's original code confused rows with columns
! when computing Delta_k, which he called 'rspace'
! (See my e-mail of 20Sept07)
!
do i = 1, 3
Delta_k_i(i) = sqrt(dot_product(recip_lattice(i, :), recip_lattice(i, :))) &
/num_points
end do
kmesh_spacing_singleinteger = maxval(Delta_k_i)
end function kmesh_spacing_singleinteger
function kmesh_spacing_mesh(mesh)
!! Same as kmesh_spacing_singleinteger, but for a kmesh with three
!! different mesh samplings along the three directions
use w90_parameters, only: recip_lattice
integer, dimension(3), intent(in) :: mesh
real(kind=dp) :: kmesh_spacing_mesh
integer :: i
real(kind=dp) :: Delta_k_i(3)
do i = 1, 3
Delta_k_i(i) = sqrt(dot_product(recip_lattice(i, :), recip_lattice(i, :))) &
/mesh(i)
end do
kmesh_spacing_mesh = maxval(Delta_k_i)
end function kmesh_spacing_mesh
!
!=========================================================!
subroutine pw90common_fourier_R_to_k(kpt, OO_R, OO, alpha)
!=========================================================!
! !
!! For alpha=0:
!! O_ij(R) --> O_ij(k) = sum_R e^{+ik.R}*O_ij(R)
!!
!! For alpha=1,2,3:
!! sum_R [cmplx_i*R_alpha*e^{+ik.R}*O_ij(R)]
!! where R_alpha is a Cartesian component of R
!! ***REMOVE EVENTUALLY*** (replace with pw90common_fourier_R_to_k_new)
! !
!=========================================================!
use w90_constants, only: dp, cmplx_0, cmplx_i, twopi
use w90_parameters, only: num_kpts, kpt_latt, num_wann, use_ws_distance
use w90_ws_distance, only: irdist_ws, crdist_ws, &
wdist_ndeg, ws_translate_dist
implicit none
! Arguments
!
real(kind=dp) :: kpt(3)
complex(kind=dp), dimension(:, :, :), intent(in) :: OO_R
complex(kind=dp), dimension(:, :), intent(out) :: OO
integer :: alpha
integer :: ir, i, j, ideg
real(kind=dp) :: rdotk
complex(kind=dp) :: phase_fac
if (use_ws_distance) CALL ws_translate_dist(nrpts, irvec)
OO(:, :) = cmplx_0
do ir = 1, nrpts
! [lp] Shift the WF to have the minimum distance IJ, see also ws_distance.F90
if (use_ws_distance) then
do j = 1, num_wann
do i = 1, num_wann
do ideg = 1, wdist_ndeg(j, i, ir)
rdotk = twopi*dot_product(kpt(:), real(irdist_ws(:, ideg, i, j, ir), dp))
!phase_fac=cmplx(cos(rdotk),sin(rdotk),dp)/real(ndegen(ir)*wdist_ndeg(i,j,ir),dp)
phase_fac = cmplx(cos(rdotk), sin(rdotk), dp)/real(ndegen(ir)*wdist_ndeg(i, j, ir), dp)
if (alpha == 0) then
OO(i, j) = OO(i, j) + phase_fac*OO_R(i, j, ir)
elseif (alpha == 1 .or. alpha == 2 .or. alpha == 3) then
OO(i, j) = OO(i, j) + cmplx_i*crdist_ws(alpha, ideg, i, j, ir)*phase_fac*OO_R(i, j, ir)
!OO(i,j)=OO(i,j)+cmplx_i*crvec(alpha,ir)*phase_fac*OO_R(i,j,ir)
else
stop 'wrong value of alpha in pw90common_fourier_R_to_k'
endif
enddo
enddo
enddo
else
! [lp] Original code, without IJ-dependent shift:
rdotk = twopi*dot_product(kpt(:), irvec(:, ir))
phase_fac = cmplx(cos(rdotk), sin(rdotk), dp)/real(ndegen(ir), dp)
if (alpha == 0) then
OO(:, :) = OO(:, :) + phase_fac*OO_R(:, :, ir)
elseif (alpha == 1 .or. alpha == 2 .or. alpha == 3) then
OO(:, :) = OO(:, :) + &
cmplx_i*crvec(alpha, ir)*phase_fac*OO_R(:, :, ir)
else
stop 'wrong value of alpha in pw90common_fourier_R_to_k'
endif
endif
enddo
end subroutine pw90common_fourier_R_to_k
! ***NEW***
!
!=========================================================!
subroutine pw90common_fourier_R_to_k_new(kpt, OO_R, OO, OO_dx, OO_dy, OO_dz)
!=======================================================!
! !
!! For OO:
!! $$O_{ij}(k) = \sum_R e^{+ik.R}.O_{ij}(R)$$
!! For $$OO_{dx,dy,dz}$$:
!! $$\sum_R [i.R_{dx,dy,dz}.e^{+ik.R}.O_{ij}(R)]$$
!! where R_{x,y,z} are the Cartesian components of R
! !
!=======================================================!
use w90_constants, only: dp, cmplx_0, cmplx_i, twopi
use w90_parameters, only: timing_level, num_kpts, kpt_latt, num_wann, use_ws_distance
use w90_ws_distance, only: irdist_ws, wdist_ndeg, ws_translate_dist
implicit none
! Arguments
!
real(kind=dp) :: kpt(3)
complex(kind=dp), dimension(:, :, :), intent(in) :: OO_R
complex(kind=dp), optional, dimension(:, :), intent(out) :: OO
complex(kind=dp), optional, dimension(:, :), intent(out) :: OO_dx
complex(kind=dp), optional, dimension(:, :), intent(out) :: OO_dy
complex(kind=dp), optional, dimension(:, :), intent(out) :: OO_dz
integer :: ir, i, j, ideg
real(kind=dp) :: rdotk
complex(kind=dp) :: phase_fac
if (use_ws_distance) CALL ws_translate_dist(nrpts, irvec)
if (present(OO)) OO = cmplx_0
if (present(OO_dx)) OO_dx = cmplx_0
if (present(OO_dy)) OO_dy = cmplx_0
if (present(OO_dz)) OO_dz = cmplx_0
do ir = 1, nrpts
! [lp] Shift the WF to have the minimum distance IJ, see also ws_distance.F90
if (use_ws_distance) then
do j = 1, num_wann
do i = 1, num_wann
do ideg = 1, wdist_ndeg(j, i, ir)
rdotk = twopi*dot_product(kpt(:), real(irdist_ws(:, ideg, i, j, ir), dp))
phase_fac = cmplx(cos(rdotk), sin(rdotk), dp)/real(ndegen(ir)*wdist_ndeg(i, j, ir), dp)
if (present(OO)) OO(i, j) = OO(i, j) + phase_fac*OO_R(i, j, ir)
if (present(OO_dx)) OO_dx(i, j) = OO_dx(i, j) + &
cmplx_i*crvec(1, ir)*phase_fac*OO_R(i, j, ir)
if (present(OO_dy)) OO_dy(i, j) = OO_dy(i, j) + &
cmplx_i*crvec(2, ir)*phase_fac*OO_R(i, j, ir)
if (present(OO_dz)) OO_dz(i, j) = OO_dz(i, j) + &
cmplx_i*crvec(3, ir)*phase_fac*OO_R(i, j, ir)
enddo
enddo
enddo
else
! [lp] Original code, without IJ-dependent shift:
rdotk = twopi*dot_product(kpt(:), irvec(:, ir))
phase_fac = cmplx(cos(rdotk), sin(rdotk), dp)/real(ndegen(ir), dp)
if (present(OO)) OO(:, :) = OO(:, :) + phase_fac*OO_R(:, :, ir)
if (present(OO_dx)) OO_dx(:, :) = OO_dx(:, :) + &
cmplx_i*crvec(1, ir)*phase_fac*OO_R(:, :, ir)
if (present(OO_dy)) OO_dy(:, :) = OO_dy(:, :) + &
cmplx_i*crvec(2, ir)*phase_fac*OO_R(:, :, ir)
if (present(OO_dz)) OO_dz(:, :) = OO_dz(:, :) + &
cmplx_i*crvec(3, ir)*phase_fac*OO_R(:, :, ir)
endif
enddo
end subroutine pw90common_fourier_R_to_k_new
!=========================================================!
subroutine pw90common_fourier_R_to_k_new_second_d(kpt, OO_R, OO, OO_da, OO_dadb)
!=======================================================!
! !
!! For OO:
!! $$O_{ij}(k) = \sum_R e^{+ik.R}.O_{ij}(R)$$
!! For $$OO_{dx,dy,dz}$$:
!! $$\sum_R [i.R_{dx,dy,dz}.e^{+ik.R}.O_{ij}(R)]$$
!! where R_{x,y,z} are the Cartesian components of R
!! For $$OO_{dx1,dy1,dz1;dx2,dy2,dz2}$$:
!! $$-\sum_R [R_{dx1,dy1,dz1}.R_{dx2,dy2,dz2}.e^{+ik.R}.O_{ij}(R)]$$
!! where R_{xi,yi,zi} are the Cartesian components of R
! !
!=======================================================!
use w90_constants, only: dp, cmplx_0, cmplx_i, twopi
use w90_parameters, only: timing_level, num_kpts, kpt_latt, num_wann, use_ws_distance
use w90_ws_distance, only: irdist_ws, wdist_ndeg, ws_translate_dist
implicit none
! Arguments
!
real(kind=dp) :: kpt(3)
complex(kind=dp), dimension(:, :, :), intent(in) :: OO_R
complex(kind=dp), optional, dimension(:, :), intent(out) :: OO
complex(kind=dp), optional, dimension(:, :, :), intent(out) :: OO_da
complex(kind=dp), optional, dimension(:, :, :, :), intent(out) :: OO_dadb
integer :: ir, i, j, ideg, a, b
real(kind=dp) :: rdotk
complex(kind=dp) :: phase_fac
if (use_ws_distance) CALL ws_translate_dist(nrpts, irvec)
if (present(OO)) OO = cmplx_0
if (present(OO_da)) OO_da = cmplx_0
if (present(OO_dadb)) OO_dadb = cmplx_0
do ir = 1, nrpts
! [lp] Shift the WF to have the minimum distance IJ, see also ws_distance.F90
if (use_ws_distance) then
do j = 1, num_wann
do i = 1, num_wann
do ideg = 1, wdist_ndeg(j, i, ir)
rdotk = twopi*dot_product(kpt(:), real(irdist_ws(:, ideg, i, j, ir), dp))
phase_fac = cmplx(cos(rdotk), sin(rdotk), dp)/real(ndegen(ir)*wdist_ndeg(i, j, ir), dp)
if (present(OO)) OO(i, j) = OO(i, j) + phase_fac*OO_R(i, j, ir)
if (present(OO_da)) then
do a = 1, 3
OO_da(i, j, a) = OO_da(i, j, a) + cmplx_i*crvec(a, ir)*phase_fac*OO_R(i, j, ir)
enddo
endif
if (present(OO_dadb)) then
do a = 1, 3
do b = 1, 3
OO_dadb(i, j, a, b) = OO_dadb(i, j, a, b) - &
crvec(a, ir)*crvec(b, ir)*phase_fac*OO_R(i, j, ir)
enddo
enddo
end if
enddo
enddo
enddo
else
! [lp] Original code, without IJ-dependent shift:
rdotk = twopi*dot_product(kpt(:), irvec(:, ir))
phase_fac = cmplx(cos(rdotk), sin(rdotk), dp)/real(ndegen(ir), dp)
if (present(OO)) OO(:, :) = OO(:, :) + phase_fac*OO_R(:, :, ir)
if (present(OO_da)) then
do a = 1, 3
OO_da(:, :, a) = OO_da(:, :, a) + cmplx_i*crvec(a, ir)*phase_fac*OO_R(:, :, ir)
enddo
endif
if (present(OO_dadb)) then
do a = 1, 3
do b = 1, 3
OO_dadb(:, :, a, b) = OO_dadb(:, :, a, b) - &
crvec(a, ir)*crvec(b, ir)*phase_fac*OO_R(:, :, ir)
enddo
enddo
end if
endif
enddo
end subroutine pw90common_fourier_R_to_k_new_second_d
subroutine pw90common_fourier_R_to_k_new_second_d_TB_conv(kpt, OO_R, oo_a_R, OO, OO_da, OO_dadb)
!=======================================================!
! modified version of pw90common_fourier_R_to_k_new_second_d, includes wannier centres in
! the exponential inside the sum (so called TB convention)
!
!! For OO:
!! $$O_{ij}(k) = \sum_R e^{+ik.(R+tau_ij)}.O_{ij}(R)$$
!! For $$OO_{dx,dy,dz}$$:
!! $$\sum_R [i.(R+tau_ij)_{dx,dy,dz}.e^{+ik.(R+tau_ij)}.O_{ij}(R)]$$
!! where R_{x,y,z} are the Cartesian components of R
!! For $$OO_{dx1,dy1,dz1;dx2,dy2,dz2}$$:
!! $$-\sum_R [(R+tau_ij)_{dx1,dy1,dz1}.(R+tau_ij)_{dx2,dy2,dz2}.e^{+ik.(R+tau_ij)}.O_{ij}(R)]$$
!! where {xi,yi,zi} denote the Cartesian components and
! tau_ij = tau_j - tau_i, being tau_i=<0i|r|0i> the individual wannier centres
!=======================================================!
use w90_constants, only: dp, cmplx_0, cmplx_i, twopi
use w90_parameters, only: timing_level, num_kpts, kpt_latt, num_wann, &
use_ws_distance, wannier_centres, recip_lattice
use w90_ws_distance, only: irdist_ws, wdist_ndeg, ws_translate_dist
use w90_utility, only: utility_cart_to_frac
implicit none
! Arguments
!
real(kind=dp) :: kpt(3)
complex(kind=dp), dimension(:, :, :), intent(in) :: OO_R
complex(kind=dp), optional, dimension(:, :), intent(out) :: OO
complex(kind=dp), optional, dimension(:, :, :), intent(out) :: OO_da
complex(kind=dp), optional, dimension(:, :, :, :), intent(out) :: OO_dadb
complex(kind=dp), dimension(:, :, :, :), intent(in) :: oo_a_R
integer :: ir, i, j, ideg, a, b
real(kind=dp) :: rdotk
complex(kind=dp) :: phase_fac
real(kind=dp) :: local_wannier_centres(3, num_wann), wannier_centres_frac(3, num_wann)
real(kind=dp) :: r_sum(3)
r_sum = 0.d0
if (use_ws_distance) CALL ws_translate_dist(nrpts, irvec)
! calculate wannier centres in cartesian
local_wannier_centres(:, :) = 0.d0
do j = 1, num_wann
do ir = 1, nrpts
if ((irvec(1, ir) .eq. 0) .and. (irvec(2, ir) .eq. 0) .and. (irvec(3, ir) .eq. 0)) then
local_wannier_centres(1, j) = real(oo_a_R(j, j, ir, 1))
local_wannier_centres(2, j) = real(oo_a_R(j, j, ir, 2))
local_wannier_centres(3, j) = real(oo_a_R(j, j, ir, 3))
endif
enddo
enddo
! rotate wannier centres from cartesian to fractional coordinates
wannier_centres_frac(:, :) = 0.d0
do ir = 1, num_wann
call utility_cart_to_frac(local_wannier_centres(:, ir), wannier_centres_frac(:, ir), recip_lattice)
enddo
if (present(OO)) OO = cmplx_0
if (present(OO_da)) OO_da = cmplx_0
if (present(OO_dadb)) OO_dadb = cmplx_0
do ir = 1, nrpts
! [lp] Shift the WF to have the minimum distance IJ, see also ws_distance.F90
if (use_ws_distance) then
do j = 1, num_wann
do i = 1, num_wann
do ideg = 1, wdist_ndeg(j, i, ir)
rdotk = twopi*dot_product(kpt(:), real(irdist_ws(:, ideg, i, j, ir) + &
wannier_centres_frac(:, j) - wannier_centres_frac(:, i), dp))
phase_fac = cmplx(cos(rdotk), sin(rdotk), dp)/real(ndegen(ir)*wdist_ndeg(i, j, ir), dp)
if (present(OO)) OO(i, j) = OO(i, j) + phase_fac*OO_R(i, j, ir)
if (present(OO_da)) then
do a = 1, 3
OO_da(i, j, a) = OO_da(i, j, a) + cmplx_i*(crvec(a, ir) + local_wannier_centres(a, j) - &
local_wannier_centres(a, i))*phase_fac*OO_R(i, j, ir)
enddo
endif
if (present(OO_dadb)) then
do a = 1, 3