wann_uHu_soc_tlo.F 16.7 KB
Newer Older
Daniel Wortmann's avatar
Daniel Wortmann committed
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
c*****************************************c
c    Spin-orbit coupling cont. to uHu     c
c  < u_{k+b1} | V_{k}^{soc} | u_{k+b2} >  c
c*****************************************c
c    Routine to set up T(lm,lmp) for all  c
c    pairs (b1,b2) and every atom n       c
c    due to local orbitals                c
c*****************************************c
c                J.-P. Hanke, Dec. 2015   c
c*****************************************c
      MODULE m_wann_uHu_soc_tlo
      CONTAINS
      SUBROUTINE wann_uHu_soc_tlo(
     >                  input,atoms,
     >                  ntypd,jmtd,lmaxd,jspd,
     >                  ntype,dx,rmsh,jri,lmax,natd,
     >                  lmd,irank,nlod,llod,
     >                  loplod,ello,llo,nlo,lo1l,l_dulo,ulo_der,
     >                  flo,flo_b,kdiff,kdiff2,nntot,nntot2,
     >                  vr,epar,jspin,jspin_b,jspins,l_spav,theta,phi,
     >                  yl1,yl2,jj1,jj2,p,p_b,q,q_b,ntyp,angso,
     >                  l_socscreen,socscreen_fac,
     <                  tuulo,tdulo,tulou,tulod,tuloulo)

      USE m_intgr, ONLY : intgr0
      USE m_sphbes
      USE m_ylm
      USE m_gaunt, ONLY: gaunt1
      USE m_sointg
      USE m_constants
      USE m_types
      

      IMPLICIT NONE
C     .. Intrinsic Functions ..
      INTRINSIC abs,cmplx,max,mod

      TYPE(t_input),INTENT(IN)   :: input
      TYPE(t_atoms),INTENT(IN)   :: atoms

C     ..
C     .. Scalar Arguments ..
      INTEGER, INTENT (IN) :: ntypd,jmtd,lmaxd,jspd,jspins
      INTEGER, INTENT (IN) :: lmd,ntype,irank,jspin,jspin_b
      INTEGER, INTENT (IN) :: nlod,llod,loplod,natd
      INTEGER, INTENT (IN) :: nntot,nntot2,ntyp
      REAL,    INTENT (IN) :: theta,phi,socscreen_fac
      LOGICAL, INTENT (IN) :: l_socscreen
C     ..
C     .. Array Arguments ..
      COMPLEX, INTENT (INOUT)::
     >        tdulo(0:lmd,nlod,-llod:llod),
     >        tuulo(0:lmd,nlod,-llod:llod),
     >        tulou(0:lmd,nlod,-llod:llod),
     >        tulod(0:lmd,nlod,-llod:llod),
     >        tuloulo(nlod,-llod:llod,nlod,-llod:llod)
      COMPLEX, INTENT (IN) :: angso(0:lmaxd,0:lmaxd)
      COMPLEX, INTENT (IN) :: yl1((lmaxd+1)**2),yl2((lmaxd+1)**2)
      INTEGER, INTENT (IN) :: llo(nlod,ntypd),nlo(ntypd)
      INTEGER, INTENT (IN) :: lo1l(0:llod,ntypd)
      INTEGER, INTENT (IN) :: jri(ntypd),lmax(ntypd)
      INTEGER, INTENT (IN) :: ulo_der(nlod,ntypd)
      REAL,    INTENT (IN) :: dx(ntypd),rmsh(jmtd,ntypd)
      REAL,    INTENT (IN) :: ello(nlod,ntypd,max(jspd,2))
      REAL,    INTENT (IN) :: epar(0:lmaxd,ntypd,max(jspd,2))
      REAL,    INTENT (IN) :: vr(jmtd,ntypd,jspd)
      REAL,    INTENT (IN) :: flo  (ntypd,jmtd,2,nlod),
     >                        flo_b(ntypd,jmtd,2,nlod)
      REAL,    INTENT (IN) :: p(jmtd,0:lmaxd,0:lmaxd),
     >                        q(jmtd,0:lmaxd,0:lmaxd),
     >                        p_b(jmtd,0:lmaxd,0:lmaxd),
     >                        q_b(jmtd,0:lmaxd,0:lmaxd)
      REAL,    INTENT (IN) :: kdiff (3,nntot)
      REAL,    INTENT (IN) :: kdiff2(3,nntot2)
      REAL,    INTENT (IN) :: jj1(0:lmaxd,jmtd),jj2(0:lmaxd,jmtd)
      LOGICAL, INTENT (IN) :: l_dulo(nlod,ntypd),l_spav
C     ..
C     .. Local Scalars ..
      REAL gc1,gc2,c,e,r0
80
      COMPLEX cil
Daniel Wortmann's avatar
Daniel Wortmann committed
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
      INTEGER i,l,l2,lh,lm,lmp,lp,lo,lop
      INTEGER lp1,lpl,m,mp,mu
      INTEGER ljj1,ljj2,mjj1,mjj2
      INTEGER ltil1,mtil1,mtil2
      INTEGER lmini,lmini2, lmaxi,lmaxi2
      INTEGER ll,llp,lljj1,lljj2,lmjj1,lmjj2
      INTEGER lltil1,lmtil1,lmtil2,sp1,sp2
      INTEGER mmin,mmax,lwn, ikpt_b,ikpt_b2
C     ..
C     .. Local Arrays ..
      REAL, ALLOCATABLE :: dvulo(:,:,:,:,:),ulovd(:,:,:,:,:)
      REAL, ALLOCATABLE :: uvulo(:,:,:,:,:),ulovu(:,:,:,:,:)
      REAL, ALLOCATABLE :: ulovulo(:,:,:,:,:)
      REAL, ALLOCATABLE :: plo(:,:,:),plo_b(:,:,:)
      REAL, ALLOCATABLE :: x(:)
      REAL, ALLOCATABLE :: v0(:),vso(:,:)
      INTEGER :: ispjsp(2)
      DATA ispjsp/1,-1/

#if (defined(CPP_MPI) && !defined(CPP_T90))
      INCLUDE 'mpif.h'
      INTEGER ierr(3)
#endif

      c = c_light(1.0)
      sp1 = ispjsp(jspin)
      sp2 = ispjsp(jspin_b)

      allocate( ulovulo(0:lmaxd,0:lmaxd,0:lmaxd,nlod,nlod) )
      allocate( ulovu(0:lmaxd,0:lmaxd,0:lmaxd,0:lmaxd,nlod) )
      allocate( ulovd(0:lmaxd,0:lmaxd,0:lmaxd,0:lmaxd,nlod) )
      allocate( uvulo(0:lmaxd,0:lmaxd,0:lmaxd,0:lmaxd,nlod) )
      allocate( dvulo(0:lmaxd,0:lmaxd,0:lmaxd,0:lmaxd,nlod) )
      allocate( x(jmtd) )
      allocate( v0(jmtd), vso(jmtd,2) )
      allocate( plo(jmtd,0:lmaxd,nlod), plo_b(jmtd,0:lmaxd,nlod) )

      lwn = lmax(ntyp)
      r0 = rmsh(1,ntyp)

       DO lo = 1, nlo(ntyp)
        DO ljj1 = 0, lwn
         DO i=1,jri(ntyp)
          plo(i,ljj1,lo) = flo(ntyp,i,1,lo)*jj1(ljj1,i)
          plo_b(i,ljj1,lo) = flo_b(ntyp,i,1,lo)*jj2(ljj1,i)
         ENDDO
        ENDDO
       ENDDO

c****************************************************c
c   compute radial integrals                         c
c     <u(l')jj1(ljj1) | vso(lh) | u(l)jj2(ljj2)>     c
c****************************************************c

! lo-lo
       DO lop = 1, nlo(ntyp)
        lp = llo(lop,ntyp)
        DO lo = 1, nlo(ntyp)
         l = llo(lop,ntyp)

        ! construct vso
        v0 = 0.
        IF (jspins.EQ.1) THEN
           v0(1:jri(ntyp)) = vr(1:jri(ntyp),ntyp,1)
           e = ( ello(lo,ntyp,1) + ello(lop,ntyp,1) )/2.
        ELSE
           DO i = 1,jri(ntyp)
             v0(i) = (vr(i,ntyp,1)+vr(i,ntyp,jspins))/2.
           END DO
           e = ( ello(lo,ntyp,1) + ello(lo,ntyp,jspins)
     >          +ello(lop,ntyp,1) + ello(lop,ntyp,jspins) )/4.
        END IF

        CALL sointg(ntyp,e,vr(:,ntyp,:),v0,atoms,input,vso)  

        ! spin-averaging
        if(l_spav)then
            DO i= 1,jmtd
             vso(i,1)= (vso(i,1)+vso(i,2))/2.
             vso(i,2)= vso(i,1)
           ENDDO        
        endif

        if(l_socscreen) vso(:,:) = vso(:,:)*socscreen_fac

        ! calculate radial integrals
         DO ljj1 = 0, lwn
          DO ljj2 = 0, lwn
           DO lh=0,lwn
              lmini = abs(lp-ljj1)
              lmini2= abs(lh-ljj2)
              lmaxi = lp+ljj1
              lmaxi2= lh+ljj2
              if((lmini.gt.lh) .or. (lmaxi.lt.lh).or.
     >           (lmini2.gt.l) .or. (lmaxi2.lt.l).or.
     >           (mod(lp+lh+ljj1,2).ne.0).or.
     >           (mod( l+lh+ljj2,2).ne.0) ) then

               ulovulo(lh,ljj2,ljj1,lo,lop) = 0.

              else

               DO i=1,jri(ntyp)
                x(i) = plo(i,ljj1,lop)*plo_b(i,ljj2,lo)*vso(i,jspin)
               ENDDO
               call intgr0(x(1:jri(ntyp)),r0,dx(ntyp),jri(ntyp),
     >                     ulovulo(lh,ljj2,ljj1,lo,lop))

              endif
           ENDDO
          ENDDO
         ENDDO
        ENDDO
       ENDDO!lh

! apw-lo
       DO lop = 1, nlo(ntyp)
        lp = llo(lop,ntyp)
        DO l = 0, lwn

        ! construct vso
        v0 = 0.
        IF (jspins.EQ.1) THEN
           v0(1:jri(ntyp)) = vr(1:jri(ntyp),ntyp,1)
           e = ( epar(l,ntyp,1) + ello(lop,ntyp,1) )/2.
        ELSE
           DO i = 1,jri(ntyp)
             v0(i) = (vr(i,ntyp,1)+vr(i,ntyp,jspins))/2.
           END DO
           e = ( epar(l,ntyp,1) + epar(l,ntyp,jspins)
     >          +ello(lop,ntyp,1) + ello(lop,ntyp,jspins) )/4.
        END IF

        CALL sointg(ntyp,e,vr(:,ntyp,:),v0,atoms,input,vso)  

        ! spin-averaging
        if(l_spav)then
            DO i= 1,jmtd
             vso(i,1)= (vso(i,1)+vso(i,2))/2.
             vso(i,2)= vso(i,1)
           ENDDO        
        endif

        ! calculate radial integrals
         DO ljj1 = 0, lwn
          DO ljj2 = 0, lwn
           DO lh=0,lwn
              lmini = abs(lp-ljj1)
              lmini2= abs(lh-ljj2)
              lmaxi = lp+ljj1
              lmaxi2= lh+ljj2
              if((lmini.gt.lh) .or. (lmaxi.lt.lh).or.
     >           (lmini2.gt.l) .or. (lmaxi2.lt.l).or.
     >           (mod(lp+lh+ljj1,2).ne.0).or.
     >           (mod( l+lh+ljj2,2).ne.0) ) then

               ulovu(lh,l,ljj2,ljj1,lop) = 0.
               ulovd(lh,l,ljj2,ljj1,lop) = 0.

              else

               DO i=1,jri(ntyp)
                x(i) = plo(i,ljj1,lop)*p_b(i,l,ljj2)*vso(i,jspin)
               ENDDO
               call intgr0(x(1:jri(ntyp)),r0,dx(ntyp),jri(ntyp),
     >                     ulovu(lh,l,ljj2,ljj1,lop))
               DO i=1,jri(ntyp)
                x(i) = plo(i,ljj1,lop)*q_b(i,l,ljj2)*vso(i,jspin)
               ENDDO
               call intgr0(x(1:jri(ntyp)),r0,dx(ntyp),jri(ntyp),
     >                     ulovd(lh,l,ljj2,ljj1,lop))

              endif

              lmini = abs(l-ljj1)
              lmini2= abs(lh-ljj2)
              lmaxi = l+ljj1
              lmaxi2= lh+ljj2
              if((lmini.gt.lh) .or. (lmaxi.lt.lh).or.
     >           (lmini2.gt.lp) .or. (lmaxi2.lt.lp).or.
     >           (mod( l+lh+ljj1,2).ne.0).or.
     >           (mod(lp+lh+ljj2,2).ne.0) ) then

               uvulo(lh,l,ljj2,ljj1,lop) = 0.
               dvulo(lh,l,ljj2,ljj1,lop) = 0.

              else

               DO i=1,jri(ntyp)
                x(i) = p(i,l,ljj1)*plo_b(i,ljj2,lop)*vso(i,jspin)
               ENDDO
               call intgr0(x(1:jri(ntyp)),r0,dx(ntyp),jri(ntyp),
     >                     uvulo(lh,l,ljj2,ljj1,lop))
               DO i=1,jri(ntyp)
                x(i) = q(i,l,ljj1)*plo_b(i,ljj2,lop)*vso(i,jspin)
               ENDDO
               call intgr0(x(1:jri(ntyp)),r0,dx(ntyp),jri(ntyp),
     >                     dvulo(lh,l,ljj2,ljj1,lop))
             
              endif
           ENDDO
          ENDDO
         ENDDO
        ENDDO
       ENDDO!lh

c***************** SOC CONTRIBUTION **********************c
c   compute product of the two Gaunt coefficients         c
c   with the radial integrals (+prefactors)               c
c*********************************************************c
c   We deal with two Gaunt coefficients:                  c
c  G1 = G( (ltil1, mtil1), (ljj1, mjj1), (lp, mp)       ) c
c  G2 = G( (l , m)       , (ljj2, mjj2), (ltil1, mtil2) ) c
c*********************************************************c
c   use Gaunt conditions to reduce number of operations.  c
c   coefficient G(L1,L2,L3) only nonzero if               c
c       a)  l1 + l2 + l3 = even                           c
c       b)  |l2-l3| <= l1 <= l2+l3                        c
c       c)  m1 = m2 + m3                                  c
c*********************************************************c

! lo-lo
       DO lop=1,nlo(ntyp)
        lp = llo(lop,ntyp)
        DO mp=-lp,lp
         DO ltil1=0,lwn
          lltil1 = ltil1*(ltil1+1)
          DO mtil1=-ltil1,ltil1
           lmtil1 = lltil1+mtil1
           mmin = max(-ltil1,mtil1-1)
           mmax = min( ltil1,mtil1+1)
           mjj1=mtil1-mp
           DO ljj1=0,lwn
            lljj1 = ljj1 * (ljj1 + 1)
            lmjj1 = lljj1 + mjj1
            lmini = abs(lp-ljj1)
            lmaxi = lp+ljj1
            ! Gaunt conditions (G1):
            if((lmini.gt.ltil1).or.(lmaxi.lt.ltil1).or.
     >         (mod(lp+ltil1+ljj1,2).ne.0).or.
     >         (abs(mjj1).gt.ljj1)) cycle

             gc1 = gaunt1(ltil1, ljj1, lp, 
     >                    mtil1, mjj1, mp, lmaxd)

        DO lo=1,nlo(ntyp)
         l = llo(lo,ntyp)
         DO m=-l,l
          DO mtil2=mmin,mmax     ! mtil1-1 <= mtil2 <= mtil1+1
           lmtil2 = lltil1+mtil2
           mjj2=m-mtil2
           DO ljj2=0,lwn
            lljj2 = ljj2 * (ljj2 + 1)
            lmjj2 = lljj2 + mjj2
            lmini2 = abs(ltil1-ljj2)
            lmaxi2 = ltil1+ljj2
             ! Gaunt conditions (G2):
             if((lmini2.gt.l).or.(lmaxi2.lt.l).or.
     >          (mod(l+ltil1+ljj2,2).ne.0).or.
     >          (abs(mjj2).gt.ljj2)) cycle

             gc2 = gaunt1( l, ljj2, ltil1,
     >                     m, mjj2, mtil2, lmaxd)

             ! set up prefactor
346
             cil =   ( ImagUnit ** (lp - l + ljj1 + ljj2) )
Daniel Wortmann's avatar
Daniel Wortmann committed
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
     +              *( (-1.0) ** ljj1 )
     +              * gc1 * gc2
     +              * conjg( yl1(lmjj1 + 1) )
     +              * conjg( yl2(lmjj2 + 1) )
     +              * conjg( angso(lmtil1,lmtil2) )

             ! compute T-coefficients
             tuloulo(lo, m, lop, mp) = tuloulo(lo, m, lop, mp)
     >         + cil * ulovulo(ltil1, ljj2, ljj1, lo, lop)

           ENDDO!ljj2
          ENDDO!m
         ENDDO!l
        ENDDO!mtil2
           ENDDO!ljj1
          ENDDO!mp
         ENDDO!lp
        ENDDO!mtil1
       ENDDO!ltil1

! apw-lo
       DO lop=1,nlo(ntyp)
        lp = llo(lop,ntyp)
        DO mp=-lp,lp
         DO ltil1=0,lwn
          lltil1 = ltil1*(ltil1+1)
          DO mtil1=-ltil1,ltil1
           lmtil1 = lltil1+mtil1
           mmin = max(-ltil1,mtil1-1)
           mmax = min( ltil1,mtil1+1)
           mjj1=mtil1-mp
           DO ljj1=0,lwn
            lljj1 = ljj1 * (ljj1 + 1)
            lmjj1 = lljj1 + mjj1
            lmini = abs(lp-ljj1)
            lmaxi = lp+ljj1
            ! Gaunt conditions (G1):
            if((lmini.gt.ltil1).or.(lmaxi.lt.ltil1).or.
     >         (mod(lp+ltil1+ljj1,2).ne.0).or.
     >         (abs(mjj1).gt.ljj1)) cycle

             gc1 = gaunt1(ltil1, ljj1, lp, 
     >                    mtil1, mjj1, mp, lmaxd)

        DO mtil2=mmin,mmax     ! mtil1-1 <= mtil2 <= mtil1+1
         lmtil2 = lltil1+mtil2
         DO l=0,lwn
          ll = l*(l+1)
          DO m=-l,l
           lm = ll+m
           mjj2=m-mtil2
           DO ljj2=0,lwn
            lljj2 = ljj2 * (ljj2 + 1)
            lmjj2 = lljj2 + mjj2
            lmini2 = abs(ltil1-ljj2)
            lmaxi2 = ltil1+ljj2
             ! Gaunt conditions (G2):
             if((lmini2.gt.l).or.(lmaxi2.lt.l).or.
     >          (mod(l+ltil1+ljj2,2).ne.0).or.
     >          (abs(mjj2).gt.ljj2)) cycle

             gc2 = gaunt1( l, ljj2, ltil1,
     >                     m, mjj2, mtil2, lmaxd)

             ! set up prefactor
412
             cil =   ( ImagUnit ** (lp - l + ljj1 + ljj2) )
Daniel Wortmann's avatar
Daniel Wortmann committed
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
     +              *( (-1.0) ** ljj1 )
     +              * gc1 * gc2
     +              * conjg( yl1(lmjj1 + 1) )
     +              * conjg( yl2(lmjj2 + 1) )
     +              * conjg( angso(lmtil1,lmtil2) )

             ! compute T-coefficients
             tulou(lm, lop, mp) = tulou(lm, lop, mp)
     >         + cil * ulovu(ltil1, l, ljj2, ljj1, lop)

             tulod(lm, lop, mp) = tulod(lm, lop, mp)
     >         + cil * ulovd(ltil1, l, ljj2, ljj1, lop)

           ENDDO!ljj2
          ENDDO!m
         ENDDO!l
        ENDDO!mtil2
           ENDDO!ljj1
          ENDDO!mp
         ENDDO!lp
        ENDDO!mtil1
       ENDDO!ltil1


       DO l=0,lwn
        ll = l*(l+1)
        DO m=-l,l
         lm = ll+m
         DO ltil1=0,lwn
          lltil1 = ltil1*(ltil1+1)
          DO mtil1=-ltil1,ltil1
           lmtil1 = lltil1+mtil1
           mmin = max(-ltil1,mtil1-1)
           mmax = min( ltil1,mtil1+1)
           mjj1=mtil1-m
           DO ljj1=0,lwn
            lljj1 = ljj1 * (ljj1 + 1)
            lmjj1 = lljj1 + mjj1
            lmini = abs(l-ljj1)
            lmaxi = l+ljj1
            ! Gaunt conditions (G1):
            if((lmini.gt.ltil1).or.(lmaxi.lt.ltil1).or.
     >         (mod(l+ltil1+ljj1,2).ne.0).or.
     >         (abs(mjj1).gt.ljj1)) cycle

             gc1 = gaunt1(ltil1, ljj1, l, 
     >                    mtil1, mjj1, m, lmaxd)

        DO mtil2=mmin,mmax     ! mtil1-1 <= mtil2 <= mtil1+1
         lmtil2 = lltil1+mtil2
         DO lop=1,nlo(ntyp)
          lp = llo(lop,ntyp)
          DO mp=-lp,lp
           mjj2=mp-mtil2
           DO ljj2=0,lwn
            lljj2 = ljj2 * (ljj2 + 1)
            lmjj2 = lljj2 + mjj2
            lmini2 = abs(ltil1-ljj2)
            lmaxi2 = ltil1+ljj2
             ! Gaunt conditions (G2):
             if((lmini2.gt.lp).or.(lmaxi2.lt.lp).or.
     >          (mod(lp+ltil1+ljj2,2).ne.0).or.
     >          (abs(mjj2).gt.ljj2)) cycle

             gc2 = gaunt1( lp, ljj2, ltil1,
     >                     mp, mjj2, mtil2, lmaxd)

             ! set up prefactor
481
             cil =   ( ImagUnit ** (l - lp + ljj1 + ljj2) )
Daniel Wortmann's avatar
Daniel Wortmann committed
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
     +              *( (-1.0) ** ljj1 )
     +              * gc1 * gc2
     +              * conjg( yl1(lmjj1 + 1) )
     +              * conjg( yl2(lmjj2 + 1) )
     +              * conjg( angso(lmtil1,lmtil2) )

             ! compute T-coefficients
             tuulo(lm, lop, mp) = tuulo(lm, lop, mp)
     >         + cil * uvulo(ltil1, l, ljj2, ljj1, lop)

             tdulo(lm, lop, mp) = tdulo(lm, lop, mp)
     >         + cil * dvulo(ltil1, l, ljj2, ljj1, lop)

           ENDDO!ljj2
          ENDDO!m
         ENDDO!l
        ENDDO!mtil2
           ENDDO!ljj1
          ENDDO!mp
         ENDDO!lp
        ENDDO!mtil1
       ENDDO!ltil1

#if (defined(CPP_MPI) && !defined(CPP_T90))
      CALL MPI_BARRIER(MPI_COMM_WORLD,ierr)
#endif

      deallocate( ulovulo )
      deallocate( uvulo, dvulo, ulovu, ulovd )
      deallocate( plo,plo_b )
      deallocate( x )
      deallocate( v0, vso )

      END SUBROUTINE wann_uHu_soc_tlo
      END MODULE m_wann_uHu_soc_tlo