wann_rmat.f 13.3 KB
Newer Older
1 2 3 4 5 6
!--------------------------------------------------------------------------------
! Copyright (c) 2016 Peter Grünberg Institut, Forschungszentrum Jülich, Germany
! This file is part of FLEUR and available as free software under the conditions
! of the MIT license as expressed in the LICENSE file in more detail.
!--------------------------------------------------------------------------------

7
      module m_wann_rmat
Daniel Wortmann's avatar
Daniel Wortmann committed
8
      USE m_juDFT
9 10 11 12 13 14 15 16 17 18 19 20 21
      contains
      subroutine wann_rmat(
     >          bmat,amat, 
     >          rvecnum,rvec,kpoints,
     >          jspins_in,nkpts,l_bzsym,film,l_onedimens,
     >          l_nocosoc,band_min,band_max,neigd,
     >          l_socmmn0,wan90version)
c***********************************************************
c     Calculate the matrix elements of the position operator
c     between wannier functions.
c
c     FF, 2009
c***********************************************************
22
      use m_constants, only:pimach, ImagUnit
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
      use m_wann_read_umatrix

      implicit none
      real, intent(in) :: bmat(:,:)
      real, intent(in) :: amat(:,:)

      integer, intent(in) :: rvecnum
      integer, intent(in) :: rvec(:,:)
      real,    intent(in) :: kpoints(:,:)
      integer, intent(in) :: jspins_in
      integer, intent(in) :: nkpts
      logical, intent(in) :: l_bzsym
      logical, intent(in) :: film
      logical, intent(in) :: l_onedimens

      logical, intent(in) :: l_nocosoc
      integer, intent(in) :: band_min(2),band_max(2),neigd

      logical, intent(in) :: l_socmmn0
      integer, intent(in) :: wan90version

      real,allocatable    :: bpunkt(:,:)
      real                :: kdiff(3)
      real,allocatable    :: wb(:)
      integer,allocatable :: gb(:,:,:),bpt(:,:)
      integer             :: ikpt,jspins,ikpt_help
      integer             :: kpts,nntot,nn
      logical             :: l_file
      integer             :: num_wann,num_kpts,num_nnmax,jspin
      integer             :: kspin,kkspin
      integer             :: num_wann2
      integer             :: i,j,k,m,info,r1,r2,r3,dummy1
      integer             :: dummy2,dummy3,dummy4
      integer             :: hopmin,hopmax,counter,m1,m2
      integer,allocatable :: iwork(:)
      real,allocatable    :: energy(:,:),ei(:)
      real,allocatable    :: eigw(:,:),rwork(:)
      complex,allocatable :: work(:),vec(:,:)
      complex,allocatable :: u_matrix(:,:,:),m_matrix(:,:,:,:)
      complex,allocatable :: posop(:,:,:,:),posop2(:,:,:,:)
      complex             :: fac,eulav,eulav1
      real                :: tmp_omi,rdotk,tpi,minenerg,maxenerg
      real, allocatable   :: minieni(:),maxieni(:)
      character           :: jobz,uplo
      integer             :: kpt,band,lee,lwork,lrwork,liwork,n,lda
      real                :: rvalue
      logical             :: um_format
      logical             :: repro_eig
      logical             :: l_chk,l_proj
      logical             :: have_disentangled
      integer,allocatable :: ndimwin(:)
      logical,allocatable :: lwindow(:,:)
      integer             :: chk_unit,nkp,ntmp,ierr
      character(len=33)   :: header
      character(len=20)   :: checkpoint
      real                :: tmp_latt(3,3), tmp_kpt_latt(3,nkpts)
      real                :: omega_invariant
      complex,allocatable :: u_matrix_opt(:,:,:)
      integer             :: num_bands,nn2
      logical             :: l_umdat
      real,allocatable    :: eigval2(:,:)
      real,allocatable    :: eigval_opt(:,:)
      real                :: scale,a,b
      character(len=2)    :: spinspin12(0:2)
      character(len=3)    :: spin12(2)
      character(len=6)    :: filename
      integer             :: jp,mp,kk
      integer             :: rvecind,rvecind_0
      real                :: realvec(3)
      complex,allocatable :: mmnk(:,:,:,:)
      logical             :: l_worksout,l_she

      data spinspin12/'  ','.1' , '.2'/
      data spin12/'WF1','WF2'/

      tpi=2*pimach()

      jspins=jspins_in
      if(l_nocosoc)jspins=1

      write(6,*)"nkpts=",nkpts

      do jspin=1,jspins  !spin loop
c*****************************************************
c        get num_bands and num_wann from the proj file
c*****************************************************
         do j=jspin,0,-1
          inquire(file=trim('proj'//spinspin12(j)),exist=l_file)
          if(l_file)then
            filename='proj'//spinspin12(j)
            exit
          endif
         enddo
         if(l_file)then
          open (203,file=trim(filename),status='old')
          rewind (203)
         else
Daniel Wortmann's avatar
Daniel Wortmann committed
120
            CALL judft_error("no proj/proj.1/proj.2",calledby
121 122 123 124 125 126 127 128 129 130 131
     +           ="wann_rmat")
         endif
         read (203,*) num_wann,num_bands
         close (203)
         write(6,*)'According to proj there are ',num_bands,' bands'
         write(6,*)"and ",num_wann," wannier functions."

c****************************************************************
c        get nntot and bk and wb from bkpts file
c****************************************************************
         inquire (file='bkpts',exist=l_file)
Daniel Wortmann's avatar
Daniel Wortmann committed
132
         if (.not.l_file)  CALL judft_error("need bkpts"
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
     +      ,calledby ="wann_rmat")
         open (202,file='bkpts',form='formatted',status='old')
         read (202,'(i4)') nntot
         allocate ( gb(3,nntot,nkpts) )
         allocate ( bpt(nntot,nkpts) )
         do ikpt=1,nkpts
          do nn=1,nntot
           read (202,'(2i6,3x,3i4)')
     &      ikpt_help,bpt(nn,ikpt),(gb(i,nn,ikpt),i=1,3)
          enddo !nn
         enddo !ikpt
         allocate( wb(nntot) )
         allocate( bpunkt(3,nntot) )
         do nn=1,nntot
            read(202,*)bpunkt(:,nn),wb(nn)
         enddo
         close (202)

         do ikpt=1,nkpts
          do nn=1,nntot
           kdiff(1)=kpoints(1,bpt(nn,ikpt))+
     +            gb(1,nn,ikpt) - kpoints(1,ikpt)
           kdiff(2)=kpoints(2,bpt(nn,ikpt))+
     +            gb(2,nn,ikpt) - kpoints(2,ikpt)
           kdiff(3)=kpoints(3,bpt(nn,ikpt))+
     +            gb(3,nn,ikpt) - kpoints(3,ikpt)
           kdiff=matmul(transpose(bmat),kdiff)
           l_worksout=.false.
           do nn2=1,nntot
             if (abs(kdiff(1)-bpunkt(1,nn2)).lt.1.e-5)then
               if (abs(kdiff(2)-bpunkt(2,nn2)).lt.1.e-5)then
                 if (abs(kdiff(3)-bpunkt(3,nn2)).lt.1.e-5)then
                     l_worksout=.true. 
                     exit
                 endif          
               endif
             endif 
           enddo !nn2
           if(.not.l_worksout)then
              write(*,*)"ikpt,nn=",ikpt,nn
              write(*,*)"kdiff(1)=",kdiff(1)
              write(*,*)"kdiff(2)=",kdiff(2)
              write(*,*)"kdiff(3)=",kdiff(3)
              stop 'worksout'
           endif   
          enddo !nn
         enddo !ikpt 

c****************************************************************
c        read in chk
c****************************************************************
         num_kpts=nkpts
         allocate( u_matrix_opt(num_bands,num_wann,nkpts) )
         allocate( u_matrix(num_wann,num_wann,nkpts) )
         allocate( lwindow(num_bands,nkpts) )
         allocate( ndimwin(nkpts) )
         allocate( m_matrix(num_wann,num_wann,nntot,num_kpts) )
         call wann_read_umatrix2(
     >            nkpts,num_wann,num_bands,
     >            um_format,jspin,wan90version,
     <            have_disentangled,
     <            lwindow,ndimwin,
     <            u_matrix_opt,u_matrix,m_matrix)
c****************************************************************
c        read in eig-file
c****************************************************************
         write(6,*)"read in eig-file"
         allocate(energy(num_bands,num_kpts))
         inquire(file=spin12(jspin)//'.eig',exist=l_umdat)
Daniel Wortmann's avatar
Daniel Wortmann committed
202
         IF(.NOT.l_umdat)  CALL judft_error
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
     +        ("Thou shall not hide your eig file",calledby
     +        ="wann_hopping")
         open(300,file=spin12(jspin)//'.eig',form='formatted')
         do i=1,num_kpts
           do j=1,num_bands
              read(300,*)band,kpt,energy(j,i)
           enddo
         enddo
         close(300)

         minenerg=minval(energy(:,:))
         maxenerg=maxval(energy(:,:))
         write(6,*)"minenerg=",minenerg
         write(6,*)"maxenerg=",maxenerg


         allocate(eigval_opt(num_bands,nkpts))
         allocate(eigval2(num_wann,nkpts))
         eigval_opt=0.0
         eigval2=0.0

         if(have_disentangled) then

           do nkp=1,num_kpts
            counter=0
            do j=1,num_bands
              if(lwindow(j,nkp)) then
                counter=counter+1
                eigval_opt(counter,nkp)=energy(j,nkp)
              end if
            end do
           end do
       
           do nkp=1,num_kpts
            do j=1,num_wann
             do m=1,ndimwin(nkp)
                eigval2(j,nkp)=eigval2(j,nkp)+eigval_opt(m,nkp)* 
     &    real(conjg(u_matrix_opt(m,j,nkp))*u_matrix_opt(m,j,nkp))
             enddo
            enddo
           enddo

         else
           eigval2 = energy
         end if                    !have_disentangled

         deallocate(eigval_opt)
         deallocate(energy)

c****************************************************************
c        Set up posop.
c****************************************************************
         write(6,*)"Set up posop."
         allocate( posop2(3,num_wann,num_wann,num_kpts) )
         posop2=cmplx(0.0,0.0)
         
         do ikpt=1,nkpts
          do nn=1,nntot
           kdiff(1)=kpoints(1,bpt(nn,ikpt))+
     +            gb(1,nn,ikpt) - kpoints(1,ikpt)
           kdiff(2)=kpoints(2,bpt(nn,ikpt))+
     +            gb(2,nn,ikpt) - kpoints(2,ikpt)
           kdiff(3)=kpoints(3,bpt(nn,ikpt))+
     +            gb(3,nn,ikpt) - kpoints(3,ikpt)
           kdiff=matmul(transpose(bmat),kdiff)          
           do i=1,num_wann
            do j=1,num_wann
             if(j.eq.i)then
              do kk=1,3
272
               posop2(kk,j,i,ikpt)=posop2(kk,j,i,ikpt)+ImagUnit*
273 274 275 276
     &             wb(nn)*kdiff(kk)*(m_matrix(j,i,nn,ikpt)-1.0)  
              enddo !kk
             else
              do kk=1,3
277
               posop2(kk,j,i,ikpt)=posop2(kk,j,i,ikpt)+ImagUnit*
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
     &             wb(nn)*kdiff(kk)*m_matrix(j,i,nn,ikpt)  
              enddo !kk
             endif
            enddo !j
           enddo !i
          enddo !nn   
         enddo !ikpt

         do ikpt=1,nkpts
          do i=1,num_wann
           do j=1,i
            do kk=1,3
              posop2(kk,j,i,ikpt)=(posop2(kk,j,i,ikpt)+
     &                        conjg(posop2(kk,i,j,ikpt)))/2.0
              posop2(kk,i,j,ikpt)=posop2(kk,j,i,ikpt)
            enddo !kk
           enddo !j
          enddo !i
         enddo !ikpt

         allocate( posop(3,num_wann,num_wann,rvecnum) )
         posop=cmplx(0.0,0.0)
         do rvecind=1,rvecnum
            do k=1,nkpts  
              rdotk=tpi*( kpoints(1,k)*rvec(1,rvecind)+
     &                    kpoints(2,k)*rvec(2,rvecind)+
     &                    kpoints(3,k)*rvec(3,rvecind) )
              fac=cmplx(cos(rdotk),-sin(rdotk))
              do i=1,num_wann
               do j=1,num_wann
                do kk=1,3
                 posop(kk,j,i,rvecind)=
     &                   posop(kk,j,i,rvecind)+
     &                   fac*posop2(kk,j,i,k)
                enddo !kk
               enddo !j
              enddo !i
            enddo !k
         enddo !rvecind
         posop=posop/cmplx(real(nkpts),0.0)
         deallocate( posop2 )

         do rvecind=1,rvecnum
          if(rvec(1,rvecind).eq.0)then
           if(rvec(2,rvecind).eq.0)then
            if(rvec(3,rvecind).eq.0)then
               rvecind_0=rvecind
               goto 123
            endif
           endif
          endif
         enddo !rvecind
         stop 'Ou est ce point-la ?'
 123     continue

         do i=1,num_wann
          do kk=1,3 
             posop(kk,i,i,rvecind_0)=0.0
          enddo !kk
         enddo !i

         do ikpt=1,nkpts
          do nn=1,nntot  
           kdiff(1)=kpoints(1,bpt(nn,ikpt))+
     +            gb(1,nn,ikpt) - kpoints(1,ikpt)
           kdiff(2)=kpoints(2,bpt(nn,ikpt))+
     +            gb(2,nn,ikpt) - kpoints(2,ikpt)
           kdiff(3)=kpoints(3,bpt(nn,ikpt))+
     +            gb(3,nn,ikpt) - kpoints(3,ikpt)
           kdiff=matmul(transpose(bmat),kdiff)/real(nkpts)             
           do i=1,num_wann
            do kk=1,3 
             posop(kk,i,i,rvecind_0)=posop(kk,i,i,rvecind_0)-
     &        wb(nn)*kdiff(kk)*aimag(log(m_matrix(i,i,nn,ikpt)))
            enddo !kk
           enddo !i
          enddo !nn
         enddo !ikpt

c********************************************
c        Print posop.
c********************************************
         open(321,file='posop'//spinspin12(jspin),form='formatted')
         do rvecind=1,rvecnum
            r3=rvec(3,rvecind)
            r2=rvec(2,rvecind)
            r1=rvec(1,rvecind)             
            do j=1,num_wann
             do i=1,num_wann
              do kk=1,3  
               write(321,'(i3,1x,i3,1x,i3,1x,i3,
     &                 1x,i3,1x,i3,1x,f20.8,1x,f20.8)')
     &                 r1,r2,r3,i,j,kk,posop(kk,i,j,rvecind) 
              enddo !kk  
             enddo !j
            enddo !i
         enddo !rvecind  
         close(321)

         deallocate( lwindow,u_matrix_opt,ndimwin )
         deallocate( u_matrix,m_matrix )
         deallocate( posop,eigval2 )
         deallocate( gb,bpt,wb,bpunkt )
      enddo !jspin

      end subroutine wann_rmat
      end module m_wann_rmat