refine Subroutine

public subroutine refine(Kd, m, q, a)

calculate the accurate characteristic value by the secant method

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: Kd
integer, intent(in) :: m

Order of Mathieu functions

real(kind=wp), intent(in) :: q

Parameter of Mathieu functions

real(kind=wp), intent(inout) :: a

Parameter of Mathieu functions


Calls

proc~~refine~~CallsGraph proc~refine specfun_module::refine proc~cvf specfun_module::cvf proc~refine->proc~cvf

Called by

proc~~refine~~CalledByGraph proc~refine specfun_module::refine proc~cva2 specfun_module::cva2 proc~cva2->proc~refine proc~mtu0 specfun_module::mtu0 proc~mtu0->proc~cva2 proc~mtu12 specfun_module::mtu12 proc~mtu12->proc~cva2

Source Code

   subroutine refine(Kd,m,q,a)

      integer,intent(in) :: Kd
      integer,intent(in) :: m !! Order of Mathieu functions
      real(wp),intent(in) :: q !! Parameter of Mathieu functions
      real(wp),intent(inout) :: a !! Parameter of Mathieu functions

      real(wp) :: ca , delta , f , f0 , f1 , x , x0 , x1
      integer :: it , mj

      real(wp),parameter :: eps = 1.0e-14_wp

      mj = 10 + m
      ca = a
      delta = 0.0_wp
      x0 = a
      call cvf(Kd,m,q,x0,mj,f0)
      x1 = 1.002_wp*a
      call cvf(Kd,m,q,x1,mj,f1)
      do it = 1 , 100
         mj = mj + 1
         x = x1 - (x1-x0)/(1.0_wp-f0/f1)
         call cvf(Kd,m,q,x,mj,f)
         if ( abs(1.0_wp-x1/x)<eps .or. f==0.0_wp ) exit
         x0 = x1
         f0 = f1
         x1 = x
         f1 = f
      enddo
      a = x

   end subroutine refine