The Elk Code
readafieldt.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2014 K. Krieger, 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 
6 subroutine readafieldt
7 use modmain
8 use modtddft
9 implicit none
10 ! local variables
11 integer ios,ntimes_,its,its_
12 real(8) times_,t1
13 ! generate the time step grid
14 call gentimes
15 ! read in the time-dependent vector potential
16 open(50,file='AFIELDT.OUT',form='FORMATTED',action='READ',status='OLD', &
17  iostat=ios)
18 if (ios /= 0) then
19  write(*,*)
20  write(*,'("Error(readafieldt): error opening AFIELDT.OUT")')
21  write(*,*)
22  stop
23 end if
24 read(50,*) ntimes_
25 if (ntimes_ < 1) then
26  write(*,*)
27  write(*,'("Error(readafieldt): ntimes < 1 : ",I8)') ntimes_
28  write(*,*)
29  stop
30 end if
31 ntimes=min(ntimes,ntimes_)
32 if (allocated(afieldt)) deallocate(afieldt)
33 allocate(afieldt(3,ntimes))
34 do its=1,ntimes
35  read(50,*) its_,times_,afieldt(:,its)
36  if (its /= its_) then
37  write(*,*)
38  write(*,'("Error(readafieldt): time step number mismatch")')
39  write(*,'(" internal : ",I8)') its
40  write(*,'(" AFIELDT.OUT : ",I8)') its_
41  write(*,*)
42  stop
43  end if
44  t1=abs(times(its)-times_)
45  if (t1 > 1.d-10) then
46  write(*,*)
47  write(*,'("Error(readafieldt): time step mismatch for step number ",&
48  &I8)') its
49  write(*,'(" internal : ",G18.10)') times(its)
50  write(*,'(" AFIELDT.OUT : ",G18.10)') times_
51  stop
52  end if
53 end do
54 close(50)
55 if (.not.tafspt) return
56 open(50,file='AFSPT.OUT',form='FORMATTED',action='READ',status='OLD',iostat=ios)
57 if (ios /= 0) then
58  write(*,*)
59  write(*,'("Error(readafieldt): error opening AFSPT.OUT")')
60  write(*,*)
61  stop
62 end if
63 read(50,*) ntimes_
64 if (ntimes /= ntimes_) then
65  write(*,*)
66  write(*,'("Error(readafieldt): differing ntimes")')
67  write(*,'(" internal : ",I8)') ntimes
68  write(*,'(" AFSPT.OUT : ",I8)') ntimes_
69  write(*,*)
70  stop
71 end if
72 if (allocated(afspt)) deallocate(afspt)
73 allocate(afspt(3,3,ntimes))
74 do its=1,ntimes
75  read(50,*) its_,times_,afspt(:,:,its)
76  if (its /= its_) then
77  write(*,*)
78  write(*,'("Error(readafieldt): time step number mismatch")')
79  write(*,'(" internal : ",I8)') its
80  write(*,'(" AFSPT.OUT : ",I8)') its_
81  write(*,*)
82  stop
83  end if
84  t1=abs(times(its)-times_)
85  if (t1 > 1.d-10) then
86  write(*,*)
87  write(*,'("Error(readafieldt): time step mismatch for step number ",&
88  &I8)') its
89  write(*,'(" internal : ",G18.10)') times(its)
90  write(*,'(" AFSPT.OUT : ",G18.10)') times_
91  stop
92  end if
93 end do
94 close(50)
95 end subroutine
96 
integer ntimes
Definition: modtddft.f90:42
subroutine readafieldt
Definition: readafieldt.f90:7
subroutine gentimes
Definition: gentimes.f90:7
logical tafspt
Definition: modtddft.f90:60
real(8), dimension(:,:), allocatable afieldt
Definition: modtddft.f90:58
real(8), dimension(:), allocatable times
Definition: modtddft.f90:48
real(8), dimension(:,:,:), allocatable afspt
Definition: modtddft.f90:62