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