Class | signfunc_mod |
In: |
ZolotarevModule/signfunc_mod.F90
|
Zolotarev coefficients for
sign(x) : (Lmin/Lmax) < |x| < 1,
is computed.
[1] T.-W.Chiu, T.-H.Hsieh, C.-H.Huang, T.-R.Huang, Phys.Rev.D66 (2002) 114502. [2] T.-W.Chiu, Phys.Rev.Lett. 90 (2003) 071601. [3] A.D.Kennedy, hep-lat/0402038.
Subroutine : | |
this : | type(zolotarev_signfunc), intent(inout) |
subroutine delete_zolotarev(this) implicit none type(zolotarev_signfunc), intent(inout) :: this if (associated(this%sn)) deallocate(this%sn) if (associated(this%cn)) deallocate(this%cn) if (associated(this%dn)) deallocate(this%dn) if (associated(this%aa)) deallocate(this%aa) if (associated(this%omg)) deallocate(this%omg) return end subroutine
Function : | |
NS : | integer |
kappa : | real(QP), intent(in) |
eps : | real(QP), intent(in) |
Search for NS satisfies err(NS) < eps using bisection method
function get_ns_exact(kappa,eps) result(NS) ! ! Search for NS satisfies err(NS) < eps ! using bisection method ! implicit none real(QP), intent(in) :: kappa,eps integer :: NS real(QP) :: tol,err0,err1,errm,err1_0,elimit integer :: i,j,iflg integer :: ns0,ns1,nsm,ns1_0 integer, parameter :: ns_max=124 elimit=EPSILON(elimit) ns0=1 ns1=10 err0 = errorbound_zol(kappa,ns0)-eps err1 = errorbound_zol(kappa,ns1)-eps iflg=0 do j=1,100 ! write(*,'(I4,2(I3,E24.15))')j,ns0,err0,ns1,err1 if (err0*err1 > 0.0_QP) then err1_0 = err1 ns1_0 = ns1 ns1 = min(ns1*2,ns_max) err1 = errorbound_zol(kappa,ns1)-eps !============================== ! no sign flip until ns=ns_max !============================== if (ns1_0 == ns1) then iflg=1 exit endif else nsm = (ns1+ns0)/2 errm = errorbound_zol(kappa,nsm)-eps if (err0*errm > 0.0_QP) then ns0 = nsm err0 = errm else ns1 = nsm err1 = errm endif endif ! write(*,*)ns0,err0,ns1,err1,elimit if (ns1-ns0 == 1) then iflg=1 exit endif enddo if (iflg==0) then write(*,*)"Error Stop" stop endif NS = ns1 return end function
Subroutine : | |
this : | type(zolotarev_signfunc), intent(inout) |
Lmax : | real(DP), intent(in) |
Lmin : | real(DP), intent(in) |
eps : | real(DP), intent(in) |
subroutine new_zolotarev(this,Lmax,Lmin,eps) implicit none type(zolotarev_signfunc), intent(inout) :: this real(DP), intent(in) :: Lmax,Lmin,eps real(QP), allocatable :: rr(:) real(QP) :: KKpm,sn,cn,dn,KK,x integer :: i,j this%kappa = Lmin/Lmax this%kappa_prime = sqrt((1.0_QP-this%kappa)*(1.0_QP+this%kappa)) this%eps = eps this%NS = get_NS_exact(this%kappa,this%eps) if (associated(this%sn)) deallocate(this%sn) if (associated(this%cn)) deallocate(this%cn) if (associated(this%dn)) deallocate(this%dn) if (associated(this%aa)) deallocate(this%aa) if (associated(this%omg)) deallocate(this%omg) allocate(this%sn(this%NS)) allocate(this%cn(this%NS)) allocate(this%dn(this%NS)) allocate(this%aa(this%NS)) allocate(this%omg(this%NS)) allocate(rr(this%NS)) if (mod(this%NS,2)==0) then this%flag=.true. ! NS is even else this%flag=.false. ! NS is odd endif !========================== ! compute pole location !========================== ! ! KKpm = K(sqrt(1-kappa^2)) = K' ! KKpm = comp_elliptic_integ_K(this%kappa) this%KKpm = KKpm do i=1,this%NS !========================================================== ! z = i*K(kappa')/NS ! sn(z,kappa'),cn(z,kappa'),dn(z,kappa'),K(kappa') !========================================================== call sncndnK2(i*KKpm/this%NS,this%kappa,sn,cn,dn,KK) this%sn(i) = sn this%cn(i) = cn this%dn(i) = dn if (cn > 0.0_QP) then this%aa(i) = -(this%kappa*sn/cn)**2 endif enddo !==================== ! scale factor MM !==================== this%MM=1.0_QP do i=1,this%NS/2-1 this%MM = this%MM*this%sn(2*i-1)/this%sn(2*i) enddo i = this%NS/2 if (this%flag) then this%MM = this%MM*this%sn(2*i-1) else this%MM = this%MM*this%sn(2*i-1)/this%sn(2*i) endif this%MM=this%MM**2 !=============================== ! compute overall factor ! via min/max error relation. ! use first and second min/max ! values !=============================== this%AC = 0.0_QP do j=1,2 ! ! min/max location x for error sqrt(x)R(x) \sim 1. ! x = this%kappa**2/(1.0_QP-(1.0_QP-this%kappa**2)*this%sn(j)**2) rr(j) = 1.0_QP do i=1,this%NS/2-1 rr(j) = rr(j) * (x-this%aa(2*i))/(x-this%aa(2*i-1)) enddo i = this%NS/2 if (this%flag) then rr(j) = rr(j) /(x-this%aa(2*i-1)) else rr(j) = rr(j) * (x-this%aa(2*i))/(x-this%aa(2*i-1)) endif rr(j) = sqrt(x)*rr(j) this%AC = this%AC+rr(j) enddo this%AC = 2.0_QP/this%AC !--------------- ! error bound ! delta !--------------- do j=1,1 x = 1.0_QP - this%AC*rr(j) if (x > 0) then this%delta = x else this%delta = -x endif enddo this%delta_ex = errorbound_zol(this%kappa,this%NS) !--------------- ! scaled modulus ! lambda !--------------- this%lambda = (1.0_QP - this%delta)/(1.0_QP + this%delta) deallocate(rr) !--------------------- ! root of ! |sign(x)-R(x)| = 0 ! omg <= x !--------------------- do i=1,this%NS x = sqrt((1.0_QP+3*this%lambda)/(1.0_QP+this%lambda)**3) sn = imcomp_elliptic_integ1(x,this%lambda) x = (-1)**(i-1)*this%MM*sn + floor(real(i,kind=KIND(1.0_QP))/2)*2*this%KKpm/this%NS call sncndnK2(x,this%kappa,sn,cn,dn,KK) this%omg(i) = dn enddo return end subroutine
Function : | |
fd : | real(DP) |
this : | type(zolotarev_signfunc), intent(in) |
xd : | real(DP), intent(in) |
approximation with Rational form
function signfunc_zolotarev(this,xd) result(fd) ! ! approximation with Rational form ! implicit none type(zolotarev_signfunc), intent(in) :: this real(DP), intent(in) :: xd real(DP) :: fd real(QP) :: f,x integer :: i x = xd f = 1.0_QP do i=1,this%NS/2-1 f = f * (x**2-this%aa(2*i))/(x**2-this%aa(2*i-1)) enddo i = this% NS/2 if (this%flag) then f = f / (x**2-this%aa(2*i-1)) else f = f * (x**2-this%aa(2*i)) / (x**2-this%aa(2*i-1)) endif f = this%AC*x*f fd = f return end function
Function : | |
fd : | real(DP) |
this : | type(zolotarev_signfunc), intent(in) |
xd : | real(DP), intent(in) |
approximation with Rational form version 2 (polar)
function signfunc_zolotarev2(this,xd) result(fd) ! ! approximation with Rational form version 2 (polar) ! implicit none type(zolotarev_signfunc), intent(in) :: this real(DP), intent(in) :: xd real(DP) :: fd real(QP) :: f,tp,tm,x integer :: i x = xd tp = 1.0_QP tm = 1.0_QP do i=1,this%NS tm = tm*(1-this%omg(i)*x/this%kappa) tp = tp*(1+this%omg(i)*x/this%kappa) enddo f = (tp - tm)/(tp + tm) fd = f return end function
Derived Type : | |
NS : | integer |
sn(:) => NULL() : | real(QP), pointer |
cn(:) => NULL() : | real(QP), pointer |
dn(:) => NULL() : | real(QP), pointer |
aa(:) => NULL() : | real(QP), pointer |
omg(:) => NULL() : | real(QP), pointer |
AC : | real(QP) |
kappa : | real(QP) |
kappa_prime : | real(QP) |
Lmax : | real(QP) |
Lmin : | real(QP) |
eps : | real(QP) |
lambda : | real(QP) |
delta : | real(QP) |
delta_ex : | real(QP) |
MM : | real(QP) |
KKpm : | real(QP) |
flag : | logical |