File: tdtemp.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 (58 lines) | stat: -rw-r--r-- 1,321 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

! Copyright (C) 2015 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 tdtemp(occsvp)
use modmain
use modtddft
implicit none
! arguments
real(8), intent(in) :: occsvp(nstsv,nkpt)
! local variables
integer, parameter :: maxit=1000
integer ik,ist,it
real(8), parameter :: eps=1.d-6
real(8) sw,dsw,sum,sp
real(8) x,t1,t2,t3
! external functions
real(8) sdelta_fd,stheta_fd
external sdelta_fd,stheta_fd
! initial smearing width
sw=1.d-6
! initial smearing width step size
dsw=1.d-6
sp=0.d0
do it=1,maxit
  t1=1.d0/sw
  sum=0.d0
  do ik=1,nkpt
    do ist=1,nstsv
      x=(efermi-evalsv(ist,ik))*t1
      t2=occmax*stheta_fd(x)
      t3=sdelta_fd(x)*x*t1
      sum=sum+(occsvp(ist,ik)-t2)*t3
    end do
  end do
  if ((sum*sp).lt.0.d0) then
    dsw=-0.5d0*dsw
    if (abs(dsw).lt.eps) goto 10
  else
    dsw=1.5d0*dsw
  end if
  sp=sum
  sw=sw+dsw
  if ((sw.lt.0.d0).or.(sw.gt.1.d6)) exit
end do
write(*,*)
write(*,'("Warning(tdtemp): could not estimate effective temperature")')
return
10 continue
! write effective temperature to file
t1=sw/kboltz
open(50,file='TDTEMP.OUT',form='FORMATTED',position='APPEND')
write(50,'(2G18.10)') times(itimes),t1
close(50)
return
end subroutine