tlmplm.F90 10.2 KB
Newer Older
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
MODULE m_tlmplm
  IMPLICIT NONE
  !*********************************************************************
  !     sets up the t(l'm',lm) matrices for each atom type. these
  !     matrices are the k-point independent quantities needed to
  !     determine the hamiltonian contributions coming from
  !     the spheres.  the results are written to direct access
  !     unit28 labelled by atom type number (this allows for
  !     multitasking).  only non=zero elements are calculated
  !     and an index array is also set up.
  !               m. weinert     1986
  !---------------------------------------------------------------------
  !     introduced lh0 and vr0 for gw-calculations (see below)
  !*********************************************************************
  CONTAINS
  SUBROUTINE tlmplm(sphhar,atoms,dimension,enpara,&
       jspin,jsp,mpi, vr,gwc,lh0,input, td,ud)

    USE m_intgr, ONLY : intgr3
    USE m_radflo
    USE m_radfun
    USE m_tlo
    USE m_gaunt, ONLY: gaunt1,gaunt2
    USE m_types
    IMPLICIT NONE
    TYPE(t_mpi),INTENT(IN)      :: mpi
    TYPE(t_dimension),INTENT(IN):: dimension
    TYPE(t_input),INTENT(IN)    :: input
    TYPE(t_sphhar),INTENT(IN)   :: sphhar
    TYPE(t_atoms),INTENT(IN)    :: atoms
    TYPE(t_enpara),INTENT(IN)   :: enpara



    !     ..
    !     .. Scalar Arguments ..
    INTEGER, INTENT (IN) :: gwc,lh0 ! normally lh0 = 1; for
    INTEGER, INTENT (IN) :: jspin,jsp !physical spin&spin index for data
    !     ..
    !     .. Array Arguments ..

Daniel Wortmann's avatar
Daniel Wortmann committed
42
    REAL,    INTENT (IN) :: vr(atoms%jmtd,0:sphhar%nlhd,atoms%ntype)   ! this is for the
43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61
    TYPE(t_tlmplm),INTENT(INOUT) :: td
    TYPE(t_usdus),INTENT(INOUT)  :: ud

    !     ..
    !     .. Local Scalars ..
    COMPLEX cil,ci
    REAL temp,wronk
    INTEGER i,l,l2,lamda,lh,lm,lmin,lmin0,lmp,lmpl,lmplm,lmx,lmxx,lp
    INTEGER lp1,lpl ,mem,mems,mp,mu,n,nh,noded,nodeu ,na,m,nsym
    LOGICAL l_write
    !     ..
    !     .. Local Arrays ..
    INTEGER indt(0:dimension%lmplmd)
    REAL dvd(0:atoms%lmaxd*(atoms%lmaxd+3)/2,lh0:sphhar%nlhd )
    REAL dvu(0:atoms%lmaxd*(atoms%lmaxd+3)/2,lh0:sphhar%nlhd )
    REAL uvd(0:atoms%lmaxd*(atoms%lmaxd+3)/2,lh0:sphhar%nlhd )
    REAL uvu(0:atoms%lmaxd*(atoms%lmaxd+3)/2,lh0:sphhar%nlhd )
    REAL f(atoms%jmtd,2,0:atoms%lmaxd),g(atoms%jmtd,2,0:atoms%lmaxd),x(atoms%jmtd)
    REAL flo(atoms%jmtd,2,atoms%nlod)
Daniel Wortmann's avatar
Daniel Wortmann committed
62 63
    REAL uuilon(atoms%nlod,atoms%ntype),duilon(atoms%nlod,atoms%ntype)
    REAL ulouilopn(atoms%nlod,atoms%nlod,atoms%ntype)
64

65
#ifdef CPP_MPI
66 67 68 69 70 71 72 73 74 75 76 77 78
     INCLUDE 'mpif.h'
    INTEGER ierr(3)
#endif
    !     ..
    ci = cmplx(0.0,1.0)
    td%tdulo(:,:,:,jsp) = cmplx(0.0,0.0)
    td%tuulo(:,:,:,jsp) = cmplx(0.0,0.0)
    td%tuloulo(:,:,:,jsp) = cmplx(0.0,0.0)

    !
    !--->    generate the wavefunctions for each l
    !
    l_write=mpi%irank==0
Daniel Wortmann's avatar
Daniel Wortmann committed
79 80 81 82 83 84 85 86
    !!$    l_write=.false.
    !!$    call gaunt2(atoms%lmaxd)
    !!$OMP PARALLEL DO DEFAULT(NONE)&
    !!$OMP PRIVATE(indt,dvd,dvu,uvd,uvu,f,g,x,flo,uuilon,duilon,ulouilopn)&
    !!$OMP PRIVATE(cil,temp,wronk,i,l,l2,lamda,lh,lm,lmin,lmin0,lmp,lmpl)&
    !!$OMP PRIVATE(lmplm,lmx,lmxx,lp,lp1,lpl,m,mem,mems,mp,mu,n,nh,noded)&
    !!$OMP PRIVATE(nodeu,nsym,na)&
    !!$OMP SHARED(dimension,atoms,gwc,lh0,jspin,jsp,sphhar,enpara,td,ud,l_write,ci,vr,mpi,input)
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
    DO  n = 1,atoms%ntype
       na=sum(atoms%neq(:n-1))+1

       IF (l_write) WRITE (6,FMT=8000) n
       IF (gwc==2) WRITE (14) atoms%rmsh(1:atoms%jri(n),n),atoms%lmax(n)
       DO l = 0,atoms%lmax(n)
          CALL radfun(l,n,jspin,enpara%el0(l,n,jspin),vr(:,0,n),atoms,&
               f(1,1,l),g(1,1,l),ud,nodeu,noded,wronk)
          IF (l_write) WRITE (6,FMT=8010) l,enpara%el0(l,n,jspin),ud%us(l,n,jspin),&
               ud%dus(l,n,jspin),nodeu,ud%uds(l,n,jspin),ud%duds(l,n,jspin),noded,ud%ddn(l,n,jspin),wronk
          IF (gwc==2) WRITE (14) f(1:atoms%jri(n),1,l),g(1:atoms%jri(n),1,l),&
               f(1:atoms%jri(n),2,l),g(1:atoms%jri(n),2,l), ud%us(l,n,jspin),&
               ud%dus(l,n,jspin),ud%uds(l,n,jspin),ud%duds(l,n,jspin), enpara%el0(l,n,jspin)
       END DO
8000   FORMAT (1x,/,/,' wavefunction parameters for atom type',i3,':',&
                      /,t32,'radial function',t79,'energy derivative',/,t3,&
                      'l',t8,'energy',t26,'value',t39,'derivative',t53,&
                      'nodes',t68,'value',t81,'derivative',t95,'nodes',t107,&
                      'norm',t119,'wronskian')
8010   FORMAT (i3,f10.5,2 (5x,1p,2e16.7,i5),1p,2e16.7)
       !
       !--->   generate the extra wavefunctions for the local orbitals,
       !--->   if there are any.
       !
       IF (atoms%nlo(n).GE.1) THEN
          CALL radflo(atoms,n,jspin,enpara%ello0(1,1,jspin), vr(:,0,n), f,g,mpi,&
               ud, uuilon,duilon,ulouilopn,flo)
          IF (gwc==2) THEN
             DO i=1,atoms%nlo(n)
                WRITE (14) flo(1:atoms%jri(n),1,i),&
                     flo(1:atoms%jri(n),2,i), ud%ulos(i,n,jspin),ud%dulos(i,n,jspin),enpara%ello0(i,n,jspin)
             ENDDO
          ENDIF
       END IF

       nsym = atoms%ntypsy(na)
       nh = sphhar%nlh(nsym)
       !
       !--->    generate the irreducible integrals (u(l'):v(lamda,nu:u(l))
       !--->    for l' .ge. l, but only thos that will contribute
       !
       DO lp = 0,atoms%lmax(n)
          lp1 = (lp* (lp+1))/2
          DO l = 0,lp
             lpl = lp1 + l
             !--->    loop over non-spherical components of the potential: must
             !--->    satisfy the triangular conditions and that l'+l+lamda even
             !--->    (conditions from the gaunt coefficient)
             DO lh = lh0, nh
                lamda = sphhar%llh(lh,nsym)
                lmin = lp - l
                lmx = lp + l
                IF ((mod(lamda+lmx,2).EQ.1) .OR. (lamda.LT.lmin) .OR. (lamda.GT.lmx)) THEN
                   uvu(lpl,lh) = 0.0
                   dvd(lpl,lh) = 0.0
                   uvd(lpl,lh) = 0.0
                   dvu(lpl,lh) = 0.0
                ELSE
                   DO i = 1,atoms%jri(n)
                      x(i) = (f(i,1,lp)*f(i,1,l)+f(i,2,lp)*f(i,2,l))* vr(i,lh,n)
                   END DO
                   CALL intgr3(x,atoms%rmsh(1,n),atoms%dx(n),atoms%jri(n),temp)
                   uvu(lpl,lh) = temp
                   DO i = 1,atoms%jri(n)
                      x(i) = (g(i,1,lp)*f(i,1,l)+g(i,2,lp)*f(i,2,l))* vr(i,lh,n)
                   END DO
                   CALL intgr3(x,atoms%rmsh(1,n),atoms%dx(n),atoms%jri(n),temp)
                   dvu(lpl,lh) = temp
                   DO i = 1,atoms%jri(n)
                      x(i) = (f(i,1,lp)*g(i,1,l)+f(i,2,lp)*g(i,2,l))* vr(i,lh,n)
                   END DO
                   CALL intgr3(x,atoms%rmsh(1,n),atoms%dx(n),atoms%jri(n),temp)
                   uvd(lpl,lh) = temp
                   DO i = 1,atoms%jri(n)
                      x(i) = (g(i,1,lp)*g(i,1,l)+g(i,2,lp)*g(i,2,l))* vr(i,lh,n)
                   END DO
                   CALL intgr3(x,atoms%rmsh(1,n),atoms%dx(n),atoms%jri(n),temp)
                   dvd(lpl,lh) = temp
                END IF
             END DO
          END DO
       END DO
       !--->    generate the various t(l'm',lm) matrices for l'm'.ge.lm
       td%tuu(:,n,jsp) = cmplx(0.0,0.0)
       td%tdd(:,n,jsp) = cmplx(0.0,0.0)
       td%tud(:,n,jsp) = cmplx(0.0,0.0)
       td%tdu(:,n,jsp) = cmplx(0.0,0.0)
       indt(:) = 0
       !--->    loop over l'm'
       DO lp = 0,atoms%lmax(n)
          lp1 = (lp* (lp+1))/2
          DO mp = -lp,lp
             lmp = lp* (lp+1) + mp
             lmpl = (lmp* (lmp+1))/2
             !--->    loop over lattice harmonics
             DO lh = lh0, nh
                lamda = sphhar%llh(lh,nsym)
                lmin0 = abs(lp-lamda)
                IF (lmin0.GT.lp) CYCLE
                !-->     ensure l+l'+lamda even
                lmxx = lp - mod(lamda,2)
                mems = sphhar%nmem(lh,nsym)
                DO mem = 1,mems
                   mu = sphhar%mlh(mem,lh,nsym)
                   m = mp - mu
                   lmin = max(lmin0,abs(m))
                   l2 = abs(lmxx-lmin)
                   lmin = lmin + mod(l2,2)
                   DO l = lmin,lmxx,2
                      lm = l* (l+1) + m
                      IF (lm.GT.lmp) CYCLE
                      lpl = lp1 + l
                      lmplm = lmpl + lm
                      cil = ((ci** (l-lp))*sphhar%clnu(mem,lh,nsym))*&
                           gaunt1(lp,lamda,l,mp,mu,m,atoms%lmaxd)
                      td%tuu(lmplm,n,jsp) = td%tuu(lmplm,n,jsp) + cil*uvu(lpl,lh)
                      td%tdd(lmplm,n,jsp) = td%tdd(lmplm,n,jsp) + cil*dvd(lpl,lh)
                      td%tud(lmplm,n,jsp) = td%tud(lmplm,n,jsp) + cil*uvd(lpl,lh)
                      td%tdu(lmplm,n,jsp) = td%tdu(lmplm,n,jsp) + cil*dvu(lpl,lh)
                      indt(lmplm) = 1
                   END DO
                END DO
             END DO
          END DO
       END DO
       !--->    include diagonal terms from muffin-tin hamiltonian
       !--->   only, if secvar = .false. If second variation is
       !--->   is switched on, the tlmplm matrix contains only the
       !--->   off-diagonal terms
       !-gu deleted: 20.11.97
       !--->    set up mapping array
       DO lp = 0,atoms%lmax(n)
          DO mp = -lp,lp
             lmp = lp* (lp+1) + mp
             DO l = 0,atoms%lmax(n)
                DO m = -l,l
                   lm = l* (l+1) + m
                   IF (lmp.GE.lm) THEN
                      lmplm = (lmp* (lmp+1))/2 + lm
                      IF (indt(lmplm).NE.0) THEN
                         td%ind(lmp,lm,n,jsp) = lmplm
                      ELSE
                         td%ind(lmp,lm,n,jsp) = -9999
                      END IF
                   ELSE
                      lmplm = (lm* (lm+1))/2 + lmp
                      IF (indt(lmplm).NE.0) THEN
                         td%ind(lmp,lm,n,jsp) = -lmplm
                      ELSE
                         td%ind(lmp,lm,n,jsp) = -9999
                      END IF
                   END IF
                END DO
             END DO
          END DO
       END DO
       !
       !--->   set up the t-matrices for the local orbitals,
       !--->   if there are any
       IF (atoms%nlo(n).GE.1) THEN
247
          CALL tlo(atoms,sphhar,jspin,jsp,n,enpara,lh0,input,vr(1,0,n),&
248
                   na,flo,f,g,ud, uuilon,duilon,ulouilopn, td)
249 250 251 252

       ENDIF

    ENDDO
Daniel Wortmann's avatar
Daniel Wortmann committed
253
    !!$OMP END PARALLEL DO
254 255 256 257


  END SUBROUTINE tlmplm
END MODULE m_tlmplm