wann_socmat.F 7.48 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 8 9
      module m_wann_socmat
      contains
      subroutine wann_socmat(
Daniel Wortmann's avatar
Daniel Wortmann committed
10
     >               mpi,enpara,input,noco,atoms,
11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27
     >               lmaxd,ntypd,nlod,natd,neigd,
     >               llod,jmtd,jspd,nlhd,soc_opt,neq,
     >               ntype,theta_in,phi_in,jspins,irank,
     >               jri,lmax,dx,rmsh,el,ello,nlo,llo,
     >               l_dulo,ulo_der,vr,
     >               acof,bcof,chelp,
     <               hsomtx)
c***********************************************************************
c     Calculate matrix elements of the spin-orbit interaction for those
c     Bloch functions out of which the Wannier functions are 
c     constructed. From these matrix elements the spin-orbit Hamiltonian
c     in the basis set of Wannier functions may be constructed.
c     
c     Frank Freimuth
c***********************************************************************
      USE m_spnorb 
      USE m_hsoham
Daniel Wortmann's avatar
Daniel Wortmann committed
28
      USE m_types
29
      implicit none
Daniel Wortmann's avatar
Daniel Wortmann committed
30 31 32 33 34 35 36 37

      TYPE(t_mpi),INTENT(IN)      :: mpi

      TYPE(t_enpara),INTENT(IN)   :: enpara
      TYPE(t_input),INTENT(IN)    :: input
      TYPE(t_noco),INTENT(IN)     :: noco
      TYPE(t_atoms),INTENT(IN)    :: atoms

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
      integer, intent(in) :: lmaxd
      integer, intent(in) :: ntypd
      integer, intent(in) :: nlod
      integer, intent(in) :: natd
      integer, intent(in) :: neigd

      integer, intent(in) :: llod
      integer, intent(in) :: jmtd
      integer, intent(in) :: jspd
      integer, intent(in) :: nlhd
      LOGICAL, INTENT (IN):: soc_opt(:) !(ntypd+2)
      integer, intent(in) :: neq(:) !ntype

      integer, intent(in) :: ntype
      real,    intent(in) :: theta_in
      real,    intent(in) :: phi_in
      integer, intent(in) :: jspins
      integer, intent(in) :: irank

      integer, intent(in) :: jri(:)  ! ntypd
      integer, intent(in) :: lmax(:) ! ntypd
      real,    intent(in) :: dx(:)   ! ntypd
      real,    intent(in) :: rmsh(:,:)  ! jmtd,ntypd
      real,    intent(in) :: el(0:,:,:) ! 0:lmaxd,ntypd,max(2,jspd)
      real,    intent(in) :: ello(:,:,:)! nlod,ntypd,max(2,jspd)
      integer, intent(in) :: nlo(:)     ! ntypd
      integer, intent(in) :: llo(:,:)   ! nlod,ntypd

      logical, intent(in) :: l_dulo(:,:)  ! l_dulo(nlod,ntypd)
      integer, intent(in) :: ulo_der(:,:) ! ulo_der(nlod,ntypd)
      real,    intent(in) :: vr(:,0:,:,:)    ! jmtd,0:nlhd,ntypd,jspd   
  
      complex, intent(in) :: acof(:,0:,:,:) !acof(noccbd,0:lmd,natd,jspd)
      complex, intent(in) :: bcof(:,0:,:,:) !bcof(noccbd,0:lmd,natd,jspd)
      complex, intent(in) :: chelp(-llod:,:,:,:,:) !chelp(-llod:llod,neigd,nlod,natd,jspd)

      complex, intent(out):: hsomtx(:,:,:,:) !(2,2,neigd,neigd)

      integer :: n,l,nwdd,nw,ispin,ie,na,ll1,m,lm,i,nsz(2)
      real    :: s(3),r2
      logical :: l_all
      CHARACTER*3 chntype
      real    :: theta,phi,pi

      real,allocatable :: ddn(:,:,:) ! 0:lmaxd,ntypd,jspd
      real,allocatable :: us(:,:,:)  ! 0:lmaxd,ntypd,jspd
      real,allocatable :: dus(:,:,:) ! 0:lmaxd,ntypd,jspd
      real,allocatable :: uds(:,:,:) ! 0:lmaxd,ntypd,jspd
      real,allocatable :: duds(:,:,:)! 0:lmaxd,ntypd,jspd

      real,allocatable :: ulos(:,:,:)  ! nlod,ntypd,jspd
      real,allocatable :: dulos(:,:,:) ! nlod,ntypd,jspd
      real,allocatable :: uulon(:,:,:) ! nlod,ntypd,jspd
      real,allocatable :: dulon(:,:,:) ! nlod,ntypd,jspd  

      COMPLEX,ALLOCATABLE :: ahelp(:,:,:,:,:)
      complex,allocatable :: bhelp(:,:,:,:,:)

      REAL,    ALLOCATABLE :: rsopdp(:,:,:,:),rsopdpd(:,:,:,:)
      REAL,    ALLOCATABLE :: rsopp(:,:,:,:),rsoppd(:,:,:,:)
      REAL,    ALLOCATABLE :: rsoplop(:,:,:,:)
      REAL,    ALLOCATABLE :: rsoplopd(:,:,:,:),rsopdplo(:,:,:,:)
      REAL,    ALLOCATABLE :: rsopplo(:,:,:,:),rsoploplop(:,:,:,:,:)
      COMPLEX, ALLOCATABLE :: soangl(:,:,:,:,:,:)

Daniel Wortmann's avatar
Daniel Wortmann committed
103 104
      TYPE(t_usdus):: usdus

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
      nwdd=1
      nw=1
      nsz=neigd

      IF (.true.) THEN
        theta= theta_In
        phi= phi_In
      ELSE
        pi= 4.*ATAN(1.)
        theta= -theta_In
        phi=   phi_In+pi
        ! now the definition of rotation matrices
        ! is equivalent to the def in the noco-routines
      ENDIF

      ALLOCATE ( ahelp(-lmaxd:lmaxd,lmaxd,natd,neigd,jspd) )
      ALLOCATE ( bhelp(-lmaxd:lmaxd,lmaxd,natd,neigd,jspd) )

      do ispin=1,jspd
       DO ie = 1, neigd
        DO na = 1, natd
         DO l = 1, lmaxd
          ll1 = l*(l+1)
          DO m = -l,l
             lm = ll1 + m
             ahelp(m,l,na,ie,ispin) = (acof(ie,lm,na,ispin))
             bhelp(m,l,na,ie,ispin) = (bcof(ie,lm,na,ispin))
          ENDDO !m
         ENDDO !l
        ENDDO !na
       ENDDO !ie
      enddo !ispin

      ALLOCATE( rsopdp(ntypd,lmaxd,2,2),rsopdpd(ntypd,lmaxd,2,2),
     +          rsopp(ntypd,lmaxd,2,2),rsoppd(ntypd,lmaxd,2,2),
     +          rsoplop(ntypd,nlod,2,2),rsoplopd(ntypd,nlod,2,2),
     +          rsopdplo(ntypd,nlod,2,2),rsopplo(ntypd,nlod,2,2),
     +          rsoploplop(ntypd,nlod,nlod,2,2),
     +          soangl(lmaxd,-lmaxd:lmaxd,2,lmaxd,-lmaxd:lmaxd,2) )

Daniel Wortmann's avatar
Daniel Wortmann committed
145 146 147 148 149 150 151 152 153 154
      ALLOCATE(usdus%us(0:atoms%lmaxd,atoms%ntype,jspd), 
     +         usdus%dus(0:atoms%lmaxd,atoms%ntype,jspd),
     +         usdus%uds(0:atoms%lmaxd,atoms%ntype,jspd),
     +         usdus%duds(0:atoms%lmaxd,atoms%ntype,jspd),
     +         usdus%ddn(0:atoms%lmaxd,atoms%ntype,jspd),
     +         usdus%ulos(atoms%nlod,atoms%ntype,jspd),
     +         usdus%dulos(atoms%nlod,atoms%ntype,jspd),
     +         usdus%uulon(atoms%nlod,atoms%ntype,jspd),
     +         usdus%dulon(atoms%nlod,atoms%ntype,jspd))

155 156 157 158 159 160 161 162
      ALLOCATE( us(0:lmaxd,ntypd,jspd), dus(0:lmaxd,ntypd,jspd),
     +          uds(0:lmaxd,ntypd,jspd),duds(0:lmaxd,ntypd,jspd),
     +          ddn(0:lmaxd,ntypd,jspd),
     +          ulos(nlod,ntypd,jspd),dulos(nlod,ntypd,jspd),
     +          uulon(nlod,ntypd,jspd),dulon(nlod,ntypd,jspd))

      soangl(:,:,:,:,:,:) = cmplx(0.0,0.0)
      CALL spnorb(
Daniel Wortmann's avatar
Daniel Wortmann committed
163 164 165
     >         atoms,noco,input,mpi, enpara,
     >         vr,rsopp,rsoppd,rsopdp,rsopdpd,usdus,
     <         rsoplop,rsoplopd,rsopdplo,rsopplo,rsoploplop,soangl)
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

      l_all=.false.
      IF (irank==0) THEN
           IF ( soc_opt(ntype+1).or.l_all ) THEN
             IF (l_all) THEN
               WRITE (6,fmt='(A)') 'Only SOC contribution of certain'
     &          //' atom types included in Hamiltonian.'
             ELSE 
               WRITE (chntype,'(i3)') ntype
               WRITE (6,fmt='(A,2x,'//chntype//'l1)') 'SOC contributi'
     &          //'on of certain atom types included in Hamiltonian:',
     &          (soc_opt(n),n=1,ntype)
             ENDIF
           ELSE
             WRITE(6,fmt='(A,1x,A)') 'SOC contribution of all atom'//
     &        ' types inculded in Hamiltonian.'
           ENDIF
           IF (soc_opt(ntype+2)) THEN
             WRITE(6,fmt='(A)')
     &        'SOC Hamiltonian is constructed by neglecting B_xc.'
           ENDIF
      ENDIF
      soangl(:,:,:,:,:,:) = conjg(soangl(:,:,:,:,:,:))

      CALL hsoham(
Daniel Wortmann's avatar
Daniel Wortmann committed
191 192
     >    atoms,noco,input,
     >    nsz,chelp,
193 194 195 196 197 198
     >    rsoplop,rsoplopd,rsopdplo,rsopplo,rsoploplop,
     >    ahelp,bhelp,rsopp,rsoppd,rsopdp,rsopdpd,soangl,
     <    hsomtx)
 
      end subroutine wann_socmat
      end module m_wann_socmat