File: writephn.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (64 lines) | stat: -rw-r--r-- 1,713 bytes parent folder | download
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
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64

! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

subroutine writephn
use modmain
use modphonon
implicit none
! local variables
integer iq,i,j,is,ia,ip
! allocatable arrays
real(8), allocatable :: w(:)
complex(8), allocatable :: dynq(:,:,:),dynr(:,:,:)
complex(8), allocatable :: dynp(:,:),ev(:,:)
! initialise universal variables
call init0
call init2
allocate(w(nbph))
allocate(dynq(nbph,nbph,nqpt))
allocate(dynr(nbph,nbph,nqptnr))
allocate(dynp(nbph,nbph))
allocate(ev(nbph,nbph))
! read in the dynamical matrices
call readdyn(dynq)
! apply the acoustic sum rule
call sumrule(dynq)
! Fourier transform the dynamical matrices to real-space
call dynqtor(dynq,dynr)
open(50,file='PHONON.OUT',form='FORMATTED')
do iq=1,nphwrt
  call dynrtoq(vqlwrt(:,iq),dynr,dynp)
  call dynev(dynp,w,ev)
  write(50,*)
  write(50,'(I6,3G18.10," : q-point, vqlwrt")') iq,vqlwrt(:,iq)
  do j=1,nbph
    write(50,*)
    write(50,'(I6,G18.10," : mode, frequency")') j,w(j)
    i=0
    do is=1,nspecies
      do ia=1,natoms(is)
        do ip=1,3
          i=i+1
          if (i.eq.1) then
            write(50,'(3I4,2G18.10," : species, atom, polarisation, &
             &eigenvector")') is,ia,ip,ev(i,j)
          else
            write(50,'(3I4,2G18.10)') is,ia,ip,ev(i,j)
          end if
        end do
      end do
    end do
  end do
  write(50,*)
end do
close(50)
write(*,*)
write(*,'("Info(writephn): phonon frequencies and eigenvectors written to &
 &PHONON.OUT")')
write(*,'(" for all q-vectors in the phwrite list")')
deallocate(w,dynq,dynr,dynp,ev)
return
end subroutine