The Elk Code
Loading...
Searching...
No Matches
rdmvaryn.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2007-2008 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3
! This file is distributed under the terms of the GNU Lesser General Public
4
! License. See the file COPYING for license details.
5
6
!BOP
7
! !ROUTINE: rdmvaryn
8
! !INTERFACE:
9
subroutine
rdmvaryn
10
! !USES:
11
use
modmain
12
use
modrdm
13
use
modmpi
14
! !DESCRIPTION:
15
! Calculates new occupation numbers from old by using the derivatives of the
16
! total energy: $n_i^{\rm new} = n_i^{\rm old}-\tau \gamma_i$, where $\tau$ is
17
! chosen such that $0 \le n_i \le n_{\rm max}$ with
18
! $$ \gamma_i=\begin{cases}
19
! g_i(n_{\rm max}-n_i) & g_i > 0 \\
20
! g_i n_i & g_i\le 0 \end{cases} $$
21
! where $g_i=\partial E/\partial n_i-\kappa$, and $\kappa$ is chosen such that
22
! $\sum_i\gamma_i=0$.
23
!
24
! !REVISION HISTORY:
25
! Created 2009 (JKD,Sharma)
26
!EOP
27
!BOC
28
implicit none
29
! local variables
30
integer
,
parameter
:: maxit=10000
31
integer
it,ik,ist
32
real
(8),
parameter
:: eps=1.d-12
33
real
(8) tau,sm,gs,gsp,dgs
34
real
(8) kapa,dkapa,t1
35
! allocatable arrays
36
real
(8),
allocatable
:: dedn(:,:),gamma(:,:)
37
! add constant to occupation numbers for charge conservation
38
sm=0.d0
39
do
ik=1,
nkpt
40
do
ist=1,
nstsv
41
sm=sm+
wkpt
(ik)*
occsv
(ist,ik)
42
end do
43
end do
44
t1=(
chgval
-sm)/dble(
nstsv
)
45
occsv
(:,:)=
occsv
(:,:)+t1
46
! redistribute charge so that occupation numbers are in the interval [0,occmax]
47
sm=0.d0
48
do
ik=1,
nkpt
49
do
ist=1,
nstsv
50
if
(
occsv
(ist,ik) >
occmax
)
then
51
sm=sm+
wkpt
(ik)*(
occsv
(ist,ik)-
occmax
)
52
occsv
(ist,ik)=
occmax
53
end if
54
if
(
occsv
(ist,ik) < 0.d0)
then
55
sm=sm+
wkpt
(ik)*
occsv
(ist,ik)
56
occsv
(ist,ik)=0.d0
57
end if
58
end do
59
end do
60
do
ist=1,
nstsv
61
do
ik=1,
nkpt
62
if
(sm > 0.d0)
then
63
t1=
wkpt
(ik)*(
occmax
-
occsv
(ist,ik))
64
t1=min(t1,sm)
65
occsv
(ist,ik)=
occsv
(ist,ik)+t1/
wkpt
(ik)
66
sm=sm-t1
67
else
68
t1=
wkpt
(ik)*
occsv
(ist,ik)
69
t1=min(t1,-sm)
70
occsv
(ist,ik)=
occsv
(ist,ik)-t1/
wkpt
(ik)
71
sm=sm+t1
72
end if
73
end do
74
end do
75
allocate
(dedn(
nstsv
,
nkpt
))
76
allocate
(gamma(
nstsv
,
nkpt
))
77
! get the derivatives
78
call
rdmdedn
(dedn)
79
! find suitable value of kapa such that sum of gamma is 0
80
gsp=0.d0
81
kapa=0.d0
82
dkapa=0.1d0
83
do
it=1,maxit
84
gs=0.d0
85
sm=0.d0
86
do
ik=1,
nkpt
87
do
ist=1,
nstsv
88
t1=dedn(ist,ik)-kapa
89
if
(t1 > 0.d0)
then
90
gamma(ist,ik)=t1*(
occmax
-
occsv
(ist,ik))
91
else
92
gamma(ist,ik)=t1*
occsv
(ist,ik)
93
end if
94
gs=gs+
wkpt
(ik)*gamma(ist,ik)
95
sm=sm+
wkpt
(ik)*gamma(ist,ik)**2
96
end do
97
end do
98
sm=sqrt(sm)
99
sm=max(sm,1.d0)
100
t1=abs(gs)/sm
101
if
(t1 < eps)
goto
10
102
if
(it >= 2)
then
103
dgs=gs-gsp
104
if
(gs*dgs > 0.d0) dkapa=-dkapa
105
if
(gs*gsp < 0.d0)
then
106
dkapa=0.5d0*dkapa
107
else
108
dkapa=1.1d0*dkapa
109
end if
110
end if
111
gsp=gs
112
kapa=kapa+dkapa
113
end do
114
write
(*,*)
115
write
(*,
'("Error(rdmvaryn): could not find offset")'
)
116
write
(*,*)
117
stop
118
10
continue
119
! write derivatives and occupation numbers to file
120
call
rdmwritededn
(dedn)
121
deallocate
(dedn)
122
! normalize gamma if sum of squares is greater than 1
123
sm=0.d0
124
do
ik=1,
nkpt
125
do
ist=1,
nstsv
126
sm=sm+
wkpt
(ik)*gamma(ist,ik)**2
127
end do
128
end do
129
if
(sm > 1.d0)
then
130
t1=1.d0/sqrt(sm)
131
gamma(:,:)=t1*gamma(:,:)
132
end if
133
! find step size which keeps occupation numbers in the interval [0,occmax]
134
tau=
taurdmn
135
20
continue
136
if
(abs(tau) < eps)
goto
30
137
do
ik=1,
nkpt
138
do
ist=1,
nstsv
139
t1=
occsv
(ist,ik)+tau*gamma(ist,ik)
140
if
(gamma(ist,ik) > 0.d0)
then
141
if
(t1 >
occmax
+eps)
then
142
tau=0.75d0*tau
143
goto
20
144
end if
145
end if
146
if
(gamma(ist,ik) < 0.d0)
then
147
if
(t1 < -eps)
then
148
tau=0.75d0*tau
149
goto
20
150
end if
151
end if
152
end do
153
end do
154
30
continue
155
! update occupation numbers and write to file
156
do
ik=1,
nkpt
157
do
ist=1,
nstsv
158
occsv
(ist,ik)=
occsv
(ist,ik)+tau*gamma(ist,ik)
159
end do
160
call
putoccsv
(
filext
,ik,
occsv
(:,ik))
161
end do
162
deallocate
(gamma)
163
end subroutine
164
!EOC
165
modmain
Definition
modmain.f90:6
modmain::wkpt
real(8), dimension(:), allocatable wkpt
Definition
modmain.f90:475
modmain::filext
character(256) filext
Definition
modmain.f90:1300
modmain::chgval
real(8) chgval
Definition
modmain.f90:722
modmain::nkpt
integer nkpt
Definition
modmain.f90:461
modmain::nstsv
integer nstsv
Definition
modmain.f90:886
modmain::occmax
real(8) occmax
Definition
modmain.f90:898
modmain::occsv
real(8), dimension(:,:), allocatable occsv
Definition
modmain.f90:902
modmpi
Definition
modmpi.f90:6
modrdm
Definition
modrdm.f90:6
modrdm::taurdmn
real(8) taurdmn
Definition
modrdm.f90:17
putoccsv
subroutine putoccsv(fext, ik, occsvp)
Definition
putoccsv.f90:7
rdmdedn
subroutine rdmdedn(dedn)
Definition
rdmdedn.f90:10
rdmvaryn
subroutine rdmvaryn
Definition
rdmvaryn.f90:10
rdmwritededn
subroutine rdmwritededn(dedn)
Definition
rdmwritededn.f90:10
rdmvaryn.f90
Generated by
1.9.8