wann_projgen.F 5.73 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
      module m_wann_projgen
      contains
      subroutine wann_projgen(
     >               ntype,neq,natd,zatom,l_nocosoc,wann)
c*****************************************
c     Generate the proj-file.
c     Frank Freimuth
c*****************************************
15

16
      use m_types, only: t_wann
17 18
      USE m_constants

19
      implicit none
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
      integer, intent(in) :: ntype
      integer, intent(in) :: neq(ntype)
      integer, intent(in) :: natd
      real,intent(in)     :: zatom(ntype)
      logical,intent(in)  :: l_nocosoc
      type(t_wann), intent(in) :: wann

      integer             :: j,nat,n,at,k,num_proj
      integer             :: num_wann,kk
      character(len=2)    :: symbol(100*ntype)
      character(len=2)    :: symb(natd)
      integer             :: lwf(100*ntype)
      integer             :: mrwf(100*ntype)
      integer             :: rwf(100*ntype)
      real                :: alpha(100*ntype),beta(100*ntype)
      real                :: gamma(100*ntype),zona(100*ntype)
      real                :: regio(100*ntype)      
      character*2         :: namat(0:103)
      integer             :: projections(-5:5)
      integer             :: projspin

      DATA namat/'va',' H','He','Li','Be',' B',' C',' N',' O',' F','Ne',
     +     'Na','Mg','Al','Si',' P',' S','Cl','Ar',' K','Ca','Sc','Ti',
     +     ' V','Cr','Mn','Fe','Co','Ni','Cu','Zn','Ga','Ge','As','Se',
     +     'Br','Kr','Rb','Sr',' Y','Zr','Nb','Mo','Tc','Ru','Rh','Pd',
     +     'Ag','Cd','In','Sn','Sb','Te',' I','Xe','Cs','Ba','La','Ce',
     +     'Pr','Nd','Pm','Sm','Eu','Gd','Tb','Dy','Ho','Er','Tm','Yb',
     +     'Lu','Hf','Ta',' W','Re','Os','Ir','Pt','Au','Hg','Tl','Pb',
     +     'Bi','Po','At','Rn','Fr','Ra','Ac','Th','Pa',' U','Np','Pu',
     +     'Am','Cm','Bk','Cf','Es','Fm','Md','No','Lw'/

      projections(:)=0
      projections(-5)=6
      projections(-4)=5
      projections(-3)=4
      projections(-2)=3
      projections(-1)=2
      projections(0)=1
      projections(1)=3
      projections(2)=5
      projections(3)=7


      nat=0
      do n=1,ntype
         do at=1,neq(n)
            nat=nat+1
            symb(nat)=namat(nint(zatom(n)))
         enddo
      enddo

      open(200,file='projgen_inp',status='old')
      j=0
      do 
         j=j+1
         read(200,*,end=299,err=299)symbol(j),lwf(j),mrwf(j),rwf(j)
         num_proj=j
      enddo   
 299  continue
      print*,"found ",num_proj,"projections"
81
      write(oUnit,*)"found ",num_proj,"projections"
82 83 84 85 86 87 88 89 90 91 92 93
      close(200)

      open(300,file='proj')
      alpha=0.0
      beta=0.0
      gamma=0.0
      zona=0.0
      regio=1.0
      num_wann=0
      do j=1,nat
         do k=1,num_proj
            print*,"symbols:",symbol(k)," ",symb(j)
94
            write(oUnit,*)"symbols:",symbol(k)," ",symb(j)
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
            if(trim(adjustl(symbol(k))).eq.
     &             trim(adjustl(symb(j))))then
               if(mrwf(k).ne.0)then
                  num_wann=num_wann+1
               else
                  num_wann=num_wann+projections(lwf(k))
               endif
            endif
         enddo
      enddo
      if(l_nocosoc)num_wann=num_wann*2
      print*,num_wann," wannier functions have been defined"
      if(l_nocosoc)then
         write(300,*) num_wann, 
     &                MAX(wann%band_max(1)-wann%band_min(1)+1,num_wann),
     &                "  t  "
      else
         write(300,*) num_wann, 
     &                MAX(wann%band_max(1)-wann%band_min(1)+1,num_wann)
      endif
      if(l_nocosoc)then
       do projspin=1,-1,-2  
        do j=1,nat
         do k=1,num_proj
            if(
     &          trim(adjustl(symbol(k)))
     &                    .eq.
     &          trim(adjustl(symb(j)))
     &                    )then
               if(mrwf(k).ne.0)then
                 write (300,'(i3,1x,i2,1x,i2,1x,i2,1x,i2)') 
     &                j,lwf(k),mrwf(k),rwf(k),projspin
                 write (300,'(2x,4f10.6,1x,f4.2)') 
     &                alpha(k),beta(k),gamma(k),
     &                zona(k),regio(k)
               else
                  do kk=1,projections(lwf(k))
                     write (300,'(i3,1x,i2,1x,i2,1x,i2,1x,i2)') 
     &                j,lwf(k),kk,rwf(k),projspin
                     write (300,'(2x,4f10.6,1x,f4.2)') 
     &                alpha(k),beta(k),gamma(k),
     &                zona(k),regio(k)
                  enddo
               endif
            endif
         enddo
        enddo
       enddo 
      else
        do j=1,nat
         do k=1,num_proj
            if(
     &          trim(adjustl(symbol(k)))
     &                    .eq.
     &          trim(adjustl(symb(j)))
     &                    )then
               if(mrwf(k).ne.0)then
                 write (300,'(i3,1x,i2,1x,i2,1x,i2)') 
     &                j,lwf(k),mrwf(k),rwf(k)
                 write (300,'(2x,4f10.6,1x,f4.2)') 
     &                alpha(k),beta(k),gamma(k),
     &                zona(k),regio(k)
               else
                  do kk=1,projections(lwf(k))
                     write (300,'(i3,1x,i2,1x,i2,1x,i2)') 
     &                j,lwf(k),kk,rwf(k)
                     write (300,'(3x,4f10.6,1x,f4.2)') 
     &                alpha(k),beta(k),gamma(k),
     &                zona(k),regio(k)
                  enddo
               endif
            endif
         enddo
        enddo
      endif
      close(300)

      end subroutine wann_projgen
      end module m_wann_projgen