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
6subroutine latvstep
7use modmain
8use modmpi
9implicit none
10! local variables
11integer i
12real(8) t1
13do 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
19 else
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)
24end do
25! compute the new unit cell volume
27! preserve the volume if required (added by Ronald Cohen)
28if (latvopt == 2) then
29 t1=(omega0/omega)**(1.d0/3.d0)
30 avec(1:3,1:3)=t1*avec(1:3,1:3)
32end if
33! each MPI process should have numerically identical lattice vectors
34call mpi_bcast(avec,9,mpi_double_precision,0,mpicom,ierror)
35end subroutine
36
subroutine latvstep
Definition latvstep.f90:7
integer nstrain
Definition modmain.f90:1012
real(8), dimension(3, 3) bvec
Definition modmain.f90:16
real(8) omega
Definition modmain.f90:20
real(8), dimension(9) stressp
Definition modmain.f90:1023
real(8), dimension(3, 3) avec
Definition modmain.f90:12
integer latvopt
Definition modmain.f90:1034
real(8) tau0latv
Definition modmain.f90:1038
real(8) omegabz
Definition modmain.f90:22
real(8), dimension(9) stress
Definition modmain.f90:1021
real(8), dimension(3, 3, 9) strain
Definition modmain.f90:1016
real(8) omega0
Definition modmain.f90:20
real(8), dimension(9) taulatv
Definition modmain.f90:1040
integer ierror
Definition modmpi.f90:19
integer mpicom
Definition modmpi.f90:11
subroutine reciplat(avec, bvec, omega, omegabz)
Definition reciplat.f90:10