File: writetddft.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (182 lines) | stat: -rw-r--r-- 4,785 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182

! Copyright (C) 2014 K. Krieger, J. K. Dewhurst, S. Sharma and E. K. U. Gross.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

subroutine writetddft
use modmain
use modtddft
use moddftu
implicit none
! local variables
integer is,ia,ias
real(8) t1
character(256) fext
! allocatable arrays
real(8), allocatable :: rvfmt(:,:,:),rvfir(:,:)
! file extension
write(fext,'("_TS",I8.8,".OUT")') itimes
! delete all files at first time-step
if (itimes.le.1) then
  open(50,file='CHARGEMT_TD.OUT')
  close(50,status='DELETE')
  open(50,file='CHARGEIR_TD.OUT')
  close(50,status='DELETE')
  open(50,file='MOMENT_TD.OUT')
  close(50,status='DELETE')
  open(50,file='MOMENTM_TD.OUT')
  close(50,status='DELETE')
  open(50,file='MOMENTMT_TD.OUT')
  close(50,status='DELETE')
  open(50,file='MOMENTIR_TD.OUT')
  close(50,status='DELETE')
  open(50,file='CURRENT_TD.OUT')
  close(50,status='DELETE')
  open(50,file='CURRENTM_TD.OUT')
  close(50,status='DELETE')
  if (tddos) then
    open(50,file='TDTEMP.OUT')
    close(50,status='DELETE')
  end if
end if
! muffin-tin charges
open(50,file='CHARGEMT_TD.OUT',form='FORMATTED',position='APPEND')
write(50,'(G18.10)') times(itimes)
do is=1,nspecies
  do ia=1,natoms(is)
    ias=idxas(ia,is)
    write(50,'(2I4,G18.10)') is,ia,chgmt(ias)
  end do
end do
write(50,*)
close(50)
! interstitial charge
open(50,file='CHARGEIR_TD.OUT',form='FORMATTED',position='APPEND')
write(50,'(2G18.10)') times(itimes),chgir
close(50)
! spin moment
open(50,file='MOMENT_TD.OUT',form='FORMATTED',position='APPEND')
write(50,'(4G18.10)') times(itimes),momtot(1:ndmag)
close(50)
open(50,file='MOMENTM_TD.OUT',form='FORMATTED',position='APPEND')
write(50,'(2G18.10)') times(itimes),momtotm
close(50)
! muffin-tin moments
open(50,file='MOMENTMT_TD.OUT',form='FORMATTED', &
 position='APPEND')
write(50,'(G18.10)') times(itimes)
do is=1,nspecies
  do ia=1,natoms(is)
    ias=idxas(ia,is)
    write(50,'(2I4,3G18.10)') is,ia,mommt(1:ndmag,ias)
  end do
end do
write(50,*)
close(50)
! interstitial moment
open(50,file='MOMENTIR_TD.OUT',form='FORMATTED', &
 position='APPEND')
write(50,'(4G18.10)') times(itimes),momir(1:ndmag)
close(50)
! total current
open(50,file='CURRENT_TD.OUT',form='FORMATTED',position='APPEND')
write(50,'(4G18.10)') times(itimes),curtot(:)
close(50)
! total current magnitude
open(50,file='CURRENTM_TD.OUT',form='FORMATTED',position='APPEND')
t1=sqrt(sum(curtot(:)**2))
write(50,'(4G18.10)') times(itimes),curtotm
close(50)
! write optional quantities
if (ntswrite.le.0) goto 10
if (mod(itimes-1,ntswrite).ne.0) goto 10
! charge density in 1D
if (tdrho1d) then
  open(50,file='RHO1D'//trim(fext),form='FORMATTED')
  open(51,file='RHOLINES'//trim(fext),form='FORMATTED')
  call plot1d(50,51,1,rhomt,rhoir)
  close(50)
  close(51)
end if
! charge density in 2D
if (tdrho2d) then
  open(50,file='RHO2D'//trim(fext),form='FORMATTED')
  call plot2d(50,1,rhomt,rhoir)
  close(50)
end if
! charge density in 3D
if (tdrho3d) then
  open(50,file='RHO3D'//trim(fext),form='FORMATTED')
  call plot3d(50,1,rhomt,rhoir)
  close(50)
end if
! magnetisation in 2D or 3D
if ((tdmag2d.or.tdmag3d).and.(spinpol)) then
  allocate(rvfmt(npmtmax,natmtot,3),rvfir(ngtot,3))
  if (ncmag) then
! non-collinear
    rvfmt(:,:,:)=magmt(:,:,:)
    rvfir(:,:)=magir(:,:)
  else
! collinear
    rvfmt(:,:,1:2)=0.d0
    rvfir(:,1:2)=0.d0
    rvfmt(:,:,3)=magmt(:,:,1)
    rvfir(:,3)=magir(:,1)
  end if
  if (tdmag2d) then
! project the magnetisation onto the plotting plane
    call proj2d(rvfmt,rvfir)
    open(50,file='MAG2D'//trim(fext),form='FORMATTED')
    call plot2d(50,3,rvfmt,rvfir)
    close(50)
  end if
  if (tdmag3d) then
    open(50,file='MAG3D'//trim(fext),form='FORMATTED')
    call plot3d(50,3,rvfmt,rvfir)
    close(50)
  end if
  deallocate(rvfmt,rvfir)
end if
! current density in 2D
if (tdcd2d) then
  allocate(rvfmt(npmtmax,natmtot,3),rvfir(ngtot,3))
  rvfmt(:,:,:)=cdmt(:,:,:)
  rvfir(:,:)=cdir(:,:)
  call proj2d(rvfmt,rvfir)
  open(50,file='CURDEN2D'//trim(fext),form='FORMATTED')
  call plot2d(50,3,rvfmt,rvfir)
  close(50)
  deallocate(rvfmt,rvfir)
end if
if (tdcd3d) then
  open(50,file='CURDEN3D'//trim(fext),form='FORMATTED')
  call plot3d(50,3,cdmt,cdir)
  close(50)
end if
! calculate and write tensor moments
if (dftu.ne.0) then
  if (tmwrite) then
    open(50,file='TMDFTU'//trim(fext),form='FORMATTED')
    if (spinorb) then
      call writetm3du(50)
    else
      call writetm2du(50)
    end if
    close(50)
  end if
end if
! write time-dependent DOS
if (tddos) then
  call writetddos(fext)
end if
! write the atomic forces
if (tforce) then
  open(50,file='FORCES'//trim(fext),form='FORMATTED')
  call writeforces(50)
  close(50)
end if
10 continue
return
end subroutine