The Elk Code
writestate.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: writestate
8 ! !INTERFACE:
9 subroutine writestate
10 ! !USES:
11 use modmain
12 use moddftu
13 ! !DESCRIPTION:
14 ! Writes the charge density, potentials and other relevant variables to the
15 ! file {\tt STATE.OUT}. Note to developers: changes to the way the variables
16 ! are written should be mirrored in {\tt readstate}.
17 !
18 ! !REVISION HISTORY:
19 ! Created May 2003 (JKD)
20 !EOP
21 !BOC
22 implicit none
23 ! local variables
24 integer idm,is,ias
25 ! allocatable arrays
26 real(8), allocatable :: rfmt(:,:,:),rvfmt(:,:,:,:),rvfcmt(:,:,:,:)
27 open(100,file='STATE'//trim(filext),form='UNFORMATTED',action='WRITE')
28 write(100) version
29 write(100) spinpol
30 write(100) nspecies
31 write(100) lmmaxo
32 write(100) nrmtmax
33 write(100) nrcmtmax
34 do is=1,nspecies
35  write(100) natoms(is)
36  write(100) nrmt(is)
37  write(100) rsp(1:nrmt(is),is)
38  write(100) nrcmt(is)
39  write(100) rcmt(1:nrcmt(is),is)
40 end do
41 write(100) ngridg
42 write(100) ngvec
43 write(100) ndmag
44 write(100) nspinor
45 write(100) fsmtype
46 write(100) ftmtype
47 write(100) dftu
48 write(100) lmmaxdm
49 write(100) xcgrad
50 write(100) efermi
51 write(100) dlefe
52 ! muffin-tin functions are unpacked to maintain backward compatibility
53 allocate(rfmt(lmmaxo,nrmtmax,natmtot))
54 if (spinpol) then
55  allocate(rvfmt(lmmaxo,nrmtmax,natmtot,ndmag))
56  allocate(rvfcmt(lmmaxo,nrcmtmax,natmtot,ndmag))
57 end if
58 ! write the density
59 do ias=1,natmtot
60  is=idxis(ias)
61  call rfmtpack(.false.,nrmt(is),nrmti(is),rhomt(:,ias),rfmt(:,:,ias))
62 end do
63 write(100) rfmt,rhoir
64 ! write the Coulomb potential
65 do ias=1,natmtot
66  is=idxis(ias)
67  call rfmtpack(.false.,nrmt(is),nrmti(is),vclmt(:,ias),rfmt(:,:,ias))
68 end do
69 write(100) rfmt,vclir
70 ! write the exchange-correlation potential
71 do ias=1,natmtot
72  is=idxis(ias)
73  call rfmtpack(.false.,nrmt(is),nrmti(is),vxcmt(:,ias),rfmt(:,:,ias))
74 end do
75 write(100) rfmt,vxcir
76 ! write the Kohn-Sham effective potential
77 do ias=1,natmtot
78  is=idxis(ias)
79  call rfmtpack(.false.,nrmt(is),nrmti(is),vsmt(:,ias),rfmt(:,:,ias))
80 end do
81 write(100) rfmt,vsir
82 if (spinpol) then
83 ! write the magnetisation, exchange-correlation and effective magnetic fields
84  do idm=1,ndmag
85  do ias=1,natmtot
86  is=idxis(ias)
87  call rfmtpack(.false.,nrmt(is),nrmti(is),magmt(:,ias,idm), &
88  rvfmt(:,:,ias,idm))
89  end do
90  end do
91  write(100) rvfmt,magir
92  do idm=1,ndmag
93  do ias=1,natmtot
94  is=idxis(ias)
95  call rfmtpack(.false.,nrmt(is),nrmti(is),bxcmt(:,ias,idm), &
96  rvfmt(:,:,ias,idm))
97  end do
98  end do
99  write(100) rvfmt,bxcir
100  do idm=1,ndmag
101  do ias=1,natmtot
102  is=idxis(ias)
103  call rfmtpack(.false.,nrcmt(is),nrcmti(is),bsmt(:,ias,idm), &
104  rvfcmt(:,:,ias,idm))
105  end do
106  end do
107  write(100) rvfcmt,bsir
108 ! write fixed spin moment magnetic fields
109  if (fsmtype /= 0) then
110  write(100) bfsmc
111  write(100) bfsmcmt
112  end if
113 end if
114 ! write the meta-GGA exchange-correlation potential
115 if (any(xcgrad == [3,4,5,6])) then
116  do ias=1,natmtot
117  is=idxis(ias)
118  call rfmtpack(.false.,nrmt(is),nrmti(is),wxcmt(:,ias),rfmt(:,:,ias))
119  end do
120  write(100) rfmt,wxcir
121 end if
122 ! write the potential matrix in each muffin-tin
123 if ((dftu /= 0).or.(ftmtype /= 0)) then
124  write(100) vmatmt
125 end if
126 ! write the fixed tensor moment potential matrix
127 if (ftmtype /= 0) then
128  write(100) vmftm
129 end if
130 close(100)
131 deallocate(rfmt)
132 if (spinpol) deallocate(rvfmt,rvfcmt)
133 end subroutine
134 !EOC
135 
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 lmmaxo
Definition: modmain.f90:203
real(8), dimension(:,:), allocatable vxcmt
Definition: modmain.f90:634
logical spinpol
Definition: modmain.f90:228
real(8), dimension(:), pointer, contiguous rhoir
Definition: modmain.f90:614
integer, parameter lmmaxdm
Definition: moddftu.f90:14
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
subroutine writestate
Definition: writestate.f90:10
real(8), dimension(3) bfsmc
Definition: modmain.f90:257
real(8), dimension(:), allocatable vsir
Definition: modmain.f90:651
real(8), dimension(:), allocatable vxcir
Definition: modmain.f90:634
integer nrcmtmax
Definition: modmain.f90:175
real(8), dimension(:,:), allocatable wxcmt
Definition: modmain.f90:676
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
integer ngvec
Definition: modmain.f90:396
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(:,:), allocatable bxcir
Definition: modmain.f90:636
complex(8), dimension(:,:,:,:,:), allocatable vmftm
Definition: moddftu.f90:82
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
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
integer nrmtmax
Definition: modmain.f90:152
real(8), dimension(:,:), pointer, contiguous vsmt
Definition: modmain.f90:649
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
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