The Elk Code
 
Loading...
Searching...
No Matches
readrhos.f90
Go to the documentation of this file.
1
2! Copyright (C) 2022 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 readrhos
7use modmain
8use modtddft
9implicit none
10integer ios
11integer natmtot_,npmtmax_,ngtot_
12! allocate static density and charge global arrays
13if (allocated(rhosmt)) deallocate(rhosmt)
14allocate(rhosmt(npmtmax,natmtot,3))
15if (allocated(rhosir)) deallocate(rhosir)
16allocate(rhosir(ngtot,3))
17if (allocated(chgsmt)) deallocate(chgsmt)
18allocate(chgsmt(natmtot,3))
19open(100,file='RHOSTAT.OUT',form='UNFORMATTED',action='READ',status='OLD', &
20 iostat=ios)
21if (ios /= 0) then
22 write(*,*)
23 write(*,'("Error(readrhos): error opening RHOSTAT.OUT")')
24 write(*,*)
25 stop
26end if
27read(100) natmtot_
28if (natmtot /= natmtot_) then
29 write(*,*)
30 write(*,'("Error(readrhos): differing natmtot")')
31 write(*,'(" current : ",I6)') natmtot
32 write(*,'(" RHOSTAT.OUT : ",I6)') natmtot_
33 write(*,*)
34 stop
35end if
36read(100) npmtmax_
37if (npmtmax /= npmtmax_) then
38 write(*,*)
39 write(*,'("Error(readrhos): differing npmtmax")')
40 write(*,'(" current : ",I6)') npmtmax
41 write(*,'(" RHOSTAT.OUT : ",I6)') npmtmax_
42 write(*,*)
43 stop
44end if
45read(100) ngtot_
46if (ngtot /= ngtot_) then
47 write(*,*)
48 write(*,'("Error(readrhos): differing ngtot")')
49 write(*,'(" current : ",I8)') ngtot
50 write(*,'(" RHOSTAT.OUT : ",I8)') ngtot_
51 write(*,*)
52 stop
53end if
54read(100) rhosmt,rhosir
55read(100) chgsmt
56read(100) chgstot
57close(100)
58end subroutine
59
integer ngtot
Definition modmain.f90:390
integer natmtot
Definition modmain.f90:40
integer npmtmax
Definition modmain.f90:216
real(8), dimension(:,:,:), allocatable rhosmt
Definition modtddft.f90:82
real(8), dimension(:,:), allocatable chgsmt
Definition modtddft.f90:86
real(8), dimension(:,:), allocatable rhosir
Definition modtddft.f90:82
real(8), dimension(3) chgstot
Definition modtddft.f90:84
subroutine readrhos
Definition readrhos.f90:7