The Elk Code
readstate.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2005 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 !BOP
7 ! !ROUTINE: readstate
8 ! !INTERFACE:
9 subroutine readstate
10 ! !USES:
11 use modmain
12 use moddftu
13 ! !DESCRIPTION:
14 ! Reads in the charge density and other relevant variables from the file
15 ! {\tt STATE.OUT}. Checks for version and parameter compatibility.
16 !
17 ! !REVISION HISTORY:
18 ! Created May 2003 (JKD)
19 !EOP
20 !BOC
21 implicit none
22 ! local variables
23 logical spinpol_
24 integer is,ia,ias,lmmax,lm,ir,jr
25 integer idm,jdm,mapidm(3),ios
26 integer i1,i2,i3,j1,j2,j3,n
27 integer version_(3)
28 integer nspecies_,natoms_,lmmaxo_
29 integer nrmt_(maxspecies),nrmtmax_
30 integer nrcmt_(maxspecies),nrcmtmax_
31 integer ngridg_(3),ngtot_,ngvec_
32 integer ndmag_,nspinor_,fsmtype_,ftmtype_
33 integer dftu_,lmmaxdm_,xcgrad_
34 real(8) t1
35 ! allocatable arrays
36 integer, allocatable :: mapir(:)
37 real(8), allocatable :: rsp_(:,:),rcmt_(:,:)
38 real(8), allocatable :: wcrmt_(:,:,:),wcrcmt_(:,:,:)
39 real(8), allocatable :: rfmt_(:,:,:),rfir_(:)
40 real(8), allocatable :: rvfmt_(:,:,:,:),rvfir_(:,:)
41 real(8), allocatable :: rvfcmt_(:,:,:,:),rfmt(:,:)
42 real(8), allocatable :: bfsmcmt_(:,:),fi(:),fo(:)
43 complex(8), allocatable :: vsig_(:)
44 complex(8), allocatable :: vmatmt_(:,:,:,:,:),vmftm_(:,:,:,:,:)
45 open(100,file='STATE'//trim(filext),form='UNFORMATTED',action='READ', &
46  status='OLD',iostat=ios)
47 if (ios /= 0) then
48  write(*,*)
49  write(*,'("Error(readstate): error opening ",A)') 'STATE'//trim(filext)
50  write(*,*)
51  stop
52 end if
53 read(100) version_
54 if (version_(1) < 2) then
55  write(*,*)
56  write(*,'("Error(readstate): unable to read STATE.OUT from versions earlier &
57  &than 2.0.0")')
58  write(*,*)
59  stop
60 end if
61 if (any(version(:) /= version_(:))) then
62  write(*,*)
63  write(*,'("Warning(readstate): different versions")')
64  write(*,'(" current : ",I0,".",I0,".",I0)') version
65  write(*,'(" STATE.OUT : ",I0,".",I0,".",I0)') version_
66 end if
67 read(100) spinpol_
68 read(100) nspecies_
69 if (nspecies /= nspecies_) then
70  write(*,*)
71  write(*,'("Error(readstate): differing nspecies")')
72  write(*,'(" current : ",I4)') nspecies
73  write(*,'(" STATE.OUT : ",I4)') nspecies_
74  write(*,*)
75  stop
76 end if
77 read(100) lmmaxo_
78 lmmax=min(lmmaxo,lmmaxo_)
79 read(100) nrmtmax_
80 read(100) nrcmtmax_
81 allocate(rsp_(nrmtmax_,nspecies))
82 allocate(rcmt_(nrcmtmax_,nspecies))
83 do is=1,nspecies
84  read(100) natoms_
85  if (natoms(is) /= natoms_) then
86  write(*,*)
87  write(*,'("Error(readstate): differing natoms for species ",I4)') is
88  write(*,'(" current : ",I4)') natoms(is)
89  write(*,'(" STATE.OUT : ",I4)') natoms_
90  write(*,*)
91  stop
92  end if
93  read(100) nrmt_(is)
94  read(100) rsp_(1:nrmt_(is),is)
95  read(100) nrcmt_(is)
96  read(100) rcmt_(1:nrcmt_(is),is)
97 end do
98 read(100) ngridg_
99 read(100) ngvec_
100 read(100) ndmag_
101 if (spinpol_.and.(ndmag_ /= 1).and.(ndmag_ /= 3)) then
102  write(*,*)
103  write(*,'("Error(readstate): invalid ndmag in STATE.OUT : ",I8)') ndmag_
104  write(*,*)
105  stop
106 end if
107 read(100) nspinor_
108 read(100) fsmtype_
109 if (version_(1) > 2) then
110  read(100) ftmtype_
111 else
112  ftmtype_=0
113 end if
114 read(100) dftu_
115 read(100) lmmaxdm_
116 if ((version_(1) > 5).or.((version_(1) == 5).and.(version_(2) > 0))) then
117  read(100) xcgrad_
118 else
119  xcgrad_=0
120 end if
121 if ((version_(1) > 9).or.((version_(1) == 9).and.(version_(2) > 5))) then
122  read(100) efermi
123 else
124  call readefm
125 end if
126 if ((version_(1) > 10).or.((version_(1) == 10).and.(version_(2) > 6))) then
127  read(100) dlefe
128 end if
129 ngtot_=ngridg_(1)*ngridg_(2)*ngridg_(3)
130 ! map from old interstitial grid to new
131 allocate(mapir(ngtot))
132 do i3=0,ngridg(3)-1
133  t1=dble(i3*ngridg_(3))/dble(ngridg(3))
134  j3=modulo(nint(t1),ngridg_(3))
135  do i2=0,ngridg(2)-1
136  t1=dble(i2*ngridg_(2))/dble(ngridg(2))
137  j2=modulo(nint(t1),ngridg_(2))
138  do i1=0,ngridg(1)-1
139  t1=dble(i1*ngridg_(1))/dble(ngridg(1))
140  j1=modulo(nint(t1),ngridg_(1))
141  jr=j3*ngridg_(2)*ngridg_(1)+j2*ngridg_(1)+j1+1
142  ir=i3*ngridg(2)*ngridg(1)+i2*ngridg(1)+i1+1
143  mapir(ir)=jr
144  end do
145  end do
146 end do
147 ! determine the spline coefficient weights on the old radial mesh
148 allocate(wcrmt_(12,nrmtmax_,nspecies))
149 allocate(wcrcmt_(12,nrcmtmax_,nspecies))
150 do is=1,nspecies
151  call wspline(nrmt_(is),rsp_(:,is),wcrmt_(:,:,is))
152  call wspline(nrcmt_(is),rcmt_(:,is),wcrcmt_(:,:,is))
153 end do
154 allocate(rfmt_(lmmaxo_,nrmtmax_,natmtot),rfir_(ngtot_))
155 allocate(rfmt(lmmaxo,nrmtmax))
156 n=max(nrmtmax,nrmtmax_)
157 allocate(fi(n),fo(n))
158 ! read the muffin-tin density
159 read(100) rfmt_,rfir_
160 ! regrid and pack the muffin-tin function
161 call rgfmt(rhomt)
162 ! regrid the interstitial function
163 rhoir(1:ngtot)=rfir_(mapir(1:ngtot))
164 ! read the Coulomb potential, regrid and pack
165 read(100) rfmt_,rfir_
166 call rgfmt(vclmt)
167 vclir(1:ngtot)=rfir_(mapir(1:ngtot))
168 ! read the exchange-correlation potential, regrid and pack
169 read(100) rfmt_,rfir_
170 call rgfmt(vxcmt)
171 vxcir(1:ngtot)=rfir_(mapir(1:ngtot))
172 ! read the Kohn-Sham effective potential, regrid and pack
173 if (version_(1) > 2) then
174  read(100) rfmt_,rfir_
175 else
176  allocate(vsig_(ngvec_))
177  read(100) rfmt_,rfir_,vsig_
178  deallocate(vsig_)
179 end if
180 call rgfmt(vsmt)
181 vsir(1:ngtot)=rfir_(mapir(1:ngtot))
182 ! compute the Kohn-Sham potential on coarse interstitial grid
183 call rfirftoc(vsir,vsirc)
184 ! read the magnetisation, exchange-correlation and effective magnetic fields
185 if (spinpol_) then
186 ! component map for spin-polarised case
187  mapidm(:)=0
188  if (ndmag == ndmag_) then
189  do idm=1,ndmag
190  mapidm(idm)=idm
191  end do
192  else
193  mapidm(ndmag)=ndmag_
194  end if
195  allocate(rvfmt_(lmmaxo_,nrmtmax_,natmtot,ndmag_))
196  allocate(rvfir_(ngtot_,ndmag_))
197  allocate(rvfcmt_(lmmaxo_,nrcmtmax_,natmtot,ndmag_))
198  read(100) rvfmt_,rvfir_
199  call rgvfmt(magmt)
200  call rgvir(magir)
201  read(100) rvfmt_,rvfir_
202  call rgvfmt(bxcmt)
203  call rgvir(bxcir)
204  read(100) rvfcmt_,rvfir_
205  call rgvfcmt(bsmt)
206  call rgvir(bsir)
207 ! compute the Kohn-Sham magnetic field on coarse interstitial grid
208  do idm=1,ndmag
209  call rfirftoc(bsir(:,idm),bsirc(:,idm))
210  end do
211  deallocate(rvfmt_,rvfir_,rvfcmt_)
212 ! read fixed spin moment effective fields
213  if (fsmtype_ /= 0) then
214  allocate(bfsmcmt_(3,natmtot))
215  read(100) bfsmc
216  read(100) bfsmcmt_
217  if (fsmtype /= 0) bfsmcmt(:,:)=bfsmcmt_(:,:)
218 ! make sure that the constraining fields are perpendicular to the fixed moments
219 ! for fixed direction calculations (Y. Kvashnin and LN)
220  if (fsmtype < 0) then
221  if (ncmag) then
222  call r3vo(momfix,bfsmc)
223  do is=1,nspecies
224  do ia=1,natoms(is)
225  ias=idxas(ia,is)
226  call r3vo(mommtfix(:,ia,is),bfsmcmt(:,ias))
227  end do
228  end do
229  else
230  bfsmc(:)=0.d0
231  bfsmcmt(:,:)=0.d0
232  end if
233  end if
234  deallocate(bfsmcmt_)
235  end if
236 end if
237 if (any(xcgrad == [3,4,5,6])) then
238  if (any(xcgrad_ == [3,4,5,6])) then
239  read(100) rfmt_,rfir_
240  call rgfmt(wxcmt)
241  wxcir(1:ngtot)=rfir_(mapir(1:ngtot))
242  else
243  wxcmt(:,:)=0.d0
244  wxcir(1:ngtot)=0.d0
245  end if
246 end if
247 deallocate(wcrmt_,wcrcmt_,rfmt_,rfir_,rfmt,fi,fo)
248 ! read DFT+U potential matrix in each muffin-tin
249 if (((dftu /= 0).and.(dftu_ /= 0)).or. &
250  ((ftmtype /= 0).and.(ftmtype_ /= 0))) then
251  allocate(vmatmt_(lmmaxdm_,nspinor_,lmmaxdm_,nspinor_,natmtot))
252  read(100) vmatmt_
253  lmmax=min(lmmaxdm,lmmaxdm_)
254  vmatmt(:,:,:,:,:)=0.d0
255  if (nspinor == nspinor_) then
256  vmatmt(1:lmmax,:,1:lmmax,:,:)=vmatmt_(1:lmmax,:,1:lmmax,:,:)
257  else if ((nspinor == 1).and.(nspinor_ == 2)) then
258  vmatmt(1:lmmax,1,1:lmmax,1,:)=0.5d0*(vmatmt_(1:lmmax,1,1:lmmax,1,:) &
259  +vmatmt_(1:lmmax,2,1:lmmax,2,:))
260  else
261  vmatmt(1:lmmax,1,1:lmmax,1,:)=vmatmt_(1:lmmax,1,1:lmmax,1,:)
262  vmatmt(1:lmmax,2,1:lmmax,2,:)=vmatmt_(1:lmmax,1,1:lmmax,1,:)
263  end if
264  deallocate(vmatmt_)
265 end if
266 ! read fixed tensor moment potential matrix elements
267 if ((ftmtype /= 0).and.(ftmtype_ /= 0)) then
268  allocate(vmftm_(lmmaxdm_,nspinor_,lmmaxdm_,nspinor_,natmtot))
269  read(100) vmftm_
270  lmmax=min(lmmaxdm,lmmaxdm_)
271  vmftm_(:,:,:,:,:)=0.d0
272  if (nspinor == nspinor_) then
273  vmftm(1:lmmax,:,1:lmmax,:,:)=vmftm_(1:lmmax,:,1:lmmax,:,:)
274  else if ((nspinor == 1).and.(nspinor_ == 2)) then
275  vmftm(1:lmmax,1,1:lmmax,1,:)=0.5d0*(vmftm_(1:lmmax,1,1:lmmax,1,:) &
276  +vmftm_(1:lmmax,2,1:lmmax,2,:))
277  else
278  vmftm(1:lmmax,1,1:lmmax,1,:)=vmftm_(1:lmmax,1,1:lmmax,1,:)
279  vmftm(1:lmmax,2,1:lmmax,2,:)=vmftm_(1:lmmax,1,1:lmmax,1,:)
280  end if
281  deallocate(vmftm_)
282 end if
283 close(100)
284 
285 contains
286 
287 subroutine rgfmt(rfmtp)
288 implicit none
289 ! arguments
290 real(8), intent(out) :: rfmtp(npmtmax,natmtot)
291 do ias=1,natmtot
292  is=idxis(ias)
293 ! regrid the muffin-tin function
294  do lm=1,lmmax
295  fi(1:nrmt_(is))=rfmt_(lm,1:nrmt_(is),ias)
296  call rfinterp(nrmt_(is),rsp_(:,is),wcrmt_(:,:,is),fi,nrmt(is),rsp(:,is),fo)
297  rfmt(lm,1:nrmt(is))=fo(1:nrmt(is))
298  end do
299  rfmt(lmmax+1:lmmaxo,1:nrmt(is))=0.d0
300 ! pack the muffin-tin function
301  call rfmtpack(.true.,nrmt(is),nrmti(is),rfmt,rfmtp(:,ias))
302 end do
303 end subroutine
304 
305 subroutine rgvfmt(rvfmt)
306 implicit none
307 ! arguments
308 real(8), intent(out) :: rvfmt(npmtmax,natmtot,ndmag)
309 do idm=1,ndmag
310  jdm=mapidm(idm)
311  if (jdm == 0) then
312  rvfmt(:,:,idm)=0.d0
313  cycle
314  end if
315  do ias=1,natmtot
316  is=idxis(ias)
317  do lm=1,lmmax
318  fi(1:nrmt_(is))=rvfmt_(lm,1:nrmt_(is),ias,jdm)
319  call rfinterp(nrmt_(is),rsp_(:,is),wcrmt_(:,:,is),fi,nrmt(is),rsp(:,is), &
320  fo)
321  rfmt(lm,1:nrmt(is))=fo(1:nrmt(is))
322  end do
323  rfmt(lmmax+1:lmmaxo,1:nrmt(is))=0.d0
324  call rfmtpack(.true.,nrmt(is),nrmti(is),rfmt,rvfmt(:,ias,idm))
325  end do
326 end do
327 end subroutine
328 
329 subroutine rgvfcmt(rvfcmt)
330 implicit none
331 ! arguments
332 real(8), intent(out) :: rvfcmt(npcmtmax,natmtot,ndmag)
333 do idm=1,ndmag
334  jdm=mapidm(idm)
335  if (jdm == 0) then
336  rvfcmt(:,:,idm)=0.d0
337  cycle
338  end if
339  do ias=1,natmtot
340  is=idxis(ias)
341  do lm=1,lmmax
342  fi(1:nrcmt_(is))=rvfcmt_(lm,1:nrcmt_(is),ias,jdm)
343  call rfinterp(nrcmt_(is),rcmt_(:,is),wcrcmt_(:,:,is),fi,nrcmt(is), &
344  rcmt(:,is),fo)
345  rfmt(lm,1:nrcmt(is))=fo(1:nrcmt(is))
346  end do
347  rfmt(lmmax+1:lmmaxo,1:nrcmt(is))=0.d0
348  call rfmtpack(.true.,nrcmt(is),nrcmti(is),rfmt,rvfcmt(:,ias,idm))
349  end do
350 end do
351 end subroutine
352 
353 subroutine rgvir(rvfir)
354 implicit none
355 ! arguments
356 real(8), intent(out) :: rvfir(ngtot,ndmag)
357 do idm=1,ndmag
358  jdm=mapidm(idm)
359  if (jdm /= 0) then
360  rvfir(1:ngtot,idm)=rvfir_(mapir(1:ngtot),jdm)
361  else
362  rvfir(1:ngtot,idm)=0.d0
363  end if
364 end do
365 end subroutine
366 
367 end subroutine
368 !EOC
369 
real(8), dimension(:,:), allocatable rcmt
Definition: modmain.f90:177
real(8), dimension(:), allocatable wxcir
Definition: modmain.f90:676
real(8) efermi
Definition: modmain.f90:907
integer, dimension(3) ngridg
Definition: modmain.f90:386
character(256) filext
Definition: modmain.f90:1301
integer ngtot
Definition: modmain.f90:390
integer lmmaxo
Definition: modmain.f90:203
real(8), dimension(:,:), allocatable vxcmt
Definition: modmain.f90:634
real(8), dimension(:), pointer, contiguous rhoir
Definition: modmain.f90:614
subroutine rgfmt(rfmtp)
Definition: readstate.f90:288
subroutine rgvfmt(rvfmt)
Definition: readstate.f90:306
integer, parameter lmmaxdm
Definition: moddftu.f90:14
integer, dimension(maxatoms, maxspecies) idxas
Definition: modmain.f90:42
integer ndmag
Definition: modmain.f90:238
real(8), dimension(:,:), allocatable vclmt
Definition: modmain.f90:624
real(8), dimension(:,:), pointer, contiguous rhomt
Definition: modmain.f90:614
real(8), dimension(3) bfsmc
Definition: modmain.f90:257
subroutine rfinterp(ni, xi, wci, fi, no, xo, fo)
Definition: rfinterp.f90:10
real(8), dimension(:), allocatable vsir
Definition: modmain.f90:651
real(8), dimension(3) momfix
Definition: modmain.f90:253
real(8), dimension(:), allocatable vxcir
Definition: modmain.f90:634
real(8), dimension(:,:), allocatable wxcmt
Definition: modmain.f90:676
subroutine wspline(n, x, wc)
Definition: wspline.f90:7
subroutine rgvir(rvfir)
Definition: readstate.f90:354
subroutine rgvfcmt(rvfcmt)
Definition: readstate.f90:330
subroutine readefm
Definition: readefm.f90:10
real(8), dimension(:,:,:), allocatable bxcmt
Definition: modmain.f90:636
real(8), dimension(:,:), allocatable bsir
Definition: modmain.f90:658
real(8), dimension(:,:,:), pointer, contiguous magmt
Definition: modmain.f90:616
integer ftmtype
Definition: moddftu.f90:70
pure subroutine rfmtpack(tpack, nr, nri, rfmt1, rfmt2)
Definition: rfmtpack.f90:7
real(8), dimension(:), allocatable vclir
Definition: modmain.f90:624
integer nspinor
Definition: modmain.f90:267
complex(8), dimension(:,:,:,:,:), allocatable vmatmt
Definition: moddftu.f90:20
real(8), dimension(3, maxatoms, maxspecies) mommtfix
Definition: modmain.f90:259
real(8), dimension(:,:), allocatable bxcir
Definition: modmain.f90:636
complex(8), dimension(:,:,:,:,:), allocatable vmftm
Definition: moddftu.f90:82
real(8), dimension(:), pointer, contiguous vsirc
Definition: modmain.f90:653
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
real(8), dimension(:,:), allocatable bfsmcmt
Definition: modmain.f90:263
integer dftu
Definition: moddftu.f90:32
subroutine readstate
Definition: readstate.f90:10
real(8) dlefe
Definition: modmain.f90:830
integer, dimension(3), parameter version
Definition: modmain.f90:1289
integer nspecies
Definition: modmain.f90:34
real(8), dimension(:,:), pointer, contiguous magir
Definition: modmain.f90:616
real(8), dimension(:,:), allocatable rsp
Definition: modmain.f90:135
integer natmtot
Definition: modmain.f90:40
integer, dimension(maxspecies) nrcmt
Definition: modmain.f90:173
integer, dimension(maxspecies) nrcmti
Definition: modmain.f90:211
pure subroutine r3vo(x, y)
Definition: r3vo.f90:7
integer nrmtmax
Definition: modmain.f90:152
logical ncmag
Definition: modmain.f90:240
real(8), dimension(:,:), pointer, contiguous bsirc
Definition: modmain.f90:660
real(8), dimension(:,:), pointer, contiguous vsmt
Definition: modmain.f90:649
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
subroutine rfirftoc(rfir, rfirc)
Definition: rfirftoc.f90:7
integer xcgrad
Definition: modmain.f90:602
real(8), dimension(:,:,:), pointer, contiguous bsmt
Definition: modmain.f90:656
integer fsmtype
Definition: modmain.f90:251
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150