wann_rw_eig.F 12.3 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 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 80 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 346 347 348 349 350 351 352 353 354 355 356 357 358 359
      MODULE m_wann_rw_eig
      use m_juDFT
#ifdef CPP_MPI
#define CPP_HYBEIG
#endif
c****************************************************************
c    write WF1.eig and determine maximum of number of bands
c      Frank Freimuth, October 2006
c****************************************************************
      CONTAINS
      SUBROUTINE wann_write_eig(
     >              ntapwf,
     >              lmaxd,ntypd,nlod,neigd,nvd,jspd,
     >              isize,jspin,nbasfcn,nlotot,
     >              l_ss,l_noco,nrec,fullnkpts,irecl,
     >              l_bzsym,l_byindex,l_bynumber,l_byenergy,
     >              irreduc,odi,band_min,band_max,numbands,
     >              e1s,e2s,ef,l_paulimag,nkpt,    
     <              nbnd,kpoints)

      use m_cdnread, only:cdn_read
      use m_od_types,only:od_inp

      IMPLICIT NONE
      integer,intent(in) :: ntapwf
      integer,intent(in) :: lmaxd,ntypd,nlod,neigd,nvd,jspd
      integer,intent(in) :: isize,jspin,nbasfcn,nlotot
      logical,intent(in) :: l_ss,l_noco
      integer,intent(in) :: nrec,fullnkpts
#ifndef CPP_MPI
      integer(8), intent (in) :: irecl
#else
      include 'mpif.h'
      INTEGER(KIND=MPI_OFFSET_KIND),intent(in) :: irecl
#endif
      logical,intent(in) :: l_byindex,l_bynumber,l_byenergy
      integer,intent(in) :: irreduc(fullnkpts)
      type (od_inp),intent(in) :: odi
      integer,intent(in) :: band_min,band_max,numbands
      logical,intent(in) :: l_bzsym
      real,intent(in)    :: e1s,e2s,ef
      logical,intent(in) :: l_paulimag
      integer,intent(in) :: nkpt

      integer,intent(out):: nbnd
      real,intent(out)   :: kpoints(fullnkpts)

      logical          :: l_eig
      character(len=3) :: spin12(2)
      data spin12/'WF1' , 'WF2'/
      integer          :: ikpt,kptibz
      integer          :: nmat,nbands,nv(jspd)
      real             :: wk, bkpt(3),eig(neigd),cp_time(9)
      integer          :: k1(nvd,jspd),k2(nvd,jspd),k3(nvd,jspd)
      integer          :: nkbnd,i
      REAL             :: ello(nlod,ntypd,jspd),evdu(2,jspd)
      REAL             :: epar(0:lmaxd,ntypd,jspd)
      INTEGER          :: kveclo(nlotot)
      integer          :: n_start,n_end,co
      integer          :: bnd1,bnd2,kpt
      complex, allocatable :: paulimat(:,:,:,:) 
      real             :: sum1,sum2,sum3
      integer          :: num_bands,err
      real,parameter   :: hartree=27.21138505

#if ( !defined(CPP_INVERSION) || defined(CPP_SOC) )
      COMPLEX  :: z(nbasfcn,neigd)
#else
      REAL     :: z(nbasfcn,neigd)
#endif

      n_start=1
      n_end=neigd

      inquire (file=spin12(jspin)//'.eig',exist=l_eig)
      if (l_eig) then
         open (306,file=spin12(jspin)//'.eig',
     &                               form='formatted',status='old')
         rewind (306)
      else
         open (306,file=spin12(jspin)//'.eig',
     &                               form='formatted',status='new')
      endif!l_eig

  
      write(*,*)"nkpt=",nkpt
      if(l_paulimag.and..false.)then
        num_bands=band_max-band_min+1
        allocate( paulimat(num_bands,num_bands,3,nkpt),stat=err )
        IF (err/=0)  CALL juDFT_error
     +       ("error allocating paulimat",calledby
     +       ="wann_rw_eig")
        open(655,file='wpmat')
        do ikpt=1,nkpt
          read(655,*)
          do co=1,3
            read(655,*)
            do bnd2=1,num_bands
               do bnd1=1,num_bands
                  read(655,*)paulimat(bnd1,bnd2,co,ikpt)
               enddo
            enddo
          enddo
        enddo      
        close(655)
      endif
      sum1=0.0
      sum2=0.0
      sum3=0.0

      nbnd=0
      do ikpt = 1,fullnkpts

        kptibz=ikpt 
        if(l_bzsym) kptibz=irreduc(ikpt)
            call wann_read_eig(
     >              ntapwf,     
     >              lmaxd,ntypd,nlod,neigd,nvd,jspd,
     >              0,isize,kptibz,jspin,nbasfcn,nlotot,
     >              l_ss,l_noco,nrec,irecl,
     <              nmat,nv,ello,evdu,epar,kveclo,
     <              k1,k2,k3,bkpt,wk,nbands,eig,z,cp_time)


        nkbnd = 0
        do i = 1,nbands
         if((eig(i).ge.e1s .and. nkbnd.lt.numbands.and.l_bynumber).or.
     &      (eig(i).ge.e1s.and.eig(i).le.e2s.and.l_byenergy ).or.
     &      (i.ge.band_min.and.i.le.band_max.and.l_byindex))  then
           nkbnd = nkbnd + 1
           if(l_bynumber.or.l_byindex)then
              write (306,'(2i12,f19.13)') nkbnd,ikpt,(eig(i)-ef)*
     &                                    hartree
           endif
           if(l_paulimag.and..false.)then
              if(eig(i).le.ef)then
                 sum1=sum1+paulimat(nkbnd,nkbnd,1,kptibz)
                 sum2=sum2+paulimat(nkbnd,nkbnd,2,kptibz)
                 sum3=sum3+paulimat(nkbnd,nkbnd,3,kptibz)
              endif
           endif   
         endif 
        enddo 
        if (nkbnd.ge.nbnd) nbnd = nkbnd

        if (odi%d1) kpoints(ikpt) = bkpt(3)

      enddo !ikpt

      if(l_paulimag)then
        write(6,*)"sum1=",sum1/fullnkpts
        write(6,*)"sum2=",sum2/fullnkpts
        write(6,*)"sum3=",sum3/fullnkpts
      endif  


      if(l_byenergy)then   !now we know the maximum of bands
         do ikpt=1,fullnkpts
            kptibz=ikpt
            if(l_bzsym)kptibz=irreduc(ikpt)
            call wann_read_eig(     
     >              ntapwf,       
     >              lmaxd,ntypd,nlod,neigd,nvd,jspd,
     >              0,isize,kptibz,jspin,nbasfcn,nlotot,
     >              l_ss,l_noco,nrec,irecl,
     <              nmat,nv,ello,evdu,epar,kveclo,
     <              k1,k2,k3,bkpt,wk,nbands,eig,z,cp_time)


              nkbnd = 0
              do i = 1,nbands
                if(eig(i).ge.e1s .and. nkbnd.lt.nbnd)then
                     nkbnd = nkbnd + 1
                  write (306,'(2i12,f19.13)')nkbnd,ikpt,(eig(i)-ef)*
     &                                       hartree
                endif
              enddo 
         enddo   
      endif
      close (306)
      end subroutine wann_write_eig
c*********************************************************************
c     each process writes its own eig file if hdf is not used for mpi
c     this routine finds the appropriate eig file for a given kpoint
c     Frank Freimuth, October 2006
c*********************************************************************
      subroutine wann_mpi_read_eig(     
     >              lmaxd,ntypd,nlod,neigd,nvd,jspd,
     >              irank,isize,kptibz,jspin,nbasfcn,nlotot,
     >              l_ss,l_noco,nrec,irecl,
     <              nmat,nv,ello,evdu,epar,kveclo,
     <              k1,k2,k3,bkpt,wk,nbands,eig,z,cp_time)

      use m_cdnread, only:cdn_read
      implicit none
      INTEGER, INTENT (IN) :: irank,isize,kptibz,nbasfcn,neigd,nlotot
      INTEGER, INTENT (IN) :: nrec,nvd,jspd,jspin
      INTEGER, INTENT (IN) :: lmaxd,ntypd,nlod
#ifndef CPP_MPI
      integer(8), intent (in) :: irecl
#else
      include 'mpif.h'
      INTEGER(KIND=MPI_OFFSET_KIND),intent(in) :: irecl
#endif
      LOGICAL, INTENT (IN) :: l_ss,l_noco
      INTEGER, INTENT (OUT) :: nbands,nmat
      REAL,    INTENT (OUT) :: wk

      INTEGER, INTENT (OUT) :: k1(nvd,jspd),k2(nvd,jspd),k3(nvd,jspd)
      INTEGER, INTENT (OUT) :: nv(jspd),kveclo(nlotot)
      REAL,    INTENT (OUT) :: bkpt(3),eig(neigd)
      REAL,    INTENT (INOUT) :: cp_time(9)
      REAL,    INTENT (OUT) :: ello(nlod,ntypd,jspd),evdu(2,jspd)
      REAL,    INTENT (OUT) :: epar(0:lmaxd,ntypd,jspd)

#if ( !defined(CPP_INVERSION) || defined(CPP_SOC) )
      COMPLEX, INTENT (OUT) :: z(nbasfcn,neigd)
#else
      REAL,    INTENT (OUT) :: z(nbasfcn,neigd)
#endif

      character(len=12) :: fname
      integer :: aoff,d1,d10
      integer :: cpu_index
      logical :: l_eig
      integer :: n_start,n_end

      n_start=1
      n_end=neigd


      do cpu_index=0,isize-1

        if (mod(kptibz-1,isize).eq.cpu_index) then



        if(cpu_index.eq.irank)then

           call cdn_read(
     >              lmaxd,ntypd,nlod,neigd,nvd,jspd,
     >              irank,isize,kptibz,jspin,nbasfcn,nlotot,
     >              l_ss,l_noco,nrec,kptibz,66,neigd,
     >              n_start,n_end,      
     <              nmat,nv,ello,evdu,epar,kveclo,
     <              k1,k2,k3,bkpt,wk,nbands,eig,z,cp_time)

        else !eigenvalues are stored in the file of another process

           d10=int((cpu_index+0.5)/10)
           d1=mod(cpu_index,10)
           aoff=iachar('1')-1
           if(d10.ge.10)aoff=iachar('7')
           fname='eig'//achar(d10+aoff)//achar(d1+iachar('1')-1)
           l_eig=.false.
           inquire(file=fname,exist=l_eig)
           IF(.NOT.l_eig)  CALL juDFT_error("one of eig files not found"
     +          ,calledby ="wann_rw_eig")
           OPEN(666,file=fname,access='direct',form='unformatted',
     &                            recl=irecl,status='old')
           call cdn_read(
     >              lmaxd,ntypd,nlod,neigd,nvd,jspd,
     >              cpu_index,isize,kptibz,jspin,nbasfcn,nlotot,
     >              l_ss,l_noco,nrec,kptibz,666,neigd,
     >              n_start,n_end,
     <              nmat,nv,ello,evdu,epar,kveclo,
     <              k1,k2,k3,bkpt,wk,nbands,eig,z,cp_time)

           close(666)
        endif   !cpu_index==irank
        endif !find process
      enddo !cpu_index


      end subroutine wann_mpi_read_eig  
c*****************************************************************
c     read in eig file
c*****************************************************************
      subroutine wann_read_eig(  
     >              ntapwf,
     >              lmaxd,ntypd,nlod,neigd,nvd,jspd,
     >              irank,isize,kptibz,jspin,nbasfcn,nlotot,
     >              l_ss,l_noco,nrec,irecl,
     <              nmat,nv,ello,evdu,epar,kveclo,
     <              k1,k2,k3,bkpt,wk,nbands,eig,z,cp_time)

      use m_cdnread, only:cdn_read
      implicit none
      integer, intent (in) :: ntapwf
      INTEGER, INTENT (IN) :: irank,isize,kptibz,nbasfcn,neigd,nlotot
      INTEGER, INTENT (IN) :: nrec,nvd,jspd,jspin
      INTEGER, INTENT (IN) :: lmaxd,ntypd,nlod
#ifndef CPP_MPI
      integer(8), intent (in) :: irecl
#else
      include 'mpif.h'
      INTEGER(KIND=MPI_OFFSET_KIND),intent(in) :: irecl
#endif
      LOGICAL, INTENT (IN) :: l_ss,l_noco
      INTEGER, INTENT (OUT) :: nbands,nmat
      REAL,    INTENT (OUT) :: wk

      INTEGER, INTENT (OUT) :: k1(nvd,jspd),k2(nvd,jspd),k3(nvd,jspd)
      INTEGER, INTENT (OUT) :: nv(jspd),kveclo(nlotot)
      REAL,    INTENT (OUT) :: bkpt(3),eig(neigd)
      REAL,    INTENT (INOUT) :: cp_time(9)
      REAL,    INTENT (OUT) :: ello(nlod,ntypd,jspd),evdu(2,jspd)
      REAL,    INTENT (OUT) :: epar(0:lmaxd,ntypd,jspd)

#if ( !defined(CPP_INVERSION) || defined(CPP_SOC) )
      COMPLEX, INTENT (OUT) :: z(nbasfcn,neigd)
#else
      REAL,    INTENT (OUT) :: z(nbasfcn,neigd)
#endif

      integer :: n_start,n_end

      n_start=1
      n_end=neigd

#if( defined(CPP_MPI) && !defined(CPP_HDF) )
! MPI without HDF
#ifdef CPP_HYBEIG
c-----if hybrid functionals are used, and mpi is defined, mpi-io is
c-----used and all processors use the same eig-file
        call cdn_read(
     >          lmaxd,ntypd,nlod,neigd,nvd,jspd,
     >          irank,isize,kptibz,jspin,nbasfcn,nlotot,
     >          l_ss,l_noco,nrec,kptibz,ntapwf,neigd,
     >          n_start,n_end,      
     <          nmat,nv,ello,evdu,epar,kveclo,
     <          k1,k2,k3,bkpt,wk,nbands,eig,z,cp_time)
#else
        call wann_mpi_read_eig(     
     >              lmaxd,ntypd,nlod,neigd,nvd,jspd,
     >              irank,isize,kptibz,jspin,nbasfcn,nlotot,
     >              l_ss,l_noco,nrec,irecl,
     <              nmat,nv,ello,evdu,epar,kveclo,
     <              k1,k2,k3,bkpt,wk,nbands,eig,z,cp_time)
#endif
#else
!MPI with HDF or no MPI
        call cdn_read(
     >              lmaxd,ntypd,nlod,neigd,nvd,jspd,
     >              0,isize,kptibz,jspin,nbasfcn,nlotot,
     >              l_ss,l_noco,nrec,kptibz,66,
     >              neigd,n_start,n_end,
     <              nmat,nv,ello,evdu,epar,kveclo,
     <              k1,k2,k3,bkpt,wk,nbands,eig,z,cp_time)

#endif
      END SUBROUTINE wann_read_eig  
      END MODULE m_wann_rw_eig