The Elk Code
tdtemp.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2015 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 
6 subroutine tdtemp(occsvp)
7 use modmain
8 use modtddft
9 implicit none
10 ! arguments
11 real(8), intent(in) :: occsvp(nstsv,nkpt)
12 ! local variables
13 integer, parameter :: maxit=1000
14 integer ik,ist,it
15 real(8), parameter :: eps=1.d-6
16 real(8) sw,dsw,sm,sp
17 real(8) x,t1,t2,t3
18 ! external functions
19 real(8), external :: sdelta_fd,stheta_fd
20 ! initial smearing width
21 sw=1.d-6
22 ! initial smearing width step size
23 dsw=1.d-6
24 sp=0.d0
25 do it=1,maxit
26  t1=1.d0/sw
27  sm=0.d0
28  do ik=1,nkpt
29  do ist=1,nstsv
30  x=(efermi-evalsv(ist,ik))*t1
31  t2=occmax*stheta_fd(x)
32  t3=sdelta_fd(x)*x*t1
33  sm=sm+(occsvp(ist,ik)-t2)*t3
34  end do
35  end do
36  if ((sm*sp) < 0.d0) then
37  dsw=-0.5d0*dsw
38  if (abs(dsw) < eps) goto 10
39  else
40  dsw=1.5d0*dsw
41  end if
42  sp=sm
43  sw=sw+dsw
44  if ((sw < 0.d0).or.(sw > 1.d6)) exit
45 end do
46 write(*,*)
47 write(*,'("Warning(tdtemp): could not estimate effective temperature")')
48 return
49 10 continue
50 ! write effective temperature to file
51 t1=sw/kboltz
52 open(50,file='TDTEMP.OUT',form='FORMATTED',position='APPEND')
53 write(50,'(2G18.10)') times(itimes),t1
54 close(50)
55 end subroutine
56 
real(8) efermi
Definition: modmain.f90:907
real(8), dimension(:,:), allocatable evalsv
Definition: modmain.f90:921
real(8), parameter kboltz
Definition: modmain.f90:1263
real(8) occmax
Definition: modmain.f90:901
subroutine tdtemp(occsvp)
Definition: tdtemp.f90:7
real(8), dimension(:), allocatable times
Definition: modtddft.f90:48
integer itimes
Definition: modtddft.f90:46