The Elk Code
bornechg.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 
6 subroutine bornechg
7 use modmain
8 use modphonon
9 use modmpi
10 use modtest
11 use moddelf
12 implicit none
13 ! local variables
14 integer ip,i
15 real(8) vc(3),pvl1(3),pvl2(3)
16 real(8) becl(3),becc(3),t1
17 ! initialise universal variables
18 call init0
19 call init1
20 ! store original parameters
21 atposl0(:,:,:)=atposl(:,:,:)
22 atposc0(:,:,:)=atposc(:,:,:)
23 ngridk0(:)=ngridk(:)
26 ! no shifting of the atomic basis
27 tshift=.false.
28 ! begin new Born effective charge task
29 10 continue
30 call bectask(80,filext)
31 ! if nothing more to do then restore original input parameters and return
32 if (isph == 0) then
33  filext='.OUT'
35  atposl(:,:,:)=atposl0(:,:,:)
36  return
37 end if
38 if (mp_mpi) then
39  write(*,'("Info(bornechg): working on ",A)') 'BEC'//trim(filext)
40 end if
41 ! dry run: just generate empty BEC files
42 if (task == 209) goto 10
43 ! apply negative atomic displacement
44 atposl(:,:,:)=atposl0(:,:,:)
45 atposc(:,:,:)=atposc0(:,:,:)
46 vc(:)=atposc(:,iaph,isph)
47 vc(ipph)=vc(ipph)-0.5d0*deltaph
48 call r3mv(ainv,vc,atposl(:,iaph,isph))
49 ! initial ground-state run should start from atomic densities
50 trdstate=.false.
51 ! run the ground-state calculation
52 call gndstate
53 ! subsequent calculations will read in the previous potential
54 trdstate=.true.
55 ! compute the first polarisation in lattice coordinates
56 call polar(pvl1)
57 ! apply positive atomic displacement
58 atposl(:,:,:)=atposl0(:,:,:)
59 atposc(:,:,:)=atposc0(:,:,:)
60 vc(:)=atposc(:,iaph,isph)
61 vc(ipph)=vc(ipph)+0.5d0*deltaph
62 call r3mv(ainv,vc,atposl(:,iaph,isph))
63 ! run the ground-state calculation again
64 call gndstate
65 ! compute the second polarisation
66 call polar(pvl2)
67 do i=1,3
68 ! add multiple of 2*pi to bring polarisation vectors into coincidence
69  pvl1(i)=modulo(pvl1(i),twopi)
70  pvl2(i)=modulo(pvl2(i),twopi)
71  t1=pvl1(i)-pvl2(i)
72  if (abs(t1-twopi) < abs(t1)) then
73  pvl1(i)=pvl1(i)-twopi
74  else if (abs(t1+twopi) < abs(t1)) then
75  pvl1(i)=pvl1(i)+twopi
76  end if
77 ! calculate the Born effective charge from the difference in polarisations
79  becl(i)=t1*(pvl2(i)-pvl1(i))
80 end do
81 ! convert from lattice to Cartesian coordinates
82 call r3mv(avec,becl,becc)
83 ! add the core and nuclear charge
84 becc(ipph)=becc(ipph)+chgcr(isph)+spzn(isph)
85 ! write Born effective charge matrix row to file
86 if (mp_mpi) then
87  do ip=1,3
88  write(80,'(G18.10," : ip = ",I4)') becc(ip),ip
89  end do
90  close(80)
91 end if
92 ! write test file if required and return
93 if (test) then
94  call writetest(208,'Born effective charge',nv=3,tol=1.d-3,rva=becc)
95  return
96 end if
97 ! delete the non-essential files
98 call delfiles(evec=.true.,eval=.true.,occ=.true.)
99 ! synchronise MPI processes
100 call mpi_barrier(mpicom,ierror)
101 goto 10
102 end subroutine
103 
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition: modtest.f90:16
subroutine gndstate
Definition: gndstate.f90:10
real(8), parameter twopi
Definition: modmain.f90:1233
character(256) filext
Definition: modmain.f90:1301
real(8), dimension(3, maxatoms, maxspecies) atposl0
Definition: modmain.f90:51
integer task
Definition: modmain.f90:1299
logical mp_mpi
Definition: modmpi.f90:17
logical tshift0
Definition: modmain.f90:352
logical tshift
Definition: modmain.f90:352
real(8), dimension(3, 3) ainv
Definition: modmain.f90:14
integer maxscl0
Definition: modmain.f90:1049
subroutine polar(pvl)
Definition: polar.f90:10
integer, dimension(3) ngridk0
Definition: modmain.f90:448
real(8), dimension(3, maxatoms, maxspecies) atposl
Definition: modmain.f90:51
real(8) occmax
Definition: modmain.f90:901
real(8), dimension(maxspecies) chgcr
Definition: modmain.f90:716
logical test
Definition: modtest.f90:11
real(8), dimension(3, 3) avec
Definition: modmain.f90:12
integer ipph
Definition: modphonon.f90:15
integer, dimension(3) ngridk
Definition: modmain.f90:448
integer isph
Definition: modphonon.f90:15
subroutine init1
Definition: init1.f90:10
subroutine delfiles(evec, devec, eval, occ, pmat, epsi)
Definition: moddelf.f90:25
Definition: modmpi.f90:6
subroutine bornechg
Definition: bornechg.f90:7
subroutine bectask(fnum, fext)
Definition: bectask.f90:7
logical trdstate
Definition: modmain.f90:682
real(8), dimension(maxspecies) spzn
Definition: modmain.f90:80
real(8) wkptnr
Definition: modmain.f90:477
subroutine init0
Definition: init0.f90:10
integer nkspolar
Definition: modmain.f90:485
integer iaph
Definition: modphonon.f90:15
pure subroutine r3mv(a, x, y)
Definition: r3mv.f90:10
real(8), dimension(3, maxatoms, maxspecies) atposc
Definition: modmain.f90:54
integer maxscl
Definition: modmain.f90:1049
real(8), dimension(3, maxatoms, maxspecies) atposc0
Definition: modmain.f90:54
integer mpicom
Definition: modmpi.f90:11
real(8) deltaph
Definition: modphonon.f90:48
integer ierror
Definition: modmpi.f90:19