wann_kptsrotate.F 5.51 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_kptsrotate
8
      use m_judft
9 10 11 12 13 14 15 16 17 18
c****************************************
c     Rotate the wave function.
c     Frank Freimuth
c****************************************
      contains
      subroutine wann_kptsrotate(
     >               natd,nlod,llod,
     >               ntypd,nlo,llo,invsat,
     >               l_noco,l_soc,
     >               ntype,neq,nlotot,
19
     >               jspin_in,
20 21 22 23
     >               oper,nop,mrot,nvd,
     >               nv,shiftkpt,
     >               tau,
     x               bkpt,k1,k2,k3,
Daniel Wortmann's avatar
Daniel Wortmann committed
24 25
     x               zMat,nsfactor)
      USE m_types
26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49
      use m_inv3
      use m_constants,only:pimach
      implicit none
      integer,intent(in)    :: natd
      integer,intent(in)    :: nlod
      integer,intent(in)    :: llod
      integer,intent(in)    :: ntypd
      INTEGER, INTENT (IN)  :: nlo(ntypd)
      integer, intent(in)   :: llo(nlod,ntypd)
      INTEGER, INTENT (IN)  :: invsat(natd)
      logical,intent(in)    :: l_noco
      logical,intent(in)    :: l_soc
      integer,intent(in)    :: ntype
      integer,intent(in)    :: neq(ntype)
      integer,intent(in)    :: nlotot
      integer,intent(in)    :: jspin_in
      integer,intent(in)    :: oper
      integer,intent(in)    :: nop
      integer,intent(in)    :: mrot(3,3,nop)
      integer,intent(in)    :: nvd,nv(:)
      integer,intent(in)    :: shiftkpt(3)
      real,intent(in)       :: tau(3,nop)
      real,intent(inout)    :: bkpt(3)
      integer,intent(inout) :: k1(:,:),k2(:,:),k3(:,:) !nvd,jspd
Daniel Wortmann's avatar
Daniel Wortmann committed
50 51 52

      TYPE(t_zmat), INTENT (INOUT) :: zMat !z(nbasfcn,noccbd) !can be real/complex

53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76
      complex,intent(out)   :: nsfactor !phase of non-symmorphic ops

      real    :: bkrot(3),tpi,arg
      integer :: rotmat(3,3),determ,ilo
      integer :: j1,j2,j3,k,j,at_ind,invoper
      real    :: shiftnonsymm(3)
      real    :: phase
      INTEGER :: kvec(2*(2*llod+1),nlod,natd)
      integer :: natom,ntyp,lm,m,l,lo,nn,jspin
      integer :: absoper,jj,jsp_start,jsp_end
      integer :: testmat(3,3)

      tpi=2.0*pimach()

      absoper=abs(oper)

      call inv3(mrot(:,:,absoper),rotmat,determ)
      shiftnonsymm(:)=matmul(rotmat,tau(:,absoper))

c      testmat=matmul(mrot(:,:,absoper),rotmat)
c      print*,testmat
c      testmat=matmul(rotmat,mrot(:,:,absoper))
c      print*,testmat

Daniel Wortmann's avatar
Daniel Wortmann committed
77 78 79 80 81 82
      IF(.NOT.zMat%l_real) THEN
         if(oper.lt.0)then
            zMat%z_c = CONJG(zMat%z_c)
            shiftnonsymm=-1.0*shiftnonsymm
         endif
      END IF
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

      if(l_noco) then
         jsp_start=1
         jsp_end=2
      else
         jsp_start=jspin_in
         jsp_end=jspin_in
      endif



      do jspin=jsp_start,jsp_end
          if(.not.(l_noco.and.(jspin.eq.2)))then
            bkrot(:)=0.0
            do k=1,3
              bkrot(:)=bkrot(:)+mrot(k,:,absoper)*bkpt(k)
            enddo
            bkpt(:)=bkrot(:)
          endif
          do j=1,nv(jspin)  !rotate reciprocal vector
               j1=mrot(1,1,absoper)*k1(j,jspin)+
     +            mrot(2,1,absoper)*k2(j,jspin)+
     +            mrot(3,1,absoper)*k3(j,jspin)
               j2=mrot(1,2,absoper)*k1(j,jspin)+
     +            mrot(2,2,absoper)*k2(j,jspin)+
     +            mrot(3,2,absoper)*k3(j,jspin)
               j3=mrot(1,3,absoper)*k1(j,jspin)+
     +            mrot(2,3,absoper)*k2(j,jspin)+
     +            mrot(3,3,absoper)*k3(j,jspin)
               k1(j,jspin)=j1
               k2(j,jspin)=j2
               k3(j,jspin)=j3
          enddo 
      enddo !jspin  

      if(oper.lt.0)then !time-inversion symmetry
         k1   = -k1
         k2   = -k2
         k3   = -k3
         bkpt = -bkpt
      endif

      do jspin=jsp_start,jsp_end
        jj=0 
        if(l_noco.and.(jspin.eq.2))then
           jj=nv(1)+nlotot
        endif
        do j=1,nv(jspin)
         phase = k1(j,jspin) * shiftnonsymm(1)
     +         + k2(j,jspin) * shiftnonsymm(2)
     +         + k3(j,jspin) * shiftnonsymm(3)
         phase = tpi*phase
         phase = cos(phase)
Daniel Wortmann's avatar
Daniel Wortmann committed
136 137 138 139 140
         IF(zMat%l_real) THEN
            zMat%z_r(j+jj,:)  = phase * zMat%z_r(j+jj,:) 
         ELSE
            zMat%z_c(j+jj,:)  = phase * zMat%z_c(j+jj,:) 
         END IF
141 142 143
        enddo    
        jj=jj+nv(jspin)
        do ilo=1,nlotot
144 145
           call judft_error("BUG in wann_kptsrotate:LOs missing")
         !j=kveclo(ilo)
146 147 148 149 150
         phase = k1(j,jspin) * shiftnonsymm(1)
     +         + k2(j,jspin) * shiftnonsymm(2)
     +         + k3(j,jspin) * shiftnonsymm(3)
         phase = tpi*phase
         phase = cos(phase)
Daniel Wortmann's avatar
Daniel Wortmann committed
151 152 153 154 155
         IF(zMat%l_real) THEN
            zMat%z_r(jj+ilo,:)  = phase * zMat%z_r(jj+ilo,:) 
         ELSE
            zMat%z_c(jj+ilo,:)  = phase * zMat%z_c(jj+ilo,:) 
         END IF
156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179
        enddo
      enddo  


      bkpt(:)   = bkpt(:)   - shiftkpt(:)
      k1(:,:)   = k1(:,:)   + shiftkpt(1)
      k2(:,:)   = k2(:,:)   + shiftkpt(2)
      k3(:,:)   = k3(:,:)   + shiftkpt(3)

      write(6,*)"in wann_kptsrotate:"
      write(6,*) "bkpt=",bkpt

      arg = tpi*(
     +        bkpt(1)*shiftnonsymm(1)+
     +        bkpt(2)*shiftnonsymm(2)+
     +        bkpt(3)*shiftnonsymm(3)  )
      
      nsfactor = cmplx(cos(arg),sin(arg))




      end subroutine wann_kptsrotate
      end module m_wann_kptsrotate