~jose-soler/siesta/unfolding

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
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
! ---
! Copyright (C) 1996-2016	The SIESTA group
!  This file is distributed under the terms of the
!  GNU General Public License: see COPYING in the top directory
!  or http://www.gnu.org/copyleft/gpl.txt .
! See Docs/Contributors.txt for a list of contributors.
! ---
module m_tbt_gf

! ----------------------------------------------------------------------
! Module that contains the routines and variables used for to obtain 
! the electrodes surface green's functions and self energies
!
! F.D.Novaes, fdnovaes@icmab.es
! ----------------------------------------------------------------------
! CONTAINS:
! 1) alloc_gf_vars
! 2) cp_gf_vars
! 3) green
! 4) sethhm2


implicit none

save

logical :: LFrstTime=.True., RFrstTime=.True., LJob, RJob, WGFFiles=.false.

complex*16, allocatable, dimension (:) :: LH00, LS00, LH01, LS01
complex*16, allocatable, dimension (:) :: RH00, RS00, RH01, RS01
complex*16, allocatable, dimension(:) :: LGS,RGS
complex*16, allocatable, dimension(:,:) :: LHAA,LSAA,LGAAq,RHAA,RSAA,RGAAq
integer, pointer, dimension (:) :: Llasto,Rlasto
!complex*16, allocatable, dimension(:) ::  Lzbulkdos,Rzbulkdos

! Saved Variables at green4.F

integer :: Lng1, Rng1,NGL2,Lngaa,NGR2,Rngaa,nuaL,nuaR

! Saved Left Variables at sethhm2
integer :: Lnuotot
integer, allocatable, dimension (:), target :: Llisth, Llisthptr, Lnumh, Lindxuo
integer, allocatable, dimension (:), target :: Lix
double precision, allocatable, dimension (:,:), target :: LH,Lxij
double precision, allocatable, dimension (:), target :: LS,Lefs
logical LGamma             ! true if Gamma
real(8), dimension(3,3) ::  Lrcell

! Saved Right Variables at sethhm2
integer :: Rnuotot
integer, allocatable, dimension (:), target :: Rlisth, Rlisthptr, Rnumh, Rindxuo
integer, allocatable, dimension (:), target :: Rix
double precision, allocatable, dimension (:,:), target :: RH, Rxij
double precision, allocatable, dimension (:), target :: RS,Refs
logical RGamma             ! true if Gamma
real(8), dimension(3,3) ::  Rrcell


contains

subroutine alloc_gf_vars(H,S,xij,indxuo,listh,listhptr,numh,efs, ix, &
        notot, nuotot, maxnh, nspin, gamma)

double precision, allocatable, dimension (:,:) :: H, xij
integer, allocatable, dimension (:) :: ix
double precision, allocatable, dimension (:) :: S,efs
integer, allocatable, dimension (:) :: listh, listhptr, numh, indxuo

! Dimesnions

integer :: notot, nuotot, maxnh, nspin
logical :: gamma

allocate(H(maxnh,nspin))
allocate(S(maxnh))
allocate(indxuo(notot))
allocate(listh(maxnh))
allocate(listhptr(nuotot))
allocate(numh(nuotot))
allocate(efs(nspin))
allocate(ix(maxnh)) 
if (.not.gamma) allocate(xij(3,maxnh))


end subroutine alloc_gf_vars

subroutine cp_gf_vars(H,S,xij,indxuo,listh,listhptr,numh,efs,ix,rcell,&
         nuotot, gamma)

use sys, only : die

double precision, dimension (:,:), pointer :: H, xij
double precision, dimension (:), pointer :: S,efs
integer, dimension (:), pointer :: listh, listhptr, numh, indxuo
integer, dimension (:), pointer :: ix
integer :: nuotot
real(8), dimension(3,3) ::  rcell
logical :: gamma


if (LJob) then

   LH=H
   if (.not. gamma) Lxij=xij
   LS=S
   Lefs=efs
   Llisth=listh
   Llisthptr=listhptr
   Lnumh=numh
   Lindxuo=indxuo
   Lix=ix
   Lnuotot=nuotot
   Lrcell=rcell
   Lgamma=gamma


else if (RJob) then
   
   RH=H
   if (.not. gamma) Rxij=xij
   RS=S
   Refs=efs
   Rlisth=listh
   Rlisthptr=listhptr
   Rnumh=numh
   Rindxuo=indxuo
   Rix=ix
   Rnuotot=nuotot
   Rrcell=rcell
   Rgamma=gamma

else

call die("In routine cp_gf_vars: It should be either LJob or RJob !!")

end if

end subroutine cp_gf_vars


! ##################################################################
! ## Driver subroutine for calculating the (ideal)                ##
! ## Left surface Greens function.                                ##          
! ##                            By                                ##
! ##              Mads Brandbyge, mbr@mic.dtu.dk                  ##
! ##################################################################

! FDN Added wq as dummy
      subroutine green(joutfile,nq,q,wq, &
          NEn,contour,wgf,efermi,zbulk,tjob,tbtjob, &
          ispin_fdf, k_fdf, zenergy_fdf,iEn_fdf)
      

! FDN
      use m_tbt_options, only : sppol, Lhsfile, Rhsfile, Lnucuse, &
                  Rnucuse
      use sys, only : die
      use fdf
! FDN


!======================================================================
      implicit none
!======================================================================

! Dimensions for input to green

      logical PRINTALOT
!      parameter(PRINTALOT=.FALSE.)
      parameter(PRINTALOT=.TRUE.)
      
!=======================================================================
! INPUT:
      integer nq                ! Number of surface q-points
! FDN q may have a z component
      real*8 q(3,nq)         ! Surface q-point
! FDN
      character*20 slabel       ! System Label (to name output files)
      character*20 llabel       ! Lead Label (to name output files)
      character*59 sname        ! System Name

      logical tjob !True if Left, False if Right
      logical tbtjob 
      integer joutfile          !unit-number of out-file

      character*33 hsfile        !name of HS-input file
      integer NEn ! no. contour points
      complex*16 contour(NEn),wgf(NEn) !energy contours and weights for GF
      real*8 efermi             ! the Fermi-energy we REQUIRE the electrode
                                ! to have (e.g. when a voltage is applied)
      integer ispin             !spin number for which to generate gf
    
!=======================================================================

!   k_|| == q-points in simple (non-repeated) cell:
      real*8, pointer:: q1(:,:),wq1(:) ! k_|| and their weights 
! FDN
      real*8, dimension(nq) :: wq
! FDN
      real*8, pointer:: q2(:,:) ! shifted q1 points
      integer na1,na2           ! Replication of unitcell
      integer nq1                ! no. q1-points (<= na1*na2 for gamma) 
      real*8 kpoint(3)          !3D k-point (q,kz)
      data kpoint /0d0, 0d0, 0d0/
!------------------------------------------------------------------     
! Hamiltonian/Overlap for given k

      integer nspin
      integer nua !no atom in uc
!      integer lasto(0:NG1) ! Index of last orbital of each atom
      integer, dimension (:), pointer:: lasto
! (we use that the number of atoms < NG1)


! ==================================================================
!     OUTPUT GAA:            


      character*33 gffile        ! name of gffile
      integer jgfu               ! output gf-file unit
!      complex*16 HAA(NGAA,nq1)  ! Real-space bulk Hamiltonian
!      complex*16 SAA(NGAA,nq1)  ! Real-space bulk Overlap
!      complex*16 GAAq(NGAA,nq1) ! Inverse ideal GF (z*SAA - HAA - SigmaAA)

!      complex*16 GS(ngaa1)
      complex*16, allocatable :: GS(:),HAA(:,:),SAA(:,:), GAAq(:,:)

      integer nucuse             !Truncate GF to *last/first* nucuse atoms
!                                !         for   *left/right*

      integer, dimension (:), allocatable :: lastou
!      integer lastou(0:NG1) ! Index of last orbital of each *used* atom
! (we use that the number of atoms < NG2)
      integer nou               !no. used orbitals     

!      complex*16 zbulkdos(NEn)
      complex*16 zbulk 

!=======================================================================
!     Helpers, workspace, tempos etc...

      character*70 gftitle      !title of gf
      character*33 paramfile    !parameter file with this "tag"

      logical tinit,tlast

!     LEFT/RIGHT sign in GF calc.: Exp(i*LRsign*Kz*Z):       
      integer LRSign

      character*5 tag
      character stag
      character*6  gfjob
      character*33 gffile_default

      logical exist1
      logical tdos,tkham
      
      complex*16, dimension (:), pointer:: H00, S00, H01, s01

!      complex*16 h00(ngaa1),s00(ngaa1)
!      complex*16 h01(ngaa1),s01(ngaa1)
!      complex*16 ab(ngaa1),ba(ngaa1)
       complex*16, allocatable :: ab(:),ba(:)


      integer i,l1,l2,ia,ia2
      integer iq1,iq
      integer ju,ng1tmp
      integer ierror

      integer ngaa,ngaa1
      integer NG1      ! Number of basis  orbitals
      integer NG2      ! Number of orbitals used

      real*8  factor
      real*8, allocatable:: eig(:)
      

      complex*16 ZEnergy
      complex*16 ZSEnergy
      complex*16 zdos
      integer iEn

! FDN 
      real*8 cell(3,3)
      integer kscell(3,3),j
      real*8 kdispl(3),tol
      complex*16 zenergy_fdf
      integer ispin_fdf,iEn_fdf
      real*8 k_fdf(3)
! FDN

      external io_assign,io_close

! FDF-stuff:
      character*33 paste,header,itemfdf     
      external paste


!=======================================================================
! BEGIN:
!=======================================================================
       tdos =.true.


       factor =fdf_convfac('Ry','eV')


      if(tjob)then
!!   left
       gfjob='LEFT  '
       LRSign=1
       tag='Left'
       stag='L'
      else
!!   right
       gfjob='RIGHT '
       LRSign=-1
       tag='Right'
       stag='R'
      endif

! output .GF files to this name

! FDN
      if ( LFrstTime .or. RFrstTime ) then 
! FDN

      if (LJob) hsfile = Lhsfile
      if (RJob) hsfile = Rhsfile 

!=======================================================================
!
! Read-in bulk parameters and do calculation:
!
!=======================================================================
        nullify(h00)
        nullify(s00)
        nullify(h01)
        nullify(s01)
        nullify(lasto)

! FDN To be changed !!!
        allocate(h00(1))
        allocate(s00(1))
        allocate(h01(1))
        allocate(s01(1))

! initialize and get nspinin         
        tinit=.true.
        tlast=.false.
        tkham=.false.
! Before read, dimension of H00,s00,h01 and s01 = 0
        ng1=0
        call sethhm2(joutfile,tinit,tkham,tlast,kpoint,ispin, &
          hsfile, nua,lasto,ng1,nspin,cell,kscell,kdispl, &
          H00,s00,h01,s01) ! ->

! FDN To be changed !!!
        deallocate(h00)
        deallocate(s00)
        deallocate(h01)
        deallocate(s01)


        if (LJob) then
          Lng1=ng1

          ngaa1=ng1*ng1
          allocate(LGS(ngaa1))

          itemfdf = paste('TS.NumUsedAtoms',tag) 
!          nucuse = fdf_integer(itemfdf,nua) !default use all atoms in uc.
          if (Lnucuse == 0) then
             nucuse=nua
          else 
             nucuse=Lnucuse
          end if

          nuaL=Lnucuse
     
!            write(*,*) "Used no. atoms:",nucuse,nua
!     Truncate lasto to used atoms: lastou:

          allocate(lastou(0:ng1))
          lastou=0

          lastou(0)=0
          ia2=0   
          if(LRsign .EQ. 1) then !Left
!     use nucuse *last* atoms            
              do ia=nua-nucuse ,nua
                 lastou(ia2)=lasto(ia)
                 ia2=ia2+1
              end do
           else                !Right
!     use nucuse *first* atoms            
              do ia=0,nucuse
                 lastou(ia2)=lasto(ia)
                 ia2=ia2+1               
              end do
           end if              !L or R            
         
!     No. used orbitals in uc:
           nou=0
           do ia=1,nucuse
              nou=nou + (lastou(ia)-lastou(ia-1))
           end do

           NGL2 = nou
           Lngaa=NGL2*NGL2
           deallocate(lastou)

!           do iEn=1,NEn
!              zbulkdos(iEn)=dcmplx(0d0,0d0)
!           end do
           

! FDN Changing the dimensions
           allocate(LHAA(LNGAA,1))  
           allocate(LSAA(LNGAA,1))
           allocate(LGAAq(LNGAA,1))

        end if
        if (RJob) then
          Rng1=ng1


          ngaa1=ng1*ng1
          allocate(RGS(ngaa1))
          itemfdf = paste('TS.NumUsedAtoms',tag) 
!          nucuse = fdf_integer(itemfdf,nua) !default use all atoms in uc.
          if (Rnucuse == 0) then
             nucuse = nua 
          else
             nucuse=Rnucuse
          end if

          nuaR=Rnucuse

!            write(*,*) "Used no. atoms:",nucuse,nua
!     Truncate lasto to used atoms: lastou:

          allocate(lastou(0:ng1))
          lastou=0

          lastou(0)=0
          ia2=0
          if(LRsign .EQ. 1) then !Left
!     use nucuse *last* atoms            
             do ia=nua-nucuse ,nua
                lastou(ia2)=lasto(ia)
                ia2=ia2+1
             end do
          else                !Right
!     use nucuse *first* atoms            
             do ia=0,nucuse
                lastou(ia2)=lasto(ia)
                ia2=ia2+1               
             end do
          end if              !L or R            
         
!     No. used orbitals in uc:
          nou=0
          do ia=1,nucuse
             nou=nou + (lastou(ia)-lastou(ia-1))
          end do

          NGR2 = nou
          Rngaa=NGR2*NGR2       
          deallocate(lastou)

!          do iEn=1,NEn
!             zbulkdos(iEn)=dcmplx(0d0,0d0)
!          end do
! FDN Changing the dimensions
          allocate(RHAA(RNGAA,1))  
          allocate(RSAA(RNGAA,1))
          allocate(RGAAq(RNGAA,1)) 

        end if 


        if(sppol.AND.nspin<2) write(joutfile,*) & 
          "WARNING/ERROR: Spin-polarized electrodes was expected"
        if(.not.sppol .AND. nspin>1)  write(joutfile,*) &
          "WARNING/ERROR: Spin-polarized electrodes was NOT expected"

      end if ! LFirstTime .or. RFirstTime

! FDN Added ....

!      if (.not. WGFFiles ) then
         iEn=iEn_fdf

         if (LJob .and. (.not. LFrstTime)) then

            if ( iEn == 1 ) then


              tinit=.false.
              tkham=.false.
              tlast=.false.
              call sethhm2(joutfile,tinit,tkham,tlast,k_fdf, &
                   ispin_fdf,hsfile,nuaL,Llasto,Lng1,nspin, &
                   cell,kscell,kdispl,LH00,Ls00,Lh01,Ls01) 
            end if ! iEn == 1
            if(iEn.eq.1) then
                  
               i=0
               do l2 = 1,NGL2
                  do l1 = 1,NGL2
                     i=i+1
! FDN iq1 ==> 1
          LHAA(i,1) = LH00(l1+(LNG1-NGL2)+LNG1*(l2+(LNG1-NGL2)-1))
          LSAA(i,1) = LS00(l1+(LNG1-NGL2)+LNG1*(l2+(LNG1-NGL2)-1))
! FDN
                  end do  ! l1
               end do     ! l2
              
               NGAA=NGL2*NGL2
               call zaxpy(NGAA,dcmplx(efermi,0.d0), &
                      LSAA(1,1),1,LHAA(1,1),1)

            end if        ! iEn.eq.1

            zsenergy = zenergy_fdf-efermi
            call calc_green(Lng1,zsenergy,Lh00,Ls00,Lh01,Ls01, &
                   LGS,zdos,joutfile,tjob,tdos)
! A Fazer .... colocar wk ...           
            if(tdos) zbulk =  zdos


            i=0
                  do l2 = 1,NGL2
                     do l1 = 1,NGL2
                        i=i+1
! FDN iq1 ==> 1
           LGAAq(i,1) = Lgs(l1+(LNG1-NGL2)+ LNG1*(l2+(LNG1-NGL2)-1) )
! FDN
                     end do     ! l1
                  end do        ! l2

         end if ! LJob


         if (RJob .and. .not.RFrstTime) then

            if ( iEn == 1 ) then



              tinit=.false.
              tkham=.false.
              tlast=.false.
              call sethhm2(joutfile,tinit,tkham,tlast,k_fdf, &
                   ispin_fdf,hsfile,nuaR,Rlasto,Rng1,nspin, &
                   cell,kscell,kdispl,RH00,Rs00,Rh01,Rs01) 
            end if ! iEn
 
            if(iEn.eq.1) then
                 
                     i=0
                     do l2 = 1,NGR2
                        do l1 = 1,NGR2
                           i=i+1
! FDN iq1 ==> 1
               RHAA(i,1) = RH00(l1+RNG1*(l2-1))
               RSAA(i,1) = RS00(l1+RNG1*(l2-1))
! FDN
                        end do  ! l1
                     end do     ! l2                  
         
!
!     Shift so efermi is energy zero
!
! FDN iq1 ==> 1
                  NGAA=NGR2*NGR2
                  call zaxpy(NGAA,dcmplx(efermi,0.d0), &
                      RSAA(1,1),1,RHAA(1,1),1)
! FDN

             end if ! iEn.eq.1 

             zsenergy = zenergy_fdf-efermi


             call calc_green(Rng1,zsenergy,Rh00,Rs00,Rh01,Rs01, &
                   RGS,zdos,joutfile,tjob,tdos) 
! A Fazer .... wk
            if(tdos) zbulk =   zdos
             i=0
                  do l2 = 1,NGR2
                     do l1 = 1,NGR2
                        i=i+1
! FDN iq1 ==> 1
                        RGAAq(i,1) = Rgs(l1+RNG1*(l2-1))
! FDN
                     end do     ! l1
                  end do        ! l2           

         end if !RJob

!      end if ! .not. WGFFiles

! FDN
      if (LFrstTime .and. tjob ) then
        LFrstTime=.false.     
! Alcocar LH00, LS00, etc ..  
      end if
      if (RFrstTime .and. (.not.tjob) ) then
        RFrstTime=.false.
! Alocar RH00, RS00, etc ...
      end if 
! FDN     


! ======================================================================
      return
      end subroutine green
! ======================================================================




      subroutine sethhm2(joutfile,tinit,tkham,tlast,kpoint,ispin, &
          hsfile,nua,lasto,nuo,nspin,cell,kscell,kdispl, &
          Hk,Sk,Hk2,Sk2)    ! ->

! ##################################################################
! ##               Setup Hamiltonian in k-space                   ##
! ##                            By                                ##
! ##              Mads Brandbyge, mbr@mic.dtu.dk                  ##
! ##################################################################
!
   

      use tsread2
! FDN
      use parallel, only : IOnode
      use m_tbt_kpts, only : reclat
! FDN
!-----------------------------------------------------------------------
      implicit none
!-----------------------------------------------------------------------

!=======================================================================
      real*8 EPS
      parameter(EPS=1.0d-8)
!=======================================================================


! INPUT

      integer joutfile          !info-out file unit
      logical tkham             ! true if full H_k hamiltonian is generated

!      integer nuo               !No. states in unitcell (expected in read-in)
      real*8 kpoint(3) 
      integer ispin             !set up H for ispin
      character*33 hsfile       !H,S parameter file


!-----------------------------------------------------------------------
! READ-IN from HS-file/OUTPUT


      integer nua               ! No. atoms in unitcell
      integer mo             ! Number of orbitals in supercell
      integer nuo            ! Number of basis  orbitals
      integer mno             ! Number of orbitals interacting
      integer nspin          ! Spin polarization (1 or 2)

!      integer numh(maxuo)        ! Number of nonzero elements of each row
!                               ! of hamiltonian matrix
!      integer listh(maxno,maxuo) ! Nonzero hamiltonian-matrix element column
!                               ! indexes for each matrix row
!      integer indxuo(maxo)      ! Index of equivalent orbital in unit cell
!                               !  Unit cell orbitals must be the first in
!                               ! orbital lists, i.e. indxuo.le.nuo, with
!                               ! nuo the number of orbitals in unit cell  

!      real*8  H(maxno,maxuo,nspin) ! Hamiltonian in sparse form
!      real*8  S(maxno,maxuo)     ! Overlap in sparse form
      real*8  qtot              ! Total number of electrons
      real*8  temp              ! Electronic temperature for Fermi smearing

      logical Gamma             ! true if Gamma

!      real*8  xij(3,maxno,maxuo) ! Vectors between orbital centers (sparse)
!                               (not read/written if only gamma point ..
!                                but must be written here!!)

      real*8 cell(3,3)          ! unit cell
! FDN
      integer kscell(3,3)
      real*8 kdispl(3)
      
! FDN

!      integer nua               ! No. atoms in unitcell
!      integer isa(maxua)         ! Species index of each atom
!      real*8 xa(3,maxua)         ! Atomic coordinates (Bohr)

    

      integer, dimension (:), pointer:: listh, listhptr, &
                             numh,indxuo,isa,lasto
      double precision, dimension (:,:), pointer:: H,xij,xa
      double precision, dimension (:), pointer:: S,efs


!-----------------------------------------------------------------------
! OUTPUT
      complex*16 Hk(:), Sk(:)
      complex*16 Hk2(:), Sk2(:)


!      integer, dimension (:), pointer:: lasto
!      integer lasto(0:maxua)   ! Index of last orbital of each atom
!-----------------------------------------------------------------------
! Helpers
      integer nb           ! No. basis orbitals on atoms (the same on all..)
!      real*8 xo(3,maxuo)        ! Atomic coordinates (Bohr)
      real*8, allocatable ::   xo(:,:) 

      integer nuotot,notot,maxnh
      integer ia
      integer i,j,io,jo,iuo,juo,j2
      real*8 k(3),kxij,rcell(3,3)
      real*8 recell(3,3)
      complex*16 cphase

!      integer ix(maxnh)
      integer, dimension (:), pointer :: ix
      integer icoi,icoa
      integer iprop,inn,it,in,ind
      real*8   xc      
      logical tinit,tlast
!-----------------------------------------------------------------------
! SAVED:
!      save numh,listh,indxuo,Gamma
!      save H,S,efs,rcell,xij,listhptr,nuotot
!      save ix

!=======================================================================
! BEGIN
      
       iprop = 3
       icoa = 0
       icoi = 0

       if (.not. LFrstTime .and. LJob) then


 
        H=>LH
        if (.not. Lgamma) xij=>Lxij
        S=>LS
        efs=>Lefs
        listh=>Llisth
        listhptr=>Llisthptr
        numh=>Lnumh
        indxuo=>Lindxuo
!        allocate(ix(size(Lix)))
        ix=>Lix
        nuotot=Lnuotot
        rcell=Lrcell
        gamma=Lgamma


       end if
       if (.not. RFrstTime .and. RJob) then


 
        H=>RH
        if (.not. Rgamma) xij=>Rxij
        S=>RS
        efs=>Refs
        listh=>Rlisth
        listhptr=>Rlisthptr
        numh=>Rnumh
        indxuo=>Rindxuo
!        allocate(ix(size(Rix)))
        ix=>Rix
        nuotot=Rnuotot
        rcell=Rrcell
        gamma=Rgamma

       end if

!-----------------------------------------------------------------
!     Read-in Hamiltonian/Overlap parameters from HS-lattice-file
!     Only read-in firsttime:
!
            
       if(tinit) then
!-----------------------------------------------------------------         
              nullify(H)
              nullify(S)
              nullify(xij)
              nullify(indxuo)
              nullify(listh)
              nullify(listhptr)
              nullify(numh)
              nullify(efs)
              nullify(xa)


! FDN kscell and kdispl added as dummys

      call TSiohs('read', &
      hsfile, gamma, nua, nuotot,notot,nspin, &
                   maxnh,numh, listhptr, listh, H, S, qtot, temp, &
                   xij, indxuo, efs, cell, isa, lasto, xa, &
                   kscell, kdispl)

! FDN

           nuo = nuotot


       
         if (LFrstTime .and. LJob) then
            allocate(LH00(nuo*nuo))
            allocate(LS00(nuo*nuo))
            allocate(LH01(nuo*nuo))
            allocate(LS01(nuo*nuo)) 
! Alocar matrizes que sao salvas
            call alloc_gf_vars(LH,LS,Lxij,Lindxuo,Llisth,Llisthptr,Lnumh &
      ,Lefs,Lix,notot,nuotot,maxnh,nspin,gamma)
         end if         

         if (RFrstTime .and. RJob) then
            allocate(RH00(nuo*nuo))
            allocate(RS00(nuo*nuo))
            allocate(RH01(nuo*nuo))
            allocate(RS01(nuo*nuo))
! Alocar matrizes que sao salvas 
            call alloc_gf_vars(RH,RS,Rxij,Rindxuo,Rlisth,Rlisthptr,Rnumh &
      ,Refs,Rix,notot,nuotot,maxnh,nspin,gamma)     
         end if
        

         if (IOnode) then
            write(joutfile,*) 'unit cell:'

            do j=1,3
               write(joutfile,'(3F8.4)') (cell(i,j),i=1,3)
            end do
         end if ! IOnode 
         
         call reclat(cell,rcell,1) !reciprocal of cell incl. 2Pi!

         call reclat(cell,recell,0)

         if(.not. Gamma) then

            
          allocate(xo(3,nuo))
          
          if (associated(ix)) nullify(ix)
          
        
          allocate (ix(maxnh))
        
!
!     Transform xij so there is no k-dep. phase within uc.
!

!     ... but first find orbital coordinate

            do ia=1,nua
               do iuo=lasto(ia-1)+1,lasto(ia)
                  xo(1,iuo)=xa(1,ia)
                  xo(2,iuo)=xa(2,ia)
                  xo(3,iuo)=xa(3,ia)
               end do           !iuo
            end do              !ia in uc


           do iuo = 1,nuo
            do j = 1,numh(iuo)
              ind = listhptr(iuo) + j
              jo = listh(ind)
              juo = indxuo(jo)
              xij(1,ind)=xij(1,ind)-(xo(1,juo)-xo(1,iuo))
              xij(2,ind)=xij(2,ind)-(xo(2,juo)-xo(2,iuo))
              xij(3,ind)=xij(3,ind)-(xo(3,juo)-xo(3,iuo))
              xc =0.0
              do  i = 1,3
                xc = xc + xij(i,ind)*recell(i,iprop)
              end do
              ix(ind) = nint(xc)
              icoa = max(icoa,nint(xc))
              icoi = min(icoi,nint(xc))
            enddo
            enddo

             deallocate( xo )
         end if                 ! not Gamma         

         
         deallocate( isa )
         deallocate( xa )
            
        if (LFrstTime .and. LJob ) then
            call cp_gf_vars(H,S,xij,indxuo,listh,listhptr,numh,efs,ix, &
      rcell,nuotot, gamma)
!          LFrstTime=.false.
        end if         

        if (RFrstTime .and. RJob ) then
            call cp_gf_vars(H,S,xij,indxuo,listh,listhptr,numh,efs,ix, &
      rcell,nuotot, gamma)
!         RFrstTime=.false.          
        end if


         tinit = .false.
         return
!-----------------------------------------------------------------
         end if                    !tinit
!-----------------------------------------------------------------
         if(tlast) then

            return
         endif
         
   


      k=kpoint
!-----------------------------------------------------------------
      if (tkham)then
!-----------------------------------------------------------------

!
! Setup H,S for this k-point:
!
      do inn = 1,nuo*nuo
        Hk(inn) = dcmplx(0.d0,0.d0)
        Sk(inn) = dcmplx(0.d0,0.d0)
      enddo

      if(.not.Gamma) then



          do iuo = 1,nuo
            do j = 1,numh(iuo)
              ind = listhptr(iuo) + j
              jo = listh(ind)
              juo = indxuo(jo)
              kxij = (k(1) * xij(1,ind) + &
                     k(2) * xij(2,ind) + &
                     k(3) * xij(3,ind) )
              cphase = cdexp(dcmplx(0d0,1d0)*kxij)
              inn = iuo+(juo-1)*nuo
              Hk(inn) = Hk(inn)+H(ind,ispin)*cphase
              Sk(inn) = Sk(inn)+S(ind)*cphase
            enddo
          enddo


      else !Gamma!

        do io = 1,nuo
          do j = 1,numh(io)
            ind = listhptr(io) + j
            jo = listh(ind)
            inn = io+(jo-1)*nuo
            Hk(inn) = Hk(inn)+H(ind,ispin)
            Sk(inn) = Sk(inn)+S(ind)
          enddo
        enddo

      end if                    !Gamma or not

!
!     Symmetrize and *make EF the energy-zero*!!!
!
      do iuo = 1,nuo
         do juo = 1,iuo-1
           it = juo+(iuo-1)*nuo
           in = iuo+(juo-1)*nuo

           Sk(it) = 0.5d0*( Sk(it) + dconjg(Sk(in)) )
           Sk(in) =  dconjg(Sk(it))

           Hk(it) = 0.5d0*( Hk(it) + dconjg(Hk(in)) ) &
                - efs(ispin)*Sk(it)
           Hk(in) =  dconjg(Hk(it))

        enddo

         in = iuo+(iuo-1)*nuo
         Sk(in) = Sk(in) - dcmplx(0d0,1d0)*dimag(Sk(in))
         Hk(in) = Hk(in) - dcmplx(0d0,1d0)*dimag(Hk(in)) &
                - efs(ispin)*Sk(in)
      enddo


!-----------------------------------------------------------------
      else                      !tkham
!-----------------------------------------------------------------

! Setup Transfer H,S for this k||-point:
!
          if(gamma)then
           write(6,*) 'Transfer matrix not possible with gamma'
           stop 'Transfer matrix not possible with gamma'
          endif


         do inn = 1,nuo*nuo
           Hk(inn) = dcmplx(0.d0,0.d0)
           Sk(inn) = dcmplx(0.d0,0.d0)
           Hk2(inn) = dcmplx(0.d0,0.d0)
           Sk2(inn) = dcmplx(0.d0,0.d0)
         enddo

          do iuo = 1,nuo
            do j = 1,numh(iuo)
              ind = listhptr(iuo) + j
              jo = listh(ind)
              juo = indxuo(jo)
               kxij = (k(1) * xij(1,ind) + &
                      k(2) * xij(2,ind) + &
                      k(3) * xij(3,ind) - &
                      k(iprop) *xij(iprop,ind) )
              cphase = cdexp(dcmplx(0d0,1d0)*kxij)

              inn = iuo+(juo-1)*nuo

              if(ix(ind).eq.0) then
                 Hk(inn) = Hk(inn)+H(ind,ispin)*cphase
                 Sk(inn) = Sk(inn)+S(ind)*cphase    
               else if(ix(ind).eq.1) then
                 Hk2(inn) = Hk2(inn)+H(ind,ispin)*cphase
                 Sk2(inn) = Sk2(inn)+S(ind)*cphase
               endif

            enddo
          enddo

!
!     Symmetrize and *make EF the energy-zero*!!!
!
      do iuo = 1,nuo
         do juo = 1,iuo-1
           it = juo+(iuo-1)*nuo
           in = iuo+(juo-1)*nuo

           Sk(it) = 0.5d0*( Sk(it) + dconjg(Sk(in)) )
           Sk(in) =  dconjg(Sk(it))

           Hk(it) = 0.5d0*( Hk(it) + dconjg(Hk(in)) ) &
                - efs(ispin)*Sk(it)
           Hk(in) =  dconjg(Hk(it))

        enddo

         in = iuo+(iuo-1)*nuo
         Sk(in)=Sk(in) - dcmplx(0d0,1d0)*dimag(Sk(in))
         Hk(in)=Hk(in) - dcmplx(0d0,1d0)*dimag(Hk(in)) &
                - efs(ispin)*Sk(in)

      enddo
      

      do iuo = 1,nuo
         do juo = 1,nuo
          in = iuo+(juo-1)*nuo
             Hk2(in)=Hk2(in) - efs(ispin)*Sk2(in)
         enddo
      enddo    
        
!-----------------------------------------------------------------
       endif                      !tkham
!-----------------------------------------------------------------
       
! ===============================================================
       return
       END subroutine sethhm2
! ===============================================================














end module m_tbt_gf