This repository has been archived by the owner on Sep 30, 2021. It is now read-only.
-
Notifications
You must be signed in to change notification settings - Fork 0
/
frprmn.f90
310 lines (274 loc) · 9 KB
/
frprmn.f90
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
subroutine frprmn(p,n,ftol,iter,nxsize,nwloc,Ns,omega,D_ev,Nbath,&
Norb,fret,comm)
implicit none
include 'mpif.h'
integer comm,nprocs,taskid,ierr,master
parameter(master=0)
integer:: iter, n, nmax, itmax,i
double precision:: fret, ftol, p(n), eps, func
parameter(nmax=60, itmax=20000, eps=1.d-10)
! USES dfunc, func, linmin
! Given a starting point p that is a vector of length n, Fletcher-Reeves-Polak-Ribiere minimization
! is performed on a function func, using its gradient as calculated by a routine dfunc.
! The convergence tolerance on the function value is input as ftol. Returned quantities are p
! (the location of the minimum), iter ( the # of iterations that were performed), and fret (the minimum value of function)
! The routine linmin is called to perform line minimizations
! parameter: nmax is the maximum anticipated value of n; itmax is the maximum allowd number of iterations;
! eps is a small number of rectify special case of converging to exactly zero function value.
integer its, j
double precision dgg, fp, gam, gg, g(nmax), h(nmax), xi(nmax),zero
integer nwloc,Ns,Nbath,Norb,nxsize
double precision omega(nwloc)
double complex D_ev(Norb,nwloc)
call mpi_comm_size(comm,nprocs,ierr)
call mpi_comm_rank(comm,taskid,ierr)
zero=0.D0
fp = func(p,Ns,nwloc,nxsize,omega,D_ev,Nbath,Norb,comm) ! initialization
if(taskid.eq.master) write(6,'("Initial difference :",e)') fp
call dfunc(p,xi,nwloc,Ns,nxsize,omega,D_ev,Nbath,Norb,comm)
if(taskid.eq.0) then
write(6,*) " x df/dx "
do i = 1, n
write(6,*) p(i),xi(i)
enddo
write(6,*)
endif
do j = 1, n
g(j) = -xi(j)
h(j) = g(j)
xi(j) = h(j)
enddo
do its = 1, itmax
iter = its
call linmin(p,xi,n,fret,Ns,nwloc,nxsize,omega,&
D_ev,Nbath,Norb,comm)
if(2.D0*abs(fret-fp).le.ftol*(abs(fret)+abs(fp)+eps)) then
if(taskid.eq.master) &
write(6,*) "After",its,"iteration converged."
return
endif
fp=fret
call dfunc(p,xi,nwloc,Ns,nxsize,omega,D_ev,Nbath,Norb,comm)
gg = 0.D0
dgg = 0.D0
do j = 1, n
gg = gg + g(j)*g(j)
! dgg = dgg + xi(j)*x(j)
dgg = dgg+(xi(j)+g(j))*xi(j)
enddo
if(gg.eq.zero) return
gam = dgg/gg
do j =1, n
g(j) = -xi(j)
h(j) = g(j)+gam*h(j)
xi(j)=h(j)
enddo
enddo
stop "frprmn maximum iterations exceeded!"
return
end
SUBROUTINE linmin(p,xi,n,fret,Ns,nwloc,nxsize,omega,D_ev,Nbath,Norb,comm)
implicit none
INTEGER n,NMAX
double precision fret,p(n),xi(n),TOL
PARAMETER (NMAX=60,TOL=1.e-8)
!U USES brent,f1dim,mnbrak
INTEGER j,ncom
double precision ax,bx,fa,fb,fx,xmin,xx,pcom(NMAX),xicom(NMAX),brent
! COMMON /f1com/ pcom,xicom,ncom
EXTERNAL f1dim
integer Ns,nwloc,nxsize,Nbath,Norb,comm
double complex D_ev(Norb,nwloc)
double precision omega(nwloc) !, ef(Norb)
do j=1,n
pcom(j)=p(j)
xicom(j)=xi(j)
enddo
ax=0.
xx=1.
call mnbrak(ax,xx,bx,fa,fx,fb,f1dim,Ns,nwloc,nxsize,omega,D_ev,&
Nbath,Norb,pcom,xicom,comm)
fret=brent(ax,xx,bx,f1dim,TOL,xmin,Ns,nwloc,nxsize,omega,D_ev,&
Nbath,Norb,pcom,xicom,comm)
do j=1,n
xi(j)=xmin*xi(j)
p(j)=p(j)+xi(j)
enddo
return
END
FUNCTION brent(ax,bx,cx,f,tol,xmin,Ns,nwloc,nxsize,omega,D_ev,&
Nbath,Norb,pcom,xicom,comm)
implicit none
INTEGER ITMAX
double precision brent,ax,bx,cx,tol,xmin,f,CGOLD,ZEPS
EXTERNAL f
PARAMETER (ITMAX=100,CGOLD=.3819660,ZEPS=1.0e-10)
INTEGER iter
double precision a,b,d,e,etemp,fu,fv,fw,fx,p,q,r,tol1,tol2,u,v,w,x,xm
double precision pcom, xicom ! cklee
integer Ns,nwloc,nxsize,Nbath,Norb,comm
double precision omega(nwloc)!,ef(Norb)
double complex D_ev(Norb,nwloc)
a=min(ax,cx)
b=max(ax,cx)
v=bx
w=v
x=v
e=0.
fx=f(x,Ns,nwloc,nxsize,omega,D_ev,Nbath,Norb,pcom,xicom,comm)
fv=fx
fw=fx
do 11 iter=1,ITMAX
xm=0.5*(a+b)
tol1=tol*abs(x)+ZEPS
tol2=2.*tol1
if(abs(x-xm).le.(tol2-.5*(b-a))) goto 3
if(abs(e).gt.tol1) then
r=(x-w)*(fx-fv)
q=(x-v)*(fx-fw)
p=(x-v)*q-(x-w)*r
q=2.*(q-r)
if(q.gt.0.) p=-p
q=abs(q)
etemp=e
e=d
if(abs(p).ge.abs(.5*q*etemp).or.p.le.q*(a-x).or.p.ge.q*(b-x)) goto 1
d=p/q
u=x+d
if(u-a.lt.tol2 .or. b-u.lt.tol2) d=sign(tol1,xm-x)
goto 2
endif
1 if(x.ge.xm) then
e=a-x
else
e=b-x
endif
d=CGOLD*e
2 if(abs(d).ge.tol1) then
u=x+d
else
u=x+sign(tol1,d)
endif
fu=f(u,Ns,nwloc,nxsize,omega,D_ev,Nbath,Norb,pcom,xicom,comm)
if(fu.le.fx) then
if(u.ge.x) then
a=x
else
b=x
endif
v=w
fv=fw
w=x
fw=fx
x=u
fx=fu
else
if(u.lt.x) then
a=u
else
b=u
endif
if(fu.le.fw .or. w.eq.x) then
v=w
fv=fw
w=u
fw=fu
else if(fu.le.fv .or. v.eq.x .or. v.eq.w) then
v=u
fv=fu
endif
endif
11 continue
stop 'brent exceed maximum iterations'
3 xmin=x
brent=fx
return
END
FUNCTION f1dim(x,Ns,nwloc,nxsize,omega,D_ev,Nbath,Norb,pcom,xicom,comm)
implicit none
INTEGER NMAX
integer Ns,nwloc,nxsize,Nbath,Norb,comm
double precision omega(nwloc)!,ef(Norb)
double complex D_ev(Norb,nwloc)
double precision f1dim,func,x
! PARAMETER (NMAX=90)
!U USES func
INTEGER j,ncom
double precision pcom(nxsize),xicom(nxsize),xt(nxsize)
! COMMON /f1com/ pcom,xicom,ncom
do j=1,nxsize
xt(j)=pcom(j)+x*xicom(j)
enddo
f1dim=func(xt,Ns,nwloc,nxsize,omega,D_ev,Nbath,Norb,comm)
return
END
SUBROUTINE mnbrak(ax,bx,cx,fa,fb,fc,func,Ns,nwloc,nxsize,&
omega,D_ev,Nbath,Norb,pcom,xicom,comm)
implicit none
integer Ns,nwloc,nxsize,Nbath,Norb,comm
double precision omega(nwloc)
double complex D_ev(Norb,nwloc)
double precision ax,bx,cx,fa,fb,fc,func,GOLD,GLIMIT,TINY
EXTERNAL func
PARAMETER (GOLD=1.618034, GLIMIT=100., TINY=1.e-20)
double precision dum,fu,q,r,u,ulim
double precision pcom, xicom
fa=func(ax,Ns,nwloc,nxsize,omega,D_ev,Nbath,Norb,pcom,xicom,comm)
fb=func(bx,Ns,nwloc,nxsize,omega,D_ev,Nbath,Norb,pcom,xicom,comm)
if(fb.gt.fa)then
dum=ax
ax=bx
bx=dum
dum=fb
fb=fa
fa=dum
endif
cx=bx+GOLD*(bx-ax)
fc=func(cx,Ns,nwloc,nxsize,omega,D_ev,Nbath,Norb,pcom,xicom,comm)
1 if(fb.ge.fc)then
r=(bx-ax)*(fb-fc)
q=(bx-cx)*(fb-fa)
u=bx-((bx-cx)*q-(bx-ax)*r)/(2.*sign(max(abs(q-r),TINY),q-r))
ulim=bx+GLIMIT*(cx-bx)
if((bx-u)*(u-cx).gt.0.)then
fu=func(u,Ns,nwloc,nxsize,omega,D_ev,Nbath,Norb,pcom,xicom,comm)
if(fu.lt.fc)then
ax=bx
fa=fb
bx=u
fb=fu
return
else if(fu.gt.fb)then
cx=u
fc=fu
return
endif
u=cx+GOLD*(cx-bx)
fu=func(u,Ns,nwloc,nxsize,omega,D_ev,Nbath,Norb,pcom,xicom,comm)
else if((cx-u)*(u-ulim).gt.0.)then
fu=func(u,Ns,nwloc,nxsize,omega,D_ev,Nbath,Norb,pcom,xicom,comm)
if(fu.lt.fc)then
bx=cx
cx=u
u=cx+GOLD*(cx-bx)
fb=fc
fc=fu
fu=func(u,Ns,nwloc,nxsize,omega,D_ev,Nbath,Norb,pcom,xicom,comm)
endif
else if((u-ulim)*(ulim-cx).ge.0.)then
u=ulim
fu=func(u,Ns,nwloc,nxsize,omega,D_ev,Nbath,Norb,pcom,xicom,comm)
else
u=cx+GOLD*(cx-bx)
fu=func(u,Ns,nwloc,nxsize,omega,D_ev,Nbath,Norb,pcom,xicom,comm)
endif
ax=bx
bx=cx
cx=u
fa=fb
fb=fc
fc=fu
goto 1
endif
return
END