cv0 Subroutine

public subroutine cv0(Kd, m, q, a0)

Compute the initial characteristic value of Mathieu functions for m ≤ 12 or q ≤ 300 or q ≥ m*m

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(out) :: a0

Characteristic value


Calls

proc~~cv0~~CallsGraph proc~cv0 specfun_module::cv0 proc~cvql specfun_module::cvql proc~cv0->proc~cvql proc~cvqm specfun_module::cvqm proc~cv0->proc~cvqm

Called by

proc~~cv0~~CalledByGraph proc~cv0 specfun_module::cv0 proc~cva2 specfun_module::cva2 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 cv0(Kd,m,q,a0)

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

      real(wp) :: q2

      q2 = q*q
      if ( m==0 ) then
         if ( q<=1.0_wp ) then
            a0 = (((.0036392_wp*q2-.0125868_wp)*q2+.0546875_wp)*q2-.5_wp)*q2
         elseif ( q<=10.0_wp ) then
            a0 = ((3.999267e-3_wp*q-9.638957e-2_wp)*q-.88297_wp)*q + .5542818_wp
         else
            call cvql(Kd,m,q,a0)
         endif
      elseif ( m==1 ) then
         if ( q<=1.0 .and. Kd==2 ) then
            a0 = (((-6.51e-4_wp*q-.015625_wp)*q-.125_wp)*q+1.0_wp)*q + 1.0_wp
         elseif ( q<=1.0 .and. Kd==3 ) then
            a0 = (((-6.51e-4_wp*q+.015625_wp)*q-.125_wp)*q-1.0_wp)*q + 1.0_wp
         elseif ( q<=10.0 .and. Kd==2 ) then
            a0 = (((-4.94603e-4_wp*q+1.92917e-2_wp)*q-.3089229_wp)*q+1.33372_wp) &
               & *q + .811752_wp
         elseif ( q<=10.0_wp .and. Kd==3 ) then
            a0 = ((1.971096e-3_wp*q-5.482465e-2_wp)*q-1.152218_wp)*q + 1.10427_wp
         else
            call cvql(Kd,m,q,a0)
         endif
      elseif ( m==2 ) then
         if ( q<=1.0 .and. Kd==1 ) then
            a0 = (((-.0036391_wp*q2+.0125888_wp)*q2-.0551939_wp)*q2+.416667_wp) &
                 *q2 + 4.0_wp
         elseif ( q<=1.0_wp .and. Kd==4 ) then
            a0 = (.0003617_wp*q2-.0833333_wp)*q2 + 4.0_wp
         elseif ( q<=15 .and. Kd==1 ) then
            a0 = (((3.200972e-4_wp*q-8.667445e-3_wp)*q-1.829032e-4_wp) &
                  *q+.9919999_wp)*q + 3.3290504_wp
         elseif ( q<=10.0_wp .and. Kd==4 ) then
            a0 = ((2.38446e-3_wp*q-.08725329_wp)*q-4.732542e-3_wp)*q + 4.00909_wp
         else
            call cvql(Kd,m,q,a0)
         endif
      elseif ( m==3 ) then
         if ( q<=1.0_wp .and. Kd==2 ) then
            a0 = ((6.348e-4_wp*q+.015625_wp)*q+.0625_wp)*q2 + 9.0_wp
         elseif ( q<=1.0_wp .and. Kd==3 ) then
            a0 = ((6.348e-4_wp*q-.015625_wp)*q+.0625_wp)*q2 + 9.0
         elseif ( q<=20.0_wp .and. Kd==2 ) then
            a0 = (((3.035731e-4_wp*q-1.453021e-2_wp)*q+.19069602_wp)*q-.1039356_wp) &
                  *q + 8.9449274_wp
         elseif ( q<=15.0_wp .and. Kd==3 ) then
            a0 = ((9.369364e-5_wp*q-.03569325_wp)*q+.2689874_wp)*q + 8.771735_wp
         else
            call cvql(Kd,m,q,a0)
         endif
      elseif ( m==4 ) then
         if ( q<=1.0 .and. Kd==1 ) then
            a0 = ((-2.1e-6_wp*q2+5.012e-4_wp)*q2+.0333333_wp)*q2 + 16.0_wp
         elseif ( q<=1.0_wp .and. Kd==4 ) then
            a0 = ((3.7e-6_wp*q2-3.669e-4_wp)*q2+.0333333_wp)*q2 + 16.0_wp
         elseif ( q<=25.0_wp .and. Kd==1 ) then
            a0 = (((1.076676e-4_wp*q-7.9684875e-3_wp)*q+.17344854_wp)*q-.5924058_wp) &
                  *q + 16.620847_wp
         elseif ( q<=20.0_wp .and. Kd==4 ) then
            a0 = ((-7.08719e-4_wp*q+3.8216144e-3_wp)*q+.1907493_wp)*q + 15.744_wp
         else
            call cvql(Kd,m,q,a0)
         endif
      elseif ( m==5 ) then
         if ( q<=1.0_wp .and. Kd==2 ) then
            a0 = ((6.8e-6_wp*q+1.42e-5_wp)*q2+.0208333_wp)*q2 + 25.0_wp
         elseif ( q<=1.0_wp .and. Kd==3 ) then
            a0 = ((-6.8e-6_wp*q+1.42e-5_wp)*q2+.0208333_wp)*q2 + 25.0_wp
         elseif ( q<=35.0_wp .and. Kd==2 ) then
            a0 = (((2.238231e-5_wp*q-2.983416e-3_wp)*q+.10706975_wp)*q-.600205_wp) &
                  *q + 25.93515_wp
         elseif ( q<=25.0_wp .and. Kd==3 ) then
            a0 = ((-7.425364e-4_wp*q+2.18225e-2_wp)*q+4.16399e-2_wp)*q + 24.897_wp
         else
            call cvql(Kd,m,q,a0)
         endif
      elseif ( m==6 ) then
         if ( q<=1.0_wp ) then
            a0 = (.4e-6_wp*q2+.0142857_wp)*q2 + 36.0_wp
         elseif ( q<=40.0_wp .and. Kd==1 ) then
            a0 = (((-1.66846e-5_wp*q+4.80263e-4_wp)*q+2.53998e-2_wp)*q-.181233_wp) &
                  *q + 36.423_wp
         elseif ( q<=35.0_wp .and. Kd==4 ) then
            a0 = ((-4.57146e-4_wp*q+2.16609e-2_wp)*q-2.349616e-2_wp)*q + 35.99251_wp
         else
            call cvql(Kd,m,q,a0)
         endif
      elseif ( m==7 ) then
         if ( q<=10.0_wp ) then
            call cvqm(m,q,a0)
         elseif ( q<=50.0_wp .and. Kd==2 ) then
            a0 = (((-1.411114e-5_wp*q+9.730514e-4_wp)*q-3.097887e-3_wp) &
               & *q+3.533597e-2_wp)*q + 49.0547_wp
         elseif ( q<=40.0_wp .and. Kd==3 ) then
            a0 = ((-3.043872e-4_wp*q+2.05511e-2_wp)*q-9.16292e-2_wp)*q + 49.19035_wp
         else
            call cvql(Kd,m,q,a0)
         endif
      elseif ( m>=8 ) then
         if ( q<=3.0_wp*m ) then
            call cvqm(m,q,a0)
         elseif ( q>m*m ) then
            call cvql(Kd,m,q,a0)
         elseif ( m==8 .and. Kd==1 ) then
            a0 = (((8.634308e-6_wp*q-2.100289e-3_wp)*q+.169072_wp)*q-4.64336_wp) &
                  *q + 109.4211_wp
         elseif ( m==8 .and. Kd==4 ) then
            a0 = ((-6.7842e-5_wp*q+2.2057e-3_wp)*q+.48296_wp)*q + 56.59_wp
         elseif ( m==9 .and. Kd==2 ) then
            a0 = (((2.906435e-6_wp*q-1.019893e-3_wp)*q+.1101965_wp)*q-3.821851_wp) &
                  *q + 127.6098_wp
         elseif ( m==9 .and. Kd==3 ) then
            a0 = ((-9.577289e-5_wp*q+.01043839_wp)*q+.06588934_wp)*q + 78.0198_wp
         elseif ( m==10 .and. Kd==1 ) then
            a0 = (((5.44927e-7_wp*q-3.926119e-4_wp)*q+.0612099_wp)*q-2.600805_wp) &
                  *q + 138.1923_wp
         elseif ( m==10 .and. Kd==4 ) then
            a0 = ((-7.660143e-5_wp*q+.01132506_wp)*q-.09746023_wp)*q + 99.29494_wp
         elseif ( m==11 .and. Kd==2 ) then
            a0 = (((-5.67615e-7_wp*q+7.152722e-6_wp)*q+.01920291_wp)*q-1.081583_wp) &
                  *q + 140.88_wp
         elseif ( m==11 .and. Kd==3 ) then
            a0 = ((-6.310551e-5_wp*q+.0119247_wp)*q-.2681195_wp)*q + 123.667_wp
         elseif ( m==12 .and. Kd==1 ) then
            a0 = (((-2.38351e-7_wp*q-2.90139e-5_wp)*q+.02023088_wp)*q-1.289_wp) &
               & *q + 171.2723_wp
         elseif ( m==12 .and. Kd==4 ) then
            a0 = (((3.08902e-7_wp*q-1.577869e-4_wp)*q+.0247911_wp)*q-1.05454_wp) &
                  *q + 161.471_wp
         endif
      endif

   end subroutine cv0