The Elk Code
 
Loading...
Searching...
No Matches
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
6subroutine tdtemp(occsvp)
7use modmain
8use modtddft
9implicit none
10! arguments
11real(8), intent(in) :: occsvp(nstsv,nkpt)
12! local variables
13integer, parameter :: maxit=1000
14integer ik,ist,it
15real(8), parameter :: eps=1.d-6
16real(8) sw,dsw,sm,sp
17real(8) x,t1,t2,t3
18! external functions
19real(8), external :: sdelta_fd,stheta_fd
20! initial smearing width
21sw=1.d-6
22! initial smearing width step size
23dsw=1.d-6
24sp=0.d0
25do 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
45end do
46write(*,*)
47write(*,'("Warning(tdtemp): could not estimate effective temperature")')
48return
4910 continue
50! write effective temperature to file
51t1=sw/kboltz
52open(50,file='TDTEMP.OUT',form='FORMATTED',position='APPEND')
53write(50,'(2G18.10)') times(itimes),t1
54close(50)
55end subroutine
56
real(8) efermi
Definition modmain.f90:904
real(8), parameter kboltz
Definition modmain.f90:1262
real(8) occmax
Definition modmain.f90:898
real(8), dimension(:,:), allocatable evalsv
Definition modmain.f90:918
real(8), dimension(:), allocatable times
Definition modtddft.f90:48
integer itimes
Definition modtddft.f90:46
elemental real(8) function sdelta_fd(x)
Definition sdelta_fd.f90:10
elemental real(8) function stheta_fd(x)
Definition stheta_fd.f90:10
subroutine tdtemp(occsvp)
Definition tdtemp.f90:7