The Elk Code
Loading...
Searching...
No Matches
genstrain.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2013 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
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
genstrain
7
use
modmain
8
implicit none
9
! local variables
10
integer
i,j,k
11
real
(8) a(3,3),b(3,3),t1
12
! set first strain equal to isotropic scaling
13
t1=norm2(
avec
(1:3,1:3))
14
strain
(1:3,1:3,1)=
avec
(1:3,1:3)/t1
15
nstrain
=1
16
do
i=1,3
17
do
j=1,3
18
! set strain tensor in lattice coordinates to delta_ij
19
a(1:3,1:3)=0.d0
20
a(i,j)=1.d0
21
! symmetrise strain tensor
22
call
symmat
(a)
23
! convert to mixed Cartesian-lattice coordinates
24
call
r3mtm
(
ainv
,a,b)
25
! orthogonalise strain tensor to previous tensors
26
do
k=1,
nstrain
27
t1=sum(b(1:3,1:3)*
strain
(1:3,1:3,k))
28
b(1:3,1:3)=b(1:3,1:3)-t1*
strain
(1:3,1:3,k)
29
end do
30
! compute the norm
31
t1=norm2(b(1:3,1:3))
32
if
(t1 <
epslat
) cycle
33
! normalise tensor and store in global array
34
nstrain
=
nstrain
+1
35
strain
(1:3,1:3,
nstrain
)=b(1:3,1:3)/t1
36
end do
37
end do
38
! zero small components
39
do
k=1,
nstrain
40
do
i=1,3
41
do
j=1,3
42
if
(abs(
strain
(i,j,k)) <
epslat
)
strain
(i,j,k)=0.d0
43
end do
44
end do
45
end do
46
! lattice optimisation case
47
if
((
task
== 2).or.(
task
== 3))
then
48
if
(
latvopt
== 2)
then
49
! remove isotropic scaling when latvopt=2
50
strain
(1:3,1:3,1)=
strain
(1:3,1:3,
nstrain
)
51
nstrain
=
nstrain
-1
52
else
if
(
latvopt
< 0)
then
53
! optimise over particular strain when latvopt < 0
54
i=abs(
latvopt
)
55
if
(i >
nstrain
)
then
56
write
(*,*)
57
write
(*,
'("Error(genstrain): |latvopt| > nstrain : ",2I8)'
) i,
nstrain
58
write
(*,*)
59
stop
60
end if
61
strain
(1:3,1:3,1)=
strain
(1:3,1:3,i)
62
nstrain
=1
63
end if
64
end if
65
end subroutine
66
genstrain
subroutine genstrain
Definition
genstrain.f90:7
modmain
Definition
modmain.f90:6
modmain::nstrain
integer nstrain
Definition
modmain.f90:1012
modmain::avec
real(8), dimension(3, 3) avec
Definition
modmain.f90:12
modmain::latvopt
integer latvopt
Definition
modmain.f90:1034
modmain::epslat
real(8) epslat
Definition
modmain.f90:24
modmain::strain
real(8), dimension(3, 3, 9) strain
Definition
modmain.f90:1016
modmain::task
integer task
Definition
modmain.f90:1298
modmain::ainv
real(8), dimension(3, 3) ainv
Definition
modmain.f90:14
r3mtm
pure subroutine r3mtm(a, b, c)
Definition
r3mtm.f90:10
symmat
subroutine symmat(al)
Definition
symmat.f90:7
genstrain.f90
Generated by
1.9.8