The Elk Code
 
Loading...
Searching...
No Matches
genstress.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 genstress
7use modmain
8use modmpi
9implicit none
10! local variables
11real(8) et0
12! store original parameters
13avec0(:,:)=avec(:,:)
15tforce=.false.
16! restore original symmetries
17call symmetry
18! generate the strain tensors
19call genstrain
20! zero the stress tensor components
21stress(:)=0.d0
22! run the ground-state calculation
23call gndstate
24! check for stop signal
25if (tstop) goto 10
26! subsequent calculations will read in the potential from STATE.OUT
27trdstate=.true.
28! store the total energy
29et0=engytot
30! loop over strain tensors
32 if (mp_mpi) then
33 write(*,'("Info(genstress): strain tensor ",I1," of ",I1)') istrain,nstrain
34 end if
35! restore the lattice vectors
36 avec(:,:)=avec0(:,:)
37! run the ground-state calculation
38 call gndstate
39! check for stop signal
40 if (tstop) goto 10
41! compute the stress tensor component
43end do
4410 continue
45! compute the maximum stress magnitude over all lattice vectors
46stressmax=maxval(abs(stress(1:nstrain)))
47! restore original parameters
48istrain=0
49avec(:,:)=avec0(:,:)
51end subroutine
52
subroutine genstrain
Definition genstrain.f90:7
subroutine genstress
Definition genstress.f90:7
subroutine gndstate
Definition gndstate.f90:10
real(8) deltast
Definition modmain.f90:1019
integer nstrain
Definition modmain.f90:1012
real(8), dimension(3, 3) avec0
Definition modmain.f90:12
real(8), dimension(3, 3) avec
Definition modmain.f90:12
logical tstop
Definition modmain.f90:1052
real(8), dimension(9) stress
Definition modmain.f90:1021
logical trdstate
Definition modmain.f90:682
integer istrain
Definition modmain.f90:1014
real(8) stressmax
Definition modmain.f90:1025
real(8) engytot
Definition modmain.f90:980
logical tforce
Definition modmain.f90:986
logical tforce0
Definition modmain.f90:986
logical mp_mpi
Definition modmpi.f90:17
subroutine symmetry
Definition symmetry.f90:7