The Elk Code
 
Loading...
Searching...
No Matches
bornecdyn.f90
Go to the documentation of this file.
1
2! Copyright (C) 2020 J. K. Dewhurst and S. Sharma.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine bornecdyn
7use modmain
8use modphonon
9use modtddft
10use modmpi
11use modtest
12implicit none
13! local variables
14integer its,iw,i
15real(8) vc(3),w1,w2,t0,t1,t2
16character(256) fext
17! allocatable arrays
18real(8), allocatable :: w(:),wt(:),jt(:,:)
19real(8), allocatable :: f1(:),f2(:)
20complex(8), allocatable :: becw(:,:)
21! initialise universal variables
22call init0
23! store original parameters
24atposl0(:,:,:)=atposl(:,:,:)
25atposc0(:,:,:)=atposc(:,:,:)
26afieldc0(:)=afieldc(:)
28! no shifting of atomic basis allowed
29tshift=.false.
30! generate the time step grid
31call gentimes
32! allocate local arrays
33allocate(w(nwplot),wt(ntimes),jt(3,ntimes))
34allocate(f1(ntimes),f2(ntimes),becw(nwplot,3))
35! generate frequency grid (always non-negative)
36w1=max(wplot(1),0.d0)
37w2=max(wplot(2),w1)
38t1=(w2-w1)/dble(nwplot)
39do iw=1,nwplot
40 w(iw)=w1+t1*dble(iw-1)
41end do
42! determine the weights for the spline integration
43call wsplint(ntimes,times,wt)
44! generate a zero A-field and write to file
45npulse=0
46nramp=0
47if (mp_mpi) call genafieldt
48! initial ground-state run should start from atomic densities
49trdstate=.false.
50! begin new Born effective charge task
5110 continue
52call bectask(80,fext)
53! if nothing more to do then restore original input parameters and return
54if (isph == 0) then
55 filext='.OUT'
56 atposl(:,:,:)=atposl0(:,:,:)
57 afieldc(:)=afieldc0(:)
59 trdatdv=.false.
60 deallocate(w,wt,jt,f1,f2,becw)
61 return
62end if
63if (mp_mpi) then
64 write(*,'("Info(bornecdyn): working on ",A)') 'BEC'//trim(fext)
65end if
66! break the crystal symmetry for the displaced atom
67atposl(:,:,:)=atposl0(:,:,:)
68atposc(:,:,:)=atposc0(:,:,:)
69vc(:)=atposc(:,iaph,isph)
70vc(ipph)=vc(ipph)-0.5d0*deltaph
71call r3mv(ainv,vc,atposl(:,iaph,isph))
72! apply a small static A-field
73afieldc(:)=0.d0
74afieldc(ipph)=1.d-4
75! run the ground-state calculation
76call gndstate
77! subsequent calculations will read in the previous potential
78trdstate=.true.
79! write zero atomic forces to file
80if (mp_mpi) call becforce
81! write displacement to file
82atdvc(:,:,:,:)=0.d0
84if (mp_mpi) call writeatdvc
85trdatdv=.true.
86! run the time evolution calculation with Ehrenfest dynamics
87task=462
88call tddft
89task=478
90! read in the total current from file
91call readjtot(jt)
92! filter the high-frequency components from the current
93do its=1,ntimes
94 t1=exp(-swidth*times(its))
95 jt(:,its)=t1*jt(:,its)
96end do
97! compute the dynamical BEC from the Fourier transformed current
98t0=1.d0/(deltaph*cos(tdphi))
99do i=1,3
100 call zftft(w,wt,3,jt(i,1),becw(:,i))
101 becw(:,i)=t0*becw(:,i)
102end do
103! static and nuclear charge
104t1=sum(chgsmt(iasph,1:3))/3.d0+spzn(isph)
105! write Born effective charge matrix row to file
106if (mp_mpi) then
107 do i=1,3
108 if (i == ipph) then
109 t2=t1
110 else
111 t2=0.d0
112 end if
113 do iw=1,nwplot
114 write(80,'(2G18.10)') w(iw),dble(becw(iw,i))+t2
115 end do
116 write(80,*)
117 do iw=1,nwplot
118 write(80,'(2G18.10)') w(iw),aimag(becw(iw,i))
119 end do
120 write(80,*)
121 end do
122 close(80)
123end if
124! synchronise MPI processes
125call mpi_barrier(mpicom,ierror)
126! write test file if required and return
127if (test) then
128 call writetest(478,'dynamical Born effective charge',nv=nwplot,tol=1.d-2, &
129 zva=becw)
130 return
131end if
132goto 10
133end subroutine
134
subroutine becforce
Definition becforce.f90:7
subroutine bectask(fnum, fext)
Definition bectask.f90:7
subroutine bornecdyn
Definition bornecdyn.f90:7
subroutine genafieldt
subroutine gentimes
Definition gentimes.f90:7
subroutine gndstate
Definition gndstate.f90:10
subroutine init0
Definition init0.f90:10
logical tshift
Definition modmain.f90:352
real(8), dimension(3) afieldc0
Definition modmain.f90:325
real(8), dimension(3, maxatoms, maxspecies) atposc
Definition modmain.f90:54
character(256) filext
Definition modmain.f90:1300
real(8), dimension(3) afieldc
Definition modmain.f90:325
integer nwplot
Definition modmain.f90:1070
logical trdatdv
Definition modmain.f90:62
real(8) swidth
Definition modmain.f90:892
logical tshift0
Definition modmain.f90:352
real(8), dimension(2) wplot
Definition modmain.f90:1076
real(8), dimension(3, maxatoms, maxspecies) atposl0
Definition modmain.f90:51
logical trdstate
Definition modmain.f90:682
integer task
Definition modmain.f90:1298
real(8), dimension(3, maxatoms, maxspecies) atposc0
Definition modmain.f90:54
real(8), dimension(3, 3) ainv
Definition modmain.f90:14
real(8), dimension(3, 0:1, maxatoms, maxspecies) atdvc
Definition modmain.f90:64
real(8), dimension(3, maxatoms, maxspecies) atposl
Definition modmain.f90:51
real(8), dimension(maxspecies) spzn
Definition modmain.f90:80
integer ierror
Definition modmpi.f90:19
integer mpicom
Definition modmpi.f90:11
logical mp_mpi
Definition modmpi.f90:17
integer iaph
Definition modphonon.f90:15
integer ipph
Definition modphonon.f90:15
integer iasph
Definition modphonon.f90:15
integer isph
Definition modphonon.f90:15
real(8) deltaph
Definition modphonon.f90:48
real(8) tdphi
Definition modtddft.f90:52
real(8), dimension(:), allocatable times
Definition modtddft.f90:48
real(8), dimension(:,:), allocatable chgsmt
Definition modtddft.f90:86
integer npulse
Definition modtddft.f90:24
integer nramp
Definition modtddft.f90:29
integer ntimes
Definition modtddft.f90:42
logical test
Definition modtest.f90:11
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition modtest.f90:16
pure subroutine r3mv(a, x, y)
Definition r3mv.f90:10
subroutine readjtot(jt)
Definition readjtot.f90:7
subroutine tddft
Definition tddft.f90:7
subroutine writeatdvc
Definition writeatdvc.f90:7
subroutine wsplint(n, x, w)
Definition wsplint.f90:7
subroutine zftft(w, wt, ld, ft, fw)
Definition zftft.f90:7