The Elk Code
Loading...
Searching...
No Matches
readforcet.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2020 J. K. Dewhurst and S. Sharma.
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
readforcet
7
use
modmain
8
use
modtddft
9
implicit none
10
! local variables
11
integer
ios,its,its_
12
integer
is,ia,ias,is_,ia_
13
real
(8) times_,t1
14
! allocate and zero the time-dependent force array
15
if
(
allocated
(
forcet
))
deallocate
(
forcet
)
16
allocate
(
forcet
(3,
natmtot
,
ntimes
))
17
forcet
(:,:,:)=0.d0
18
! read in the time-dependent total atomic forces
19
open
(50,file=
'FORCETOT_TD.OUT'
,form=
'FORMATTED'
,action=
'READ'
,status=
'OLD'
, &
20
iostat=ios)
21
if
(ios /= 0)
then
22
write
(*,*)
23
write
(*,
'("Error(readforcet): error opening FORCETOT_TD.OUT")'
)
24
write
(*,*)
25
stop
26
end if
27
do
its=1,
ntimes
-1,
ntsforce
28
read
(50,*) its_,times_
29
if
(its /= its_)
then
30
write
(*,*)
31
write
(*,
'("Error(readforcet): time step number mismatch")'
)
32
write
(*,
'(" internal : ",I8)'
) its
33
write
(*,
'(" FORCETOT_TD.OUT : ",I8)'
) its_
34
write
(*,*)
35
stop
36
end if
37
t1=abs(
times
(its)-times_)
38
if
(t1 > 1.d-10)
then
39
write
(*,*)
40
write
(*,
'("Error(readforcet): time step mismatch for step number ",I8)'
) its
41
write
(*,
'(" internal : ",G18.10)'
)
times
(its)
42
write
(*,
'(" FORCETOT_TD.OUT : ",G18.10)'
) times_
43
stop
44
end if
45
do
is=1,
nspecies
46
do
ia=1,
natoms
(is)
47
ias=
idxas
(ia,is)
48
read
(50,*) is_,ia_,
forcet
(:,ias,its)
49
if
((is /= is_).or.(ia /= ia_))
then
50
write
(*,*)
51
write
(*,
'("Error(readforcet): species or atom number mismatch for time &
52
&step number ",I8)'
) its
53
write
(*,
'(" internal : ",2I4)'
) is,ia
54
write
(*,
'(" FORCETOT_TD.OUT : ",2I4)'
) is_,ia_
55
write
(*,*)
56
stop
57
end if
58
end do
59
end do
60
end do
61
close
(50)
62
end subroutine
63
modmain
Definition
modmain.f90:6
modmain::natoms
integer, dimension(maxspecies) natoms
Definition
modmain.f90:36
modmain::idxas
integer, dimension(maxatoms, maxspecies) idxas
Definition
modmain.f90:42
modmain::natmtot
integer natmtot
Definition
modmain.f90:40
modmain::nspecies
integer nspecies
Definition
modmain.f90:34
modtddft
Definition
modtddft.f90:6
modtddft::times
real(8), dimension(:), allocatable times
Definition
modtddft.f90:48
modtddft::ntsforce
integer ntsforce
Definition
modtddft.f90:98
modtddft::forcet
real(8), dimension(:,:,:), allocatable forcet
Definition
modtddft.f90:102
modtddft::ntimes
integer ntimes
Definition
modtddft.f90:42
readforcet
subroutine readforcet
Definition
readforcet.f90:7
readforcet.f90
Generated by
1.9.8