The Elk Code
acgwse.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 acgwse(ist,jst,sem,wr,ser)
7 use modmain
8 use modgw
9 implicit none
10 ! arguments
11 integer ist,jst
12 complex(8), intent(in) :: sem(nstsv,nstsv,0:nwfm)
13 real(8), intent(in) :: wr(nwplot)
14 complex(8), intent(out) :: ser(nstsv,nstsv,nwplot)
15 ! automatic arrays
16 complex(8) zm(0:nwfm),zwr(nwplot),zr(nwplot)
17 zm(0:nwfm)=sem(ist,jst,0:nwfm)
18 zwr(1:nwplot)=wr(1:nwplot)
19 select case(actype)
20 case(1)
21 ! fit a multipole model
22  call acpole(zm,zwr,zr)
23 case(10)
24 ! stabilised Pade approximant
25  call pades(nspade,swidth,nwfm+1,wfm,zm,nwplot,zwr,zr)
26 case default
27  write(*,*)
28  write(*,'("Error(acgwse): actype not defined : ",I8)') actype
29  write(*,*)
30  stop
31 end select
32 ser(ist,jst,1:nwplot)=zr(1:nwplot)
33 end subroutine
34 
integer actype
Definition: modgw.f90:30
real(8) swidth
Definition: modmain.f90:895
subroutine acpole(zm, zwr, zr)
Definition: acpole.f90:7
Definition: modgw.f90:6
integer nspade
Definition: modgw.f90:35
subroutine pades(ns, r, ni, zi, ui, no, zo, uo)
Definition: pades.f90:7
subroutine acgwse(ist, jst, sem, wr, ser)
Definition: acgwse.f90:7
complex(8), dimension(:), allocatable wfm
Definition: modgw.f90:25