The Elk Code
 
Loading...
Searching...
No Matches
writetddft.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
6subroutine writetddft
7use modmain
8use modtddft
9use moddftu
10use moddelf
11implicit none
12! local variables
13integer is,ia,ias
14character(32) fext
15! allocatable arrays
16real(8), allocatable :: rvfmt(:,:,:),rvfir(:,:)
17! file extension
18write(fext,'("_TS",I8.8,".OUT")') itimes
19! delete files at first time step
20if (itimes <= 1) then
21 call delfile('TOTENERGY_TD.OUT')
22 call delfile('CHARGEMT_TD.OUT')
23 call delfile('CHARGEIR_TD.OUT')
24 if (spinpol) then
25 call delfile('MOMENT_TD.OUT')
26 call delfile('MOMENTM_TD.OUT')
27 call delfile('MOMENTMT_TD.OUT')
28 call delfile('MOMENTIR_TD.OUT')
29 end if
30 call delfile('JTOT_TD.OUT')
31 call delfile('JTOTM_TD.OUT')
32 if (tddos) call delfile('TDTEMP.OUT')
33 if (tdlsj) call delfile('TDLSJ.OUT')
34 if (tforce) call delfile('FORCETOT_TD.OUT')
35 if (tatdisp) then
36 call delfile('ATDISPL_TD.OUT')
37 call delfile('ATDISPC_TD.OUT')
38 end if
39 if (tafindt) call delfile('AFIND_TD.OUT')
40end if
41! total energy
42call writetdengy
43! muffin-tin charges
44open(50,file='CHARGEMT_TD.OUT',form='FORMATTED',position='APPEND')
45write(50,'(G18.10)') times(itimes)
46do is=1,nspecies
47 do ia=1,natoms(is)
48 ias=idxas(ia,is)
49 write(50,'(2I4,G18.10)') is,ia,chgmt(ias)
50 end do
51end do
52write(50,*)
53close(50)
54! interstitial charge
55open(50,file='CHARGEIR_TD.OUT',form='FORMATTED',position='APPEND')
56write(50,'(2G18.10)') times(itimes),chgir
57close(50)
58! write spin moments to file if required
59if (spinpol) call writemomtd
60! total current
61open(50,file='JTOT_TD.OUT',form='FORMATTED',position='APPEND')
62write(50,'(4G18.10)') times(itimes),jtot(:)
63close(50)
64! total current magnitude
65open(50,file='JTOTM_TD.OUT',form='FORMATTED',position='APPEND')
66write(50,'(4G18.10)') times(itimes),jtotm
67close(50)
68! write the time-dependent atomic forces
69if (tforce.and.ttsforce) then
70 call writetdforces
71 open(50,file='FORCES'//trim(fext),form='FORMATTED',action='WRITE')
72 call writeforces(50)
73 close(50)
74end if
75! write the time-dependent atomic displacements
76if (tatdisp.and.ttsforce) call writeatdisp
77! write the time-dependent induced A-field
78if (tafindt) then
79 open(50,file='AFIND_TD.OUT',form='FORMATTED',position='APPEND')
80 write(50,'(4G18.10)') times(itimes),afindt(:,0)
81 close(50)
82end if
83! write optional quantities
84if (.not.ttswrite) return
85! charge density in 1D
86if (tdrho1d) then
87 open(50,file='RHO1D'//trim(fext),form='FORMATTED',action='WRITE')
88 open(51,file='RHOLINES.OUT',form='FORMATTED',action='WRITE')
89 call plot1d(50,51,1,rhomt,rhoir)
90 close(50)
91 close(51)
92end if
93! charge density in 2D
94if (tdrho2d) then
95 open(50,file='RHO2D'//trim(fext),form='FORMATTED',action='WRITE')
96 call plot2d(.false.,50,1,rhomt,rhoir)
97 close(50)
98end if
99! charge density in 3D
100if (tdrho3d) then
101 open(50,file='RHO3D'//trim(fext),form='FORMATTED',action='WRITE')
102 call plot3d(50,1,rhomt,rhoir)
103 close(50)
104end if
105! magnetisation in 1D, 2D or 3D
106if ((tdmag1d.or.tdmag2d.or.tdmag3d).and.spinpol) then
107 allocate(rvfmt(npmtmax,natmtot,3),rvfir(ngtot,3))
108 if (ncmag) then
109! non-collinear
110 rvfmt(:,:,:)=magmt(:,:,:)
111 rvfir(:,:)=magir(:,:)
112 else
113! collinear
114 rvfmt(:,:,1:2)=0.d0
115 rvfir(:,1:2)=0.d0
116 rvfmt(:,:,3)=magmt(:,:,1)
117 rvfir(:,3)=magir(:,1)
118 end if
119 if (tdmag1d) then
120 open(50,file='MAG1D'//trim(fext),form='FORMATTED',action='WRITE')
121 open(51,file='MAGLINES.OUT',form='FORMATTED',action='WRITE')
122 call plot1d(50,51,3,rvfmt,rvfir)
123 close(50)
124 close(51)
125 end if
126 if (tdmag2d) then
127 open(50,file='MAG2D'//trim(fext),form='FORMATTED',action='WRITE')
128 call plot2d(.true.,50,3,rvfmt,rvfir)
129 close(50)
130 end if
131 if (tdmag3d) then
132 open(50,file='MAG3D'//trim(fext),form='FORMATTED',action='WRITE')
133 call plot3d(50,3,rvfmt,rvfir)
134 close(50)
135 end if
136 deallocate(rvfmt,rvfir)
137end if
138! gauge-invariant current density in 1D
139if (tdjr1d) then
140 open(50,file='JR1D'//trim(fext),form='FORMATTED',action='WRITE')
141 open(51,file='JRLINES.OUT',form='FORMATTED',action='WRITE')
142 call plot1d(50,51,3,jrmt,jrir)
143 close(50)
144 close(51)
145end if
146! gauge-invariant current density in 2D
147if (tdjr2d) then
148 open(50,file='JR2D'//trim(fext),form='FORMATTED',action='WRITE')
149 call plot2d(.true.,50,3,jrmt,jrir)
150 close(50)
151end if
152! gauge-invariant current density in 3D
153if (tdjr3d) then
154 open(50,file='JR3D'//trim(fext),form='FORMATTED',action='WRITE')
155 call plot3d(50,3,jrmt,jrir)
156 close(50)
157end if
158! calculate and write tensor moments
159if (dftu /= 0) then
160 if (tmwrite) call writetdtm3
161end if
162end subroutine
163
subroutine delfile(fname)
Definition moddelf.f90:15
logical tmwrite
Definition moddftu.f90:66
integer dftu
Definition moddftu.f90:32
real(8), dimension(3) jtot
Definition modmain.f90:748
integer ngtot
Definition modmain.f90:390
logical ncmag
Definition modmain.f90:240
real(8), dimension(:,:,:), pointer, contiguous magmt
Definition modmain.f90:616
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
real(8), dimension(:), pointer, contiguous rhoir
Definition modmain.f90:614
integer, dimension(maxatoms, maxspecies) idxas
Definition modmain.f90:42
logical spinpol
Definition modmain.f90:228
real(8) jtotm
Definition modmain.f90:748
integer natmtot
Definition modmain.f90:40
real(8), dimension(:,:), pointer, contiguous magir
Definition modmain.f90:616
logical tatdisp
Definition modmain.f90:59
real(8), dimension(:), allocatable chgmt
Definition modmain.f90:732
integer npmtmax
Definition modmain.f90:216
real(8), dimension(:,:,:), allocatable jrmt
Definition modmain.f90:622
logical tforce
Definition modmain.f90:986
real(8) chgir
Definition modmain.f90:730
real(8), dimension(:,:), pointer, contiguous rhomt
Definition modmain.f90:614
integer nspecies
Definition modmain.f90:34
real(8), dimension(:,:), allocatable jrir
Definition modmain.f90:622
logical tdjr2d
Definition modtddft.f90:93
logical tdrho1d
Definition modtddft.f90:91
logical tdrho3d
Definition modtddft.f90:91
real(8), dimension(3, 0:1) afindt
Definition modtddft.f90:67
logical ttswrite
Definition modtddft.f90:80
real(8), dimension(:), allocatable times
Definition modtddft.f90:48
logical tdlsj
Definition modtddft.f90:94
logical tdmag2d
Definition modtddft.f90:92
logical tdrho2d
Definition modtddft.f90:91
logical tddos
Definition modtddft.f90:94
logical tafindt
Definition modtddft.f90:72
logical tdjr3d
Definition modtddft.f90:93
logical ttsforce
Definition modtddft.f90:100
logical tdjr1d
Definition modtddft.f90:93
integer itimes
Definition modtddft.f90:46
logical tdmag1d
Definition modtddft.f90:92
logical tdmag3d
Definition modtddft.f90:92
subroutine plot1d(fnum1, fnum2, nf, rfmt, rfir)
Definition plot1d.f90:10
subroutine plot2d(tproj, fnum, nf, rfmt, rfir)
Definition plot2d.f90:10
subroutine plot3d(fnum, nf, rfmt, rfir)
Definition plot3d.f90:10
subroutine writeatdisp
subroutine writeforces(fnum)
subroutine writemomtd
Definition writemomtd.f90:7
subroutine writetddft
Definition writetddft.f90:7
subroutine writetdengy
subroutine writetdforces
subroutine writetdtm3
Definition writetdtm3.f90:7