wann_updown.F 39.8 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 10
      MODULE m_wann_updown
      use m_juDFT
      CONTAINS
      SUBROUTINE wann_updown(
11
     >      mpi,input,sym,atoms,stars,vacuum,sphhar,oneD,noco,cell,vTot,
12
     >      enpara,
13
     >      eig_id,l_real,mpi_comm,l_dulo,l_noco,l_ss,lmaxd,ntypd,
14
     >      neigd,natd,nop,nvd,jspd,nbasfcn,llod,nlod,ntype,
15
     >      omtil,nlo,llo,lapw_l,invtab,mrot,ngopr,neq,lmax,
16 17 18 19 20 21
     >      invsat,invsatnr,nkpt,taual,rmt,amat,bmat,bbmat,alph,
     >      beta,qss,sk2,phi2,odi,ods,irank,isize,n3d,nmzxyd,nmzd,
     >      jmtd,nlhd,nq3,nvac,invs,invs2,film,nlh,jri,ntypsd,
     >      ntypsy,jspins,nkptd,dx,n2d,rmsh,e1s,e2s,ulo_der,
     >      ustep,ig,k1d,k2d,k3d,rgphs,slice,kk,nnne,
     >      z1,nv2d,nmzxy,nmz,delz,ig2,area,tau,zatom,nq2,nop2,
22
     >      volint,symor,pos,ef,l_soc,
23
     >      memd,lnonsph,clnu,lmplmd,mlh,nmem,llh,lo1l,
24
     >      theta,phi)
25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43
c****************************************************************************
c     Calculate spin-up-spin-down matrix elements.
c     Useful in the case of spin-orbit coupling and noncollinear magnetism.
c     Implemented matrix elements:
c         (i)   Overlaps => Pauli matrix elements.
c         (ii)  Momentum operator => Spin current.
c         (iii) H_{SOC} => Include SOC in the basis set of Wannier functions.
c         (iv)  Spin-current into the MT-spheres => Torque.
c         (v)   Matrix elements of exchange field => Torque.
c     Frank Freimuth
c****************************************************************************
      use m_types
      use m_juDFT
      use m_wann_rw_eig
      use m_abcof
      use m_radfun
      use m_radflo
      use m_cdnread, only : cdn_read0, cdn_read
      use m_constants, only : pimach
44
!      use m_wann_mmk0_od_vac
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
      use m_wann_mmk0_vac
      use m_wann_mmk0_updown_sph
      use m_wann_abinv
      use m_wann_kptsrotate
      use m_wann_read_inp
      use m_wann_radovlp_integrals
      use m_matmul,only : matmul3,matmul3r
      use m_wann_socmat
      use m_wann_write_hsomtx
      use m_wann_write_nabla
      use m_wann_write_amn
      use m_wann_write_nabla
      use m_wann_mmkb_int
#ifdef CPP_TOPO
      use m_wann_perpmag
      use m_wann_nabla_updown
      use m_wann_surfcurr_updown
      use m_wann_surfcurr_int_updown
      use m_wann_gabeffgagaunt
#endif
      IMPLICIT NONE
#include "cpp_double.h"
#ifdef CPP_MPI
      include 'mpif.h'
      integer ierr(3)
      integer cpu_index
      integer stt(MPI_STATUS_SIZE)

#endif
74 75 76 77 78 79 80 81 82 83 84

      TYPE(t_mpi),INTENT(IN)    :: mpi
      TYPE(t_input),INTENT(IN)  :: input
      TYPE(t_sym),INTENT(IN)    :: sym
      TYPE(t_atoms),INTENT(IN)  :: atoms
      TYPE(t_stars),INTENT(IN)  :: stars
      TYPE(t_vacuum),INTENT(IN) :: vacuum
      TYPE(t_sphhar),INTENT(IN) :: sphhar
      TYPE(t_oneD),INTENT(IN)   :: oneD
      TYPE(t_noco),INTENT(IN)   :: noco
      TYPE(t_cell),INTENT(IN)   :: cell
85
      TYPE(t_potden),INTENT(IN) :: vTot
86 87
      TYPE(t_enpara),INTENT(IN) :: enpara

88
      logical, intent (in) :: invs,invs2,film,slice,symor,l_real
89
      integer, intent (in) :: lmaxd,ntypd,neigd,nkptd,kk,nnne,mpi_comm
90
      integer, intent (in) :: natd,nop,nvd,jspd,nbasfcn,nq2,nop2,eig_id
91
      integer, intent (in) :: llod,nlod,ntype,n3d,n2d
92 93 94 95
      integer, intent (in) :: nmzxyd,nmzd,jmtd,nlhd,nq3,nvac
      integer, intent (in) :: ntypsd,jspins,k1d,k2d,k3d
      real,    intent (in) :: omtil,e1s,e2s,delz,area,z1,volint
      integer, intent (in) :: ig(-k1d:k1d,-k2d:k2d,-k3d:k3d)
96
      complex, intent (in) :: rgphs(-k1d:k1d,-k2d:k2d,-k3d:k3d)
97 98 99 100
      integer, intent (in) :: nlh(ntypsd),jri(ntypd),ntypsy(natd)
      integer, intent (in) :: nlo(ntypd),llo(nlod,ntypd),lapw_l(ntypd)
      integer, intent (in) :: invtab(nop),mrot(3,3,nop),ngopr(natd)
      integer, intent (in) :: neq(ntypd),lmax(ntypd)
101
      integer, intent (in) :: invsat(natd),invsatnr(natd),nkpt
102 103 104 105 106
      integer, intent (in) :: irank,isize,nv2d,nmzxy,nmz
      integer, intent (in) :: ulo_der(nlod,ntypd),ig2(n3d)
      real,    intent (in) :: taual(3,natd),rmt(ntypd),dx(ntypd)
      real,    intent (in) :: amat(3,3),bmat(3,3),bbmat(3,3)
      real,    intent (in) :: rmsh(jmtd,ntypd),tau(3,nop),zatom(ntype)
107
      real,    intent (in) :: alph(ntypd),beta(ntypd),qss(3)
108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124
      real,    intent (in) :: pos(3,natd)
      real, intent(in)     :: ef
      complex, intent (in) :: ustep(n3d)
      logical, intent (in) :: l_dulo(nlod,ntypd),l_noco,l_ss
      logical, intent (in) :: l_soc
      integer, intent (in) :: memd
      integer, intent (in) :: lnonsph(ntypd)
      complex, intent (in) :: clnu(memd,0:nlhd,ntypsd)
      integer, intent (in) :: lmplmd
      integer, intent (in) :: mlh(memd,0:nlhd,ntypsd)
      integer, intent (in) :: nmem(0:nlhd,ntypsd)
      integer, intent (in) :: llh(0:nlhd,ntypsd)
      integer, intent (in) :: lo1l(0:llod,ntypd)
      real,    intent (in) :: sk2(n2d),phi2(n2d)
      type (od_inp), intent (in) :: odi
      type (od_sym), intent (in) :: ods
      real, intent(in) :: theta,phi
125
     
126 127 128
cccccccccccccccccc   local variables   cccccccccccccccccccc
      type(t_wann) :: wann
      integer :: lmd,nlotot,n,nmat,iter,ikpt,ikpt_b,nmat_b
129
      integer :: addnoco,addnoco2,funbas,loplod
130 131 132 133 134 135 136
      integer :: noccbd,noccbd_b,nn,nkpts,i,jspin,j,l,i_rec,m,nwf,nwfp
      integer :: jsp_start,jsp_end,nrec,nrec1,nbands,nbands_b
      integer :: nodeu,noded,n_size,na,n_rank,nbnd,numbands
      integer :: i1,i2,i3,in,lda
      integer :: n_bands(0:neigd),nslibd,nslibd_b
      character(len=8) :: dop,iop,name(10)
      real    :: bkpt(3),bkpt_b(3),sfp,tpi,wronk,wk,wk_b,phase
137
      complex :: c_phase, zzConjgTemp
138 139 140
      real    :: eig(neigd),cp_time(9),efermi
      logical :: l_p0,l_bkpts,l_proj,l_amn,l_mmn
!!! energy window boundaries
141 142
      integer, allocatable :: nv(:)
      integer, allocatable :: nv_b(:)
143 144
      integer, allocatable :: k1(:,:),k2(:,:),k3(:,:)
      integer, allocatable :: k1_b(:,:),k2_b(:,:),k3_b(:,:)
145

146 147 148 149 150 151 152 153 154
      real,    allocatable :: eigg(:,:)
      real kpoints(nkptd)
!!! a and b coeff. constructed for each k-point
      complex, allocatable :: acof(:,:,:,:)
      complex, allocatable :: bcof(:,:,:,:)
      complex, allocatable :: ccof(:,:,:,:,:)
!!! the parameters for the number of wfs
      integer :: nwfs
!!! the potential in the spheres and the vacuum
155
      real, allocatable :: vr(:,:,:),vz(:,:,:)
156 157 158 159 160 161
!!! bkpts data
      integer nntot,ikpt_help
      integer, allocatable :: gb(:,:,:),bpt(:,:)
!!! radial wavefunctions in the muffin-tins and more ...
      real,    allocatable :: flo(:,:,:,:,:)
      real,    allocatable :: ff(:,:,:,:,:),gg(:,:,:,:,:)
162

163 164 165 166 167 168 169
      real     :: uuilon(nlod,ntypd),duilon(nlod,ntypd)
      real     :: ulouilopn(nlod,nlod,ntypd)
!!! energy parameters
      real    :: ello(nlod,ntypd,max(2,jspd)),evac(2,max(2,jspd))
      real    :: epar(0:lmaxd,ntypd,max(2,jspd)),evdu(2,max(jspd,2))
!!! the Mmn matrices
      complex, allocatable :: surfcurr(:,:,:,:)
170
      complex, allocatable :: hsomtx(:,:,:,:,:),hsomtx_l(:,:,:,:)
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
      complex, allocatable :: mmnk(:,:,:,:),mmn(:,:,:)
      complex, allocatable :: perpmag(:,:,:)
      complex, allocatable :: amn(:,:,:),nablamat(:,:,:,:)
      complex, allocatable :: socamn(:,:,:)
      complex, allocatable :: socmmn(:,:,:)
      complex, allocatable :: socmmnk(:,:,:,:)
      complex, allocatable :: omn(:,:,:),vec(:,:),a(:)
      complex, allocatable :: smn(:,:,:)
      complex, allocatable :: psiw(:,:,:)
c..wf-hamiltonian in reciprocal space Hwf(k)
      complex, allocatable :: hwfk(:,:,:)
c..egenvalues of the wf-hamiltonian Hwf(k)
      real, allocatable :: eigw(:,:)
c..wf-hamiltonian in real space (hopping in the same unit cell)
      complex, allocatable :: hwfr(:,:)
      real, allocatable :: ei(:)
      complex, allocatable :: work(:)
      real,allocatable::centers(:,:,:)
      integer, allocatable :: iwork(:)
      integer :: info,lwork,lrwork,liwork,lee
      character :: jobz, uplo
      logical :: l_file
      logical :: l_amn2, l_conjugate
      character(len=3) :: spin12(2)
      data   spin12/'WF1' , 'WF2'/
      character(len=30)  :: task
      integer,allocatable::irreduc(:)
      integer,allocatable::mapkoper(:)
      integer :: fullnkpts,kpt,kptibz,kptibz_b,j1,j2,j3,oper,oper_b,k
      real :: bkrot(3)
      integer :: ios,kplot,kplotoper,plotoper,gfcut
      integer :: funbas_start,funbas_stop
      complex :: phasust
      integer,allocatable::pair_to_do(:,:)
      integer :: ngopr1(natd)
      integer,allocatable::maptopair(:,:,:)
      integer :: ldu,ldv,wannierspin,jspin2
      complex, allocatable :: amn_copy(:,:)
      complex, allocatable :: sunit(:,:)
      real,    allocatable :: sdiag(:)
      complex, allocatable :: umn(:,:),vmn(:,:)
      real, allocatable :: rwork(:)
      character :: jobu*1,jobv*1
      real,allocatable::kdiff(:,:)
      integer,allocatable :: shiftkpt(:,:)
      integer :: unigrid(6),gfthick
      complex,allocatable::ujug(:,:,:,:),ujdg(:,:,:,:)
      complex,allocatable::djug(:,:,:,:),djdg(:,:,:,:)
      complex,allocatable::ujulog(:,:,:,:)
      complex,allocatable::djulog(:,:,:,:)
      complex,allocatable::ulojug(:,:,:,:)
      complex,allocatable::ulojdg(:,:,:,:)
      complex,allocatable::ulojulog(:,:,:,:)
224 225 226 227
      real delta,delta1,time_interstitial,time_mmn
      real time_total,delta2,delta3
      real time_lapw_expand,time_rw,time_symm,time_film
      real time_lapw_plot,time_ujugaunt,time_abcof
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
      integer :: n_start,n_end,mlotot,mlolotot,err
      integer :: mlot_d,mlolot_d,ilo,dir,l1,l2,ii,jj
      character(len=2)    :: spin012(0:2)
      data spin012/'  ', '.1', '.2'/
      character(len=6)    :: filename
      real                :: arg,hescale
      complex             :: nsfactor,nsfactor_b,value
      real                :: b1(3),b2(3),b3(3)
      real,allocatable    :: radial1_ff(:,:,:,:,:)
      real,allocatable    :: radial1_gg(:,:,:,:,:)
      real,allocatable    :: radial1_fg(:,:,:,:,:)
      real,allocatable    :: radial1_gf(:,:,:,:,:)
      real,allocatable    :: radial1_flo(:,:,:,:,:)
      real,allocatable    :: radial1_glo(:,:,:,:,:)
      real,allocatable    :: radial1_lof(:,:,:,:,:)
      real,allocatable    :: radial1_log(:,:,:,:,:)
      real,allocatable    :: radial1_lolo(:,:,:,:,:)

      real,parameter      :: bohrtocm=0.529177e-8
      real,parameter      :: condquant=7.7480917e-5
      real                :: dirfacs(3)
      complex,allocatable :: ubeffug(:,:,:) 
      complex,allocatable :: ubeffdg(:,:,:) 
      complex,allocatable :: dbeffug(:,:,:) 
      complex,allocatable :: dbeffdg(:,:,:) 
      complex,allocatable :: ubeffulog(:,:,:) 
      complex,allocatable :: dbeffulog(:,:,:) 
      complex,allocatable :: ulobeffug(:,:,:) 
      complex,allocatable :: ulobeffdg(:,:,:) 
      complex,allocatable :: ulobeffulog(:,:,:) 
      integer :: spinloop1,spinloop2,ilop

260
      TYPE(t_usdus) :: usdus
261
      TYPE(t_mat)   :: zMat(2),zzMat(2)
262
      TYPE(t_lapw)  :: lapw
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

      IF(l_ss) CALL juDFT_error("spin-spiral not yet in this version"
     +     ,calledby ="wann_updown")

c-----initializations
      time_interstitial=0.0
      time_mmn=0.0
      time_total=0.0
      time_ujugaunt=0.0
      time_abcof=0.0
      time_rw=0.0
      time_symm=0.0
      time_film=0.0
      ngopr1(:)=1

      call cpu_time(time_total)

      l_p0 = .false.
      if (irank.eq.0) l_p0 = .true.

      sfp = 2* sqrt( pimach() )
      tpi = 2* pimach()
      lmd = lmaxd*(lmaxd+2)

!!!   should be changed in case the windows are really used
288
      nkpts = nkpt
289 290

c-----read the input file to determine what to do
291
      call wann_read_inp(input,l_p0,wann)
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

      IF(wann%l_byenergy.AND.wann%l_byindex) CALL juDFT_error
     +     ("byenergy.and.byindex",calledby ="wann_updown")
      IF(wann%l_byenergy.AND.wann%l_bynumber) CALL juDFT_error
     +     ("byenergy.and.bynumber",calledby ="wann_updown")
      IF(wann%l_bynumber.AND.wann%l_byindex) CALL juDFT_error
     +     ("bynumber.and.byindex",calledby ="wann_updown")
      IF(.NOT.(wann%l_bynumber.OR.wann%l_byindex.OR.wann%l_byenergy))
     &     CALL juDFT_error("no rule to sort bands",calledby
     +     ="wann_updown")

      efermi=ef
      if(.not.wann%l_fermi)efermi=0.0

#ifdef CPP_MPI
         call MPI_BARRIER(mpi_comm,ierr)
#endif

c**************************************************************
c   for bzsym=.true.: determine mapping between kpts and w90kpts
c**************************************************************
      if (wann%l_bzsym) then
         l_file=.false.
         inquire(file='w90kpts',exist=l_file)
         IF(.NOT.l_file)  CALL juDFT_error
     +        ("w90kpts not found, needed if bzsym",calledby
     +        ="wann_updown")
         open(412,file='w90kpts',form='formatted')
         read(412,*)fullnkpts
         close(412)
         if(l_p0)print*,"fullnkpts=",fullnkpts
         IF(fullnkpts<nkpts) CALL juDFT_error("fullnkpts.lt.nkpts"
     +        ,calledby ="wann_updown")
         allocate(irreduc(fullnkpts),mapkoper(fullnkpts))
         allocate(shiftkpt(3,fullnkpts))
         l_file=.false.
         inquire(file='kptsmap',exist=l_file)
         IF(.NOT.l_file)  CALL juDFT_error
     +        ("kptsmap not found, needed if bzsym",calledby
     +        ="wann_updown")
         open(713,file='kptsmap')
         do i=1,fullnkpts
            read(713,*)kpt,irreduc(i),mapkoper(i),shiftkpt(:,i)
335 336
            IF(kpt/=i) CALL juDFT_error("kpt.ne.i",
     +                                  calledby ="wann_updown")
337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362
            if(l_p0)print*,i,irreduc(i),mapkoper(i)
         enddo
         close(713)
         IF(MAXVAL(irreduc(:))/=nkpts) CALL juDFT_error
     +        ("max(irreduc(:))/=nkpts",calledby ="wann_updown")
      else
         fullnkpts=nkpts
      endif

      nrec = 0
      if(l_p0)then
      write (*,*) 'fermi energy:',efermi
      write (*,*) 'emin,emax=',e1s,e2s
      write (*,*) 'nbasfcn =',nbasfcn
      endif
      nlotot = 0
      do n = 1, ntype
        do l = 1,nlo(n)
          nlotot = nlotot + neq(n) * ( 2*llo(l,n) + 1 )
        enddo
      enddo

c*********************************************************
cccccccccccccccc   initialize the potential   cccccccccccc
c*********************************************************

363 364
      allocate (vz(nmzd,2,jspd))
      allocate (vr(jmtd,ntypd,jspd))
365

366
      vz = vTot%vacz
367 368 369 370

      do jspin = 1,jspins
        do n = 1, ntype
          do j = 1,jri(n)
371
            vr(j,n,jspin) = vTot%mt(j,0,n,jspin)
372 373 374 375 376 377 378 379
          enddo
        enddo
      enddo
    
cccccccccccccccc   end of the potential part  ccccccccccc
      wannierspin=jspd
      if(l_soc) wannierspin=2

380
      allocate ( nv(wannierspin) )
381 382 383 384
      allocate ( k1(nvd,wannierspin),k2(nvd,wannierspin),
     &           k3(nvd,wannierspin) )
      allocate ( ff(ntypd,jmtd,2,0:lmaxd,2) )
      allocate ( gg(ntypd,jmtd,2,0:lmaxd,2) )
385 386 387 388 389 390 391 392 393 394
      allocate ( usdus%us(0:lmaxd,ntypd,2) )
      allocate ( usdus%uds(0:lmaxd,ntypd,2) )
      allocate ( usdus%dus(0:lmaxd,ntypd,2) )
      allocate ( usdus%duds(0:lmaxd,ntypd,2) )
      allocate ( usdus%ddn(0:lmaxd,ntypd,2) )
      allocate ( usdus%ulos(nlod,ntypd,2) )
      allocate ( usdus%dulos(nlod,ntypd,2) )
      allocate ( usdus%uulon(nlod,ntypd,2) )
      allocate ( usdus%dulon(nlod,ntypd,2) )
      allocate ( usdus%uloulopn(nlod,nlod,ntypd,2) )
395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418
c****************************************************
c     jspin=1 
c****************************************************
      jspin=1

c...read number of bands and wannier functions from file proj

c..reading the proj.1 / proj.2 / proj file
       l_proj=.false.  
       do j=jspin,0,-1
         inquire(file=trim('proj'//spin012(j)),exist=l_proj)
         if(l_proj)then
            filename='proj'//spin012(j)
            exit
         endif
       enddo

       if(l_proj)then
         open (203,file=trim(filename),status='old')
         rewind (203)
         read (203,*) nwfs,numbands
         rewind (203)
         close (203)
       else
419 420
          CALL juDFT_error("no proj/proj.1/proj.2",
     +                     calledby ="wann_updown")
421 422 423 424 425 426 427 428 429 430 431
       endif  


       jspin2=jspin
       if(l_soc .and. jspins.eq.1)jspin2=1
       jsp_start = jspin ; jsp_end = jspin

cccccccccccc   read in the eigenvalues and vectors   cccccc


       do jspin=1,jspins
432 433 434
          CALL judft_error("TODO:wann_updown")
!         call cdn_read0(eig_id,irank,isize,jspin,jspd,l_noco,
 !    <                  ello,evac,epar,wk,n_bands,n_size)
435 436 437 438 439 440 441 442 443 444 445 446 447 448
       enddo

       jspin=1

c.. now we want to define the maximum number of the bands by all kpts

      nbnd = 0

      i_rec = 0 ; n_rank = 0
c*************************************************************
c..writing down the eig.1 and/or eig.2 files
c*************************************************************
      if(l_p0)then
         call wann_write_eig(
449
     >            eig_id,l_real,
450 451
     >            lmaxd,ntypd,nlod,neigd,nvd,wannierspin,
     >            isize,jspin,nbasfcn,nlotot,
452
     >            l_ss,l_noco,nrec,fullnkpts,
453 454 455 456 457 458
     >            wann%l_bzsym,wann%l_byindex,wann%l_bynumber,
     >            wann%l_byenergy,
     >            irreduc,odi,wann%band_min(jspin),
     >            wann%band_max(jspin),
     >            numbands,
     >            e1s,e2s,efermi,wann%l_pauli,nkpts,
459
     <            nbnd,kpoints,.false.,-1)
460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495
         if(odi%d1)then
            kpoints(:)=kpoints(:)*bmat(3,3)
         endif
      endif!l_p0

! nbnd is calculated for process zero and is sent here to the others
#ifdef CPP_MPI
      if(l_p0)then
         do cpu_index=1,isize-1
      call MPI_SEND(nbnd,1,MPI_INTEGER,cpu_index,1,mpi_comm,ierr)
         enddo
      else
       call MPI_RECV(nbnd,1,MPI_INTEGER,0,1,mpi_comm,stt,ierr)
      endif
#endif

       print*,"process: ",irank," nbnd= ",nbnd
c##################################################################

      allocate ( mmn(nbnd,nbnd,fullnkpts) )
      mmn(:,:,:) = cmplx(0.,0.)
      if((l_soc.or.l_noco) .and. (jspin.eq.1))
     &   allocate( socmmn(nbnd,nbnd,fullnkpts) )

      if(wann%l_nabla)then
         allocate ( nablamat(3,nbnd,nbnd,fullnkpts) )
         nablamat = cmplx(0.,0.)
      endif

      if(wann%l_surfcurr)then
         allocate ( surfcurr(3,nbnd,nbnd,fullnkpts) )
         surfcurr = cmplx(0.,0.)
      endif

      if(wann%l_socmat)then
         allocate ( hsomtx(2,2,nbnd,nbnd,fullnkpts) )
496
         allocate ( hsomtx_l(nbnd,nbnd,2,2) )
497 498 499 500 501 502 503 504 505 506 507 508 509 510 511
      endif

      write (*,*) 'nwfs=',nwfs

      allocate ( flo(ntypd,jmtd,2,nlod,2) )

c***********************************************************
c     Calculate the radial wave functions for the two spins.
c***********************************************************
      do jspin=1,wannierspin
        jspin2=jspin
        if(l_soc .and. jspins.eq.1)jspin2=1
        do n = 1,ntype
          do l = 0,lmax(n)
             call radfun(
512 513 514
     >              l,n,jspin,epar(l,n,jspin2),vr(1,n,jspin2),atoms,
     <              ff(n,:,:,l,jspin),gg(n,:,:,l,jspin),usdus,
     <              nodeu,noded,wronk)
515 516
          enddo !l
          do ilo = 1, nlo(n)
517

518
             call radflo(
519 520 521 522
     >             atoms,n,jspin,ello(:,:,jspin2),vr(1,n,jspin2),
     >             ff(n,1:,1:,0:,jspin),gg(n,1:,1:,0:,jspin),mpi,
     <             usdus,uuilon,duilon,ulouilopn,flo(n,:,:,:,jspin))

523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701
          enddo !ilo
        enddo !n
      enddo !jspin

c****************************************************************
c     Calculate the radial overlap integrals.
c****************************************************************
      if(wann%l_mmn0)then
       allocate( radial1_ff(2,2,0:lmaxd,0:lmaxd,ntype) )
       allocate( radial1_fg(2,2,0:lmaxd,0:lmaxd,ntype) )
       allocate( radial1_gf(2,2,0:lmaxd,0:lmaxd,ntype) )
       allocate( radial1_gg(2,2,0:lmaxd,0:lmaxd,ntype) )

       allocate( radial1_flo(2,2,0:lmaxd,1:nlod,ntype) )
       allocate( radial1_glo(2,2,0:lmaxd,1:nlod,ntype) )
       allocate( radial1_lof(2,2,1:nlod,0:lmaxd,ntype) )
       allocate( radial1_log(2,2,1:nlod,0:lmaxd,ntype) )
       allocate( radial1_lolo(2,2,1:nlod,1:nlod,ntype) )

       do n=1,ntype
        do l1=0,lmaxd
         do l2=0,lmaxd
          do spinloop1=1,2
           do spinloop2=1,2
            call wann_radovlp_integrals(
     >         rmsh(:,n),jri(n),dx(n),
     >         ff(n,:,:,l1,spinloop2),
     >         ff(n,:,:,l2,spinloop1),
     <         radial1_ff(spinloop2,spinloop1,l1,l2,n))
            call wann_radovlp_integrals(
     >         rmsh(:,n),jri(n),dx(n),
     >         ff(n,:,:,l1,spinloop2),
     >         gg(n,:,:,l2,spinloop1),
     <         radial1_fg(spinloop2,spinloop1,l1,l2,n))
            call wann_radovlp_integrals(
     >         rmsh(:,n),jri(n),dx(n),
     >         gg(n,:,:,l1,spinloop2),
     >         ff(n,:,:,l2,spinloop1),
     <         radial1_gf(spinloop2,spinloop1,l1,l2,n))
            call wann_radovlp_integrals(
     >         rmsh(:,n),jri(n),dx(n),
     >         gg(n,:,:,l1,spinloop2),
     >         gg(n,:,:,l2,spinloop1),
     <         radial1_gg(spinloop2,spinloop1,l1,l2,n))
           enddo !spinloop2
          enddo !spinloop1
         enddo !l2
        enddo !l1
       enddo !n

       do n=1,ntype
        do ilo=1,nlo(n)
         do ilop=1,nlo(n)
          do spinloop1=1,2
           do spinloop2=1,2
             call wann_radovlp_integrals(
     > rmsh(:,n),jri(n),dx(n),
     > flo(n,:,:,ilo,spinloop2),
     > flo(n,:,:,ilop,spinloop1),
     < radial1_lolo(spinloop2,spinloop1,ilo,ilop,n))
           enddo !spinloop2
          enddo !spinloop1
         enddo!ilop
        enddo !ilo
       enddo !n

       do n=1,ntype
        do ilo=1,nlo(n)
         do l2=0,lmaxd
          do spinloop1=1,2
           do spinloop2=1,2
             call wann_radovlp_integrals(
     > rmsh(:,n),jri(n),dx(n),
     > flo(n,:,:,ilo,spinloop2),
     > ff(n,:,:,l2,spinloop1),
     < radial1_lof(spinloop2,spinloop1,ilo,l2,n))
           enddo !spinloop2
          enddo !spinloop1
         enddo!l2
        enddo !ilo
       enddo !n

       do n=1,ntype
        do ilo=1,nlo(n)
         do l2=0,lmaxd
          do spinloop1=1,2
           do spinloop2=1,2
             call wann_radovlp_integrals(
     > rmsh(:,n),jri(n),dx(n),
     > flo(n,:,:,ilo,spinloop2),
     > gg(n,:,:,l2,spinloop1),
     < radial1_log(spinloop2,spinloop1,ilo,l2,n))
           enddo !spinloop2
          enddo !spinloop1
         enddo!l2
        enddo !ilo
       enddo !n

       do n=1,ntype
        do ilo=1,nlo(n)
         do l2=0,lmaxd
          do spinloop1=1,2
           do spinloop2=1,2
             call wann_radovlp_integrals(
     > rmsh(:,n),jri(n),dx(n),
     > ff(n,:,:,l2,spinloop2),
     > flo(n,:,:,ilo,spinloop1),
     < radial1_flo(spinloop2,spinloop1,l2,ilo,n))
           enddo !spinloop2
          enddo !spinloop1
         enddo!l2
        enddo !ilo
       enddo !n

       do n=1,ntype
        do ilo=1,nlo(n)
         do l2=0,lmaxd
          do spinloop1=1,2
           do spinloop2=1,2
             call wann_radovlp_integrals(
     > rmsh(:,n),jri(n),dx(n),
     > gg(n,:,:,l2,spinloop2),
     > flo(n,:,:,ilo,spinloop1),
     < radial1_glo(spinloop2,spinloop1,l2,ilo,n))
           enddo !spinloop2
          enddo !spinloop1
         enddo!l2
        enddo !ilo
       enddo !n

      endif


      jspin=1

c****************************************************************
c     compute matrix elements of exchange field
c****************************************************************
#ifdef CPP_TOPO
      if(wann%l_perpmag)then
         allocate( perpmag(nbnd,nbnd,fullnkpts) )
         perpmag=0.0
         allocate( ubeffug(0:lmd,0:lmd,1:ntype) )
         allocate( ubeffdg(0:lmd,0:lmd,1:ntype) )
         allocate( dbeffug(0:lmd,0:lmd,1:ntype) )
         allocate( dbeffdg(0:lmd,0:lmd,1:ntype) )
         allocate( ubeffulog(0:lmd,0:15,1:ntype) )
         allocate( dbeffulog(0:lmd,0:15,1:ntype) )
         allocate( ulobeffug(0:lmd,0:15,1:ntype) )
         allocate( ulobeffdg(0:lmd,0:15,1:ntype) )
         allocate( ulobeffulog(0:15,0:15,1:ntype) )
         call wann_gabeffgagaunt(
     >            neq,mlh,nlhd,nlh,ntypsy,llh,lmax,
     >            nmem,ntype,ntypd,bbmat,bmat,
     >            nlod,nlo,llo,flo,ff,gg,jri,rmsh,dx,jmtd,
     >            lmaxd,lmd,clnu,
     <            ubeffug,ubeffdg,dbeffug,dbeffdg,ubeffulog,
     <            dbeffulog,
     <            ulobeffug,ulobeffdg,ulobeffulog)
      endif   
#endif      


      i_rec = 0 ; n_rank = 0

c****************************************************************
c..   loop by kpoints starts!      each may be a separate task
c****************************************************************
      do 10 ikpt = wann%ikptstart,fullnkpts  ! loop by k-points starts
        kptibz=ikpt
        if(wann%l_bzsym) kptibz=irreduc(ikpt)
        if(wann%l_bzsym) oper=mapkoper(ikpt)

        i_rec = i_rec + 1
      if (mod(i_rec-1,isize).eq.irank) then

c****************************************************************
c     Obtain eigenvectors for the two spin channels.
c****************************************************************
702
      
703 704 705 706 707 708 709 710
      allocate ( eigg(neigd,2) )
      call cpu_time(delta)
      do jspin=1,wannierspin
         if(jspin.eq.1)nrec=0
         if(jspin.eq.2)nrec=nkpts
         if(l_noco)nrec=0
         n_start=1
         n_end=neigd
711

712
      zzMat(jspin)%l_real = l_real
713 714
      zzMat(jspin)%matsize1 = nbasfcn
      zzMat(jspin)%matsize2 = neigd
715
      IF(l_real) THEN
716 717
         ALLOCATE (zzMat(jspin)%data_r(zzMat(jspin)%matsize1,
     +                                 zzMat(jspin)%matsize2))
718
      ELSE
719 720
         ALLOCATE (zzMat(jspin)%data_c(zzMat(jspin)%matsize1,
     +                                 zzMat(jspin)%matsize2))
721 722
      END IF

723 724 725 726 727 728 729
!         CALL cdn_read(
!     >              eig_id,
!     >              nvd,jspd,irank,isize,kptibz,jspin,nbasfcn,
!     >              l_ss,l_noco,neigd,n_start,n_end,
!     <              nmat,nv,ello,evdu,epar,
!     <              k1,k2,k3,bkpt,wk,nbands,eigg(:,jspin),zzMat(jspin))
      call judft_error("cdn_read in wann_updown")
730 731 732 733 734 735 736 737 738 739
         nrec=0
      enddo !jspin
      call cpu_time(delta1)
      time_rw=time_rw+delta1-delta
      nslibd = 0
      
      jspin=1

c...we work only within the energy window

740 741
      DO jspin = 1, wannierspin
         zMat(jspin)%l_real = zzMat(jspin)%l_real
742 743
         zMat(jspin)%matsize1 = zzMat(jspin)%matsize1
         zMat(jspin)%matsize2 = zzMat(jspin)%matsize2
744
         IF (zzMat(jspin)%l_real) THEN
745 746 747 748
	       IF(.not.allocated(zmat(jspin)%data_r))
     +      ALLOCATE (zMat(jspin)%data_r(zMat(jspin)%matsize1,
     +                                zMat(jspin)%matsize2))
            zMat(jspin)%data_r = 0.0
749
         ELSE
750 751 752 753
	       IF(.not.allocated(zMat(jspin)%data_c))
     +      ALLOCATE (zMat(jspin)%data_c(zMat(jspin)%matsize1,
     +                                zMat(jspin)%matsize2))
            zMat(jspin)%data_c = CMPLX(0.0,0.0)
754 755 756
         END IF
      END DO

757 758 759 760
      eig(:) = 0.

      do jspin=1,wannierspin
       nslibd=0
761
       !print*,"bands used:"
762 763 764 765 766 767 768 769 770
       do i = 1,nbands
        if ((eigg(i,1).ge.e1s.and.nslibd.lt.numbands.
     &       and.wann%l_bynumber)
     &    .or.(eigg(i,1).ge.e1s.and.eigg(i,1).le.e2s.
     &       and.wann%l_byenergy)
     &    .or.(i.ge.wann%band_min(jspin).and.
     &        (i.le.wann%band_max(jspin))
     &       .and.wann%l_byindex))then

771
           !print*,i
772 773 774 775 776 777 778 779
           nslibd = nslibd + 1
           eig(nslibd) = eigg(i,1)
           if(l_noco)then
             funbas=nv(1)+nlotot
             funbas=funbas+nv(2)+nlotot
           else
             funbas=nv(jspin)+nlotot
           endif
780 781
           IF (zzMat(jspin)%l_real) THEN
              do j = 1,funbas
782
                 zMat(jspin)%data_r(j,nslibd) = zzMat(jspin)%data_r(j,i)
783 784 785
              enddo
           ELSE
              do j = 1,funbas
786
                 zMat(jspin)%data_c(j,nslibd) = zzMat(jspin)%data_c(j,i)
787 788
              enddo
           END IF
789 790 791 792 793 794 795 796 797 798 799 800 801 802 803
        endif
       enddo
      enddo
      jspin=1

c***********************************************************
c              rotate the wavefunction
c***********************************************************
      do jspin=1,wannierspin
        if (wann%l_bzsym.and.oper.ne.1) then  !rotate bkpt
          call wann_kptsrotate(
     >            natd,nlod,llod,
     >            ntypd,nlo,llo,invsat,
     >            l_noco,l_soc,
     >            ntype,neq,nlotot,
804
     >            jspin,
805 806 807 808
     >            oper,nop,mrot,nvd,nv,
     >            shiftkpt(:,ikpt),
     >            tau,
     x            bkpt,k1(:,:),k2(:,:),k3(:,:),
809
     x            zMat(jspin),nsfactor)
810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831
        else
          nsfactor=cmplx(1.0,0.0)
        endif
      enddo !jspin 
c      print*,"bkpt1=",bkpt

c******************************************************************

c...the overlap matrix Mmn which is computed for each k- and b-point

      noccbd = nslibd

      allocate ( acof(noccbd,0:lmd,natd,wannierspin),
     &           bcof(noccbd,0:lmd,natd,wannierspin),
     &           ccof(-llod:llod,noccbd,nlod,natd,wannierspin))

      acof(:,:,:,:) = cmplx(0.,0.) ; bcof(:,:,:,:) = cmplx(0.,0.)
      ccof(:,:,:,:,:) = cmplx(0.,0.)

c...generation the A,B,C coefficients in the spheres
c...for the lapws and local orbitals, summed by the basis functions

832 833 834 835 836 837 838 839 840 841
      ALLOCATE(lapw%k1(SIZE(k1,1),SIZE(k1,2)))
      ALLOCATE(lapw%k2(SIZE(k1,1),SIZE(k1,2)))
      ALLOCATE(lapw%k3(SIZE(k1,1),SIZE(k1,2)))
      lapw%nv = nv
      lapw%nmat = nmat
      lapw%k1 = k1
      lapw%k2 = k2
      lapw%k3 = k3
      ! I think the other variables of lapw are not needed here.

842 843
      call cpu_time(delta)
      do jspin=1,wannierspin
844

845 846
         CALL abcof(input,atoms,sym,cell,lapw,noccbd,usdus,
     +              noco,jspin,oneD,acof(1:,0:,1:,jspin),
847 848 849
     +              bcof(1:,0:,1:,jspin),ccof(-llod:,1:,1:,1:,jspin),
     +              zMat(jspin))

850 851 852
         call wann_abinv(atoms,
     X        acof(:,0:,:,jspin),bcof(:,0:,:,jspin),
     <       ccof(-llod:,:,:,:,jspin))
853 854
      enddo !jspin

855
      DEALLOCATE(lapw%k1,lapw%k2,lapw%k3)
856 857 858 859 860 861 862

      jspin=1
      call cpu_time(delta1)
      time_abcof=time_abcof+delta1-delta

      if(wann%l_socmat)then
         call wann_socmat(
863
     >       mpi,enpara,input,noco,atoms,
864
     >       lmaxd,ntypd,nlod,natd,noccbd,
865
     >       llod,jmtd,jspd,nlhd,neq,
866 867
     >       ntype,theta,phi,jspins,irank,
     >       jri,lmax,dx,rmsh,epar,ello,nlo,llo,
868
     >       l_dulo,ulo_der,vTot%mt,
869
     >       acof,bcof,ccof,
870 871 872 873 874 875 876 877 878 879
     <       hsomtx_l)
        do i1 = 1,2
         do i2 = 1,2
          do i = 1,nbnd
           do j = 1,nbnd
             hsomtx(i1,i2,i,j,ikpt) = hsomtx_l(i,j,i1,i2)
           enddo
          enddo
         enddo
        enddo
880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906
      endif

#ifdef CPP_TOPO

      if(wann%l_perpmag)then
         call wann_perpmag(
     >            nslibd,neq,mlh,nlhd,nlh,ntypsy,llh,lmax,
     >            nmem,ntype,ntypd,bbmat,bmat,
     >            nlod,nlo,llo,flo,ff,gg,jri,rmsh,dx,jmtd,
     >            lmaxd,lmd,clnu,
     >            ubeffug,ubeffdg,dbeffug,dbeffdg,ubeffulog,
     >            dbeffulog,
     >            ulobeffug,ulobeffdg,ulobeffulog,
     >            acof,bcof,
     <            perpmag(:,:,ikpt))
      endif   

      if(wann%l_surfcurr)then
         call wann_surfcurr_int_updown(
     >       nv2d,jspin,odi,ods,n3d,nmzxyd,n2d,ntypsd,
     >       ntype,lmaxd,jmtd,ntypd,natd,nmzd,neq,nq3,nvac,
     >       nmz,nmzxy,nq2,nop,nop2,volint,film,slice,symor,
     >       invs,invs2,z1,delz,ngopr,ntypsy,jri,
     >       pos,taual,zatom,rmt,
     >       lmax,mrot,tau,rmsh,invtab,amat,bmat,bbmat,ikpt,
     >       nvd,nv(jspin),bkpt,omtil,
     >       k1(:,jspin),k2(:,jspin),k3(:,jspin),
907
     >       nslibd,nbasfcn,neigd,zMat,
908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961
     >       dirfacs,
     >       surfcurr(:,:,:,ikpt))
         call wann_surfcurr_updown(
     >       dirfacs,amat,
     >       ntype,lmaxd,lmax,natd,
     >       neq,noccbd,lmd,natd,llod,nlod,
     >       nlo,llo, 
     >       acof,bcof,ccof,
     >       us,dus,duds,uds,
     >       ulos,dulos,
     >       rmt,pos, 
     &       surfcurr(:,:,:,ikpt))
        write(6,*)"dirfacs=",dirfacs
      endif

      if(wann%l_nabla)then
         call wann_nabla_updown(
     >       jspin,ntype,lmaxd,lmax,natd,
     >       jmtd,jri,rmsh,dx,neq,
     >       noccbd,lmd,natd,llod,nlod, 
     >       ff(:,:,:,:,1),gg(:,:,:,:,1),
     >       ff(:,:,:,:,2),gg(:,:,:,:,2),
     >       acof(:,:,:,1),bcof(:,:,:,1),
     >       ccof(:,:,:,:,1),
     >       acof(:,:,:,2),bcof(:,:,:,2),
     >       ccof(:,:,:,:,2),
     &       nablamat(:,:,:,ikpt))

         addnoco=0
         do 41 i = n_rank+1,nv(jspin),n_size
           b1(1)=bkpt(1)+k1(i,jspin)
           b1(2)=bkpt(2)+k2(i,jspin)
           b1(3)=bkpt(3)+k3(i,jspin)
           b2(1)=b1(1)*bmat(1,1)+b1(2)*bmat(2,1)+b1(3)*bmat(3,1)
           b2(2)=b1(1)*bmat(1,2)+b1(2)*bmat(2,2)+b1(3)*bmat(3,2)
           b2(3)=b1(1)*bmat(1,3)+b1(2)*bmat(2,3)+b1(3)*bmat(3,3)
           do 42 j = n_rank+1,nv(jspin),n_size
c            b1(1)=bkpt(1)+k1(j,jspin)
c            b1(2)=bkpt(2)+k2(j,jspin)
c            b1(3)=bkpt(3)+k3(j,jspin)
c            b3(1)=b1(1)*bmat(1,1)+b1(2)*bmat(2,1)+b1(3)*bmat(3,1)
c            b3(2)=b1(1)*bmat(1,2)+b1(2)*bmat(2,2)+b1(3)*bmat(3,2)
c            b3(3)=b1(1)*bmat(1,3)+b1(2)*bmat(2,3)+b1(3)*bmat(3,3)
c-->     determine index and phase factor
            i1 = k1(j,jspin) - k1(i,jspin)
            i2 = k2(j,jspin) - k2(i,jspin)
            i3 = k3(j,jspin) - k3(i,jspin)
            in = ig(i1,i2,i3)
            if (in.eq.0) goto 42
            phase   = rgphs(i1,i2,i3)
            phasust = cmplx(phase,0.0)*ustep(in)

            do m = 1,nslibd
             do n = 1,nslibd
962 963
              do dir=1,3
               IF (zMat%l_real) THEN
964 965
                  zzConjgTemp = cmplx(zMat(1)%data_r(i+addnoco,m) *
     +                                zMat(2)%data_r(j+addnoco,n), 0.0)
966
               ELSE
967 968
                  zzConjgTemp =       zMat(1)%data_c(i+addnoco,m) *
     +                          CONJG(zMat(2)%data_c(j+addnoco,n))
969 970 971 972
               END IF
               value=phasust*zzConjgTemp
               nablamat(dir,m,n,ikpt) = nablamat(dir,m,n,ikpt) -
     +                                  value*b2(dir)
973 974 975 976 977 978 979 980 981 982 983 984 985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020
              enddo 
             enddo
            enddo

 42        continue
 41      continue
      endif   
#endif

c***************************************************
c     Calculate matrix elements of the Pauli matrix.
c***************************************************
      if(wann%l_mmn0)then

c------> interstitial contribution to updown.mmn0-matrix
       addnoco=0
       if(l_noco)then
          addnoco=nv(1)+nlotot
       endif

c$$$       do i = n_rank+1,nv(jspin),n_size
c$$$         do 22 j = n_rank+1,nv(jspin),n_size
c$$$
c$$$c-->     determine index and phase factor
c$$$            i1 = k1(j,jspin) - k1(i,jspin)
c$$$            i2 = k2(j,jspin) - k2(i,jspin)
c$$$            i3 = k3(j,jspin) - k3(i,jspin)
c$$$            in = ig(i1,i2,i3)
c$$$            if (in.eq.0) goto 22
c$$$            phase   = rgphs(i1,i2,i3)
c$$$            phasust = cmplx(phase,0.0)*ustep(in)
c$$$            do m = 1,nslibd
c$$$             do n = 1,nslibd
c$$$#if ( !defined(CPP_INVERSION) || defined(CPP_SOC) )
c$$$              mmn(m,n,ikpt) =
c$$$     =        mmn(m,n,ikpt) + 
c$$$     +          phasust*z(i,m,1)*conjg(z(j+addnoco,n,2))
c$$$#else
c$$$              mmn(m,n,ikpt) =
c$$$     =        mmn(m,n,ikpt) + 
c$$$     +          phasust*cmplx(z(i,m,1)*z(j+addnoco,n,2),0.0)
c$$$#endif
c$$$             enddo !n
c$$$            enddo !m
c$$$
c$$$  22     continue
c$$$       enddo !i  

1021
      addnoco2 = addnoco      ! TODO: correct for addnoco2
1022
      call wann_mmkb_int(
1023
     >         cmplx(1.,0.),addnoco,addnoco2,
1024 1025
     >         nvd,k1d,k2d,k3d,
     >         n3d,k1(:,jspin),k2(:,jspin),k3(:,jspin),
1026
     >         nv(jspin),neigd,nbasfcn,zMat(1),nslibd,
1027
     >         k1(:,jspin),k2(:,jspin),k3(:,jspin),
1028
     >         nv(jspin),zMat(2),nslibd,
1029 1030 1031 1032 1033 1034 1035 1036
     >         nbnd,
     >         rgphs,ustep,ig,(/ 0,0,0 /),
     <         mmn(:,:,ikpt))


c---> spherical contribution to updown.mmn0-matrix
       call wann_mmk0_updown_sph(
     >          l_noco,alph,beta,
1037
     >          llod,noccbd,nlod,natd,ntypd,lmaxd,lmax,lmd,
1038 1039 1040 1041 1042 1043 1044
     >          ntype,neq,nlo,llo,
     >          radial1_ff,radial1_gg,radial1_fg,radial1_gf,
     >          radial1_flo,radial1_glo,
     >          radial1_lof,radial1_log,
     >          radial1_lolo,
     >          acof(1:noccbd,:,:,:),
     >          bcof(1:noccbd,:,:,:),ccof(:,1:noccbd,:,:,:),
1045
     >          usdus%ddn,usdus%uulon,usdus%dulon,usdus%uloulopn,
1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073
     =          mmn(1:noccbd,1:noccbd,ikpt))

c---> vacuum contribution to updown.mmn0-matrix

       if (film .and. .not.odi%d1) then

c         call wann_mmk0_vac(
c     >           z1,nmzd,nv2d,k1d,k2d,k3d,n3d,nvac,
c     >           ig,nmz,delz,ig2,area,bmat,
c     >           bbmat,evac(:,jspin),bkpt,vz(:,:,jspin2),
c     >           nslibd,jspin,k1,k2,k3,wannierspin,nvd,
c     >           nbasfcn,neigd,z,nv,omtil,
c     <           mmn(1:noccbd,1:noccbd,ikpt))

       elseif (odi%d1) then

c          call wann_mmk0_od_vac(
c     >           z1,nmzxyd,nmzd,nv2d,k1d,k2d,k3d,n2d,n3d,
c     >           ig,nmzxy,nmz,delz,ig2,odi%n2d,
c     >           bbmat,evac(1,jspin),bkpt,odi%M,odi%mb,
c     >           vz(:,1,jspin2),odi,
c     >           nslibd,jspin,k1,k2,k3,wannierspin,nvd,area,
c     >           nbasfcn,neigd,z,nv,sk2,phi2,omtil,qss,
c     <           mmn(1:noccbd,1:noccbd,ikpt))

       endif
      endif !l_mmn0

1074
      deallocate ( acof,bcof,ccof,eigg )
1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116

      write (*,*) 'nslibd=',nslibd


      endif   ! loop by processors

10    continue ! end of cycle by the k-points

#ifdef CPP_MPI
      call MPI_BARRIER(mpi_comm,ierr)
#endif
c******************************************************
c     Write down the projections.
c******************************************************
  5   continue

      if(wann%l_nabla)then
        hescale=tpi*condquant/bohrtocm/omtil
        hescale=sqrt(hescale)
        nablamat=nablamat*hescale
        call wann_write_nabla(
     >      mpi_comm,l_p0,'nabl.updown',
     >      'Matrix elements of nabla operator',
     >      nbnd,fullnkpts,nbnd,
     >      irank,isize,
     <      nablamat)
      endif

      if(wann%l_socmat)then
        call wann_write_hsomtx(
     >      mpi_comm,l_p0,'WF1.hsomtx',
     >      'Matrix elements of spin-orbit coupling',
     >      2,2,nbnd,nbnd,fullnkpts,
     >      irank,isize,
     <      hsomtx)
      endif

      if(wann%l_mmn0)then
        call wann_write_amn(
     >      mpi_comm,l_p0,'updown.mmn0',
     >      'Overlaps of the wavefunct. at the same kpoint',
     >      nbnd,fullnkpts,nbnd,
1117 1118
     >      irank,isize,.false.,
     <      mmn,wann%l_unformatted)
1119 1120 1121 1122 1123 1124 1125
      endif !l_soc and l_mmn0  

      if(wann%l_perpmag)then
        call wann_write_amn(
     >      mpi_comm,l_p0,'updown.perpmag',
     >      'Overlaps of the wavefunct. with Beff',
     >      nbnd,fullnkpts,nbnd,
1126 1127
     >      irank,isize,.true.,
     <      perpmag,wann%l_unformatted)
1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147
      endif !l_soc and l_mmn0  

      if(wann%l_surfcurr)then
        surfcurr=conjg(surfcurr/cmplx(0.0,2.0))
        hescale=tpi*condquant/bohrtocm/omtil
        hescale=sqrt(hescale)
        surfcurr=surfcurr*hescale
        call wann_write_nabla(
     >      mpi_comm,l_p0,'surfcurr.updown',
     >      'Surface currents',
     >      nbnd,fullnkpts,nbnd,
     >      irank,isize,
     <      surfcurr)
      endif 

      if( allocated (nablamat) ) deallocate( nablamat )
      deallocate ( mmn )

      deallocate(flo)
 
1148
      deallocate ( vr,vz,nv,k1,k2,k3 )
1149
      deallocate ( ff,gg )
1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173
      if (wann%l_bzsym)  deallocate(irreduc,mapkoper)

      call cpu_time(delta)
      time_total=delta-time_total
      write(6,*)"time_total=",time_total
      write(6,*)"time_mmn=",time_mmn
      write(6,*)"time_interstitial=",time_interstitial
      write(6,*)"time_ujugaunt=",time_ujugaunt
      write(6,*)"time_abcof=",time_abcof
      write(6,*)"time_symm=",time_symm
      write(6,*)"time_rw=",time_rw
      write(6,*)"time_film=",time_film

#ifdef CPP_MPI
      call MPI_BARRIER(mpi_comm,ierr)
#endif

#ifdef CPP_MPI
      call MPI_BARRIER(mpi_comm,ierr)
#endif


      END SUBROUTINE wann_updown
      END MODULE m_wann_updown