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
6subroutine readforcet
7use modmain
8use modtddft
9implicit none
10! local variables
11integer ios,its,its_
12integer is,ia,ias,is_,ia_
13real(8) times_,t1
14! allocate and zero the time-dependent force array
15if (allocated(forcet)) deallocate(forcet)
16allocate(forcet(3,natmtot,ntimes))
17forcet(:,:,:)=0.d0
18! read in the time-dependent total atomic forces
19open(50,file='FORCETOT_TD.OUT',form='FORMATTED',action='READ',status='OLD', &
20 iostat=ios)
21if (ios /= 0) then
22 write(*,*)
23 write(*,'("Error(readforcet): error opening FORCETOT_TD.OUT")')
24 write(*,*)
25 stop
26end if
27do 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
60end do
61close(50)
62end subroutine
63
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
integer, dimension(maxatoms, maxspecies) idxas
Definition modmain.f90:42
integer natmtot
Definition modmain.f90:40
integer nspecies
Definition modmain.f90:34
real(8), dimension(:), allocatable times
Definition modtddft.f90:48
integer ntsforce
Definition modtddft.f90:98
real(8), dimension(:,:,:), allocatable forcet
Definition modtddft.f90:102
integer ntimes
Definition modtddft.f90:42
subroutine readforcet
Definition readforcet.f90:7