The Elk Code
readdvs.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2008 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
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 readdvs(iq,is,ia,ip,dvsmt,dvsir)
7 use modmain
8 implicit none
9 ! arguments
10 integer, intent(in) :: iq,is,ia,ip
11 complex(8), intent(out) :: dvsmt(npmtmax,natmtot),dvsir(ngtot)
12 ! local variables
13 integer js,jas,ios
14 integer version_(3),nspecies_
15 integer lmmaxo_,natoms_
16 integer nrmt_,ngridg_(3)
17 character(256) fext,fname
18 ! allocatable arrays
19 complex(8), allocatable :: zfmt(:,:,:)
20 call dynfext(iq,is,ia,ip,fext)
21 fname='DVS'//trim(fext)
22 open(150,file=trim(fname),form='UNFORMATTED',action='READ',status='OLD', &
23  iostat=ios)
24 if (ios /= 0) then
25  write(*,*)
26  write(*,'("Error(readdvs): error opening ",A)') trim(fname)
27  write(*,*)
28  stop
29 end if
30 read(150) version_
31 if ((version(1) /= version_(1)).or.(version(2) /= version_(2)) &
32  .or.(version(3) /= version_(3))) then
33  write(*,*)
34  write(*,'("Warning(readdvs): different versions")')
35  write(*,'(" current : ",I0,".",I0,".",I0)') version
36  write(*,'(" file : ",I0,".",I0,".",I0)') version_
37  write(*,'(" in file ",A)') trim(fname)
38 end if
39 read(150) nspecies_
40 if (nspecies /= nspecies_) then
41  write(*,*)
42  write(*,'("Error(readdvs): differing nspecies")')
43  write(*,'(" current : ",I4)') nspecies
44  write(*,'(" file : ",I4)') nspecies_
45  write(*,'(" in file ",A)') trim(fname)
46  write(*,*)
47  stop
48 end if
49 read(150) lmmaxo_
50 if (lmmaxo /= lmmaxo_) then
51  write(*,*)
52  write(*,'("Error(readdvs): differing lmmaxo")')
53  write(*,'(" current : ",I4)') lmmaxo
54  write(*,'(" file : ",I4)') lmmaxo_
55  write(*,'(" in file ",A)') trim(fname)
56  write(*,*)
57  stop
58 end if
59 do js=1,nspecies
60  read(150) natoms_
61  if (natoms(js) /= natoms_) then
62  write(*,*)
63  write(*,'("Error(readdvs): differing natoms for species ",I4)') js
64  write(*,'(" current : ",I4)') natoms(js)
65  write(*,'(" file : ",I4)') natoms_
66  write(*,'(" in file ",A)') trim(fname)
67  write(*,*)
68  stop
69  end if
70  read(150) nrmt_
71  if (nrmt(js) /= nrmt_) then
72  write(*,*)
73  write(*,'("Error(readdvs): differing nrmt for species ",I4)') js
74  write(*,'(" current : ",I6)') nrmt(js)
75  write(*,'(" file : ",I6)') nrmt_
76  write(*,'(" in file ",A)') trim(fname)
77  write(*,*)
78  stop
79  end if
80 end do
81 read(150) ngridg_
82 if ((ngridg(1) /= ngridg_(1)).or.(ngridg(2) /= ngridg_(2)).or. &
83  (ngridg(3) /= ngridg_(3))) then
84  write(*,*)
85  write(*,'("Error(readdvs): differing ngridg")')
86  write(*,'(" current : ",3I6)') ngridg
87  write(*,'(" file : ",3I6)') ngridg_
88  write(*,'(" in file ",A)') trim(fname)
89  write(*,*)
90  stop
91 end if
92 allocate(zfmt(lmmaxo,nrmtmax,natmtot))
93 read(150) zfmt,dvsir
94 do jas=1,natmtot
95  js=idxis(jas)
96  call zfmtpack(.true.,nrmt(js),nrmti(js),zfmt(:,:,jas),dvsmt(:,jas))
97 end do
98 close(150)
99 deallocate(zfmt)
100 end subroutine
101 
integer, dimension(3) ngridg
Definition: modmain.f90:386
integer lmmaxo
Definition: modmain.f90:203
subroutine dynfext(iq, is, ia, ip, fext)
Definition: dynfext.f90:7
subroutine readdvs(iq, is, ia, ip, dvsmt, dvsir)
Definition: readdvs.f90:7
pure subroutine zfmtpack(tpack, nr, nri, zfmt1, zfmt2)
Definition: zfmtpack.f90:7
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
integer, dimension(3), parameter version
Definition: modmain.f90:1289
integer nspecies
Definition: modmain.f90:34
integer nrmtmax
Definition: modmain.f90:152
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150