The Elk Code
acpole.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2017 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 acpole(zm,zwr,zr)
7 use modmain
8 use modgw
9 implicit none
10 ! arguments
11 complex(8), intent(in) :: zm(0:nwfm),zwr(nwplot)
12 complex(8), intent(out) :: zr(nwplot)
13 ! local variables
14 integer, parameter :: maxit=1000
15 integer n2,iter,iw,i,j
16 real(8), parameter :: eps=1.d-5
17 ! external functions
18 complex(8), external :: zfpole
19 n2=2*(2*npole+1)
20 block
21 real(8) x(n2,n2+1)
22 ! intialise simplex
23 x(:,:)=0.d0
24 ! fit the constant
25 x(1,2)=0.5d0
26 x(2,3)=0.5d0
27 call minf_nm(1,zm,n2,x,maxit,iter,eps)
28 ! fit the constant and the first pole
29 x(3,1)=1.d0
30 do i=1,6
31  x(i,i+1)=x(i,1)+0.1d0
32 end do
33 call minf_nm(1,zm,n2,x,maxit,iter,eps)
34 ! fit the remaining poles one-by-one
35 i=7
36 do j=2,npole
37  x(i,1)=1.d0
38  x(i,i+1)=x(i,1)+0.1d0
39  i=i+1
40  x(i,i+1)=0.1d0
41  i=i+1
42  x(i,i+1)=0.1d0
43  i=i+1
44  x(i,i+1)=0.1d0
45  i=i+1
46  call minf_nm(1,zm,n2,x,maxit,iter,eps)
47 end do
48 ! fit the constant and the first pole again
49 do i=1,6
50  x(i,i+1)=x(i,1)+0.1d0
51 end do
52 call minf_nm(1,zm,n2,x,maxit,iter,eps)
53 ! fit everything together
54 if (npole > 1) then
55  do i=1,n2
56  x(i,i+1)=x(i,1)+0.1d0
57  end do
58  call minf_nm(1,zm,n2,x,maxit,iter,eps)
59 end if
60 do iw=1,nwplot
61  zr(iw)=zfpole(npole,x(:,1),zwr(iw))
62 end do
63 end block
64 end subroutine
65 
subroutine minf_nm(id, rd, n, x, maxit, iter, eps)
Definition: minf_nm.f90:7
subroutine acpole(zm, zwr, zr)
Definition: acpole.f90:7
integer npole
Definition: modgw.f90:32
Definition: modgw.f90:6