bandstr1.F 14.4 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
      MODULE m_bandstr1
      use m_juDFT
9
      USE m_types
10 11 12 13
!----------------------------------------------------------------------
!----------------------------------------------------------------------
      CONTAINS
      SUBROUTINE bandstr1(
14 15
     >                    idsyst,idtype,bmat,kpts,input,l_fillArrays)

16 17
      IMPLICIT NONE

18
      TYPE(t_input),INTENT(IN) :: input
19 20 21
      TYPE(t_kpts),INTENT(INOUT)  :: kpts

      INTEGER, INTENT (IN)  :: idsyst,idtype
22
      REAL,    INTENT (IN)  :: bmat(3,3)
23
      LOGICAL, INTENT (IN)  :: l_fillArrays
24 25 26 27

      REAL, POINTER             :: syp(:,:) 
      CHARACTER(len=1), POINTER :: ssy(:)
      REAL,    ALLOCATABLE      :: rsyp(:,:),del(:),d(:)
28
      REAL,    ALLOCATABLE      :: bkTemp(:,:)
29 30 31 32 33
      INTEGER, ALLOCATABLE      :: nk(:)
      INTEGER nosyp,i,j,n,ntot
      REAL dk(3),syp1(3),rsyp1(3)
      
      CALL get_points(
34
     >                idsyst,idtype,bmat,kpts,
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
     <                nosyp,syp,ssy)
!
! --> transform to cartesian units
!
      ALLOCATE ( rsyp(3,nosyp),del(nosyp),d(nosyp),nk(nosyp) ) 
      DO i = 1,nosyp
        syp1(:) = syp(:,i)
        !CALL cotra3(syp1,rsyp1,bmat)
        rsyp1=matmul(syp1,bmat)
        rsyp(:,i) = rsyp1(:)
      ENDDO
!
! --> calculate length between points
!
      d(1) = 0.0
      DO i = 2,nosyp
        del(i) = ( rsyp(1,i) - rsyp(1,i-1) )**2 +
     +           ( rsyp(2,i) - rsyp(2,i-1) )**2 +
     +           ( rsyp(3,i) - rsyp(3,i-1) )**2
        del(i) = sqrt(del(i))
        d(i) = d(i-1) + del(i)
      ENDDO
!
! --> distibute points evenly
!
      ntot = nosyp
      DO i = 2,nosyp
62
        nk(i) = NINT ( (kpts%nkpt-nosyp)*( del(i) / d(nosyp) ) )
63
        ntot = ntot + nk(i)
64 65 66 67
      ENDDO
      kpts%nkpt = ntot

      ALLOCATE (bkTemp(3,kpts%nkpt))
68 69 70 71 72 73 74 75
!
! --> generate k-points mesh 
!
      n = 1
      DO i = 2,nosyp
         dk(1) = ( syp(1,i) - syp(1,i-1) ) / (nk(i)+1)
         dk(2) = ( syp(2,i) - syp(2,i-1) ) / (nk(i)+1)
         dk(3) = ( syp(3,i) - syp(3,i-1) ) / (nk(i)+1)
76
         bkTemp(:,n) = syp(:,i-1)
77 78
         n = n + 1
         DO j = 1, nk(i)
79
            bkTemp(:,n) = bkTemp(:,n-1) + dk(:)
80 81 82
            n = n + 1
         ENDDO
      ENDDO
83
      bkTemp(:,n) = syp(:,nosyp)
Daniel Wortmann's avatar
Daniel Wortmann committed
84
      kpts%nkpt = kpts%nkpt
85
      kpts%posScale = 1.0
86

87 88 89 90
      IF(l_fillArrays) THEN
         IF(ALLOCATED(kpts%bk)) THEN
            DEALLOCATE(kpts%bk)
         END IF
Daniel Wortmann's avatar
Daniel Wortmann committed
91 92
         IF(ALLOCATED(kpts%wtkpt)) THEN
            DEALLOCATE(kpts%wtkpt)
93
         END IF
Daniel Wortmann's avatar
Daniel Wortmann committed
94
         ALLOCATE (kpts%bk(3,kpts%nkpt), kpts%wtkpt(kpts%nkpt))
95
         kpts%bk(:,:) = bkTemp(:,:)
Daniel Wortmann's avatar
Daniel Wortmann committed
96
         kpts%wtkpt = 1.0
97
      ELSE
98 99 100 101 102 103 104 105 106 107 108 109 110
         OPEN (41,file='kpts',form='formatted',status='new')
         IF (.NOT.input%film) THEN
            WRITE(41,'(i5,f20.10)') kpts%nkpt,kpts%posScale
            DO n = 1, kpts%nkpt
               WRITE (41,'(4f10.5)') bkTemp(:,n),1.0
            ENDDO
         ELSE
           WRITE(41,'(i5,f20.10,3x,l1)') kpts%nkpt,kpts%posScale,.false.
           DO n = 1, kpts%nkpt
              WRITE (41,'(3f10.5)') bkTemp(1:2,n),1.0
           ENDDO
         END IF
         CLOSE (41)
111 112 113
      END IF

      CALL write_gnu(
114
     >               nosyp,d,ssy,input)
115

116
      DEALLOCATE ( rsyp,syp,del,nk,ssy,d,bkTemp )
117 118 119 120 121 122 123

      END SUBROUTINE bandstr1
!----------------------------------------------------------------------
! once the file "bands.1" and "bands.2" are created, activate with:
! gnuplot < band.gnu > band.ps
!----------------------------------------------------------------------
      SUBROUTINE write_gnu(
124
     >                     nosyp,d,ssy,input)
125 126
!
      IMPLICIT NONE
127 128 129

      TYPE(t_input),INTENT(IN) :: input
      INTEGER, INTENT (IN) :: nosyp
130 131 132 133
      REAL,    INTENT (IN) :: d(nosyp)
      CHARACTER(len=1), INTENT (IN) :: ssy(nosyp)
      
      INTEGER n,aoff,adel
134
      CHARACTER(LEN=200) tempTitle
135 136
      aoff = iachar('a')-1
      adel = iachar('a')-iachar('A')
Daniel Wortmann's avatar
Daniel Wortmann committed
137
      !write(*,*) aoff,adel 
138 139 140 141 142 143

      OPEN (27,file='band.gnu',status='unknown')
      WRITE (27,900)
      WRITE (27,901)
      WRITE (27,902)
      WRITE (27,903)
144 145 146 147 148
      WRITE(tempTitle,'(10a)') input%comment
      IF(TRIM(ADJUSTL(tempTitle)).EQ.'') THEN
         tempTitle = "Fleur Bandstructure"
      END IF
      WRITE (27,904) TRIM(ADJUSTL(tempTitle))
149 150 151 152 153 154 155 156
      DO n = 1, nosyp
        WRITE (27,905) d(n),d(n)
      ENDDO
      WRITE (27,906) d(1),d(nosyp)
!
! nomal labels
!
      IF (iachar(ssy(1)) < aoff ) THEN
157
        WRITE (27,907) ssy(1),d(1),achar(92)
158
      ELSE
159
        WRITE (27,907) " ",d(1),achar(92)
160 161 162
      ENDIF
      DO n = 2, nosyp-1
        IF (iachar(ssy(n)) < aoff ) THEN 
163
          WRITE (27,908) ssy(n),d(n),achar(92)
164
        ELSE
165
          WRITE (27,908) " ",d(n),achar(92)
166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184
        ENDIF
      ENDDO
      IF (iachar(ssy(nosyp)) < aoff ) THEN
        WRITE (27,909) ssy(nosyp),d(nosyp)
      ELSE
        WRITE (27,909) " ",d(nosyp)
      ENDIF
!
! greek labels
!
      DO n = 1, nosyp
        IF (iachar(ssy(n)) > aoff ) THEN
          WRITE (27,914) achar(iachar(ssy(n))-adel),d(n)
        ENDIF
      ENDDO
!
! now write the rest
!
      WRITE (27,910) 
185 186
      WRITE (27,911) d(nosyp)+0.00001,achar(92)
      IF (input%jspins == 2) WRITE (27,912) achar(92)
187 188 189 190 191 192 193
      WRITE (27,913)
      CLOSE (27)

 900  FORMAT ('set terminal postscript enhanced "Times-Roman" 20')
 901  FORMAT ('set xlabel ""')
 902  FORMAT ('set ylabel "E - E_F (eV)"')
 903  FORMAT ('set nokey')
194
 904  FORMAT ('set title "',a,'"')
195 196
 905  FORMAT ('set arrow from',f9.5,', -9.0 to',f9.5,',  5.0 nohead')
 906  FORMAT ('set arrow from',f9.5,', 0.0 to',f9.5,', 0.0 nohead lt 3')
197 198
 907  FORMAT ('set xtics ("',a1,'"',f9.5,', ',a)
 908  FORMAT ('           "',a1,'"',f9.5,', ',a)
199 200
 909  FORMAT ('           "',a1,'"',f9.5,'  )')
 910  FORMAT ('set ytics -8,2,4')
201 202
 911  FORMAT ('plot [0:',f9.5,'] [-9:5] ',a)
 912  FORMAT ('"bands.2" using 1:($2+0.00) w p pt 12 ps 0.5,',a)
203 204 205 206 207 208 209 210 211 212
 913  FORMAT ('"bands.1" using 1:($2+0.00)  w p pt  7 ps 0.5')
 914  FORMAT ('set label "',a1,'" at ',f9.5,
     +        ', -9.65 center font "Symbol,20"')

      END SUBROUTINE write_gnu
!----------------------------------------------------------------------
! given a bravais-lattice, determine <nosyp> symmetry points (syp)
! and their names (lowercase = greek)
!----------------------------------------------------------------------
      SUBROUTINE get_points(
213
     >                      idsyst,idtype,bmat,kpts,
214 215 216
     <                      nosyp,syp,ssy)

      IMPLICIT NONE
217 218 219 220 221 222

      TYPE(t_kpts),INTENT(IN)   :: kpts
      INTEGER,     INTENT (IN)  :: idsyst,idtype
      REAL,        INTENT (IN)  :: bmat(3,3)
      INTEGER,     INTENT (OUT) :: nosyp
      REAL,             POINTER :: syp(:,:)  ! actually intent out
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
      CHARACTER(len=1), POINTER :: ssy(:)

      LOGICAL         :: band_inp_file
      INTEGER         :: n
      REAL, PARAMETER :: f12 = 1./2., f14 = 1./4., zro = 0.0
      REAL, PARAMETER :: f34 = 3./4., f38 = 3./8., one = 1.0
      REAL, PARAMETER :: f13 = 1./3., f23 = 2./3.

      !Check if band_inp file exists for lines
      INQUIRE(file ="band_inp",exist= band_inp_file)
      IF (band_inp_file) THEN
         OPEN(99,file ="band_inp",status ="old")
         nosyp=0
         !count the lines
         DO
            READ(99,*,END = 100)
            nosyp=nosyp+1
         ENDDO
 100     REWIND 99
         ALLOCATE(syp(3,nosyp),ssy(nosyp) )
         DO n = 1,nosyp
            READ(99,*,err = 110,END=110) ssy(n),syp(:,n)
         ENDDO
         CLOSE(99)
         RETURN
 110     WRITE(*,*) "Error reading band_inp file"
          CALL juDFT_error("Bandstr1",calledby="bandstr1")
      ENDIF

252 253 254 255 256 257 258 259 260 261
      IF(kpts%numSpecialPoints.GE.2) THEN
         nosyp = kpts%numSpecialPoints
         ALLOCATE(syp(3,nosyp),ssy(nosyp) )
         DO n = 1, nosyp
            ssy(n) = TRIM(ADJUSTL(kpts%specialPointNames(n)))
            syp(:,n) = kpts%specialPoints(:,n)
         END DO
         RETURN
      END IF

262
      !No band_inp file -> determine default lines for bandstructure
Daniel Wortmann's avatar
Daniel Wortmann committed
263
      !write(*,*) idsyst,idtype
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 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
      nosyp = -1
      IF ( (idsyst == 1).AND.(idtype ==  3) ) THEN       ! fcc
        nosyp = 7
        ALLOCATE ( syp(3,nosyp),ssy(nosyp) )
        syp(:,1) = (/f12,f12,one/)  ; ssy(1) = "X"
        syp(:,2) = (/f38,f38,f34/)  ; ssy(2) = "K"
        syp(:,3) = (/zro,zro,zro/)  ; ssy(3) = "g"
        syp(:,4) = (/f12,f12,f12/)  ; ssy(4) = "L"
        syp(:,5) = (/f12,f14,f34/)  ; ssy(5) = "W"
        syp(:,6) = (/f12,zro,f12/)  ; ssy(6) = "X"
        syp(:,7) = (/zro,zro,zro/)  ; ssy(7) = "g"
      ENDIF
      IF ( (idsyst == 5).AND.(idtype ==  1) ) THEN       ! rhombohedric (trigonal)
        nosyp = 8
        ALLOCATE ( syp(3,nosyp),ssy(nosyp) )
        syp(:,1) = (/f12,f12, f12/)  ; ssy(1) = "Z"
        syp(:,2) = (/zro,zro, zro/)  ; ssy(2) = "g"
        syp(:,3) = (/f14,f14,-f14/)  ; ssy(3) = "K"
        syp(:,4) = (/f12,f12,-f12/)  ; ssy(4) = "Z"
        syp(:,5) = (/f14,f12,-f14/)  ; ssy(5) = "W"
        syp(:,6) = (/zro,f12, zro/)  ; ssy(6) = "L"
        syp(:,7) = (/zro,zro, zro/)  ; ssy(7) = "g"
        syp(:,8) = (/f12,f12, zro/)  ; ssy(8) = "F"
      ENDIF
      IF ( (idsyst == 4).AND.(idtype ==  1) ) THEN       ! hexagonal
        nosyp = 8
        ALLOCATE ( syp(3,nosyp),ssy(nosyp) )
        IF (bmat(1,1)*bmat(2,1)+bmat(1,2)*bmat(2,2) > 0.0) THEN
          syp(:,1) = (/zro,zro, zro/)  ; ssy(1) = "g"
          syp(:,2) = (/zro,f12, zro/)  ; ssy(2) = "M"
          syp(:,3) = (/f13,f13, zro/)  ; ssy(3) = "K"
          syp(:,4) = (/zro,zro, zro/)  ; ssy(4) = "g"
          syp(:,5) = (/zro,zro, f12/)  ; ssy(5) = "A"
          syp(:,6) = (/zro,f12, f12/)  ; ssy(6) = "L"
          syp(:,7) = (/f13,f13, f12/)  ; ssy(7) = "H"
          syp(:,8) = (/zro,zro, f12/)  ; ssy(8) = "A"
        ELSE                                             ! hexagonal (angle = 60)
          syp(:,1) = (/zro,zro, zro/)  ; ssy(1) = "g"
          syp(:,2) = (/f12,f12, zro/)  ; ssy(2) = "M"
          syp(:,3) = (/f13,f23, zro/)  ; ssy(3) = "K"
          syp(:,4) = (/zro,zro, zro/)  ; ssy(4) = "g"
          syp(:,5) = (/zro,zro, f12/)  ; ssy(5) = "A"
          syp(:,6) = (/f12,f12, f12/)  ; ssy(6) = "L"
          syp(:,7) = (/f13,f23, f12/)  ; ssy(7) = "H"
          syp(:,8) = (/zro,zro, f12/)  ; ssy(8) = "A"
        ENDIF
      ENDIF
      IF ( (idsyst == 1).AND.(idtype ==  1) ) THEN       ! simple cubic
        nosyp = 8
        ALLOCATE ( syp(3,nosyp),ssy(nosyp) )
        syp(:,1) = (/f12,f12, zro/)  ; ssy(1) = "M"
        syp(:,2) = (/zro,zro, zro/)  ; ssy(2) = "g"
        syp(:,3) = (/f12,zro, zro/)  ; ssy(3) = "X"
        syp(:,4) = (/f12,f12, zro/)  ; ssy(4) = "M"
        syp(:,5) = (/f12,f12, f12/)  ; ssy(5) = "R"
        syp(:,6) = (/f12,zro, zro/)  ; ssy(6) = "X"
        syp(:,7) = (/zro,zro, zro/)  ; ssy(7) = "g"
        syp(:,8) = (/f12,f12, f12/)  ; ssy(8) = "R"
      ENDIF
      IF ( (idsyst == 1).AND.(idtype ==  2) ) THEN       ! body centered cubic
        nosyp = 6
        ALLOCATE ( syp(3,nosyp),ssy(nosyp) )
        syp(:,1) = (/zro,zro, zro/)  ; ssy(1) = "g"
        syp(:,2) = (/f12,-f12,f12/)  ; ssy(2) = "H"
        syp(:,3) = (/zro,zro, f12/)  ; ssy(3) = "N"
        syp(:,4) = (/f14,f14, f14/)  ; ssy(4) = "P"
        syp(:,5) = (/zro,zro, zro/)  ; ssy(5) = "g"
        syp(:,6) = (/zro,zro, f12/)  ; ssy(6) = "N"
      ENDIF
      IF ( (idsyst == 2).AND.(idtype ==  2) ) THEN       ! body centered tetragonal (a > c)
        nosyp = 8
        ALLOCATE ( syp(3,nosyp),ssy(nosyp) )
        syp(:,1) = (/f12,f12,-f12/)  ; ssy(1) = "Z"    ! via Lambda and V
        syp(:,2) = (/zro,zro, zro/)  ; ssy(2) = "g"    ! via Sigma
        syp(:,3) = (/-f12,f12,f12/)  ; ssy(3) = "Z"    ! via Y
        syp(:,4) = (/zro,zro, f12/)  ; ssy(4) = "X"    ! via Delta
        syp(:,5) = (/zro,zro, zro/)  ; ssy(5) = "g"    
        syp(:,6) = (/zro,f12, zro/)  ; ssy(6) = "N"    ! via Q
        syp(:,7) = (/f14,f14, f14/)  ; ssy(7) = "P"    ! via W
        syp(:,8) = (/zro,zro, f12/)  ; ssy(8) = "X"
      ENDIF
      IF ( (idsyst == 2).AND.(idtype ==  2) ) THEN       ! body centered tetragonal (a < c)
        nosyp = 9
        ALLOCATE ( syp(3,nosyp),ssy(nosyp) )
        syp(:,1) = (/-f12,f12,f12/)  ; ssy(1) = "Z"    ! via F and Sigma
        syp(:,2) = (/zro,zro, zro/)  ; ssy(2) = "g"    ! via Delta
        syp(:,3) = (/zro,zro, f12/)  ; ssy(3) = "X"    ! via W
        syp(:,4) = (/f14,f14, f14/)  ; ssy(4) = "P"    ! via Q
        syp(:,5) = (/zro,f12, zro/)  ; ssy(5) = "N"     
        syp(:,6) = (/zro,zro, zro/)  ; ssy(6) = "g"    ! via Lambda
        syp(:,7) = (/f12,f12,-f12/)  ; ssy(7) = "Z"    ! via U and Y
        syp(:,8) = (/f12,f12, zro/)  ; ssy(8) = "X"
        syp(:,9) = (/f14,f14, f14/)  ; ssy(9) = "P"
      ENDIF
      IF ( (idsyst == 2).AND.(idtype ==  1) ) THEN       ! primitive tetragonal (a < c)
        nosyp = 8
        ALLOCATE ( syp(3,nosyp),ssy(nosyp) )
        syp(:,1) = (/zro,zro, zro/)  ; ssy(1) = "g"    ! via Delta
        syp(:,2) = (/f12,zro, zro/)  ; ssy(2) = "X"    ! via Y
        syp(:,3) = (/f12,f12, zro/)  ; ssy(3) = "M"    ! via Sigma
        syp(:,4) = (/zro,zro, zro/)  ; ssy(4) = "g"    ! via Lambda
        syp(:,5) = (/zro,zro, f12/)  ; ssy(5) = "Z"    ! via U
        syp(:,6) = (/f12,zro, f12/)  ; ssy(6) = "R"    ! via T
        syp(:,7) = (/f12,f12, f12/)  ; ssy(7) = "A"    ! via S
        syp(:,8) = (/zro,zro, f12/)  ; ssy(8) = "Z"
      ENDIF
      IF ( (idsyst == 3).AND.(idtype ==  1) ) THEN       ! primitive tetragonal (a < c)
        nosyp = 10
        ALLOCATE ( syp(3,nosyp),ssy(nosyp) )
        syp(:,1) = (/zro,zro, zro/)  ; ssy(1) = "g"    ! via Sigma
        syp(:,2) = (/f12,zro, zro/)  ; ssy(2) = "X"    ! via D
        syp(:,3) = (/f12,f12, zro/)  ; ssy(3) = "S"    ! via C
        syp(:,4) = (/zro,f12, zro/)  ; ssy(4) = "Y"    ! via Delta
        syp(:,5) = (/zro,zro, zro/)  ; ssy(5) = "g"    ! via Lambda
        syp(:,6) = (/zro,zro, f12/)  ; ssy(6) = "Z"    ! via A
        syp(:,7) = (/f12,zro, f12/)  ; ssy(7) = "U"    ! via P
        syp(:,8) = (/f12,f12, f12/)  ; ssy(8) = "R"    ! via E
        syp(:,9) = (/zro,f12, f12/)  ; ssy(8) = "T"    ! via B
        syp(:,10)= (/zro,zro, f12/)  ; ssy(8) = "Z"
      ENDIF

      IF ( nosyp == -1 )  CALL juDFT_error
     +     ("lattice not implemented",calledby ="bandstr1")

      END SUBROUTINE get_points
!----------------------------------------------------------------------
      END MODULE m_bandstr1