The Elk Code
Loading...
Searching...
No Matches
latvstep.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
latvstep
7
use
modmain
8
use
modmpi
9
implicit none
10
! local variables
11
integer
i
12
real
(8) t1
13
do
i=1,
nstrain
14
! compute product of current and previous stress tensor components
15
t1=
stress
(i)*
stressp
(i)
16
! if component is in the same direction then increase step size parameter
17
if
(t1 > 0.d0)
then
18
taulatv
(i)=
taulatv
(i)+
tau0latv
19
else
20
taulatv
(i)=
tau0latv
21
end if
22
t1=
taulatv
(i)*(
stress
(i)+
stressp
(i))
23
avec
(1:3,1:3)=
avec
(1:3,1:3)-t1*
strain
(1:3,1:3,i)
24
end do
25
! compute the new unit cell volume
26
call
reciplat
(
avec
,
bvec
,
omega
,
omegabz
)
27
! preserve the volume if required (added by Ronald Cohen)
28
if
(
latvopt
== 2)
then
29
t1=(
omega0
/
omega
)**(1.d0/3.d0)
30
avec
(1:3,1:3)=t1*
avec
(1:3,1:3)
31
omega
=
omega0
32
end if
33
! each MPI process should have numerically identical lattice vectors
34
call
mpi_bcast(
avec
,9,mpi_double_precision,0,
mpicom
,
ierror
)
35
end subroutine
36
latvstep
subroutine latvstep
Definition
latvstep.f90:7
modmain
Definition
modmain.f90:6
modmain::nstrain
integer nstrain
Definition
modmain.f90:1012
modmain::bvec
real(8), dimension(3, 3) bvec
Definition
modmain.f90:16
modmain::omega
real(8) omega
Definition
modmain.f90:20
modmain::stressp
real(8), dimension(9) stressp
Definition
modmain.f90:1023
modmain::avec
real(8), dimension(3, 3) avec
Definition
modmain.f90:12
modmain::latvopt
integer latvopt
Definition
modmain.f90:1034
modmain::tau0latv
real(8) tau0latv
Definition
modmain.f90:1038
modmain::omegabz
real(8) omegabz
Definition
modmain.f90:22
modmain::stress
real(8), dimension(9) stress
Definition
modmain.f90:1021
modmain::strain
real(8), dimension(3, 3, 9) strain
Definition
modmain.f90:1016
modmain::omega0
real(8) omega0
Definition
modmain.f90:20
modmain::taulatv
real(8), dimension(9) taulatv
Definition
modmain.f90:1040
modmpi
Definition
modmpi.f90:6
modmpi::ierror
integer ierror
Definition
modmpi.f90:19
modmpi::mpicom
integer mpicom
Definition
modmpi.f90:11
reciplat
subroutine reciplat(avec, bvec, omega, omegabz)
Definition
reciplat.f90:10
latvstep.f90
Generated by
1.9.8