The Elk Code
 
Loading...
Searching...
No Matches
phononsc.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2008 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine phononsc
7use modmain
8use modphonon
9use modmpi
10use moddelf
11implicit none
12! local variables
13integer is,ia,ja,ias,jas
14integer ip,nph,p,n,i
15real(8) a,b,t1
16real(8) ft(3,maxatoms*maxspecies)
17complex(8) z1,z2
18! allocatable arrays
19real(8), allocatable :: vsmt0(:,:),vsir0(:)
20complex(8), allocatable :: dyn(:,:)
21! store original parameters
22avec0(:,:)=avec(:,:)
23natoms0(:)=natoms(:)
24atposl0(:,:,:)=atposl(:,:,:)
25bfcmt00(:,:,:)=bfcmt0(:,:,:)
26mommtfix0(:,:,:)=mommtfix(:,:,:)
31ngridk0(:)=ngridk(:)
32! no shifting of atomic basis allowed
33tshift=.false.
34! require forces
35tforce=.true.
36! determine k-point grid size from radkpt
37autokpt=.true.
38! no reduction to primitive cell
39primcell=.false.
40! initialise universal variables
41call init0
42! initialise q-point dependent variables
43call init2
44! store original parameters
47bvec0(:,:)=bvec(:,:)
48binv0(:,:)=binv(:,:)
49atposc0(:,:,:)=atposc(:,:,:)
50ngridg0(:)=ngridg(:)
52if (allocated(ivg0)) deallocate(ivg0)
53allocate(ivg0(3,ngtot0))
54ivg0(:,:)=ivg(:,:)
55if (allocated(igfft0)) deallocate(igfft0)
56allocate(igfft0(ngtot0))
57igfft0(:)=igfft(:)
58! allocate the Kohn-Sham potential derivative arrays
59if (allocated(dvsbs)) deallocate(dvsbs)
61allocate(dvsbs(n))
62dvsmt(1:npmtmax,1:natmtot) => dvsbs(1:)
63if (allocated(dvsir)) deallocate(dvsir)
64allocate(dvsir(ngtot))
65! allocate supercell offset vector array
66if (allocated(vscph)) deallocate(vscph)
67allocate(vscph(3,nqptnr))
68! allocate dynamical matrix column
69allocate(dyn(3,natmtot))
70! begin new phonon task
7110 continue
72natoms(:)=natoms0(:)
73! find a dynamical matrix to calculate
74call dyntask(80,filext)
75! if nothing more to do then restore original input parameters and return
76if (iqph == 0) then
77 filext='.OUT'
78 natoms(:)=natoms0(:)
79 avec(:,:)=avec0(:,:)
80 atposl(:,:,:)=atposl0(:,:,:)
81 bfcmt0(:,:,:)=bfcmt00(:,:,:)
82 mommtfix(:,:,:)=mommtfix0(:,:,:)
87 ngridk(:)=ngridk0(:)
88 deallocate(ivg0,igfft0)
89 return
90end if
91if (mp_mpi) write(*,'("Info(phononsc): working on ",A)') 'DYN'//trim(filext)
92! dry run: just generate empty DYN files
93if (task == 202) goto 10
94! zero the dynamical matrix row
95dyn(:,:)=0.d0
96! zero the Kohn-Sham potential derivative
97dvsmt(:,:)=0.d0
98dvsir(:)=0.d0
99! check to see if mass is considered infinite
100if (spmass(isph) <= 0.d0) goto 20
101! loop over phases: 0 = cos and 1 = sin displacements
102if ((ivq(1,iqph) == 0).and.(ivq(2,iqph) == 0).and.(ivq(3,iqph) == 0)) then
103 nph=0
104else
105 nph=1
106end if
107! initialise or read the charge density and potentials from file
108trdstate=(task == 201)
109! loop over cos and sin displacements
110do p=0,nph
111! generate the supercell with negative displacement
112 call genscph(p,-0.5d0*deltaph)
113! run the ground-state calculation
114 call gndstate
115! subsequent calculations will read in this supercell potential
116 trdstate=.true.
117! store the total force for this displacement
118 do ias=1,natmtot
119 ft(:,ias)=forcetot(:,ias)
120 end do
121! store the Kohn-Sham potential for this displacement
122 allocate(vsmt0(npmtmax,natmtot),vsir0(ngtot))
123 vsmt0(:,:)=vsmt(:,:)
124 vsir0(:)=vsir(:)
125! generate the supercell again with positive displacement
126 call genscph(p,0.5d0*deltaph)
127! run the ground-state calculation again
128 call gndstate
129! compute the complex Kohn-Sham potential derivative with implicit q-phase
130 call phscdvs(p,vsmt0,vsir0)
131 deallocate(vsmt0,vsir0)
132! Fourier transform the force differences to obtain the dynamical matrix
133 z1=1.d0/(dble(nscph)*deltaph)
134! multiply by i for sin-like displacement
135 if (p == 1) z1=z1*zi
136 ias=0
137 jas=0
138 do is=1,nspecies
139 ja=0
140 do ia=1,natoms0(is)
141 ias=ias+1
142 do i=1,nscph
143 ja=ja+1
144 jas=jas+1
145 t1=-dot_product(vqc(:,iqph),vscph(:,i))
146 z2=z1*cmplx(cos(t1),sin(t1),8)
147 do ip=1,3
148 t1=-(forcetot(ip,jas)-ft(ip,jas))
149 dyn(ip,ias)=dyn(ip,ias)+z2*t1
150 end do
151 end do
152 end do
153 end do
154end do
15520 continue
156! write dynamical matrix row to file
157if (mp_mpi) then
158 ias=0
159 do is=1,nspecies
160 do ia=1,natoms0(is)
161 ias=ias+1
162 do ip=1,3
163 a=dble(dyn(ip,ias))
164 b=aimag(dyn(ip,ias))
165 if (abs(a) < 1.d-12) a=0.d0
166 if (abs(b) < 1.d-12) b=0.d0
167 write(80,'(2G18.10," : is = ",I4,", ia = ",I4,", ip = ",I4)') a,b,is, &
168 ia,ip
169 end do
170 end do
171 end do
172 close(80)
173! write the complex Kohn-Sham potential derivative to file
174 natoms(:)=natoms0(:)
177 ngridg(:)=ngridg0(:)
178 call writedvs(filext)
179end if
180! delete the non-essential files
181call delfiles(evec=.true.,eval=.true.,occ=.true.)
182! synchronise MPI processes
183call mpi_barrier(mpicom,ierror)
184goto 10
185end subroutine
186
subroutine dyntask(fnum, fext)
Definition dyntask.f90:7
subroutine genscph(p, dph)
Definition genscph.f90:7
subroutine gndstate
Definition gndstate.f90:10
subroutine init0
Definition init0.f90:10
subroutine init2
Definition init2.f90:7
subroutine delfiles(evec, devec, eval, occ, pmat, epsi)
Definition moddelf.f90:25
integer, dimension(:,:), allocatable ivg0
Definition modmain.f90:400
integer natmtot0
Definition modmain.f90:40
real(8), dimension(3, 3) avec0
Definition modmain.f90:12
integer ngtot
Definition modmain.f90:390
integer, dimension(3) ngridg
Definition modmain.f90:386
logical autokpt0
Definition modmain.f90:444
logical tshift
Definition modmain.f90:352
integer, dimension(3) ngridk0
Definition modmain.f90:448
real(8), dimension(3, maxatoms, maxspecies) atposc
Definition modmain.f90:54
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
character(256) filext
Definition modmain.f90:1300
real(8), dimension(3, 3) bvec
Definition modmain.f90:16
integer, dimension(:,:), allocatable ivq
Definition modmain.f90:529
complex(8), parameter zi
Definition modmain.f90:1239
real(8), dimension(3, maxatoms, maxspecies) mommtfix0
Definition modmain.f90:259
real(8), dimension(3, 3) binv
Definition modmain.f90:18
integer nqptnr
Definition modmain.f90:527
real(8), dimension(3, 3) bvec0
Definition modmain.f90:16
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
real(8), dimension(3, 3) avec
Definition modmain.f90:12
real(8), dimension(:,:), allocatable vqc
Definition modmain.f90:547
logical tshift0
Definition modmain.f90:352
logical autokpt
Definition modmain.f90:444
integer natmtot
Definition modmain.f90:40
real(8), dimension(3, maxatoms, maxspecies) bfcmt00
Definition modmain.f90:275
integer ngtot0
Definition modmain.f90:390
integer, dimension(3) ngridk
Definition modmain.f90:448
real(8), dimension(3, maxatoms, maxspecies) atposl0
Definition modmain.f90:51
integer, dimension(:), allocatable igfft
Definition modmain.f90:406
logical trdstate
Definition modmain.f90:682
real(8), dimension(3, maxatoms, maxspecies) mommtfix
Definition modmain.f90:259
real(8), dimension(3, 3) binv0
Definition modmain.f90:18
integer task
Definition modmain.f90:1298
real(8), dimension(3, maxatoms, maxspecies) atposc0
Definition modmain.f90:54
integer npmtmax
Definition modmain.f90:216
real(8), dimension(:,:), allocatable forcetot
Definition modmain.f90:990
integer, dimension(maxatoms *maxspecies) idxis0
Definition modmain.f90:44
real(8), dimension(:,:), pointer, contiguous vsmt
Definition modmain.f90:649
integer, dimension(3) ngridg0
Definition modmain.f90:386
integer, dimension(:,:), allocatable ivg
Definition modmain.f90:400
integer, dimension(maxspecies) natoms0
Definition modmain.f90:36
logical tforce
Definition modmain.f90:986
real(8), dimension(:), allocatable vsir
Definition modmain.f90:651
integer nspecies
Definition modmain.f90:34
logical primcell0
Definition modmain.f90:49
logical tforce0
Definition modmain.f90:986
real(8), dimension(3, maxatoms, maxspecies) atposl
Definition modmain.f90:51
real(8), dimension(3, maxatoms, maxspecies) bfcmt0
Definition modmain.f90:275
real(8), dimension(maxspecies) spmass
Definition modmain.f90:101
integer, dimension(:), allocatable igfft0
Definition modmain.f90:406
logical primcell
Definition modmain.f90:49
integer ierror
Definition modmpi.f90:19
integer mpicom
Definition modmpi.f90:11
logical mp_mpi
Definition modmpi.f90:17
complex(8), dimension(:), allocatable dvsir
complex(8), dimension(:,:), pointer, contiguous dvsmt
integer iqph
Definition modphonon.f90:15
real(8), dimension(:,:), allocatable vscph
Definition modphonon.f90:46
integer isph
Definition modphonon.f90:15
integer nscph
Definition modphonon.f90:44
complex(8), dimension(:), allocatable, target dvsbs
real(8) deltaph
Definition modphonon.f90:48
subroutine phononsc
Definition phononsc.f90:7
subroutine phscdvs(p, vsmt0, vsir0)
Definition phscdvs.f90:7
subroutine writedvs(fext)
Definition writedvs.f90:7