-
Notifications
You must be signed in to change notification settings - Fork 15
/
find_nn.f90
868 lines (795 loc) · 47.3 KB
/
find_nn.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
#include "alias.inc"
subroutine set_nn_table(NN_TABLE, PINPT, PPRAM, PGEOM, flag_print_nn_table, flag_replot)
use parameters, only : hopping, incar, params, poscar
use mpi_setup
use print_io
implicit none
type (poscar) :: PGEOM
type (incar) :: PINPT
type (params) :: PPRAM
type (hopping) :: NN_TABLE
integer*4 mpierr
logical flag_print_nn_table, flag_replot
!Set neighbor hopping index
if(.not. flag_replot) then
allocate(NN_TABLE%local_moment(3,PGEOM%neig))
allocate(NN_TABLE%local_moment_rot(3,PGEOM%neig))
NN_TABLE%local_moment = PGEOM%local_moment
NN_TABLE%local_moment_rot = PGEOM%local_moment_rot
allocate(NN_TABLE%local_charge(PGEOM%neig))
NN_TABLE%local_charge = PGEOM%local_charge
allocate(NN_TABLE%site_cindex(PGEOM%n_atom))
allocate(NN_TABLE%flag_site_cindex(PGEOM%n_atom))
NN_TABLE%site_cindex = PGEOM%site_cindex
NN_TABLE%flag_site_cindex = PGEOM%flag_site_cindex
call find_nn(PINPT,PPRAM,PGEOM, NN_TABLE)
if(flag_print_nn_table .and. myid .eq. 0) then
call print_nn_table(NN_TABLE,PINPT,PPRAM)
endif
if(PINPT%flag_load_nntable .and. .not. PINPT%flag_tbfit .and. .not. PPRAM%flag_use_overlap) then
! LOAD from FILE
call load_nn_table(NN_TABLE, PINPT, PPRAM)
elseif(PINPT%flag_load_nntable .and. (PINPT%flag_tbfit .or. PPRAM%flag_use_overlap) ) then
write(message,'(A)')' !WARN! Reading hopping file cannot be combined with parameter fitting procedure or ' ; write_msgi
write(message,'(A)')' with overlap matrix constructions, i.e., using overlap integrals.' ; write_msgi
write(message,'(A)')' Please turn off LOAD_HOP or TBFIT option or do not use overlap integrals in' ; write_msgi
write(message,'(A)')' your PARAM_FIT.dat file. Exit..' ; write_msgi
kill_job
endif
endif
return
endsubroutine
subroutine find_nn(PINPT,PPRAM,PGEOM,NN_TABLE)
use parameters
use mpi_setup
use time
use memory
use print_io
implicit none
type (poscar) :: PGEOM
type (incar) :: PINPT
type (params) :: PPRAM
type (hopping) :: NN_TABLE
type (hopping) :: NN_TABLE_dummy
character(*), parameter :: func = 'find_nn'
integer*4, parameter :: max_neighbor = 30
integer*4 nn, i, j, iorb, jorb, ix, iy, iz, imatrix, jmatrix, max_nn, ii
integer*4 max_nn_pair
integer*4 add_overlap
integer*4 index_sigma,index_pi,index_delta !sk param index
integer*4 index_sigma_scale,index_pi_scale,index_delta_scale ! sk scale param index
integer*4 index_stoner !stoner param index
integer*4 index_custom
integer*4 index_custom_soc
integer*4 stoner_I_param_index ! stoner I param index
integer*4 local_U_param_index ! local U param index (for example, staggered potential)
integer*4 plus_U_param_index ! plus U param index
integer*4 l_onsite_param_index ! scaling parameter for local atomic density evaluation (valid if SK_SCALE_TYPE >=11)
integer*4 index_lambda !soc param index
integer*4 max_x, max_y, max_z
integer*4 size_NN_TABLE
real*8 max_nn_dist
real*8 Rij_(3),pos_i(3), pos_j(3), Dij_, Dij0_, Dijc_
real*8 R_(3)
real*8 enorm
real*8 tij_sk
real*8 tij_cc
real*8 onsite_tol, a1(3), a2(3), a3(3)
integer*4 onsite_param_index
character*2 param_class
character*20 soc_type
integer*4 nn_class
external enorm
external tij_sk, tij_cc
logical flag_use_site_cindex
#ifdef MPI
integer*4 nn_temp, nn_mpi(nprocs), nn_mpi_disp(0:nprocs-1), mpierr
nn_mpi = 0
#endif
if(PPRAM%flag_use_overlap) then
add_overlap = 6
else
add_overlap = 0
endif
max_nn= PGEOM%n_atom * max_neighbor * PGEOM%max_orb * PGEOM%max_orb * 10
if(PPRAM%slater_koster_type .gt. 10) max_nn_pair = 2000 ! temporary
allocate( NN_TABLE_dummy%i_atom(max_nn) )
allocate( NN_TABLE_dummy%j_atom(max_nn) )
allocate( NN_TABLE_dummy%i_coord(3,max_nn) )
allocate( NN_TABLE_dummy%j_coord(3,max_nn) )
allocate( NN_TABLE_dummy%Rij(3,max_nn) )
allocate( NN_TABLE_dummy%R (3,max_nn) )
allocate( NN_TABLE_dummy%Dij(max_nn) )
allocate( NN_TABLE_dummy%Dij0(max_nn) )
allocate( NN_TABLE_dummy%Dijc(max_nn) )
allocate( NN_TABLE_dummy%i_matrix(max_nn) )
allocate( NN_TABLE_dummy%ci_orb(max_nn) )
allocate( NN_TABLE_dummy%i_sign(max_nn) )
allocate( NN_TABLE_dummy%j_matrix(max_nn) )
allocate( NN_TABLE_dummy%cj_orb(max_nn) )
allocate( NN_TABLE_dummy%j_sign(max_nn) )
allocate( NN_TABLE_dummy%p_class(max_nn) )
allocate( NN_TABLE_dummy%n_class(max_nn) )
if(PPRAM%slater_koster_type .gt. 10) then
allocate( NN_TABLE%n_nn( PGEOM%n_atom) ) ; NN_TABLE%n_nn = 0 ! initialize
allocate( NN_TABLE_dummy%n_nn( PGEOM%n_atom) ) ; NN_TABLE_dummy%n_nn = 0
allocate( NN_TABLE_dummy%R_nn(max_nn_pair,PGEOM%n_atom) ) ; NN_TABLE_dummy%R_nn = 0d0
allocate( NN_TABLE_dummy%R0_nn(max_nn_pair,PGEOM%n_atom) ) ; NN_TABLE_dummy%R0_nn = 0d0
allocate( NN_TABLE_dummy%j_nn(max_nn_pair,PGEOM%n_atom) ) ; NN_TABLE_dummy%j_nn = 0
else
allocate( NN_TABLE%n_nn( PGEOM%n_atom) ) ; NN_TABLE%n_nn = 0
allocate( NN_TABLE_dummy%n_nn( PGEOM%n_atom) ) ; NN_TABLE_dummy%n_nn = 0
endif
if( PPRAM%flag_slater_koster) allocate( NN_TABLE_dummy%sk_index_set(0:6+add_overlap,max_nn) )
if( PPRAM%flag_slater_koster) NN_TABLE_dummy%sk_index_set(0,1:max_nn) = -9999d0
if(.not.PPRAM%flag_slater_koster) allocate( NN_TABLE_dummy%cc_index_set(0:3,max_nn) )
if(.not.PPRAM%flag_slater_koster) NN_TABLE_dummy%cc_index_set(0,1:max_nn) = -9999d0
allocate( NN_TABLE_dummy%tij(max_nn) ) ; NN_TABLE_dummy%tij = 0d0
if( PPRAM%flag_use_overlap) allocate( NN_TABLE_dummy%sij(max_nn) )
if(allocated(NN_TABLE_dummy%sij))NN_TABLE_dummy%sij = 0d0
allocate( NN_TABLE_dummy%soc_param_index(max_nn) )
allocate( NN_TABLE_dummy%site_cindex(PGEOM%n_atom) ) ! this argument will be used only in this routine 'find_nn'
NN_TABLE_dummy%site_cindex = NN_TABLE%site_cindex
NN_TABLE_dummy%soc_param_index = 0
allocate( NN_TABLE%stoner_I_param_index(PGEOM%neig) )
allocate( NN_TABLE%local_U_param_index(PGEOM%neig) )
allocate( NN_TABLE%plus_U_param_index(PGEOM%neig) )
NN_TABLE%stoner_I_param_index = 0
NN_TABLE%local_U_param_index = 0
NN_TABLE%plus_U_param_index = 0
if(PPRAM%slater_koster_type .gt. 10) then
allocate( NN_TABLE_dummy%l_onsite_param_index(PGEOM%n_atom) )
NN_TABLE_dummy%l_onsite_param_index = 0
endif
#ifdef MPI
allocate( NN_TABLE_dummy%stoner_I_param_index(PGEOM%neig) )
allocate( NN_TABLE_dummy%local_U_param_index(PGEOM%neig) )
allocate( NN_TABLE_dummy%plus_U_param_index(PGEOM%neig) )
NN_TABLE_dummy%stoner_I_param_index = 0
NN_TABLE_dummy%local_U_param_index = 0
NN_TABLE_dummy%plus_U_param_index = 0
#endif
onsite_tol = NN_TABLE%onsite_tolerance
a1=PGEOM%a_latt(1:3,1)
a2=PGEOM%a_latt(1:3,2)
a3=PGEOM%a_latt(1:3,3)
write(message,*)' ' ; write_msgi
write(message,*)'#- SETUP NEIGHBOR ATOM/ORBITAL PAIR & HOPPING CLASS' ; write_msgi
write(message,'(A)')' => Hopping pair (NN_PAIR) with following conditions will be searched::' ; write_msgi
do i = 1, PGEOM%n_nn_type
write(message,'(A,A12,2(A,1F8.4))')' NN_PAIR: ',PGEOM%nn_pair(i),' R0_max:',PGEOM%nn_dist(i),' R0:',PGEOM%nn_r0(i) ; write_msgi
enddo
call time_check(t1,t0,'init')
max_x = PINPT%nn_max(1)
max_y = PINPT%nn_max(2)
max_z = PINPT%nn_max(3)
max_nn_dist = maxval(PGEOM%nn_dist(:))
nn=0;
loop_i:do i=1+myid, PGEOM%n_atom, nprocs
if(PGEOM%n_orbital(i) .eq. 0) cycle loop_i
pos_i= PGEOM%a_coord(1,i)*a1(:) + &
PGEOM%a_coord(2,i)*a2(:) + &
PGEOM%a_coord(3,i)*a3(:)
if(PPRAM%slater_koster_type .gt. 10) then
call get_l_onsite_param_index(l_onsite_param_index, PPRAM, PGEOM%c_spec(PGEOM%spec(i)) )
NN_TABLE_dummy%l_onsite_param_index(i) = l_onsite_param_index
endif
do iorb=1,PGEOM%n_orbital(i)
imatrix= sum( PGEOM%n_orbital(1:i) ) - PGEOM%n_orbital(i) + iorb
do ix=-max_x,max_x
do iy=-max_y,max_y
do iz=-max_z,max_z
loop_j:do j=1,PGEOM%n_atom
if(PGEOM%n_orbital(j) .eq. 0) cycle loop_j
pos_j= (PGEOM%a_coord(1,j) + ix)*a1(:) + &
(PGEOM%a_coord(2,j) + iy)*a2(:) + &
(PGEOM%a_coord(3,j) + iz)*a3(:)
R_ = ix * a1(:) + iy * a2(:) + iz * a3(:)
Rij_= pos_j - pos_i
Dij_=enorm(3,Rij_)
!if(i .eq. 1 .and. j .eq. 5) write(6,'(A,5I3,6F10.5)')"ZZZZ ", i, j, ix,iy,iz,PGEOM%a_coord(:,i),PGEOM%a_coord(:,j)
if(Dij_ .gt. max_nn_dist) cycle loop_j
!call get_nn_class(PGEOM, i,j, Dij_, onsite_tol, nn_class, Dij0_) !, Dij0_cut)
call get_nn_class(PGEOM, i,j, Dij_, onsite_tol, nn_class, Dij0_, Dijc_)
if(nn_class .ne. -9999) then
if(nn_class .gt. 0 .and. PPRAM%slater_koster_type .gt. 10) then
NN_TABLE_dummy%n_nn(i) = NN_TABLE_dummy%n_nn(i) + 1
NN_TABLE_dummy%j_nn(NN_TABLE_dummy%n_nn(i),i) = j
NN_TABLE_dummy%R_nn(NN_TABLE_dummy%n_nn(i),i) = Dij_
NN_TABLE_dummy%R0_nn(NN_TABLE_dummy%n_nn(i),i)= Dij0_ ! I think it should be Dijc_... need to check later on: Jun. 08. 2022, KHJ
endif
do jorb=1,PGEOM%n_orbital(j)
jmatrix= sum( PGEOM%n_orbital(1:j) ) - PGEOM%n_orbital(j) + jorb
if(jmatrix .ge. imatrix) then
call get_param_class(PGEOM,iorb,jorb,i,j,param_class)
nn = nn + 1
NN_TABLE_dummy%i_atom(nn) = i
NN_TABLE_dummy%j_atom(nn) = j
NN_TABLE_dummy%i_coord(:,nn)= pos_i
NN_TABLE_dummy%j_coord(:,nn)= pos_j
NN_TABLE_dummy%Rij(:,nn) = Rij_(:)
NN_TABLE_dummy%R (:,nn) = R_ (:)
NN_TABLE_dummy%Dij(nn) = Dij_
NN_TABLE_dummy%Dij0(nn) = Dij0_
NN_TABLE_dummy%Dijc(nn) = Dijc_
NN_TABLE_dummy%i_matrix(nn) = imatrix
NN_TABLE_dummy%ci_orb(nn) = PGEOM%c_orbital(iorb,i)
NN_TABLE_dummy%i_sign(nn) = PGEOM%orb_sign(iorb,i)
NN_TABLE_dummy%j_matrix(nn) = jmatrix
NN_TABLE_dummy%cj_orb(nn) = PGEOM%c_orbital(jorb,j)
NN_TABLE_dummy%j_sign(nn) = PGEOM%orb_sign(jorb,j)
NN_TABLE_dummy%p_class(nn) = param_class
NN_TABLE_dummy%n_class(nn) = nn_class
if( PPRAM%flag_slater_koster) NN_TABLE_dummy%sk_index_set(0:6+add_overlap,nn)= 0 ! initialize
if(.not. PPRAM%flag_slater_koster) NN_TABLE_dummy%cc_index_set(0:3,nn)= 0 ! initialize
if(nn_class .eq. 0) then
!SET ONSITE energies
call get_onsite_param_index(onsite_param_index, PPRAM, &
PGEOM%c_orbital(iorb,i), &
PGEOM%c_orbital(jorb,j), &
PGEOM%c_spec(PGEOM%spec(i)))
if( PPRAM%flag_slater_koster) then
NN_TABLE_dummy%sk_index_set(0,nn) = onsite_param_index
if(PPRAM%slater_koster_type .le. 10) then
NN_TABLE_dummy%tij(nn) = tij_sk(NN_TABLE_dummy,nn,PPRAM,onsite_tol,.false. )
endif
elseif(.not.PPRAM%flag_slater_koster) then
if(param_class .eq. 'cc') then
NN_TABLE_dummy%cc_index_set(0,nn) = onsite_param_index
NN_TABLE_dummy%tij(nn) = tij_cc(NN_TABLE_dummy,nn,PPRAM,onsite_tol)
else
NN_TABLE_dummy%tij(nn) = 0d0
endif
endif
!SET local potentail
if(PINPT%flag_local_charge) then
call get_local_U_param_index(local_U_param_index, PPRAM, nn_class, param_class, PGEOM%c_spec(PGEOM%spec(i)) )
NN_TABLE%local_U_param_index(imatrix)=local_U_param_index
endif
!SET Hubbard type +U
if(PINPT%flag_plus_U) then
call get_plus_U_param_index(plus_U_param_index, PPRAM, nn_class, param_class, PGEOM%c_spec(PGEOM%spec(i)) )
NN_TABLE%plus_U_param_index(imatrix)= plus_U_param_index
endif
!SET local magnetic moment
if(PINPT%flag_collinear .or. PINPT%flag_noncollinear) then
call get_stoner_I_param_index(stoner_I_param_index, PPRAM, nn_class, param_class, PGEOM%c_spec(PGEOM%spec(i)) )
NN_TABLE%stoner_I_param_index(imatrix)=stoner_I_param_index
endif
!SET SOC parameter
if(PINPT%flag_soc) then
if(param_class .eq. 'pp' .or. param_class .eq. 'dd' .or. param_class .eq. 'xx') then
call get_soc_param_index(index_lambda, PGEOM%c_orbital(iorb,i), PGEOM%c_orbital(jorb,j), &
PGEOM%c_spec(PGEOM%spec(i)), PPRAM, param_class )
NN_TABLE_dummy%soc_param_index(nn) = index_lambda
elseif(param_class .eq. 'cc') then
! in-plane SOC in the lattice model
soc_type = 'lsoc'
call get_soc_cc_param_index(index_custom_soc, NN_TABLE_dummy, nn, PPRAM, &
PGEOM%c_spec(PGEOM%spec(i)), PGEOM%c_spec(PGEOM%spec(j)), soc_type )
NN_TABLE_dummy%cc_index_set(2,nn) = index_custom_soc ! SOC due to in-plane symmetry breaking
! out-of-plane SOC in the lattice model => rashba SOC due to out-of-plane symmetry breaking or E_field(z)
soc_type = 'lrashba'
call get_soc_cc_param_index(index_custom_soc, NN_TABLE_dummy, nn, PPRAM, &
PGEOM%c_spec(PGEOM%spec(i)), PGEOM%c_spec(PGEOM%spec(j)), soc_type )
if(index_custom_soc .eq. 0) then
soc_type = 'lR'
call get_soc_cc_param_index(index_custom_soc, NN_TABLE_dummy, nn, PPRAM, &
PGEOM%c_spec(PGEOM%spec(i)), PGEOM%c_spec(PGEOM%spec(j)), soc_type )
endif
NN_TABLE_dummy%cc_index_set(3,nn) = index_custom_soc
endif
endif
elseif( nn_class .ge. 1) then
!SET HOPPING energies
if(PPRAM%flag_slater_koster) then
! CASE: SLATER_KOSTER TYPE HOPPING
flag_use_site_cindex = .false. ! BE CAREFUL ! This is for TaS2 system only...
call get_sk_index_set(index_sigma,index_pi,index_delta, &
index_sigma_scale,index_pi_scale,index_delta_scale, &
PPRAM, param_class, nn_class, &
PGEOM, i, j, &
NN_TABLE%site_cindex(i), NN_TABLE%site_cindex(j), flag_use_site_cindex, .false. )
NN_TABLE_dummy%sk_index_set(1,nn) = index_sigma
NN_TABLE_dummy%sk_index_set(2,nn) = index_pi
NN_TABLE_dummy%sk_index_set(3,nn) = index_delta
NN_TABLE_dummy%sk_index_set(4,nn) = index_sigma_scale
NN_TABLE_dummy%sk_index_set(5,nn) = index_pi_scale
NN_TABLE_dummy%sk_index_set(6,nn) = index_delta_scale
if(PPRAM%slater_koster_type .le. 10) then
NN_TABLE_dummy%tij(nn) = tij_sk(NN_TABLE_dummy,nn,PPRAM,onsite_tol,.false.)
endif
if(PPRAM%flag_use_overlap) then
call get_sk_index_set(index_sigma,index_pi,index_delta, &
index_sigma_scale,index_pi_scale,index_delta_scale, &
PPRAM, param_class, nn_class, &
PGEOM, i, j, &
NN_TABLE%site_cindex(i), NN_TABLE%site_cindex(j), flag_use_site_cindex, .true. )
NN_TABLE_dummy%sk_index_set(1+add_overlap,nn) = index_sigma
NN_TABLE_dummy%sk_index_set(2+add_overlap,nn) = index_pi
NN_TABLE_dummy%sk_index_set(3+add_overlap,nn) = index_delta
NN_TABLE_dummy%sk_index_set(4+add_overlap,nn) = index_sigma_scale
NN_TABLE_dummy%sk_index_set(5+add_overlap,nn) = index_pi_scale
NN_TABLE_dummy%sk_index_set(6+add_overlap,nn) = index_delta_scale
if(PPRAM%slater_koster_type .le. 10) then
NN_TABLE_dummy%sij(nn) = tij_sk(NN_TABLE_dummy,nn,PPRAM,onsite_tol,.true.)
endif
endif
elseif(.not. PPRAM%flag_slater_koster) then ! cc index : custum orbital hopping index
! CASE: USER DEFINED CUSTOM HOPPING
if(param_class .eq. 'cc') then
! SET T_IJ
call get_cc_index_set(index_custom, NN_TABLE_dummy, nn, PPRAM, &
PGEOM%c_spec(PGEOM%spec(i)), PGEOM%c_spec(PGEOM%spec(j)) )
NN_TABLE_dummy%cc_index_set(1,nn) = index_custom
NN_TABLE_dummy%tij(nn) = tij_cc(NN_TABLE_dummy,nn,PPRAM,onsite_tol)
if(PINPT%flag_soc) then
! in-plane SOC in the lattice model
soc_type = 'lsoc'
call get_soc_cc_param_index(index_custom_soc, NN_TABLE_dummy, nn, PPRAM, &
PGEOM%c_spec(PGEOM%spec(i)), PGEOM%c_spec(PGEOM%spec(j)), soc_type )
NN_TABLE_dummy%cc_index_set(2,nn) = index_custom_soc
! out-of-plane SOC in the lattice model => rashba SOC due to out-of-plane symmetry breaking or E_field(z)
soc_type = 'lrashba'
call get_soc_cc_param_index(index_custom_soc, NN_TABLE_dummy, nn, PPRAM, &
PGEOM%c_spec(PGEOM%spec(i)), PGEOM%c_spec(PGEOM%spec(j)), soc_type )
if(index_custom_soc .eq. 0) then
soc_type = 'lR'
call get_soc_cc_param_index(index_custom_soc, NN_TABLE_dummy, nn, PPRAM, &
PGEOM%c_spec(PGEOM%spec(i)), PGEOM%c_spec(PGEOM%spec(j)), soc_type )
endif
NN_TABLE_dummy%cc_index_set(3,nn) = index_custom_soc
endif
else
NN_TABLE_dummy%tij(nn) = 0d0
endif
endif
endif ! nn_class
endif
enddo !jorb
endif
enddo loop_j
enddo
enddo
enddo !cell
enddo !iorb
enddo loop_i
#ifdef MPI
nn_mpi(myid+1) = nn
call MPI_ALLGATHER(nn,1, MPI_INTEGER4, nn_mpi,1,MPI_INTEGER4, mpi_comm_earth, mpierr)
nn_mpi_disp(0) = 0
do ii = 1, nprocs-1
nn_mpi_disp(ii)= nn_mpi_disp(ii - 1) + nn_mpi(ii)
enddo
call MPI_ALLREDUCE(nn, nn_temp, 1, MPI_INTEGER4, MPI_SUM, mpi_comm_earth, mpierr)
nn = nn_temp
if(PPRAM%slater_koster_type .gt. 10) then
call MPI_ALLREDUCE(NN_TABLE_dummy%n_nn, NN_TABLE%n_nn, PGEOM%n_atom, MPI_INTEGER4, MPI_SUM, mpi_comm_earth, mpierr)
NN_TABLE%max_nn_pair = maxval(NN_TABLE%n_nn) ! update max_nn_pair
max_nn_pair = NN_TABLE%max_nn_pair
endif
#else
if(PPRAM%slater_koster_type .gt. 10) then
NN_TABLE%max_nn_pair = maxval(NN_TABLE%n_nn) ! update max_nn_pair
max_nn_pair = NN_TABLE%max_nn_pair
endif
#endif
NN_TABLE%n_neighbor = nn
if (nn .gt. max_nn) then
write(message,'(A,I8,A,A)')' !WARN! Total number of Neighbor pair is exeed MAX_NN=100*N_ATOM*MAX_ORB=',max_nn, &
' Exit... Please recompile with larger MAX_NN', func ; write_msgi
kill_job
endif
allocate( NN_TABLE%i_atom(nn) )
allocate( NN_TABLE%j_atom(nn) )
allocate( NN_TABLE%i_coord(3,nn))
allocate( NN_TABLE%j_coord(3,nn))
allocate( NN_TABLE%Rij(3,nn) )
allocate( NN_TABLE%R (3,nn) )
allocate( NN_TABLE%Dij(nn) )
allocate( NN_TABLE%Dij0(nn) )
allocate( NN_TABLE%Dijc(nn) )
allocate( NN_TABLE%i_matrix(nn) )
allocate( NN_TABLE%ci_orb(nn) )
allocate( NN_TABLE%i_sign(nn) )
allocate( NN_TABLE%j_matrix(nn) )
allocate( NN_TABLE%cj_orb(nn) )
allocate( NN_TABLE%j_sign(nn) )
allocate( NN_TABLE%p_class(nn) )
allocate( NN_TABLE%n_class(nn) )
if(PPRAM%slater_koster_type .gt. 10) then
if(max_nn_pair .ge. 1) then
if(allocated(NN_TABLE_dummy%R_nn ))allocate( NN_TABLE%R_nn(max_nn_pair,PGEOM%n_atom) )
if(allocated(NN_TABLE_dummy%R0_nn))allocate( NN_TABLE%R0_nn(max_nn_pair,PGEOM%n_atom))
if(allocated(NN_TABLE_dummy%j_nn ))allocate( NN_TABLE%j_nn(max_nn_pair,PGEOM%n_atom) )
else
if(allocated(NN_TABLE_dummy%R_nn ))allocate( NN_TABLE%R_nn(1,PGEOM%n_atom) )
if(allocated(NN_TABLE_dummy%R0_nn))allocate( NN_TABLE%R0_nn(1,PGEOM%n_atom) )
if(allocated(NN_TABLE_dummy%j_nn ))allocate( NN_TABLE%j_nn(1,PGEOM%n_atom) )
endif
allocate( NN_TABLE%l_onsite_param_index(PGEOM%n_atom) )
endif
if( PPRAM%flag_slater_koster) allocate( NN_TABLE%sk_index_set(0:6+add_overlap,nn) )
if(.not.PPRAM%flag_slater_koster) allocate( NN_TABLE%cc_index_set(0:3,nn) )
allocate( NN_TABLE%tij(nn) )
if(allocated(NN_TABLE_dummy%sij)) allocate( NN_TABLE%sij(nn) )
if( PINPT%flag_load_nntable ) allocate( NN_TABLE%tij_file(nn) )
if( PINPT%flag_load_nntable .and. allocated(NN_TABLE%sij) ) allocate( NN_TABLE%sij_file(nn) )
allocate( NN_TABLE%soc_param_index(nn) )
#ifdef MPI
call MPI_ALLGATHERV(NN_TABLE_dummy%i_atom(1:nn_mpi(myid+1)), nn_mpi(myid+1), MPI_INTEGER4, NN_TABLE%i_atom, &
nn_mpi, nn_mpi_disp, MPI_INTEGER4, mpi_comm_earth, mpierr)
call MPI_ALLGATHERV(NN_TABLE_dummy%j_atom(1:nn_mpi(myid+1)), nn_mpi(myid+1), MPI_INTEGER4, NN_TABLE%j_atom, &
nn_mpi, nn_mpi_disp, MPI_INTEGER4, mpi_comm_earth, mpierr)
call MPI_ALLGATHERV(NN_TABLE_dummy%Dij(1:nn_mpi(myid+1)), nn_mpi(myid+1), MPI_REAL8, NN_TABLE%Dij, &
nn_mpi, nn_mpi_disp, MPI_REAL8, mpi_comm_earth, mpierr)
call MPI_ALLGATHERV(NN_TABLE_dummy%Dij0(1:nn_mpi(myid+1)), nn_mpi(myid+1), MPI_REAL8, NN_TABLE%Dij0, &
nn_mpi, nn_mpi_disp, MPI_REAL8, mpi_comm_earth, mpierr)
call MPI_ALLGATHERV(NN_TABLE_dummy%Dijc(1:nn_mpi(myid+1)), nn_mpi(myid+1), MPI_REAL8, NN_TABLE%Dijc, &
nn_mpi, nn_mpi_disp, MPI_REAL8, mpi_comm_earth, mpierr)
call MPI_ALLGATHERV(NN_TABLE_dummy%i_matrix(1:nn_mpi(myid+1)), nn_mpi(myid+1), MPI_INTEGER4, NN_TABLE%i_matrix, &
nn_mpi, nn_mpi_disp, MPI_INTEGER4, mpi_comm_earth, mpierr)
call MPI_ALLGATHERV(NN_TABLE_dummy%j_matrix(1:nn_mpi(myid+1)), nn_mpi(myid+1), MPI_INTEGER4, NN_TABLE%j_matrix, &
nn_mpi, nn_mpi_disp, MPI_INTEGER4, mpi_comm_earth, mpierr)
call MPI_ALLGATHERV(NN_TABLE_dummy%n_class(1:nn_mpi(myid+1)), nn_mpi(myid+1), MPI_INTEGER4, NN_TABLE%n_class, &
nn_mpi, nn_mpi_disp, MPI_INTEGER4, mpi_comm_earth, mpierr)
call MPI_ALLGATHERV(NN_TABLE_dummy%tij(1:nn_mpi(myid+1)), nn_mpi(myid+1), MPI_REAL8, NN_TABLE%tij, &
nn_mpi, nn_mpi_disp, MPI_REAL8, mpi_comm_earth, mpierr)
if(allocated(NN_TABLE_dummy%sij)) then
call MPI_ALLGATHERV(NN_TABLE_dummy%sij(1:nn_mpi(myid+1)), nn_mpi(myid+1), MPI_REAL8, NN_TABLE%sij, &
nn_mpi, nn_mpi_disp, MPI_REAL8, mpi_comm_earth, mpierr)
endif
call MPI_ALLGATHERV(NN_TABLE_dummy%soc_param_index(1:nn_mpi(myid+1)), nn_mpi(myid+1), MPI_INTEGER4, NN_TABLE%soc_param_index, &
nn_mpi, nn_mpi_disp, MPI_INTEGER4, mpi_comm_earth, mpierr)
do ii = 1, 3
call MPI_ALLGATHERV(NN_TABLE_dummy%i_coord(ii,1:nn_mpi(myid+1)), size(NN_TABLE_dummy%i_coord(ii,1:nn_mpi(myid+1))), MPI_REAL8, &
NN_TABLE%i_coord(ii,:), nn_mpi, nn_mpi_disp, MPI_REAL8, mpi_comm_earth, mpierr)
call MPI_ALLGATHERV(NN_TABLE_dummy%j_coord(ii,1:nn_mpi(myid+1)), size(NN_TABLE_dummy%j_coord(ii,1:nn_mpi(myid+1))), MPI_REAL8, &
NN_TABLE%j_coord(ii,:), nn_mpi, nn_mpi_disp, MPI_REAL8, mpi_comm_earth, mpierr)
call MPI_ALLGATHERV(NN_TABLE_dummy%Rij(ii,1:nn_mpi(myid+1)), size(NN_TABLE_dummy%Rij(ii,1:nn_mpi(myid+1))), MPI_REAL8, &
NN_TABLE%Rij(ii,:), nn_mpi, nn_mpi_disp, MPI_REAL8, mpi_comm_earth, mpierr)
call MPI_ALLGATHERV(NN_TABLE_dummy%R(ii,1:nn_mpi(myid+1)), size(NN_TABLE_dummy%R(ii,1:nn_mpi(myid+1))), MPI_REAL8, &
NN_TABLE%R(ii,:), nn_mpi, nn_mpi_disp, MPI_REAL8, mpi_comm_earth, mpierr)
enddo
if(PPRAM%flag_slater_koster) then
do ii = 0, 6+add_overlap
call MPI_ALLGATHERV(NN_TABLE_dummy%sk_index_set(ii,1:nn_mpi(myid+1)), size(NN_TABLE_dummy%sk_index_set(ii,1:nn_mpi(myid+1))), MPI_INTEGER4, &
NN_TABLE%sk_index_set(ii,:), nn_mpi, nn_mpi_disp, MPI_INTEGER4, mpi_comm_earth, mpierr)
enddo
elseif(.not. PPRAM%flag_slater_koster) then
do ii = 0, 3
call MPI_ALLGATHERV(NN_TABLE_dummy%cc_index_set(ii,1:nn_mpi(myid+1)), size(NN_TABLE_dummy%cc_index_set(ii,1:nn_mpi(myid+1))), MPI_INTEGER4, &
NN_TABLE%cc_index_set(ii,:), nn_mpi, nn_mpi_disp, MPI_INTEGER4, mpi_comm_earth, mpierr)
enddo
endif
if(PPRAM%slater_koster_type .gt. 10 .and. max_nn_pair .gt. 0) then
do ii = 1, PGEOM%n_atom
call MPI_ALLREDUCE(NN_TABLE_dummy%R_nn(1:max_nn_pair,ii),NN_TABLE%R_nn(1:max_nn_pair,ii), max_nn_pair, MPI_REAL8 , MPI_SUM, mpi_comm_earth, mpierr)
call MPI_ALLREDUCE(NN_TABLE_dummy%R0_nn(1:max_nn_pair,ii),NN_TABLE%R0_nn(1:max_nn_pair,ii),max_nn_pair, MPI_REAL8 , MPI_SUM, mpi_comm_earth, mpierr)
call MPI_ALLREDUCE(NN_TABLE_dummy%j_nn(1:max_nn_pair,ii),NN_TABLE%j_nn(1:max_nn_pair,ii), max_nn_pair, MPI_INTEGER4, MPI_SUM, mpi_comm_earth, mpierr)
enddo
call MPI_ALLREDUCE(NN_TABLE_dummy%l_onsite_param_index, NN_TABLE%l_onsite_param_index, PGEOM%n_atom, MPI_INTEGER4, MPI_SUM, mpi_comm_earth, mpierr)
endif
call MPI_ALLREDUCE(NN_TABLE%stoner_I_param_index, NN_TABLE_dummy%stoner_I_param_index, PGEOM%neig, MPI_INTEGER4, MPI_SUM, mpi_comm_earth, mpierr)
call MPI_ALLREDUCE(NN_TABLE%local_U_param_index, NN_TABLE_dummy%local_U_param_index, PGEOM%neig, MPI_INTEGER4, MPI_SUM, mpi_comm_earth, mpierr)
call MPI_ALLREDUCE(NN_TABLE%plus_U_param_index, NN_TABLE_dummy%plus_U_param_index, PGEOM%neig, MPI_INTEGER4, MPI_SUM, mpi_comm_earth, mpierr)
NN_TABLE%stoner_I_param_index(1:PGEOM%neig) = NN_TABLE_dummy%stoner_I_param_index(1:PGEOM%neig)
NN_TABLE%local_U_param_index(1:PGEOM%neig) = NN_TABLE_dummy%local_U_param_index(1:PGEOM%neig)
NN_TABLE%plus_U_param_index(1:PGEOM%neig) = NN_TABLE_dummy%plus_U_param_index(1:PGEOM%neig)
call MPI_ALLGATHERV(NN_TABLE_dummy%ci_orb(1:nn_mpi(myid+1)), nn_mpi(myid+1)*8, MPI_CHAR, NN_TABLE%ci_orb, &
nn_mpi*8, nn_mpi_disp*8, MPI_CHAR, mpi_comm_earth, mpierr)
call MPI_ALLGATHERV(NN_TABLE_dummy%cj_orb(1:nn_mpi(myid+1)), nn_mpi(myid+1)*8, MPI_CHAR, NN_TABLE%cj_orb, &
nn_mpi*8, nn_mpi_disp*8, MPI_CHAR, mpi_comm_earth, mpierr)
call MPI_ALLGATHERV(NN_TABLE_dummy%i_sign(1:nn_mpi(myid+1)), size(NN_TABLE_dummy%i_sign(1:nn_mpi(myid+1))), MPI_REAL8, NN_TABLE%i_sign, &
nn_mpi, nn_mpi_disp, MPI_REAL8, mpi_comm_earth, mpierr)
call MPI_ALLGATHERV(NN_TABLE_dummy%j_sign(1:nn_mpi(myid+1)), size(NN_TABLE_dummy%j_sign(1:nn_mpi(myid+1))), MPI_REAL8, NN_TABLE%j_sign, &
nn_mpi, nn_mpi_disp, MPI_REAL8, mpi_comm_earth, mpierr)
call MPI_ALLGATHERV(NN_TABLE_dummy%p_class(1:nn_mpi(myid+1)), nn_mpi(myid+1)*2, MPI_CHAR, NN_TABLE%p_class, &
nn_mpi*2, nn_mpi_disp*2, MPI_CHAR, mpi_comm_earth, mpierr)
#else
NN_TABLE%i_atom(1:nn) = NN_TABLE_dummy%i_atom(1:nn)
NN_TABLE%j_atom(1:nn) = NN_TABLE_dummy%j_atom(1:nn)
NN_TABLE%i_coord(:,1:nn) = NN_TABLE_dummy%i_coord(:,1:nn)
NN_TABLE%j_coord(:,1:nn) = NN_TABLE_dummy%j_coord(:,1:nn)
NN_TABLE%Rij(:,1:nn) = NN_TABLE_dummy%Rij(:,1:nn)
NN_TABLE%R (:,1:nn) = NN_TABLE_dummy%R (:,1:nn)
NN_TABLE%Dij(1:nn) = NN_TABLE_dummy%Dij(1:nn)
NN_TABLE%Dij0(1:nn) = NN_TABLE_dummy%Dij0(1:nn)
NN_TABLE%Dijc(1:nn) = NN_TABLE_dummy%Dijc(1:nn)
NN_TABLE%i_matrix(1:nn) = NN_TABLE_dummy%i_matrix(1:nn)
NN_TABLE%ci_orb(1:nn) = NN_TABLE_dummy%ci_orb(1:nn)
NN_TABLE%i_sign(1:nn) = NN_TABLE_dummy%i_sign(1:nn)
NN_TABLE%j_matrix(1:nn) = NN_TABLE_dummy%j_matrix(1:nn)
NN_TABLE%cj_orb(1:nn) = NN_TABLE_dummy%cj_orb(1:nn)
NN_TABLE%j_sign(1:nn) = NN_TABLE_dummy%j_sign(1:nn)
NN_TABLE%p_class(1:nn) = NN_TABLE_dummy%p_class(1:nn)
NN_TABLE%n_class(1:nn) = NN_TABLE_dummy%n_class(1:nn)
if(PPRAM%slater_koster_type .gt. 10) then
if(max_nn_pair .gt. 0) then
NN_TABLE%R_nn(1:max_nn_pair,1:PGEOM%n_atom) = NN_TABLE_dummy%R_nn(1:max_nn_pair,1:PGEOM%n_atom)
NN_TABLE%R0_nn(1:max_nn_pair,1:PGEOM%n_atom) = NN_TABLE_dummy%R0_nn(1:max_nn_pair,1:PGEOM%n_atom)
NN_TABLE%j_nn(1:max_nn_pair,1:PGEOM%n_atom) = NN_TABLE_dummy%j_nn(1:max_nn_pair,1:PGEOM%n_atom)
else ! find no nearest neighbor pair
NN_TABLE%R_nn(1,1:PGEOM%n_atom) = 0d0
NN_TABLE%R0_nn(1,1:PGEOM%n_atom) = 0d0
NN_TABLE%j_nn(1,1:PGEOM%n_atom) = 0
endif
NN_TABLE%l_onsite_param_index(1:PGEOM%n_atom) = NN_TABLE_dummy%l_onsite_param_index(1:PGEOM%n_atom)
endif
if( PPRAM%flag_slater_koster) NN_TABLE%sk_index_set(0:6+add_overlap,1:nn) = NN_TABLE_dummy%sk_index_set(0:6+add_overlap,1:nn)
if(.not.PPRAM%flag_slater_koster) NN_TABLE%cc_index_set(0:3,1:nn) = NN_TABLE_dummy%cc_index_set(0:3,1:nn)
NN_TABLE%tij(1:nn) = NN_TABLE_dummy%tij(1:nn)
if(allocated(NN_TABLE%sij)) NN_TABLE%sij(1:nn) = NN_TABLE_dummy%sij(1:nn)
NN_TABLE%soc_param_index(1:nn) = NN_TABLE_dummy%soc_param_index(1:nn)
#endif
deallocate( NN_TABLE_dummy%i_atom )
deallocate( NN_TABLE_dummy%j_atom )
deallocate( NN_TABLE_dummy%i_coord )
deallocate( NN_TABLE_dummy%j_coord )
deallocate( NN_TABLE_dummy%Rij )
deallocate( NN_TABLE_dummy%R )
deallocate( NN_TABLE_dummy%Dij )
deallocate( NN_TABLE_dummy%Dij0 )
deallocate( NN_TABLE_dummy%Dijc )
deallocate( NN_TABLE_dummy%i_matrix )
deallocate( NN_TABLE_dummy%ci_orb )
deallocate( NN_TABLE_dummy%i_sign )
deallocate( NN_TABLE_dummy%j_matrix )
deallocate( NN_TABLE_dummy%cj_orb )
deallocate( NN_TABLE_dummy%j_sign )
deallocate( NN_TABLE_dummy%p_class )
deallocate( NN_TABLE_dummy%n_class )
if( PPRAM%flag_slater_koster) deallocate( NN_TABLE_dummy%sk_index_set )
if(.not.PPRAM%flag_slater_koster) deallocate( NN_TABLE_dummy%cc_index_set )
deallocate( NN_TABLE_dummy%tij )
if(allocated(NN_TABLE_dummy%sij)) deallocate( NN_TABLE_dummy%sij )
deallocate( NN_TABLE_dummy%soc_param_index )
if(allocated(NN_TABLE_dummy%n_nn)) deallocate( NN_TABLE_dummy%n_nn )
if(allocated(NN_TABLE_dummy%R_nn)) deallocate( NN_TABLE_dummy%R_nn )
if(allocated(NN_TABLE_dummy%R0_nn)) deallocate( NN_TABLE_dummy%R0_nn )
if(allocated(NN_TABLE_dummy%j_nn)) deallocate( NN_TABLE_dummy%j_nn )
if(allocated(NN_TABLE_dummy%l_onsite_param_index)) deallocate( NN_TABLE_dummy%l_onsite_param_index )
#ifdef MPI
deallocate( NN_TABLE_dummy%stoner_I_param_index )
deallocate( NN_TABLE_dummy%local_U_param_index )
deallocate( NN_TABLE_dummy%plus_U_param_index )
#endif
write(message,*)' ' ; write_msgi
write(message,'(A,I0,A)')' N_NEIGH: found ',NN_TABLE%n_neighbor, ' hopping pair' ; write_msgi
write(message,'(A,A )')' The hopping information will be written in ', 'hopping'//trim(PINPT%title(NN_TABLE%mysystem))//'.dat'; write_msgi
size_NN_TABLE=sizeof(NN_TABLE)
!if_main call report_memory(int8(size_NN_TABLE), 1, 'NN_TABLE')
call time_check(t1,t0,'end')
write(message,'(A,F12.6)')" TIME for FINDING NEIGHBOR PAIRS (s)", t1 ; write_msgi
write(message,*)' ' ; write_msgi
write(message,*)'#- END SETUP NEIGHBOR ATOM PAIR & HOPPING CLASS' ; write_msgi
write(message,*)' ' ; write_msgi
return
endsubroutine
subroutine load_nn_table(NN_TABLE, PINPT, PPRAM)
use parameters, only : hopping, incar, pid_nntable, params
use mpi_setup
implicit none
type (hopping) :: NN_TABLE
type (incar ) :: PINPT
type (params ) :: PPRAM
integer*4 ii
integer*4 iatom,jatom,mi,mj
integer*4 i_onsite,i_sig,i_pi,i_del
integer*4 n_class, i_sigs,i_pis,i_dels
integer*4 i_lambda,i_stoner,i_localU
integer*4 i_tij, i_soc, i_rashba
real*8 tij
real*8 R(3),D,D0
character*5 ci,cj,ptype
logical flag_soc, flag_slater_koster, flag_local_charge, flag_plus_U, flag_collinear
flag_soc = PINPT%flag_soc
flag_slater_koster = PPRAM%flag_slater_koster
flag_local_charge = PINPT%flag_local_charge
flag_plus_U = PINPT%flag_plus_U
flag_collinear = PINPT%flag_collinear
open(pid_nntable, file=NN_TABLE%nnfilenm, status='old')
read(pid_nntable,*) ! ignore fist line
! IMPORTANT!!!
! In the current version, if USE_OVERALP = .TRUE., the overlap matrix construction does not supported (need to write some routines for reading 'overlap.dat' file as 'hopping.dat')
! This is future work. Need to be updated... (HJ Kim. 2019. Sep.)
! some routines for reading "NN_TABLE%sij_file" is needed.....
do ii=1, NN_TABLE%n_neighbor
if(flag_soc) then
if(flag_slater_koster) then
if(flag_local_charge) then
read(pid_nntable,*)iatom,jatom,R(:),D,D0,mi,ci,mj,cj,ptype,i_onsite, i_sig, i_pi, i_del, i_sigs, i_pis, i_dels, n_class, NN_TABLE%tij_file(ii), i_lambda, i_stoner, i_localU
else
read(pid_nntable,*)iatom,jatom,R(:),D,D0,mi,ci,mj,cj,ptype,i_onsite, i_sig, i_pi, i_del, i_sigs, i_pis, i_dels, n_class, NN_TABLE%tij_file(ii), i_lambda, i_stoner
endif
else
if(flag_local_charge) then
read(pid_nntable,*)iatom,jatom,R(:),D,D0,mi,ci,mj,cj,ptype,i_onsite, i_tij, n_class, NN_TABLE%tij_file(ii), i_soc, i_rashba, i_stoner, i_localU
else
read(pid_nntable,*)iatom,jatom,R(:),D,D0,mi,ci,mj,cj,ptype,i_onsite, i_tij, n_class, NN_TABLE%tij_file(ii), i_soc, i_rashba, i_stoner
endif
endif
elseif(.not. flag_soc) then
if(flag_slater_koster) then
if(flag_collinear) then
if(flag_local_charge) then
read(pid_nntable,*)iatom,jatom,R(:),D,D0,mi,ci,mj,cj,ptype,i_onsite, i_sig, i_pi, i_del, i_sigs, i_pis, i_dels, n_class, NN_TABLE%tij_file(ii), i_stoner, i_localU
else
read(pid_nntable,*)iatom,jatom,R(:),D,D0,mi,ci,mj,cj,ptype,i_onsite, i_sig, i_pi, i_del, i_sigs, i_pis, i_dels, n_class, NN_TABLE%tij_file(ii), i_stoner
endif
else
if(flag_local_charge) then
read(pid_nntable,*)iatom,jatom,R(:),D,D0,mi,ci,mj,cj,ptype,i_onsite, i_sig, i_pi, i_del, i_sigs, i_pis, i_dels, n_class, NN_TABLE%tij_file(ii), i_localU
else
read(pid_nntable,*)iatom,jatom,R(:),D,D0,mi,ci,mj,cj,ptype,i_onsite, i_sig, i_pi, i_del, i_sigs, i_pis, i_dels, n_class, NN_TABLE%tij_file(ii)
endif
endif
else
if(flag_collinear) then
if(flag_local_charge) then
read(pid_nntable,*)iatom,jatom,R(:),D,D0,mi,ci,mj,cj,ptype,i_onsite, i_tij, n_class, NN_TABLE%tij_file(ii), i_stoner, i_localU
else
read(pid_nntable,*)iatom,jatom,R(:),D,D0,mi,ci,mj,cj,ptype,i_onsite, i_tij, n_class, NN_TABLE%tij_file(ii), i_stoner
endif
else
if(flag_local_charge) then
read(pid_nntable,*)iatom,jatom,R(:),D,D0,mi,ci,mj,cj,ptype,i_onsite, i_tij, n_class, NN_TABLE%tij_file(ii), i_localU
else
read(pid_nntable,*)iatom,jatom,R(:),D,D0,mi,ci,mj,cj,ptype,i_onsite, i_tij, n_class, NN_TABLE%tij_file(ii)
endif
endif
endif
endif
enddo
close(pid_nntable)
return
endsubroutine
subroutine print_nn_table(NN_TABLE, PINPT, PPRAM)
use parameters, only : hopping, incar, params, pid_nntable
use mpi_setup
use print_io
implicit none
type (hopping) :: NN_TABLE
type (incar ) :: PINPT
type (params ) :: PPRAM
integer*4 ii, i, i_check, mpierr
logical flag_soc, flag_slater_koster, flag_local_charge, flag_plus_U, flag_collinear
real*8 tij_sk
external tij_sk
flag_soc = PINPT%flag_soc
flag_slater_koster = PPRAM%flag_slater_koster
flag_local_charge = PINPT%flag_local_charge
flag_plus_U = PINPT%flag_plus_U
flag_collinear = PINPT%flag_collinear
if(PINPT%nsystem .gt. 1) then
NN_TABLE%nnfilenmo = 'hopping'//trim(PINPT%title(NN_TABLE%mysystem))//'.dat'
endif
open(pid_nntable, file=trim(NN_TABLE%nnfilenmo), status='unknown')
if(flag_soc) then
if(flag_slater_koster) then
write(pid_nntable,'(A,A)',ADVANCE='no')'# Iatom Jatom RIJ(x, y, z) |RIJ| |RIJ0|(ang)',&
' M_I "ORB_I" M_J "ORB_J" param_type e_o sig pi del sig_s pi_s del_s nn_class t_IJ(eV) lambda_i stoner_i'
elseif(.not. flag_slater_koster) then
write(pid_nntable,'(A,A)',ADVANCE='no')'# Iatom Jatom RIJ(x, y, z) |RIJ| |RIJ0|(ang)',&
' M_I "ORB_I" M_J "ORB_J" param_type e_o t_IJ nn_class t_IJ(eV) lsoc_i lrashba_i stoner_i'
endif
if(flag_local_charge) then
write(pid_nntable,'(A)',ADVANCE='no')' local_U'
endif
if(flag_plus_U) then
write(pid_nntable,'(A)',ADVANCE='no')' plus_U'
endif
write(pid_nntable,'(A)',ADVANCE='yes')' '
else
if(flag_slater_koster) then
write(pid_nntable,'(A,A)',ADVANCE='no')'# Iatom Jatom RIJ(x, y, z) |RIJ| |RIJ0|(ang)',&
' M_I "ORB_I" M_J "ORB_J" param_type e_o sig pi del sig_s pi_s del_s nn_class t_IJ(eV)'
elseif(.not. flag_slater_koster) then
write(pid_nntable,'(A,A)',ADVANCE='no')'# Iatom Jatom RIJ(x, y, z) |RIJ| |RIJ0|(ang)',&
' M_I "ORB_I" M_J "ORB_J" param_type e_o t_IJ nn_class t_IJ(eV)'
endif
if(flag_collinear) then
write(pid_nntable,'(A)',ADVANCE='no')' stoner_i'
endif
if(flag_local_charge) then
write(pid_nntable,'(A)',ADVANCE='no')' local_U'
endif
if(flag_plus_U) then
write(pid_nntable,'(A)',ADVANCE='no')' plus_U'
endif
write(pid_nntable,'(A)',ADVANCE='yes')' '
endif
do ii = 1, NN_TABLE%n_neighbor
if(flag_soc) then
if(flag_slater_koster) then
if(PPRAM%slater_koster_type .gt. 10) NN_TABLE%tij(ii) = tij_sk(NN_TABLE,ii,PPRAM,NN_TABLE%onsite_tolerance,.false.) ! calculate@here if SK_SCALE_TYPE > 10
write(pid_nntable,98,ADVANCE='no')NN_TABLE%i_atom(ii) , NN_TABLE%j_atom(ii), NN_TABLE%Rij(1:3,ii), NN_TABLE%Dij(ii), NN_TABLE%Dij0(ii), &
NN_TABLE%i_matrix(ii), NN_TABLE%ci_orb(ii), NN_TABLE%j_matrix(ii), NN_TABLE%cj_orb(ii),&
NN_TABLE%p_class(ii), NN_TABLE%sk_index_set(0:6,ii), NN_TABLE%n_class(ii), NN_TABLE%tij(ii), &
NN_TABLE%soc_param_index(ii)
elseif(.not. flag_slater_koster) then
write(pid_nntable,96,ADVANCE='no')NN_TABLE%i_atom(ii) , NN_TABLE%j_atom(ii), NN_TABLE%Rij(1:3,ii), NN_TABLE%Dij(ii), NN_TABLE%Dij0(ii), &
NN_TABLE%i_matrix(ii), NN_TABLE%ci_orb(ii), NN_TABLE%j_matrix(ii), NN_TABLE%cj_orb(ii),&
NN_TABLE%p_class(ii), NN_TABLE%cc_index_set(0:1,ii), NN_TABLE%n_class(ii), NN_TABLE%tij(ii), &
NN_TABLE%cc_index_set(2:3,ii)
endif
if( (NN_TABLE%i_matrix(ii) .eq. NN_TABLE%j_matrix(ii)) .and. (NN_TABLE%Dij(ii) .lt. NN_TABLE%onsite_tolerance) ) then
write(pid_nntable,'(8x,I3)',ADVANCE='no')NN_TABLE%stoner_I_param_index(NN_TABLE%i_matrix(ii))
else
write(pid_nntable,'(8x,I3)',ADVANCE='no')0
endif
if(flag_local_charge) then
if( (NN_TABLE%i_matrix(ii) .eq. NN_TABLE%j_matrix(ii)) .and. (NN_TABLE%Dij(ii) .lt. NN_TABLE%onsite_tolerance) ) then
write(pid_nntable,'(8x,I3)',ADVANCE='no')NN_TABLE%local_U_param_index(NN_TABLE%i_matrix(ii))
else
write(pid_nntable,'(8x,I3)',ADVANCE='no')0
endif
endif
if(flag_plus_U) then
if( (NN_TABLE%i_matrix(ii) .eq. NN_TABLE%j_matrix(ii)) .and. (NN_TABLE%Dij(ii) .lt. NN_TABLE%onsite_tolerance) ) then
write(pid_nntable,'(8x,I3)',ADVANCE='no')NN_TABLE%plus_U_param_index(NN_TABLE%i_matrix(ii))
else
write(pid_nntable,'(8x,I3)',ADVANCE='no')0
endif
endif
write(pid_nntable,'(A)',ADVANCE='yes')' '
elseif(.not. flag_soc) then
if(flag_slater_koster) then
if(PPRAM%slater_koster_type .gt. 10) NN_TABLE%tij(ii) = tij_sk(NN_TABLE,ii,PPRAM,NN_TABLE%onsite_tolerance,.false.) ! calculate@here if SK_SCALE_TYPE > 10
write(pid_nntable,99,ADVANCE='no')NN_TABLE%i_atom(ii) , NN_TABLE%j_atom(ii), NN_TABLE%Rij(1:3,ii), NN_TABLE%Dij(ii), NN_TABLE%Dij0(ii), &
NN_TABLE%i_matrix(ii), NN_TABLE%ci_orb(ii), NN_TABLE%j_matrix(ii), NN_TABLE%cj_orb(ii),&
NN_TABLE%p_class(ii), NN_TABLE%sk_index_set(0:6,ii), NN_TABLE%n_class(ii), NN_TABLE%tij(ii)
elseif(.not. flag_slater_koster) then
write(pid_nntable,97,ADVANCE='no')NN_TABLE%i_atom(ii) , NN_TABLE%j_atom(ii), NN_TABLE%Rij(1:3,ii), NN_TABLE%Dij(ii), NN_TABLE%Dij0(ii), &
NN_TABLE%i_matrix(ii), NN_TABLE%ci_orb(ii), NN_TABLE%j_matrix(ii), NN_TABLE%cj_orb(ii),&
NN_TABLE%p_class(ii), NN_TABLE%cc_index_set(0:1,ii), NN_TABLE%n_class(ii), NN_TABLE%tij(ii)
endif
if(flag_collinear) then
if( (NN_TABLE%i_matrix(ii) .eq. NN_TABLE%j_matrix(ii)) .and. (NN_TABLE%Dij(ii) .lt. NN_TABLE%onsite_tolerance) ) then
write(pid_nntable,'(8x,I3)',ADVANCE='no')NN_TABLE%stoner_I_param_index(NN_TABLE%i_matrix(ii))
else
write(pid_nntable,'(8x,I3)',ADVANCE='no')0
endif
endif
if(flag_local_charge) then
if( (NN_TABLE%i_matrix(ii) .eq. NN_TABLE%j_matrix(ii)) .and. (NN_TABLE%Dij(ii) .lt. NN_TABLE%onsite_tolerance) ) then
write(pid_nntable,'(8x,I3)',ADVANCE='no')NN_TABLE%local_U_param_index(NN_TABLE%i_matrix(ii))
else
write(pid_nntable,'(8x,I3)',ADVANCE='no')0
endif
endif
if(flag_plus_U) then
if( (NN_TABLE%i_matrix(ii) .eq. NN_TABLE%j_matrix(ii)) .and. (NN_TABLE%Dij(ii) .lt. NN_TABLE%onsite_tolerance) ) then
write(pid_nntable,'(8x,I3)',ADVANCE='no')NN_TABLE%plus_U_param_index(NN_TABLE%i_matrix(ii))
else
write(pid_nntable,'(8x,I3)',ADVANCE='no')0
endif
endif
write(pid_nntable,'(A)',ADVANCE='yes')' '
endif
i_check = 0
if(flag_slater_koster) then
do i = 0, 6
if(NN_TABLE%sk_index_set(i,ii) .eq. 0) i_check = i_check + 1
enddo
if (i_check .eq. 0) then
write(message,'(A)')' !WARNING! SK-parameter is not set properly! p_class=',NN_TABLE%p_class(ii),' n_class=',NN_TABLE%n_class(ii) ; write_msgi
kill_job
endif
elseif(.not. flag_slater_koster) then
do i = 0, 3
if(NN_TABLE%cc_index_set(i,ii) .eq. 0) i_check = i_check + 1
enddo
if (i_check .eq. 0) then
write(message,'(A)')' !WARNING! CC-parameter is not set properly! p_class=',NN_TABLE%p_class(ii),' n_class=',NN_TABLE%n_class(ii) ; write_msgi
kill_job
endif
endif
enddo
99 format( 1x,I6,I6,3F10.5,2F10.5,I6,3x,A5,I6,3x,A5,6X, A4,2X, 7I5, I6, 3x, F12.5)
97 format( 1x,I6,I6,3F10.5,2F10.5,I6,3x,A5,I6,3x,A5,6X, A4,2X, 2I5, I6, 3x, F12.5)
98 format( 1x,I6,I6,3F10.5,2F10.5,I6,3x,A5,I6,3x,A5,6X, A4,2X, 7I5, I6, 3x, F12.5,3x,I3)
96 format( 1x,I6,I6,3F10.5,2F10.5,I6,3x,A5,I6,3x,A5,6X, A4,2X, 2I5, I6, 3x, F12.5,3x,I3,6x,I3)
close(pid_nntable)
return
endsubroutine