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
6subroutine genstrain
7use modmain
8implicit none
9! local variables
10integer i,j,k
11real(8) a(3,3),b(3,3),t1
12! set first strain equal to isotropic scaling
13t1=norm2(avec(1:3,1:3))
14strain(1:3,1:3,1)=avec(1:3,1:3)/t1
15nstrain=1
16do 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
35 strain(1:3,1:3,nstrain)=b(1:3,1:3)/t1
36 end do
37end do
38! zero small components
39do 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
45end do
46! lattice optimisation case
47if ((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)
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
64end if
65end subroutine
66
subroutine genstrain
Definition genstrain.f90:7
integer nstrain
Definition modmain.f90:1012
real(8), dimension(3, 3) avec
Definition modmain.f90:12
integer latvopt
Definition modmain.f90:1034
real(8) epslat
Definition modmain.f90:24
real(8), dimension(3, 3, 9) strain
Definition modmain.f90:1016
integer task
Definition modmain.f90:1298
real(8), dimension(3, 3) ainv
Definition modmain.f90:14
pure subroutine r3mtm(a, b, c)
Definition r3mtm.f90:10
subroutine symmat(al)
Definition symmat.f90:7