cva2 Subroutine

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

Calculate a specific characteristic value of Mathieu functions

Arguments

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

Case code:

  • KD=1 for cem(x,q) ( m = 0,2,4,...)
  • KD=2 for cem(x,q) ( m = 1,3,5,...)
  • KD=3 for sem(x,q) ( m = 1,3,5,...)
  • KD=4 for sem(x,q) ( m = 2,4,6,...)
integer, intent(in) :: m

Order of Mathieu functions

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

Parameter of Mathieu functions

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

Characteristic value


Calls

proc~~cva2~~CallsGraph proc~cva2 specfun_module::cva2 proc~cv0 specfun_module::cv0 proc~cva2->proc~cv0 proc~cvql specfun_module::cvql proc~cva2->proc~cvql proc~cvqm specfun_module::cvqm proc~cva2->proc~cvqm proc~refine specfun_module::refine proc~cva2->proc~refine proc~cv0->proc~cvql proc~cv0->proc~cvqm proc~cvf specfun_module::cvf proc~refine->proc~cvf

Called by

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

Source Code

   subroutine cva2(Kd,m,q,a)

      integer,intent(in) :: m !! Order of Mathieu functions
      real(wp),intent(in) :: q !! Parameter of Mathieu functions
      integer,intent(in) :: Kd !! Case code:
                               !!
                               !! * KD=1 for cem(x,q)  ( m = 0,2,4,...)
                               !! * KD=2 for cem(x,q)  ( m = 1,3,5,...)
                               !! * KD=3 for sem(x,q)  ( m = 1,3,5,...)
                               !! * KD=4 for sem(x,q)  ( m = 2,4,6,...)
      real(wp),intent(out) :: a !! Characteristic value

      real(wp) :: a1 , a2 , delta , q1 , q2 , qq
      integer :: i , iflag , ndiv , nn

      if ( m<=12 .or. q<=3.0_wp*m .or. q>m*m ) then
         call cv0(Kd,m,q,a)
         if ( q/=0.0_wp .and. m/=2 ) call refine(Kd,m,q,a)
         if ( q>2.0e-3_wp .and. m==2 ) call refine(Kd,m,q,a)
      else
         ndiv = 10
         delta = (m-3.0_wp)*m/ndiv
         if ( (q-3.0_wp*m)<=(m*m-q) ) then
            do
               nn = int((q-3.0_wp*m)/delta) + 1
               delta = (q-3.0_wp*m)/nn
               q1 = 2.0_wp*m
               call cvqm(m,q1,a1)
               q2 = 3.0_wp*m
               call cvqm(m,q2,a2)
               qq = 3.0_wp*m
               do i = 1 , nn
                  qq = qq + delta
                  a = (a1*q2-a2*q1+(a2-a1)*qq)/(q2-q1)
                  iflag = 1
                  if ( i==nn ) iflag = -1
                  call refine(Kd,m,qq,a)
                  q1 = q2
                  q2 = qq
                  a1 = a2
                  a2 = a
               enddo
               if ( iflag/=-10 ) exit
               ndiv = ndiv*2
               delta = (m-3.0_wp)*m/ndiv
            end do
         else
            do
              nn = int((m*m-q)/delta) + 1
               delta = (m*m-q)/nn
               q1 = m*(m-1.0_wp)
               call cvql(Kd,m,q1,a1)
               q2 = m*m
               call cvql(Kd,m,q2,a2)
               qq = m*m
               do i = 1 , nn
                  qq = qq - delta
                  a = (a1*q2-a2*q1+(a2-a1)*qq)/(q2-q1)
                  iflag = 1
                  if ( i==nn ) iflag = -1
                  call refine(Kd,m,qq,a)
                  q1 = q2
                  q2 = qq
                  a1 = a2
                  a2 = a
               enddo
               if ( iflag/=-10 ) exit
               ndiv = ndiv*2
               delta = (m-3.0_wp)*m/ndiv
            end do
         endif

      endif

      end subroutine cva2