wann_socmat.F 4.92 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
     >               lmaxd,natd,neigd,
     >               llod,jspd,
     >               theta_in,phi_in,jspins,irank,
     >               vr,
15 16 17 18 19 20 21 22 23 24
     >               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***********************************************************************
25
      USE m_types
26 27 28
      USE m_spnorb 
      USE m_hsoham
      implicit none
Daniel Wortmann's avatar
Daniel Wortmann committed
29 30 31 32 33 34 35 36

      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

37 38 39 40 41 42 43
      integer, intent(in) :: lmaxd
      integer, intent(in) :: natd
      integer, intent(in) :: neigd

      integer, intent(in) :: llod
      integer, intent(in) :: jspd

44

45 46 47 48 49 50 51 52 53 54 55 56
      real,    intent(in) :: theta_in
      real,    intent(in) :: phi_in
      integer, intent(in) :: jspins
      integer, intent(in) :: irank


      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)

Frank Freimuth's avatar
Frank Freimuth committed
57
      complex, intent(out):: hsomtx(:,:,:,:) !(neigd,neigd,2,2)
58 59 60 61 62 63 64

      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

Frank Freimuth's avatar
Frank Freimuth committed
65 66 67 68 69
!      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
70

Frank Freimuth's avatar
Frank Freimuth committed
71 72 73 74
!      real,allocatable :: ulos(:,:,:)  ! nlod,ntypd,jspd
!      real,allocatable :: dulos(:,:,:) ! nlod,ntypd,jspd
!      real,allocatable :: uulon(:,:,:) ! nlod,ntypd,jspd
!      real,allocatable :: dulon(:,:,:) ! nlod,ntypd,jspd  
75

76 77
      COMPLEX,ALLOCATABLE :: ahelp(:,:,:,:)
      complex,allocatable :: bhelp(:,:,:,:)
78

79
      TYPE(t_rsoc)::rsoc
Daniel Wortmann's avatar
Daniel Wortmann committed
80 81
      TYPE(t_usdus):: usdus

82 83 84 85 86 87 88 89 90 91 92 93 94 95 96
      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

97 98
      ALLOCATE ( ahelp(lmaxd*(lmaxd+2),natd,neigd,input%jspins) )
      ALLOCATE ( bhelp(lmaxd*(lmaxd+2),natd,neigd,input%jspins) )
99

100
      do ispin=1,input%jspins
101 102 103 104 105 106
       DO ie = 1, neigd
        DO na = 1, natd
         DO l = 1, lmaxd
          ll1 = l*(l+1)
          DO m = -l,l
             lm = ll1 + m
107 108
             ahelp(lm,na,ie,ispin) = (acof(ie,lm,na,ispin))
             bhelp(lm,na,ie,ispin) = (bcof(ie,lm,na,ispin))
109 110 111 112 113 114
          ENDDO !m
         ENDDO !l
        ENDDO !na
       ENDDO !ie
      enddo !ispin

Daniel Wortmann's avatar
Daniel Wortmann committed
115 116 117 118 119 120 121 122 123 124
      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))

Frank Freimuth's avatar
Frank Freimuth committed
125 126 127 128 129
!      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))
130 131

      CALL spnorb(
Daniel Wortmann's avatar
Daniel Wortmann committed
132
     >         atoms,noco,input,mpi, enpara,
133 134 135 136
     >         vr,usdus,rsoc,.true.)
      
      rsoc%soangl= conjg(rsoc%soangl)

137 138
      CALL hsoham(atoms,noco,input,nsz,neigd,chelp,rsoc,ahelp,
     >            bhelp,1,natd,mpi%n_rank,mpi%n_size,mpi%SUB_COMM,
139
     <            hsomtx)
140 141 142
 
      end subroutine wann_socmat
      end module m_wann_socmat