-
Notifications
You must be signed in to change notification settings - Fork 15
/
get_eig.f90
1276 lines (1157 loc) · 53.7 KB
/
get_eig.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
#include "alias.inc"
subroutine get_eig(NN_TABLE, kp, nkp, PINPT, PPRAM, E, V, SV, neig, iband, nband, flag_vector, flag_sparse, flag_stat, flag_phase)
use parameters, only: hopping, incar, energy, spmat, params
use mpi_setup
use time
use memory
use reorder_band
use print_io
implicit none
type (hopping) :: NN_TABLE
type (incar ) :: PINPT
type (params ) :: PPRAM
type (energy ) :: EE
type (spmat ) :: SHm, SHs
integer*4 mpierr, iadd, ii
integer*4 ik, my_ik, neig
integer*4 iband,nband ! if sparse: iband = 1, nband = feast_nemax
integer*4 nkp, is, ie,fe, im, fm
real*8 percent, kp(3,nkp)
real*8 E(nband*PINPT%nspin,nkp)
complex*16 V(neig*PINPT%ispin,nband*PINPT%nspin,nkp)
complex*16 SV(neig*PINPT%ispin,nband*PINPT%nspin,nkp) ! overlap integeral multiplied S * V
complex*16 Hm(neig*PINPT%ispinor,neig*PINPT%ispinor) ! collinear magnetism hamiltonian (k-independent)
complex*16 Hs(neig*PINPT%ispinor,neig*PINPT%ispinor) ! 1st-order SO coupling hamiltonian (k-dependent if .not. SK)
logical, intent(in) :: flag_vector, flag_stat
real*8 t0, t1
character*4 timer
logical flag_phase, flag_init, flag_sparse
logical flag_sparse_SHm, flag_sparse_SHs ! if .false. sparse Hamiltonian SHm(collinear magnetic) and SHs(SOC) will not
! be added up and constructed along with the k_loop due to the total number
! of non-zero element is zero. This is determined in the get_ham_mag(soc)_sparse
integer*4 feast_ne(PINPT%nspin, nkp)
integer*4 ncpu, id
integer*4, allocatable :: ourjob(:), ourjob_disp(:), ourgroup(:)
integer*4 mygroup(0:nprocs-1),groupid
if(flag_stat) then
write(message,'( A)') ' ' ; write_msg
write(message,'( A)') ' #--- START: BAND STRUCTURE EVALUATION -----------' ; write_msg
write(message,'(2A)') ' # SYSTEM: ', trim(PINPT%ifilenm(NN_TABLE%mysystem)) ; write_msg
endif
timer = 'init'
if(COMM_KOREA%flag_split) then
ncpu = COMM_KOREA%nprocs
id = COMM_KOREA%myid
else
ncpu = nprocs
id = myid
endif
#ifdef PSPARSE
if(flag_sparse) then
#ifdef MPI
call MPI_BARRIER(mpi_comm_earth, mpierr)
call get_npar(trim(PINPT%ifilenm(1)))
allocate(ourgroup(npar))
allocate(ourjob(npar))
allocate(ourjob_disp(0:npar-1))
call mpi_job_distribution_group(npar, nkp, ourgroup, mygroup, ourjob, ourjob_disp)
call mpi_comm_anmeldung(COMM_JUELICH, npar, mygroup, COMM_JUELICH_RATHAUS)
groupid = COMM_JUELICH%color
id = groupid
call MPI_BARRIER(mpi_comm_earth, mpierr)
#else
write(message,'(A)')' DPSARSE and DMPI option should be activated together. Currently, only DPSPARSE is activated.' ; write_msg
write(message,'(A)')' Please set DPSPARSE option in your makefile OPTIONS and recompile the code. Exit...' ; write_msg
kill_job
#endif
else
allocate(ourjob(ncpu))
allocate(ourjob_disp(0:ncpu-1))
endif
#else
allocate(ourjob(ncpu))
allocate(ourjob_disp(0:ncpu-1))
#endif
#ifdef PSPARSE
if(flag_sparse) then
call report_job_distribution_group(flag_stat, ourjob, ourgroup)
else
call mpi_job_distribution_chain(nkp, ncpu, ourjob, ourjob_disp)
if(.not. COMM_KOREA%flag_split) call report_job_distribution(flag_stat, ourjob)
endif
#else
call mpi_job_distribution_chain(nkp, ncpu, ourjob, ourjob_disp)
if(.not. COMM_KOREA%flag_split) call report_job_distribution(flag_stat, ourjob)
#endif
call initialize_all (EE, neig, nband, nkp, ourjob(id+1), &
PINPT, flag_vector, PPRAM%flag_use_overlap, flag_sparse, flag_stat, &
ii, iadd, t1, t0, flag_init)
if_main call report_memory_total(PINPT%ispinor, PINPT%ispin, PINPT%nspin, neig, nband, nkp, &
flag_stat, flag_sparse, ncpu)
if(flag_stat) then
write(message,'(A)' ) ' ' ; write_msg
write(message,'(A)') ' STATUS: '
call write_log(trim(message), 1, myid)
call write_log(trim(message),22, myid)
endif
k_loop:do ik= sum(ourjob(1:id))+1, sum(ourjob(1:id+1))
my_ik = ik - sum(ourjob(1:id))
if(flag_sparse) then
#ifdef MKL_SPARSE
!NOTE: SHm is k-independent -> keep unchankged from first call
! SHs is k-independent if flag_slater_koster
! k-dependent if .not. slater_koster
call cal_eig_Hk_sparse(SHm, SHs, EE%E(:,ik), EE%V(:,:,my_ik), PINPT, PPRAM, NN_TABLE, kp(:,ik), &
neig, nband, flag_vector, flag_init, flag_phase, &
PINPT%feast_ne(1:PINPT%nspin,ik),&
flag_sparse_SHm, flag_sparse_SHs, timer)
#else
write(message,'(A)')' !WARN! The EWINDOW tag is only available if you have put -DMKL_SPARSE option' ; write_msg
write(message,'(A)')' in your make file. Please find the details in the instruction. Exit program...' ; write_msg
kill_job
#endif
elseif(.not.flag_sparse) then
call cal_eig_Hk_dense ( Hm, Hs, EE%E(:,ik), EE%V(:,:,my_ik), EE%SV(:,:,my_ik), PINPT, PPRAM, NN_TABLE, kp(:,ik), &
neig, iband, nband, flag_vector, flag_init, flag_phase,ik)
endif
#ifdef PSPARSE
if(flag_sparse) then
if(flag_stat .and. myid .eq. 0) call print_eig_status(ik, ii, iadd, ourjob)
elseif(.not. flag_sparse) then
if(flag_stat .and. id .eq. 0) call print_eig_status(ik, ii, iadd, ourjob)
endif
#else
if(flag_stat .and. id .eq. 0) call print_eig_status(ik, ii, iadd, ourjob)
#endif
enddo k_loop
#ifdef MPI
if(flag_stat .and. myid .eq. 0) then
write(message,'(A)')' ' ; write_msg
write(message,'(A)')' Gathering all results to main node 0 ...'
call write_log(trim(message), 23, myid)
endif
if(.not. COMM_KOREA%flag_split) then
if(.not. flag_sparse) then
call MPI_ALLREDUCE(EE%E, E, size(E,kind=4), MPI_REAL8, MPI_SUM, mpi_comm_earth, mpierr) ! share all results
if(flag_vector .and. nkp .ne. 1) then
call MPI_GATHERV(EE%V,size(EE%V, kind=4), MPI_COMPLEX16, V, &
ourjob *neig*PINPT%ispin*nband*PINPT%nspin, &
ourjob_disp *neig*PINPT%ispin*nband*PINPT%nspin, &
MPI_COMPLEX16, 0, mpi_comm_earth, mpierr) ! only main node keep wave vector information
if(PPRAM%flag_use_overlap) then
call MPI_GATHERV(EE%SV,size(EE%SV,kind=4), MPI_COMPLEX16, SV, &
ourjob * neig*PINPT%ispin*nband*PINPT%nspin , &
ourjob_disp * neig*PINPT%ispin*nband*PINPT%nspin , &
MPI_COMPLEX16, 0, mpi_comm_earth, mpierr) ! only main node keep overlap matrix information
endif
endif
elseif(flag_sparse) then
#ifdef PSPARSE
if(COMM_JUELICH_RATHAUS%mpi_comm .ne. MPI_COMM_NULL) then
call MPI_ALLREDUCE(EE%E, E, size(E,kind=4), MPI_REAL8, MPI_SUM, COMM_JUELICH_RATHAUS%mpi_comm, mpierr) ! share all results
endif
if(flag_vector .and. nkp .ne. 1 .and. COMM_JUELICH_RATHAUS%mpi_comm .ne. MPI_COMM_NULL) then
call MPI_GATHERV(EE%V,size(EE%V, kind=4), MPI_COMPLEX16, V, &
ourjob *neig*PINPT%ispin*nband*PINPT%nspin, &
ourjob_disp *neig*PINPT%ispin*nband*PINPT%nspin, &
MPI_COMPLEX16, 0, COMM_JUELICH_RATHAUS%mpi_comm, mpierr) ! only main node keep wave vector information
if(PPRAM%flag_use_overlap) then
call MPI_GATHERV(EE%SV,size(EE%SV,kind=4), MPI_COMPLEX16, SV, &
ourjob * neig*PINPT%ispin*nband*PINPT%nspin , &
ourjob_disp * neig*PINPT%ispin*nband*PINPT%nspin , &
MPI_COMPLEX16, 0, COMM_JUELICH_RATHAUS%mpi_comm, mpierr) ! only main node keep overlap matrix information
endif
endif
#else
call MPI_ALLREDUCE(EE%E, E, size(E,kind=4), MPI_REAL8, MPI_SUM, mpi_comm_earth, mpierr) ! share all results
if(flag_vector .and. nkp .ne. 1) then
call MPI_GATHERV(EE%V,size(EE%V, kind=4), MPI_COMPLEX16, V, &
ourjob *neig*PINPT%ispin*nband*PINPT%nspin, &
ourjob_disp *neig*PINPT%ispin*nband*PINPT%nspin, &
MPI_COMPLEX16, 0, mpi_comm_earth, mpierr) ! only main node keep wave vector information
if(PPRAM%flag_use_overlap) then
call MPI_GATHERV(EE%SV,size(EE%SV,kind=4), MPI_COMPLEX16, SV, &
ourjob * neig*PINPT%ispin*nband*PINPT%nspin , &
ourjob_disp * neig*PINPT%ispin*nband*PINPT%nspin , &
MPI_COMPLEX16, 0, mpi_comm_earth, mpierr) ! only main node keep overlap matrix information
endif
endif
#endif
endif
elseif(COMM_KOREA%flag_split) then
call MPI_ALLREDUCE(EE%E, E, size(E,kind=4), MPI_REAL8, MPI_SUM, COMM_KOREA%mpi_comm, mpierr) ! share all results
if(flag_vector .and. nkp .ne. 1) then
call MPI_GATHERV(EE%V,size(EE%V,kind=4), MPI_COMPLEX16, V, &
ourjob *neig*PINPT%ispin*nband*PINPT%nspin, &
ourjob_disp *neig*PINPT%ispin*nband*PINPT%nspin, &
MPI_COMPLEX16, 0, COMM_KOREA%mpi_comm, mpierr) ! only main node keep wave vector information
if(PPRAM%flag_use_overlap) then
call MPI_GATHERV(EE%SV,size(EE%SV,kind=4), MPI_COMPLEX16, SV, &
ourjob *neig*PINPT%ispin*nband*PINPT%nspin, &
ourjob_disp *neig*PINPT%ispin*nband*PINPT%nspin, &
MPI_COMPLEX16, 0, COMM_KOREA%mpi_comm, mpierr) ! only main node keep overlap matrix information
endif
endif
endif
if(flag_stat) then
write(message,'(A)')' done!' ; write_msg
write(message,'(A )')' ' ; write_msg
endif
#ifdef MKL_SPARSE
if(flag_sparse) then
if(.not.COMM_KOREA%flag_split) then
#ifdef PSPARSE
call MPI_ALLREDUCE(PINPT%feast_ne, feast_ne, size(feast_ne), MPI_INTEGER4, MPI_SUM, COMM_JUELICH%mpi_comm, mpierr)
if(COMM_JUELICH_RATHAUS%mpi_comm .ne. MPI_COMM_NULL) then
call MPI_ALLREDUCE(PINPT%feast_ne, feast_ne, size(feast_ne), MPI_INTEGER4, MPI_SUM, COMM_JUELICH_RATHAUS%mpi_comm, mpierr)
endif
#else
call MPI_ALLREDUCE(PINPT%feast_ne, feast_ne, size(feast_ne), MPI_INTEGER4, MPI_SUM, mpi_comm_earth, mpierr)
#endif
elseif(COMM_KOREA%flag_split) then
call MPI_ALLREDUCE(PINPT%feast_ne, feast_ne, size(feast_ne), MPI_INTEGER4, MPI_SUM, COMM_KOREA%mpi_comm, mpierr)
endif
PINPT%feast_ne = feast_ne
if(flag_stat) then
write(message,'(A,I0)')' MAX_NE_FOUND (NE_MAX): ',maxval(PINPT%feast_ne) ; write_msg
endif
endif
#endif
#else
E = EE%E ; if(flag_vector) V = EE%V
#ifdef MKL_SPARSE
if(flag_stat) then
write(message,'(A,I0)')' MAX_NE_FOUND (NE_MAX): ',maxval(PINPT%feast_ne) ; write_msg
endif
#endif
#endif
call finalize_all(EE, SHm, SHs, t1, t0, PINPT, flag_stat, flag_vector, flag_sparse)
if(flag_stat) then
write(message,'(A)')' #--- END: BAND STRUCTURE EVALUATION -----------' ; write_msg
endif
return
endsubroutine
subroutine get_eig_sepk(NN_TABLE, kp, nkp, PINPT, PPRAM, E, neig, iband, nband, flag_vector, flag_sparse, flag_stat, flag_phase)
use parameters, only: hopping, incar, energy, spmat, params
use mpi_setup
use time
use memory
use reorder_band
use print_io
implicit none
type (hopping) :: NN_TABLE
type (incar ) :: PINPT
type (params ) :: PPRAM
type (energy ) :: EE
type (spmat ) :: SHm, SHs
integer*4 mpierr, iadd, ii
integer*4 ik, neig, ijob
integer*4 iband,nband ! if sparse: iband = 1, nband = feast_nemax
integer*4 nkp, is, ie,fe, im, fm
real*8 percent, kp(3,nkp)
real*8 E(nband*PINPT%nspin, nkp)
complex*16 V(neig*PINPT%ispin,nband*PINPT%nspin)
complex*16 SV(neig*PINPT%ispin,nband*PINPT%nspin) ! overlap integeral multiplied S * V
complex*16 Hm(neig*PINPT%ispinor,neig*PINPT%ispinor) ! collinear magnetism hamiltonian (k-independent)
complex*16 Hs(neig*PINPT%ispinor,neig*PINPT%ispinor) ! 1st-order SO coupling hamiltonian (k-dependent if .not. SK)
logical, intent(in) :: flag_vector, flag_stat
real*8 t0, t1
character*4 timer
character*80 fname_header, fname
character*20,external :: int2str
logical flag_file_exist
logical flag_phase, flag_init, flag_sparse
logical flag_sparse_SHm, flag_sparse_SHs ! if .false. sparse Hamiltonian SHm(collinear magnetic) and SHs(SOC) will not
! be added up and constructed along with the k_loop due to the total number
! of non-zero element is zero. This is determined in the get_ham_mag(soc)_sparse
integer*4 feast_ne(PINPT%nspin, nkp)
integer*4 ncpu, id
integer*4, allocatable :: ourjob(:), ourjob_disp(:), ourgroup(:)
integer*4, allocatable :: non_skip_list(:)
integer*4 mygroup(0:nprocs-1), groupid
integer*4 n_non_skip
if(flag_stat) then
write(message,'(A)') ' ' ; write_msg
write(message,'(A)') ' #--- START: BAND STRUCTURE EVALUATION -----------' ; write_msg
endif
timer = 'init'
EE%E = 0d0
if(flag_vector) V = (0.d0,0.d0)
if(flag_vector .and. PPRAM%flag_use_overlap) SV = (0.d0,0.d0)
if(COMM_KOREA%flag_split) then
ncpu = COMM_KOREA%nprocs
id = COMM_KOREA%myid
else
ncpu = nprocs
id = myid
endif
! check non-skip list
allocate(non_skip_list(nkp)) ; non_skip_list = 0 ; n_non_skip = 0
do ik = 1, nkp
! check if the calculation is already done.. check restart
fname_header = 'band_structure_TBA'//trim(PINPT%title(NN_TABLE%mysystem))//'.kp_'//trim(ADJUSTL(int2str(ik)))
if(.not. PINPT%flag_write_unformatted) then
call get_fname(fname_header, fname, 1, PINPT%flag_collinear, PINPT%flag_noncollinear)
elseif(PINPT%flag_write_unformatted) then
call get_fname_bin(fname_header, fname, 1, PINPT%flag_collinear, PINPT%flag_noncollinear)
endif
inquire(file=trim(fname),exist=flag_file_exist)
if(flag_file_exist) then
write(message,'(A,A,A,I0)')' Band structure file:', trim(fname),' is already exist => skip to calculate KPT: ',ik ; write_msg
else
n_non_skip = n_non_skip + 1
non_skip_list(n_non_skip) = ik
cycle
endif
if(PINPT%nspin .eq. 2 .and. flag_file_exist) then
if(.not. PINPT%flag_write_unformatted) then
call get_fname(fname_header, fname, 2, PINPT%flag_collinear, PINPT%flag_noncollinear)
elseif(PINPT%flag_write_unformatted) then
call get_fname_bin(fname_header, fname, 2, PINPT%flag_collinear, PINPT%flag_noncollinear)
endif
inquire(file=trim(fname),exist=flag_file_exist)
if(flag_file_exist) then
write(message,'(A,A,A,I0)')' Band structure file:', trim(fname),' is already exist => skip to calculate KPT: ',ik ; write_msg
else
n_non_skip = n_non_skip + 1
non_skip_list(n_non_skip) = ik
cycle
endif
endif
enddo
#ifdef PSPARSE
if(flag_sparse) then
#ifdef MPI
call MPI_BARRIER(mpi_comm_earth, mpierr)
call get_npar(trim(PINPT%ifilenm(1)))
allocate(ourgroup(npar))
allocate(ourjob(npar))
allocate(ourjob_disp(0:npar-1))
call mpi_job_distribution_group(npar, n_non_skip, ourgroup, mygroup, ourjob, ourjob_disp)
call mpi_comm_anmeldung(COMM_JUELICH, npar, mygroup, COMM_JUELICH_RATHAUS)
groupid = COMM_JUELICH%color
id = groupid
call MPI_BARRIER(mpi_comm_earth, mpierr)
#else
write(message,'(A)')' DPSARSE and DMPI option should be activated together. Currently, only DPSPARSE is activated.' ; write_msg
write(message,'(A)')' Please set DPSPARSE option in your makefile OPTIONS and recompile the code. Exit...' ; write_msg
kill_job
#endif
else
allocate(ourjob(ncpu))
allocate(ourjob_disp(0:ncpu-1))
endif
#else
allocate(ourjob(ncpu))
allocate(ourjob_disp(0:ncpu-1))
#endif
#ifdef PSPARSE
if(flag_sparse) then
call report_job_distribution_group(flag_stat, ourjob, ourgroup)
else
call mpi_job_distribution_chain(n_non_skip, ncpu, ourjob, ourjob_disp)
if(.not. COMM_KOREA%flag_split) call report_job_distribution(flag_stat, ourjob)
endif
#else
call mpi_job_distribution_chain(n_non_skip, ncpu, ourjob, ourjob_disp)
if(.not. COMM_KOREA%flag_split) call report_job_distribution(flag_stat, ourjob)
#endif
call initialize_all (EE, neig, nband, nkp, ourjob(id+1), PINPT, flag_vector, PPRAM%flag_use_overlap, flag_sparse, flag_stat, &
ii, iadd, t1, t0, flag_init)
if_main call report_memory_singlek(PINPT%ispinor, PINPT%ispin, PINPT%nspin, neig, nband, nkp, &
flag_stat, flag_sparse, ncpu)
if(flag_stat) then
write(message,'(A)' ) ' ' ; write_msg
write(message,'(A)') ' STATUS: '
call write_log(trim(message), 1, myid)
call write_log(trim(message),22, myid)
endif
k_loop:do ijob= sum(ourjob(1:id))+1, sum(ourjob(1:id+1))
ik = non_skip_list(ijob)
if(flag_sparse) then
#ifdef MKL_SPARSE
!NOTE: SHm is k-independent -> keep unchankged from first call
! SHs is k-independent if flag_slater_koster
! k-dependent if .not. slater_koster
call cal_eig_Hk_sparse(SHm, SHs, EE%E(:,ik), V, PINPT, PPRAM, NN_TABLE, kp(:,ik), &
neig, nband, flag_vector, flag_init, flag_phase, &
PINPT%feast_ne(1:PINPT%nspin,ik), &
flag_sparse_SHm, flag_sparse_SHs, timer)
#else
write(message,'(A)')' !WARN! The EWINDOW tag is only available if you have put -DMKL_SPARSE option' ; write_msg
write(message,'(A)')' in your make file. Please find the details in the instruction. Exit program...' ; write_msg
kill_job
#endif
elseif(.not.flag_sparse) then
call cal_eig_Hk_dense ( Hm, Hs, EE%E(:,ik), V, SV, PINPT, PPRAM, NN_TABLE, kp(:,ik), &
neig, iband, nband, flag_vector, flag_init, flag_phase,ik)
endif
if(flag_stat .and. id .eq. 0) call print_eig_status(ijob, ii, iadd, ourjob)
call print_energy_singlek(EE%E(:,ik), V, SV, neig, iband, nband, ik, kp(:,ik), &
PINPT, PPRAM%flag_use_overlap, NN_TABLE%mysystem)
enddo k_loop
#ifdef MPI
! call MPI_BARRIER(mpi_comm_earth, mpierr)
if(flag_stat .and. myid .eq. 0) then
write(message,'(A)')' ' ; write_msg
write(message,'(A)')' Gathering all results to main node 0 ...'
call write_log(trim(message), 23, myid)
endif
if(.not. COMM_KOREA%flag_split) then
if(.not. flag_sparse) then
call MPI_ALLREDUCE(EE%E, E, size(E), MPI_REAL8, MPI_SUM, mpi_comm_earth, mpierr) ! share all results
elseif(flag_sparse) then
#ifdef PSPARSE
if(COMM_JUELICH_RATHAUS%mpi_comm .ne. MPI_COMM_NULL) then
call MPI_ALLREDUCE(EE%E, E, size(E), MPI_REAL8, MPI_SUM, COMM_JUELICH_RATHAUS%mpi_comm, mpierr)
endif
#else
call MPI_ALLREDUCE(EE%E, E, size(E), MPI_REAL8, MPI_SUM, mpi_comm_earth, mpierr) ! share all results
#endif
endif
elseif(COMM_KOREA%flag_split) then
call MPI_ALLREDUCE(EE%E, E, size(E), MPI_REAL8, MPI_SUM, COMM_KOREA%mpi_comm, mpierr) ! share all results
endif
if(flag_stat) then
write(message,'(A)')' done!' ; write_msg
write(message,'(A )')' ' ; write_msg
endif
#ifdef MKL_SPARSE
if(flag_sparse) then
if(.not.COMM_KOREA%flag_split) then
#ifdef PSPARSE
if(COMM_JUELICH_RATHAUS%mpi_comm .ne. MPI_COMM_NULL) then
call MPI_ALLREDUCE(PINPT%feast_ne, feast_ne, size(feast_ne), MPI_INTEGER4, MPI_SUM, COMM_JUELICH_RATHAUS%mpi_comm, mpierr)
endif
#else
call MPI_ALLREDUCE(PINPT%feast_ne, feast_ne, size(feast_ne), MPI_INTEGER4, MPI_SUM, mpi_comm_earth, mpierr)
#endif
elseif(COMM_KOREA%flag_split) then
call MPI_ALLREDUCE(PINPT%feast_ne, feast_ne, size(feast_ne), MPI_INTEGER4, MPI_SUM, COMM_KOREA%mpi_comm, mpierr)
endif
PINPT%feast_ne = feast_ne
if(flag_stat) then
write(message,'(A,I0)')' MAX_NE_FOUND (NE_MAX): ',maxval(PINPT%feast_ne) ; write_msg
endif
endif
#endif
#else
E = EE%E
#ifdef MKL_SPARSE
if(flag_stat) then
write(message,'(A,I0)')' MAX_NE_FOUND (NE_MAX): ',maxval(PINPT%feast_ne) ; write_msg
endif
#endif
#endif
call finalize_all(EE, SHm, SHs, t1, t0, PINPT, flag_stat, flag_vector, flag_sparse)
if(flag_stat) then
write(message,'(A)')' #--- END: BAND STRUCTURE EVALUATION -----------' ; write_msg
endif
! NEED TO BE UPDATED HERE !!! HJ KIM 21.Oct 2020
if(PINPT%flag_print_energy_singlek) then
write(message,'(A)')' !WARN! The calculation finished sucessfully. ' ; write_msg
write(message,'(A)')' Band structure information is printed into separate file "band_structure_TBA.kp_*'//trim(PINPT%title(NN_TABLE%mysystem))//'.dat" by request.' ; write_msg
write(message,'(A)')' However, due to some technical things, program will stop at this point. In the future release it will be updated.' ; write_msg
write(message,'(A)')' ' ; write_msg
write(message,'(A)')' Program stops...' ; write_msg
kill_job
endif
return
endsubroutine
#ifdef MKL_SPARSE
subroutine cal_eig_Hk_sparse(SHm, SHs, E, V, PINPT, PPRAM, NN_TABLE, kp, neig, &
nemax, flag_vector, flag_init, flag_phase, ne_found,&
flag_sparse_SHm, flag_sparse_SHs, timer)
use parameters, only : incar, hopping, spmat, params
use mpi_setup
use sparse_tool
use time
use print_io
implicit none
type(hopping) :: NN_TABLE
type(incar ) :: PINPT
type(params ) :: PPRAM
type(spmat ) :: SHk, SSk, SH0, SS0, SHm, SHs
integer*4 mpierr
integer*4 neig
integer*4 nemax !nemax <= neig * nspin. Choosing optimal value is critical for performance.
integer*4 ne_found(PINPT%nspin)
integer*4 ik
integer*4 ie, fe, im, fm, is
integer*4 feast_info
real*8 kp(3)
real*8 emin, emax
real*8 E(nemax*PINPT%nspin) ! will store all the energy eigenvalue for each spin
complex*16 V(neig*PINPT%ispin,nemax*PINPT%nspin) ! will store all the spin block at once in the first dimension
logical flag_vector, flag_init, flag_phase, flag_sparse_SHm, flag_sparse_SHs
real*8 t1, t0
character*4 timer
emin = PINPT%feast_emin ; emax = PINPT%feast_emax
#ifdef MPI
#ifdef PSPARSE
call MPI_BARRIER(COMM_JUELICH%mpi_comm, mpierr)
#endif
#endif
sp:do is = 1, PINPT%nspin
! if(flag_init) SHm, SHs will be kept during ik-run. (SHs will be modified if flag_slater_koster=.false.)
call time_check(t1,t0,timer)
if(.not. PPRAM%flag_use_overlap) then
call get_hamk_sparse(SHk, SH0, SHm, SHs, is, kp, &
PINPT, PPRAM, neig, NN_TABLE, flag_init, flag_phase, flag_sparse_SHm, flag_sparse_SHs)
elseif( PPRAM%flag_use_overlap) then
call get_hamk_sparse_overlap(SHk, SSk, SH0, SS0, SHm, SHs, is, kp, &
PINPT, PPRAM, neig, NN_TABLE, flag_init, flag_phase, flag_sparse_SHm, flag_sparse_SHs)
endif
call time_check(t1,t0)
if(timer .eq. 'init' .and. myid .eq. 0) then
write(message,'(A,F10.4,A)')' TIME for SPARSE MATRIX CONSTRUCTION: ',t1, ' (sec)' ; write_msg
timer = 'off'
else
timer = 'off'
endif
call get_matrix_index(ie, fe, im, fm, is, nemax, neig, PINPT%ispinor)
if(.not. PPRAM%flag_use_overlap) then
#ifdef PSPARSE
call cal_eig_hermitianx_psparse(SHk, emin, emax, nemax, ne_found(is), &
PINPT%feast_neguess, E(ie:fe), V(im:fm,ie:fe), flag_vector,&
PINPT%feast_fpm, feast_info)
#else
call cal_eig_hermitianx_sparse( SHk, emin, emax, nemax, ne_found(is), &
PINPT%feast_neguess, E(ie:fe), V(im:fm,ie:fe), flag_vector,&
PINPT%feast_fpm, feast_info)
#endif
elseif( PPRAM%flag_use_overlap) then
#ifdef PSPARSE
call cal_gen_eig_hermitianx_psparse(SHk, SSk, emin, emax, nemax, ne_found(is), &
PINPT%feast_neguess, E(ie:fe),V(im:fm,ie:fe), flag_vector, &
PINPT%feast_fpm, feast_info)
#else
call cal_gen_eig_hermitianx_sparse( SHk, SSk, emin, emax, nemax, ne_found(is), &
PINPT%feast_neguess, E(ie:fe),V(im:fm,ie:fe), flag_vector, &
PINPT%feast_fpm, feast_info)
#endif
endif
call adjust_ne_guess(feast_info, is, ne_found(is), kp, neig, nemax, PINPT)
if(allocated(SHk%H)) deallocate(SHk%H)
if(allocated(SHk%I)) deallocate(SHk%I)
if(allocated(SHk%J)) deallocate(SHk%J) ! SHk should be deallocated for each run
if(allocated(SSk%H)) deallocate(SSk%H)
if(allocated(SSk%I)) deallocate(SSk%I)
if(allocated(SSk%J)) deallocate(SSk%J)
enddo sp
if(PINPT%flag_soc .and. .not. PPRAM%flag_slater_koster) then
if(allocated(SHs%H)) deallocate(SHs%H)
if(allocated(SHs%I)) deallocate(SHs%I)
if(allocated(SHs%J)) deallocate(SHs%J)
endif
deallocate(SH0%H)
deallocate(SH0%I)
deallocate(SH0%J)
if(allocated(SS0%H)) deallocate(SS0%H)
if(allocated(SS0%I)) deallocate(SS0%I)
if(allocated(SS0%J)) deallocate(SS0%J)
#ifdef MPI
#ifdef PSPARSE
call MPI_BARRIER(COMM_JUELICH%mpi_comm, mpierr)
#endif
#endif
return
endsubroutine
#endif
subroutine cal_eig_Hk_dense(Hm, Hs, E, V, SV, PINPT, PPRAM, NN_TABLE, kp, neig, iband, nband, flag_vector, flag_init, flag_phase,ik)
use parameters, only : incar, hopping, pauli_0, params
use kronecker_prod, only: kproduct
use mpi_setup
use print_matrix
use time
use do_math
use phase_factor
type (hopping) :: NN_TABLE
type (incar ) :: PINPT
type (params ) :: PPRAM
integer*4 neig, iband, nband
integer*4 nkpoint, is, ie,fe, im, fm
integer*4 ik
real*8 kp(3)
real*8 E(nband*PINPT%nspin) ! will store all the energy eigenvalue for each spin
complex*16 V(neig*PINPT%ispin,nband*PINPT%nspin) ! will store all the spin block at once in the first dimension
complex*16 SV(neig*PINPT%ispin,nband*PINPT%nspin) ! will store all the spin block at once in the first dimension
complex*16 H0(neig,neig),S0(neig,neig) ! slater-koster hopping (k-dependent) hamiltonian and overlap matrix
complex*16 Hk(neig*PINPT%ispinor,neig*PINPT%ispinor) ! total hamiltonian (k-dependent)
complex*16 Sk(neig*PINPT%ispinor,neig*PINPT%ispinor) ! total overlap matrix (k-dependent)
complex*16 Sk_(neig*PINPT%ispinor,neig*PINPT%ispinor) ! temp. overlap matrix
complex*16 Hm(neig*PINPT%ispinor,neig*PINPT%ispinor) ! collinear magnetism hamiltonian (k-independent)
complex*16 Hs(neig*PINPT%ispinor,neig*PINPT%ispinor) ! 1st-order SO coupling hamiltonian (k-dependent if .not. SK)
logical flag_vector, flag_init, flag_phase, flag_load_nntable
real*8 t1, t0
real*8 E_(4)
character*20,external :: int2str
! This routine calculates all the eigenvalues within [iband:iband+nband-1] using Hamiltonian Hk with dense matrix format.
E_ = 0d0
flag_load_nntable = PINPT%flag_load_nntable
sp:do is = 1, PINPT%nspin
! if(flag_init) Hm, Hs will be kept during ik-run. (Hs will be modified if flag_slater_koster=.false.)
call get_hamk_dense(Hk, H0, Hm, Hs, is, kp, PINPT, PPRAM, neig, NN_TABLE, flag_init, flag_phase)
if(PINPT%flag_print_hamk) then
call print_matrix_c(Hk,neig*PINPT%ispinor, neig*PINPT%ispinor, &
'Hk_K'//trim(ADJUSTL(int2str(ik)))//'_SP'//trim(ADJUSTL(int2str(is))),1, 'F15.10')
endif
if(PPRAM%flag_use_overlap .and. is .eq. 1) then
call set_ham0(S0, kp, PPRAM, neig, NN_TABLE, F_IJ, flag_phase, .true., flag_load_nntable) ! .true. : flag_set_overlap
if(PINPT%flag_noncollinear) then
Sk = kproduct(pauli_0, S0, 2, 2, neig, neig)
else
Sk = S0
endif
if(PINPT%flag_print_hamk) then
call print_matrix_c(Sk,neig*PINPT%ispinor, neig*PINPT%ispinor, &
'Sk_K'//trim(ADJUSTL(int2str(ik)))//'_SP'//trim(ADJUSTL(int2str(is))),1, 'F15.10')
endif
endif
call get_matrix_index(ie, fe, im, fm, is, nband, neig, PINPT%ispinor)
if(.not. PPRAM%flag_use_overlap) then
call cal_eig(Hk, neig, PINPT%ispinor, PINPT%ispin, iband, nband, E(ie:fe), V(im:fm,ie:fe), flag_vector)
elseif( PPRAM%flag_use_overlap) then
if(flag_vector) then
Sk_ = Sk ! save temp (since Sk is modified after calling ZHEGE or ZHEGEX)
endif
call cal_gen_eig(Hk, Sk, neig, PINPT%ispinor, PINPT%ispin, iband, nband, E(ie:fe), V(im:fm,ie:fe), flag_vector)
if(flag_vector) then ! calculate V' = Sk|V> which will be used to calculate Mulliken charge rho = <V|Sk|V>
SV = matprod(fm-im+1, fm-im+1, 'N', Sk_, fm-im+1, fe-ie+1, 'N', V(im:fm,ie:fe))
endif
endif
enddo sp
return
endsubroutine
subroutine initialize_all(EE, neig, nband, nkp, my_nkp, PINPT, flag_vector, flag_overlap, flag_sparse, flag_stat, &
ii, iadd, t1, t0, flag_init)
use parameters, only: incar, energy
use mpi_setup
use time
implicit none
type(incar) :: PINPT
type(energy):: EE
integer*4 neig, nband, nkp, my_nkp
integer*4 iadd, ii
logical flag_vector, flag_stat, flag_init, flag_sparse, flag_overlap
real*8 t1, t0
integer*4 nL3
flag_init = .true.
call time_check(t1,t0,'init')
#ifdef MKL_SPARSE
if(flag_sparse) then
#ifdef PSPARSE
nL3 = 1
#ifdef MPI
call pfeastinit(PINPT%feast_fpm, COMM_JUELICH%mpi_comm, nL3)
#else
call feastinit(PINPT%feast_fpm)
#endif
#else
call feastinit(PINPT%feast_fpm)
#endif
if(allocated(PINPT%feast_ne)) deallocate(PINPT%feast_ne)
allocate(PINPT%feast_ne(PINPT%nspin, nkp))
PINPT%feast_ne = 0 !initialize to zero
PINPT%feast_neguess = PINPT%feast_nemax ! initialize to nemax. During calculations ne_guess will be adjusted using the ne_found in
! the previous step.
PINPT%feast_fpm(1) = 0 ! Specifies whether Extended Eigensolver routines print runtime status (0:F, 1:T)
PINPT%feast_fpm(2) = 4 ! The number of contour points N_e (see the description of FEAST algorithm
! Ref: E. Polizzi, Phys. Rev. B 79, 115112 (2009)
PINPT%feast_fpm(3) = 12 ! Error trace double precisiion stopping criteria e ( e = 10^-feast_fpm(3) )
PINPT%feast_fpm(4) = 50 ! Maximum number of Extended Eigensolver refinement loops allowed.
! If no convergence is reached within fpm(4) refinement loops,
! Extended Eigensolver routines return info=2.
PINPT%feast_fpm(5) = 0 ! User initial subspace. If fpm(5)=0 then Extended Eigensolver routines generate
! initial subspace, if fpm(5)=1 the user supplied initial subspace is used.
PINPT%feast_fpm(6) = 1 ! Extended Eigensolver stopping test.
! fpm(6)=0 : Extended Eigensolvers are stopped if the residual stopping test is satisfied.
! fpm(6)=1 : Extended Eigensolvers are stopped if this trace stopping test is satisfied.
PINPT%feast_fpm(14)= 0 ! If 1, return the computed eigenvectors subspace after one single contour integration.
endif
#endif
if(flag_stat) then
call initialize_eig_status(ii, iadd, nkp)
endif
allocate(EE%E(nband*PINPT%nspin, nkp))
EE%E = 0d0
if(.not. PINPT%flag_print_energy_singlek) then
if(allocated(EE%V)) deallocate(EE%V)
if(allocated(EE%SV)) deallocate(EE%SV)
allocate(EE%V(neig*PINPT%ispin, nband*PINPT%nspin, my_nkp))
allocate(EE%SV(neig*PINPT%ispin, nband*PINPT%nspin, my_nkp))
if(flag_vector) EE%V = (0.d0,0.d0)
if(flag_vector .and. flag_overlap) EE%SV = (0.d0,0.d0)
endif
return
endsubroutine
subroutine finalize_all(EE, SHm, SHs, t1, t0, PINPT, flag_stat, flag_vector, flag_sparse)
use parameters, only : energy, spmat, incar
use mpi_setup
use time
use print_io
implicit none
type(energy) :: EE
type(spmat ) :: SHm, SHs
type(incar ) :: PINPT
logical flag_vector, flag_stat, flag_sparse
real*8 t1, t0
call time_check(t1, t0)
if(flag_stat) then
write(message,*)' ' ; write_msg
write(message,'(A,F20.6)')" TIME for EIGENVALUE SOLVE (s)", t1 ; write_msg
endif
if(allocated(EE%E)) deallocate(EE%E)
if(allocated(EE%V)) deallocate(EE%V)
if(allocated(EE%SV)) deallocate(EE%SV)
if(flag_sparse) then
if(allocated(SHm%H)) deallocate(SHm%H)
if(allocated(SHm%I)) deallocate(SHm%I)
if(allocated(SHm%J)) deallocate(SHm%J)
if(allocated(SHs%H)) deallocate(SHs%H)
if(allocated(SHs%I)) deallocate(SHs%I)
if(allocated(SHs%J)) deallocate(SHs%J)
endif
return
endsubroutine
subroutine get_hamk_dense(Hk, H0, Hm, Hs, is, kpoint, PINPT, PPRAM, neig, NN_TABLE, flag_init, flag_phase)
use parameters, only: hopping, incar, pauli_0, pauli_x, pauli_y, pauli_z, params
use kronecker_prod, only: kproduct
use mpi_setup
use phase_factor
use do_math
use print_matrix
implicit none
type (hopping) :: NN_TABLE
type (incar ) :: PINPT
type (params ) :: PPRAM
logical flag_init, flag_phase, flag_load_nntable
integer*4 neig, is
real*8 kpoint(3)
complex*16 H0(neig,neig) ! slater-koster hopping (k-dependent)
complex*16 Hm(neig*PINPT%ispinor,neig*PINPT%ispinor) ! collinear magnetism hamiltonian (k-independent)
complex*16 Hs(neig*PINPT%ispinor,neig*PINPT%ispinor) ! 1st-order SO coupling hamiltonian (k-dependent if .not. SK)
complex*16 Hk(neig*PINPT%ispinor,neig*PINPT%ispinor) ! total hamiltonian (k-dependent)
logical flag_set_overlap
flag_load_nntable = PINPT%flag_load_nntable
flag_set_overlap = .false. ! in the first run set H
if(is .eq. 1) call set_ham0(H0, kpoint, PPRAM, neig, NN_TABLE, F_IJ, flag_phase, flag_set_overlap, flag_load_nntable)
if(flag_init) then
if(PINPT%flag_collinear) then
call set_ham_mag(Hm, NN_TABLE, PPRAM, neig, PINPT%ispinor, PINPT%flag_collinear, PINPT%flag_noncollinear)
elseif(PINPT%flag_noncollinear) then
call set_ham_mag(Hm, NN_TABLE, PPRAM, neig, PINPT%ispinor, PINPT%flag_collinear, PINPT%flag_noncollinear)
if(PINPT%flag_soc .and. PPRAM%flag_slater_koster) &
call set_ham_soc(Hs, 0d0, PPRAM, neig, NN_TABLE, F_IJ, flag_phase)
endif
flag_init = .false.
endif
if(PINPT%flag_collinear) then
Hk = H0 + ((-1d0)**(is+1)) * Hm
elseif(PINPT%flag_noncollinear) then
if(PINPT%flag_soc) then
!set up k-dependent SOC in the case of 'cc' orbitals
if(.not. PPRAM%flag_slater_koster) &
call set_ham_soc(Hs, kpoint, PPRAM, neig, NN_TABLE, F_IJ, flag_phase)
Hk = kproduct(pauli_0, H0, 2, 2, neig, neig) + Hm + Hs
else
Hk = kproduct(pauli_0, H0, 2, 2, neig, neig) + Hm
endif
elseif(.not. PINPT%flag_collinear .and. .not. PINPT%flag_noncollinear) then
Hk = H0
endif
return
endsubroutine
subroutine get_matrix_index(ie, fe, im, fm, is, nband, neig, ispinor)
implicit none
integer*4, intent(out):: ie, fe, im, fm
integer*4, intent(in) :: is, nband, neig, ispinor
! initial and final index for row index of E(ie:fe) and column index of V(:,ie:fe)
ie = 1 + (is-1)*nband
fe = nband + (is-1)*nband
! initial and final index for column index of V(im:fm,:)
im = 1 + (is-1)*neig
fm = neig*ispinor + (is-1)*neig
return
endsubroutine
subroutine set_ham0(H, kpoint, PPRAM, neig, NN_TABLE, FIJ, flag_phase, flag_set_overlap, flag_load_nntable)
use parameters, only : zi, hopping, params
use phase_factor
implicit none
interface
function FIJ(k,R)
complex*16 FIJ
real*8, intent(in) :: k(3)
real*8, intent(in) :: R(3)
endfunction
end interface
type (hopping) :: NN_TABLE
type (params ) :: PPRAM
integer*4 neig , i, j, nn
integer*4 ivel_axis
real*8 kpoint(3), tol, tij_sk, tij_cc
complex*16 H(neig ,neig)
complex*16 Eij
external tij_sk, tij_cc
logical flag_phase
logical flag_set_overlap
logical flag_load_nntable
H=0.0d0
tol=NN_TABLE%onsite_tolerance
do nn=1,NN_TABLE%n_neighbor
i=NN_TABLE%i_matrix(nn)
j=NN_TABLE%j_matrix(nn)
call get_hopping_integral(Eij, NN_TABLE, nn, PPRAM, tol, kpoint, FIJ, flag_phase, flag_set_overlap, flag_load_nntable)
if(.not. flag_set_overlap) then
if(i .eq. j .and. NN_TABLE%Dij(nn) <= tol) then
if(PPRAM%slater_koster_type .gt. 10) then
if(nint(PPRAM%param_const_nrl(4,1,NN_TABLE%local_U_param_index(i))) .ge. 1) then
H(i,j) = H(i,j) + Eij + NN_TABLE%local_charge(i)*PPRAM%param_const_nrl(5,1,(NN_TABLE%local_U_param_index(i)))
else
H(i,j) = H(i,j) + Eij + NN_TABLE%local_charge(i)*PPRAM%param_nrl(1,(NN_TABLE%local_U_param_index(i)))
endif
else
if(nint(PPRAM%param_const(4,NN_TABLE%local_U_param_index(i))) .ge. 1) then
H(i,j) = H(i,j) + Eij + NN_TABLE%local_charge(i)*PPRAM%param_const(5,(NN_TABLE%local_U_param_index(i)))
else
H(i,j) = H(i,j) + Eij + NN_TABLE%local_charge(i)*PPRAM%param((NN_TABLE%local_U_param_index(i)))
endif
endif
elseif(i .eq. j .and. NN_TABLE%Dij(nn) > tol) then
H(i,j) = H(i,j) + Eij
else
H(i,j) = H(i,j) + Eij
H(j,i) = H(j,i) + conjg(Eij)
endif
elseif( flag_set_overlap) then
if(i .eq. j .and. NN_TABLE%Dij(nn) <= tol) then
H(i,j) = H(i,j) + 1d0
elseif(i .eq. j .and. NN_TABLE%Dij(nn) > tol) then
H(i,j) = H(i,j) + Eij
else
H(i,j) = H(i,j) + Eij
H(j,i) = H(j,i) + conjg(Eij)
endif
endif
enddo
return
endsubroutine
subroutine set_ham_mag(H, NN_TABLE, PPRAM, neig, ispinor, flag_collinear, flag_noncollinear)
use parameters, only: zi, hopping, incar, params
implicit none
type (hopping) :: NN_TABLE
type (incar ) :: PINPT
type (params ) :: PPRAM
integer*4 neig, ispinor
integer*4 i, ii
complex*16 H(neig*ispinor,neig*ispinor)
logical flag_collinear, flag_noncollinear
H=0d0
if(flag_collinear) then
do i = 1, neig
if(NN_TABLE%stoner_I_param_index(i) .gt. 0) then ! if stoner parameter has been set...
if(PPRAM%slater_koster_type .gt. 10) then
if(nint(PPRAM%param_const_nrl(4,1,NN_TABLE%stoner_I_param_index(i))) .eq. 1) then ! if i-th basis has constraint .true.
H(i,i) = -0.5d0 * NN_TABLE%local_moment(1,i) * PPRAM%param_const_nrl(5,1,NN_TABLE%stoner_I_param_index(i))
else
H(i,i) = -0.5d0 * NN_TABLE%local_moment(1,i) * PPRAM%param_nrl(1,NN_TABLE%stoner_I_param_index(i))
endif
else
if(nint(PPRAM%param_const(4,NN_TABLE%stoner_I_param_index(i))) .eq. 1) then ! if i-th basis has constraint .true.
H(i,i) = -0.5d0 * NN_TABLE%local_moment(1,i) * PPRAM%param_const(5,NN_TABLE%stoner_I_param_index(i))
else
H(i,i) = -0.5d0 * NN_TABLE%local_moment(1,i) * PPRAM%param(NN_TABLE%stoner_I_param_index(i))
endif
endif
endif
enddo
elseif(flag_noncollinear) then
do i = 1, neig ! Hx
if(NN_TABLE%stoner_I_param_index(i) .gt. 0) then ! if stoner parameter has been set...
if(PPRAM%slater_koster_type .gt. 10) then
if(nint(PPRAM%param_const_nrl(4,1,NN_TABLE%stoner_I_param_index(i))) .eq. 1) then ! if i-th basis has constraint .true.
H(i,i+neig) = H(i,i+neig) - 0.5d0 * NN_TABLE%local_moment_rot(1,i) &
* PPRAM%param_const_nrl(5,1,NN_TABLE%stoner_I_param_index(i))
H(i+neig,i) = H(i+neig,i) - 0.5d0 * NN_TABLE%local_moment_rot(1,i) &
* PPRAM%param_const_nrl(5,1,NN_TABLE%stoner_I_param_index(i))
else
H(i,i+neig) = H(i,i+neig) - 0.5d0 * NN_TABLE%local_moment_rot(1,i) &
* PPRAM%param_nrl(1,NN_TABLE%stoner_I_param_index(i))
H(i+neig,i) = H(i+neig,i) - 0.5d0 * NN_TABLE%local_moment_rot(1,i) &
* PPRAM%param_nrl(1,NN_TABLE%stoner_I_param_index(i))
endif
else
if(nint(PPRAM%param_const(4,NN_TABLE%stoner_I_param_index(i))) .eq. 1) then ! if i-th basis has constraint .true.
H(i,i+neig) = H(i,i+neig) - 0.5d0 * NN_TABLE%local_moment_rot(1,i) &
* PPRAM%param_const(5,NN_TABLE%stoner_I_param_index(i))
H(i+neig,i) = H(i+neig,i) - 0.5d0 * NN_TABLE%local_moment_rot(1,i) &
* PPRAM%param_const(5,NN_TABLE%stoner_I_param_index(i))
else
H(i,i+neig) = H(i,i+neig) - 0.5d0 * NN_TABLE%local_moment_rot(1,i) &
* PPRAM%param(NN_TABLE%stoner_I_param_index(i))
H(i+neig,i) = H(i+neig,i) - 0.5d0 * NN_TABLE%local_moment_rot(1,i) &
* PPRAM%param(NN_TABLE%stoner_I_param_index(i))
endif
endif
endif
enddo
do i = 1, neig ! Hy
if(NN_TABLE%stoner_I_param_index(i) .gt. 0) then ! if stoner parameter has been set...
if(PPRAM%slater_koster_type .gt. 10) then
if(nint(PPRAM%param_const_nrl(4,1,NN_TABLE%stoner_I_param_index(i))) .eq. 1) then ! if i-th basis has constraint .true.
H(i,i+neig) = H(i,i+neig) + 0.5d0 * NN_TABLE%local_moment_rot(2,i) &
* PPRAM%param_const_nrl(5,1,NN_TABLE%stoner_I_param_index(i)) * zi
H(i+neig,i) = H(i+neig,i) - 0.5d0 * NN_TABLE%local_moment_rot(2,i) &
* PPRAM%param_const_nrl(5,1,NN_TABLE%stoner_I_param_index(i)) * zi
else
H(i,i+neig) = H(i,i+neig) + 0.5d0 * NN_TABLE%local_moment_rot(2,i) &
* PPRAM%param_nrl(1,NN_TABLE%stoner_I_param_index(i)) * zi
H(i,i+neig) = H(i,i+neig) - 0.5d0 * NN_TABLE%local_moment_rot(2,i) &
* PPRAM%param_nrl(1,NN_TABLE%stoner_I_param_index(i)) * zi
endif
else