The Elk Code
minf_nm.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2007 J. K. Dewhurst and D. W. H. Rankin.
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 minf_nm(id,rd,n,x,maxit,iter,eps)
7 implicit none
8 ! arguments
9 integer, intent(in) :: id(*)
10 real(8), intent(in) :: rd(*)
11 integer, intent(in) :: n
12 real(8), intent(inout) :: x(n,n+1)
13 integer, intent(in) :: maxit
14 integer, intent(out) :: iter
15 real(8), intent(in) :: eps
16 ! local variables
17 integer i,j,il,iu
18 ! Nelder-Mead parmeters
19 real(8), parameter :: alpha=1.d0,gamma=2.d0
20 real(8), parameter :: beta=0.5d0,sigma=0.5d0
21 real(8) fr,fe,fc,sm,t1
22 ! automatic arrays
23 real(8) f(n+1),xm(n),xr(n),xe(n),xc(n)
24 ! external functions
25 real(8), external :: fmin_nm
26 if (n < 1) then
27  write(*,*)
28  write(*,'("Error(minf_nm): n < 1 : ",I8)') n
29  write(*,*)
30  stop
31 end if
32 ! evaluate the function at each vertex
33 do i=1,n+1
34  f(i)=fmin_nm(id,rd,x(:,i))
35 end do
36 iter=0
37 10 continue
38 iter=iter+1
39 if (iter >= maxit) return
40 ! find the lowest and highest vertex
41 il=1
42 iu=1
43 do i=2,n+1
44  if (f(i) < f(il)) il=i
45  if (f(i) > f(iu)) iu=i
46 end do
47 ! check for convergence
48 if ((f(iu)-f(il)) < eps) return
49 ! compute the mean of the n lowest vertices
50 t1=1.d0/dble(n)
51 do i=1,n
52  sm=0.d0
53  do j=1,iu-1
54  sm=sm+x(i,j)
55  end do
56  do j=iu+1,n+1
57  sm=sm+x(i,j)
58  end do
59  xm(i)=t1*sm
60 end do
61 xr(:)=xm(:)+alpha*(xm(:)-x(:,iu))
62 fr=fmin_nm(id,rd,xr)
63 if (f(il) > fr) goto 30
64 if ((f(il) <= fr).and.(fr < f(iu))) then
65 ! reflection
66  x(:,iu)=xr(:)
67  f(iu)=fr
68  goto 10
69 else
70  goto 40
71 end if
72 30 continue
73 xe(:)=xm(:)+gamma*(xr(:)-xm(:))
74 fe=fmin_nm(id,rd,xe)
75 if (fr > fe) then
76 ! expansion
77  x(:,iu)=xe(:)
78  f(iu)=fe
79 else
80 ! reflection
81  x(:,iu)=xr(:)
82  f(iu)=fr
83 end if
84 goto 10
85 40 continue
86 xc(:)=xm(:)+beta*(x(:,iu)-xm(:))
87 fc=fmin_nm(id,rd,xc)
88 if (fc < f(iu)) then
89 ! contraction
90  x(:,iu)=xc(:)
91  f(iu)=fc
92  goto 10
93 end if
94 ! shrinkage
95 do j=1,il-1
96  x(:,j)=x(:,il)+sigma*(x(:,j)-x(:,il))
97  f(j)=fmin_nm(id,rd,x(1,j))
98 end do
99 do j=il+1,n+1
100  x(:,j)=x(:,il)+sigma*(x(:,j)-x(:,il))
101  f(j)=fmin_nm(id,rd,x(1,j))
102 end do
103 goto 10
104 end subroutine
105 
subroutine minf_nm(id, rd, n, x, maxit, iter, eps)
Definition: minf_nm.f90:7