File: writestress.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (44 lines) | stat: -rw-r--r-- 1,193 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44

! Copyright (C) 2013 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

subroutine writestress
use modmain
use modmpi
implicit none
! local variables
integer i,j,k
! initialise universal variables
call init0
! start from the atomic densities
trdstate=.false.
! generate the stress tensors
call genstress
! write the stress tensors to file
if (mp_mpi) then
  open(50,file='STRESS.OUT',form='FORMATTED')
  write(50,*)
  write(50,'("Lattice vector matrix, A, changed by")')
  write(50,*)
  write(50,'("     A --> A + e_i dt,")')
  write(50,*)
  write(50,'("where dt is an infinitesimal scalar and e_i is a strain tensor")')
  write(50,*)
  write(50,'("Stress is given by the derivative of the total energy dE/dt")')
  do k=1,nstrain
    write(50,*)
    write(50,'("Strain tensor : ",I1)') k
    do j=1,3
      write(50,'(3G18.10)') (strain(i,j,k),i=1,3)
    end do
    write(50,'("Stress : ",G18.10)') stress(k)
  end do
  close(50)
  write(*,*)
  write(*,'("Info(writestress):")')
  write(*,'(" Stress tensor components written to STRESS.OUT")')
end if
return
end subroutine