The Elk Code
Loading...
Searching...
No Matches
piezoelt.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2023 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
piezoelt
7
use
modmain
8
use
modmpi
9
use
modtest
10
implicit none
11
! local variables
12
integer
i,j,k
13
real
(8) pvl1(3),pvl2(3)
14
real
(8) vc(3),t1
15
! allocatable arrays
16
real
(8),
allocatable
:: pelt(:,:)
17
! initialise universal variables
18
call
init0
19
call
init1
20
! store original parameters
21
avec0
(:,:)=
avec
(:,:)
22
ngridk0
(:)=
ngridk
(:)
23
maxscl0
=
maxscl
24
tshift0
=
tshift
25
! no shifting of the atomic basis
26
tshift
=.false.
27
! generate the strain tensors
28
call
genstrain
29
! allocate the piezoelectric tensor array
30
allocate
(pelt(3,
nstrain
))
31
! initial ground-state run should start from atomic densities
32
trdstate
=.false.
33
! run the ground-state calculation
34
call
gndstate
35
! subsequent calculations will read in the previous potential
36
trdstate
=.true.
37
! compute the first polarisation in lattice coordinates
38
call
polar
(pvl1)
39
! loop over strain tensors
40
do
istrain
=1,
nstrain
41
if
(
mp_mpi
)
then
42
write
(*,
'("Info(piezoelt): working on strain tensor ",I1," of ",I1)'
) &
43
istrain
,
nstrain
44
end if
45
! restore the lattice vectors
46
avec
(:,:)=
avec0
(:,:)
47
! run the ground-state calculation
48
call
gndstate
49
! compute the second polarisation
50
call
polar
(pvl2)
51
do
i=1,3
52
! add multiple of 2*pi to bring polarisation vectors into coincidence
53
pvl1(i)=modulo(pvl1(i),
twopi
)
54
pvl2(i)=modulo(pvl2(i),
twopi
)
55
t1=pvl1(i)-pvl2(i)
56
if
(abs(t1-
twopi
) < abs(t1))
then
57
pvl1(i)=pvl1(i)-
twopi
58
else
if
(abs(t1+
twopi
) < abs(t1))
then
59
pvl1(i)=pvl1(i)+
twopi
60
end if
61
! calculate the piezoelectric tensor element from difference in polarisations
62
t1=
wkptnr
*
occmax
*dble(
nkspolar
*
ngridk
(i))/(
twopi
*
deltast
)
63
pelt(i,
istrain
)=t1*(pvl2(i)-pvl1(i))
64
end do
65
end do
66
if
(
mp_mpi
)
then
67
open
(50,file=
'PIEZOELT.OUT'
,form=
'FORMATTED'
)
68
write
(50,*)
69
write
(50,
'("Lattice vector matrix, A, changed by")'
)
70
write
(50,*)
71
write
(50,→ₖ
'(" A A + e dt,")'
)
72
write
(50,*)
73
write
(50,ₖ
'("where dt is an infinitesimal scalar and e is a strain tensor")'
)
74
write
(50,*)
75
write
(50,
'("The piezoelectric tensor is the derivative of the polarisation &
76
ᵢ
&vector dP/dt, i=1...3")'
)
77
do
k=1,
nstrain
78
write
(50,*)
79
write
(50,
'("Strain tensor k : ",I1)'
) k
80
do
j=1,3
81
write
(50,
'(3G18.10)'
) (
strain
(i,j,k),i=1,3)
82
end do
83
write
(50,ᵢ
'("Piezoelectric tensor components dP/dt, i=1...3 :")'
)
84
write
(50,
'(" lattice coordinates : ",3G18.10)'
) pelt(:,k)
85
call
r3mv
(
avec
,pelt(:,k),vc)
86
write
(50,
'(" Cartesian coordinates : ",3G18.10)'
) vc
87
write
(50,
'(" length : ",G18.10)'
) norm2(vc(1:3))
88
end do
89
close
(50)
90
write
(*,*)
91
write
(*,
'("Info(piezoelt):")'
)
92
write
(*,
'(" Piezoelectric tensor written to PIEZOELT.OUT")'
)
93
end if
94
! write test file if required
95
call
writetest
(380,
'Piezoelectric tensor'
,nv=3*
nstrain
,tol=5.d-5,rva=pelt)
96
deallocate
(pelt)
97
! restore original parameters
98
istrain
=0
99
avec
(:,:)=
avec0
(:,:)
100
tshift
=
tshift0
101
end subroutine
102
genstrain
subroutine genstrain
Definition
genstrain.f90:7
gndstate
subroutine gndstate
Definition
gndstate.f90:10
init0
subroutine init0
Definition
init0.f90:10
init1
subroutine init1
Definition
init1.f90:10
modmain
Definition
modmain.f90:6
modmain::deltast
real(8) deltast
Definition
modmain.f90:1019
modmain::nstrain
integer nstrain
Definition
modmain.f90:1012
modmain::avec0
real(8), dimension(3, 3) avec0
Definition
modmain.f90:12
modmain::tshift
logical tshift
Definition
modmain.f90:352
modmain::ngridk0
integer, dimension(3) ngridk0
Definition
modmain.f90:448
modmain::twopi
real(8), parameter twopi
Definition
modmain.f90:1230
modmain::avec
real(8), dimension(3, 3) avec
Definition
modmain.f90:12
modmain::tshift0
logical tshift0
Definition
modmain.f90:352
modmain::nkspolar
integer nkspolar
Definition
modmain.f90:485
modmain::strain
real(8), dimension(3, 3, 9) strain
Definition
modmain.f90:1016
modmain::ngridk
integer, dimension(3) ngridk
Definition
modmain.f90:448
modmain::maxscl0
integer maxscl0
Definition
modmain.f90:1046
modmain::trdstate
logical trdstate
Definition
modmain.f90:682
modmain::istrain
integer istrain
Definition
modmain.f90:1014
modmain::maxscl
integer maxscl
Definition
modmain.f90:1046
modmain::occmax
real(8) occmax
Definition
modmain.f90:898
modmain::wkptnr
real(8) wkptnr
Definition
modmain.f90:477
modmpi
Definition
modmpi.f90:6
modmpi::mp_mpi
logical mp_mpi
Definition
modmpi.f90:17
modtest
Definition
modtest.f90:6
modtest::writetest
subroutine writetest(id, descr, nv, iv, iva, tol, rv, rva, zv, zva)
Definition
modtest.f90:16
piezoelt
subroutine piezoelt
Definition
piezoelt.f90:7
polar
subroutine polar(pvl)
Definition
polar.f90:10
r3mv
pure subroutine r3mv(a, x, y)
Definition
r3mv.f90:10
piezoelt.f90
Generated by
1.9.8