File: tddftsplr.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (325 lines) | stat: -rw-r--r-- 8,772 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
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325

! Copyright (C) 2013 S. Sharma, J. K. Dewhurst 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 tddftsplr
use modmain
use modmpi
use modomp
implicit none
! local variables
integer ik,isym,iq,iw
integer ig,jg,i,j,n
integer info,nthd
real(8) v(3)
complex(8) a(4,4),b(4,4),z1
character(256) fname
! allocatable arrays
integer, allocatable :: ipiv(:)
integer(8), allocatable :: lock(:)
real(8), allocatable :: vgqc(:,:),gqc(:),gclgq(:),jlgqr(:,:,:)
complex(8), allocatable :: ylmgq(:,:),sfacgq(:,:)
complex(8), allocatable :: chi(:,:,:,:,:),chit(:),fxc(:,:,:,:)
complex(8), allocatable :: c(:,:),d(:,:,:,:),work(:)
if (.not.spinpol) then
  write(*,*)
  write(*,'("Error(tddftsplr): spin-unpolarised calculation")')
  write(*,*)
  stop
end if
! initialise global variables
call init0
call init1
call init2
call init3
! check q-vector is commensurate with k-point grid
v(:)=dble(ngridk(:))*vecql(:)
v(:)=abs(v(:)-nint(v(:)))
if ((v(1).gt.epslat).or.(v(2).gt.epslat).or.(v(3).gt.epslat)) then
  write(*,*)
  write(*,'("Error(tddftsplr): q-vector incommensurate with k-point grid")')
  write(*,'(" ngridk : ",3I6)') ngridk
  write(*,'(" vecql : ",3G18.10)') vecql
  write(*,*)
  stop
end if
! find the equivalent reduced q-point
call findqpt(vecql,isym,iq)
! read density and potentials from file
call readstate
! read Fermi energy from a file
call readfermi
! find the new linearisation energies
call linengy
! generate the APW radial functions
call genapwfr
! generate the local-orbital radial functions
call genlofr
! get the eigenvalues and occupancies from file
do ik=1,nkpt
  call getevalsv(filext,ik,vkl(:,ik),evalsv(:,ik))
  call getoccsv(filext,ik,vkl(:,ik),occsv(:,ik))
end do
! generate the G+q-vectors and related quantities
allocate(vgqc(3,ngrf),gqc(ngrf),jlgqr(njcmax,nspecies,ngrf))
allocate(ylmgq(lmmaxo,ngrf),sfacgq(ngrf,natmtot))
call gengqrf(vecqc,vgqc,gqc,jlgqr,ylmgq,sfacgq)
deallocate(vgqc)
! initialise the OpenMP locks
allocate(lock(nwrf))
do iw=1,nwrf
  call omp_init_lock(lock(iw))
end do
! compute chi0
allocate(chi(ngrf,4,ngrf,4,nwrf))
chi(:,:,:,:,:)=0.d0
call omp_hold(nkptnr/np_mpi,nthd)
!$OMP PARALLEL DEFAULT(SHARED) &
!$OMP NUM_THREADS(nthd)
!$OMP DO
do ik=1,nkptnr
! distribute among MPI processes
  if (mod(ik-1,np_mpi).ne.lp_mpi) cycle
!$OMP CRITICAL(tddftsplr_)
  write(*,'("Info(tddftsplr): ",I6," of ",I6," k-points")') ik,nkptnr
!$OMP END CRITICAL(tddftsplr_)
  call genspchi0(ik,lock,scissor,vecql,jlgqr,ylmgq,sfacgq,chi)
end do
!$OMP END DO
!$OMP END PARALLEL
call omp_free(nthd)
! destroy the OpenMP locks
do iw=1,nwrf
  call omp_destroy_lock(lock(iw))
end do
deallocate(lock)
! add chi0 from each process and redistribute
if (np_mpi.gt.1) then
  n=ngrf*4*ngrf*4*nwrf
  call mpi_allreduce(mpi_in_place,chi,n,mpi_double_complex,mpi_sum,mpicom, &
   ierror)
end if
! transform chi0 from 2x2 to 1x3 basis
do iw=1,nwrf
  do ig=1,ngrf
    do jg=1,ngrf
      a(:,:)=chi(ig,:,jg,:,iw)
      call tfm2213(a,b)
      chi(ig,:,jg,:,iw)=b(:,:)
    end do
  end do
end do
! generate transverse chi0 for the collinear case
if (.not.ncmag) then
  allocate(chit(nwrf))
  do iw=1,nwrf
    a(:,:)=chi(1,:,1,:,iw)
    call tfm13t(a,b)
    chit(iw)=b(2,2)
  end do
end if
! write chi0 to file
if (mp_mpi) then
! write chi0 to file in 1x3 basis
  do i=1,4
    do j=1,4
      write(fname,'("CHI0_",2I1,".OUT")') i-1,j-1
      open(50,file=trim(fname),form='FORMATTED')
      do iw=1,nwrf
        write(50,'(2G18.10)') dble(wrf(iw)),dble(chi(1,i,1,j,iw))
      end do
      write(50,*)
      do iw=1,nwrf
        write(50,'(2G18.10)') dble(wrf(iw)),aimag(chi(1,i,1,j,iw))
      end do
      close(50)
    end do
  end do
! write transverse chi0 for collinear case
  if (.not.ncmag) then
    open(50,file='CHI0_T.OUT',form='FORMATTED')
    do iw=1,nwrf
      write(50,'(2G18.10)') dble(wrf(iw)),dble(chit(iw))
    end do
    write(50,*)
    do iw=1,nwrf
      write(50,'(2G18.10)') dble(wrf(iw)),aimag(chit(iw))
      end do
    close(50)
  end if
end if
! compute f_xc in G-space
allocate(fxc(ngrf,4,ngrf,4))
call genspfxcg(fxc)
! generate the Coulomb Green's function in G+q-space regularised for q=0
allocate(gclgq(ngrf))
call gengclgq(.true.,iq,ngrf,gqc,gclgq)
! add the regularised Coulomb interaction to f_xc to give f_Hxc
do ig=1,ngrf
  fxc(ig,1,ig,1)=fxc(ig,1,ig,1)+gclgq(ig)
end do
deallocate(gclgq)
! matrix size
n=4*ngrf
allocate(ipiv(n),work(n))
allocate(c(n,n),d(ngrf,4,ngrf,4))
! loop over frequencies
do iw=1,nwrf
! multiply f_Hxc by -chi0 from the left
  z1=-1.d0
  call zgemm('N','N',n,n,n,z1,chi(:,:,:,:,iw),n,fxc,n,zzero,c,n)
! add the identity
  do i=1,n
    c(i,i)=c(i,i)+1.d0
  end do
! invert the matrix
  call zgetrf(n,n,c,n,ipiv,info)
  if (info.eq.0) call zgetri(n,c,n,ipiv,work,n,info)
  if (info.ne.0) then
    write(*,*)
    write(*,'("Error(tddftsplr): unable to invert matrix")')
    write(*,'(" for frequency ",I6)') iw
    write(*,*)
    stop
  end if
! multiply by chi0 on the right and store in chi
  call zgemm('N','N',n,n,n,zone,c,n,chi(:,:,:,:,iw),n,zzero,d,n)
  chi(:,:,:,:,iw)=d(:,:,:,:)
end do
deallocate(ipiv,work,c,d)
! generate transverse chi for the collinear case
if (.not.ncmag) then
  do iw=1,nwrf
    a(:,:)=chi(1,:,1,:,iw)
    call tfm13t(a,b)
    chit(iw)=b(2,2)
  end do
end if
if (mp_mpi) then
! write the complete chi matrix if required
  if (task.eq.331) then
    open(50,file='CHI.OUT',form='UNFORMATTED')
    write(50) chi
    close(50)
  end if
! write chi for G = G' = 0 in the 1x3 basis
  do i=1,4
    do j=1,4
      write(fname,'("CHI_",2I1,".OUT")') i-1,j-1
      open(50,file=trim(fname),form='FORMATTED')
      do iw=1,nwrf
        write(50,'(2G18.10)') dble(wrf(iw)),dble(chi(1,i,1,j,iw))
      end do
      write(50,*)
      do iw=1,nwrf
        write(50,'(2G18.10)') dble(wrf(iw)),aimag(chi(1,i,1,j,iw))
      end do
      close(50)
    end do
  end do
! write transverse chi for collinear case
  if (.not.ncmag) then
    open(50,file='CHI_T.OUT',form='FORMATTED')
    do iw=1,nwrf
      write(50,'(2G18.10)') dble(wrf(iw)),dble(chit(iw))
    end do
    write(50,*)
    do iw=1,nwrf
      write(50,'(2G18.10)') dble(wrf(iw)),aimag(chit(iw))
    end do
    close(50)
  end if
  write(*,*)
  write(*,'("Info(tddftsplr):")')
  write(*,'(" Spin-dependent response function chi_ij(G,G'',q,w) written to &
   &CHI_ij.OUT")')
  write(*,'(" for i,j = 0-3; G = G'' = 0; and all wplot frequencies")')
  write(*,'(" q-vector (lattice coordinates) : ")')
  write(*,'(3G18.10)') vecql
  write(*,'(" q-vector length : ",G18.10)') gqc(1)
  write(*,*)
  write(*,'(" The elements of chi labeled by (i,j) form the 4x4 matrix :")')
  write(*,*)
  write(*,'("                  (_|_ _ _)")')
  write(*,'("  chi(G,G'',q,w) = ( |     )")')
  write(*,'("                  ( |     )")')
  write(*,'("                  ( |     )")')
  write(*,*)
  write(*,'(" (0,0) is the charge-charge response drho/dv")')
  write(*,'(" (0,1-3) is the charge-magnetisation response drho/dB")')
  write(*,'(" (1-3,0) is the magnetisation-charge response dm/v")')
  write(*,'(" (1-3,1-3) is the magnetisation-magnetisation response dm/dB")')
  write(*,*)
  write(*,'(" Non-interacting Kohn-Sham reponse function written to &
   &CHI0_ij.OUT")')
  if (.not.ncmag) then
    write(*,*)
    write(*,'(" Transverse components corresponding to m_+- = m_x +- im_y")')
    write(*,'(" written to CHI_T.OUT and CHI0_T.OUT")')
  end if
  if (task.eq.331) then
    write(*,*)
    write(*,'(" Complete response function for all G, G'' written to binary &
     &file CHI.OUT")')
    write(*,'(" (array index ordering changed from version 4.5.16 onwards)")')
  end if
end if
deallocate(gqc,ylmgq,sfacgq,chi,fxc)
if (.not.ncmag) deallocate(chit)
return

contains

subroutine tfm2213(a,b)
implicit none
! arguments
complex(8), intent(in) :: a(4,4)
complex(8), intent(out) :: b(4,4)
! local variables
integer i,j
complex(8) c(4,4),z1
do i=1,4
  c(i,1)=a(i,1)+a(i,4)
  c(i,2)=a(i,2)+a(i,3)
  z1=a(i,2)-a(i,3)
  c(i,3)=cmplx(aimag(z1),-dble(z1),8)
  c(i,4)=a(i,1)-a(i,4)
end do
do j=1,4
  b(1,j)=c(1,j)+c(4,j)
  b(2,j)=c(2,j)+c(3,j)
  z1=c(2,j)-c(3,j)
  b(3,j)=cmplx(-aimag(z1),dble(z1),8)
  b(4,j)=c(1,j)-c(4,j)
end do
return
end subroutine

subroutine tfm13t(a,b)
implicit none
! arguments
complex(8), intent(in) :: a(4,4)
complex(8), intent(out) :: b(4,4)
! local variables
integer i,j
complex(8) c(4,4),z1
do i=1,4
  c(i,1)=a(i,1)
  z1=a(i,3)
  c(i,2)=a(i,2)+cmplx(aimag(z1),-dble(z1),8)
  c(i,3)=a(i,2)+cmplx(-aimag(z1),dble(z1),8)
  c(i,4)=a(i,4)
end do
do j=1,4
  b(1,j)=c(1,j)
  z1=c(3,j)
  b(2,j)=c(2,j)+cmplx(-aimag(z1),dble(z1),8)
  b(3,j)=c(2,j)+cmplx(aimag(z1),-dble(z1),8)
  b(4,j)=c(4,j)
end do
return
end subroutine

end subroutine