cvqm Subroutine

public subroutine cvqm(m, q, a0)

Compute the characteristic value of Mathieu functions for q ≤ m*m

Arguments

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

Order of Mathieu functions

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

Parameter of Mathieu functions

real(kind=wp), intent(out) :: a0

Initial characteristic value


Called by

proc~~cvqm~~CalledByGraph proc~cvqm specfun_module::cvqm proc~cv0 specfun_module::cv0 proc~cv0->proc~cvqm proc~cva2 specfun_module::cva2 proc~cva2->proc~cvqm proc~cva2->proc~cv0 proc~mtu0 specfun_module::mtu0 proc~mtu0->proc~cva2 proc~mtu12 specfun_module::mtu12 proc~mtu12->proc~cva2

Source Code

   subroutine cvqm(m,q,a0)

      integer,intent(in) :: m !! Order of Mathieu functions
      real(wp),intent(in) :: q !! Parameter of Mathieu functions
      real(wp),intent(out) :: a0 !! Initial characteristic value

      real(wp) :: hm1 , hm3 , hm5

      hm1 = 0.5_wp*q/(m*m-1.0_wp)
      hm3 = 0.25_wp*hm1**3/(m*m-4.0_wp)
      hm5 = hm1*hm3*q/((m*m-1.0_wp)*(m*m-9.0_wp))
      a0 = m*m + q*(hm1+(5.0_wp*m*m+7.0_wp)*hm3+(9.0_wp*m**4+58.0_wp*m*m+29.0_wp)*hm5)

   end subroutine cvqm