rswfp Subroutine

public subroutine rswfp(m, n, c, x, cv, Kf, R1f, R1d, R2f, R2d)

Compute prolate spheriodal radial functions of the first and second kinds, and their derivatives

Arguments

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

Mode parameter, m = 0,1,2,...

integer, intent(in) :: n

Mode parameter, n = m,m+1,m+2,...

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

Spheroidal parameter

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

Argument of radial function ( x > 1.0 )

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

Characteristic value

integer, intent(in) :: Kf

Function code: * KF=1 for the first kind * KF=2 for the second kind * KF=3 for both the first and second kinds

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

Radial function of the first kind

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

Derivative of the radial function of the first kind

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

Radial function of the second kind

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

Derivative of the radial function of the second kind


Calls

proc~~rswfp~~CallsGraph proc~rswfp specfun_module::rswfp proc~rmn1 specfun_module::rmn1 proc~rswfp->proc~rmn1 proc~rmn2l specfun_module::rmn2l proc~rswfp->proc~rmn2l proc~rmn2sp specfun_module::rmn2sp proc~rswfp->proc~rmn2sp proc~sdmn specfun_module::sdmn proc~rswfp->proc~sdmn proc~sckb specfun_module::sckb proc~rmn1->proc~sckb proc~sphj specfun_module::sphj proc~rmn1->proc~sphj proc~sphy specfun_module::sphy proc~rmn2l->proc~sphy proc~kmn specfun_module::kmn proc~rmn2sp->proc~kmn proc~lpmns specfun_module::lpmns proc~rmn2sp->proc~lpmns proc~lqmns specfun_module::lqmns proc~rmn2sp->proc~lqmns proc~msta1 specfun_module::msta1 proc~sphj->proc~msta1 proc~msta2 specfun_module::msta2 proc~sphj->proc~msta2 proc~envj specfun_module::envj proc~msta1->proc~envj proc~msta2->proc~envj

Source Code

   subroutine rswfp(m,n,c,x,Cv,Kf,R1f,R1d,R2f,R2d)

      integer,intent(in) :: m  !! Mode parameter, `m = 0,1,2,...`
      integer,intent(in) :: n  !! Mode parameter, `n = m,m+1,m+2,...`
      real(wp),intent(in) :: c  !! Spheroidal parameter
      real(wp),intent(in) :: x  !! Argument of radial function ( `x > 1.0` )
      real(wp),intent(in) :: cv !! Characteristic value
      integer,intent(in) :: Kf !! Function code:
                               !!  * `KF=1` for the first kind
                               !!  * `KF=2` for the second kind
                               !!  * `KF=3` for both the first and second kinds
      real(wp),intent(out) :: R1f !! Radial function of the first kind
      real(wp),intent(out) :: R1d !! Derivative of the radial function of
                                  !! the first kind
      real(wp),intent(out) :: R2f !! Radial function of the second kind
      real(wp),intent(out) :: R2d !! Derivative of the radial function of
                                  !! the second kind

      real(wp) :: df(200)
      integer :: id , kd

      kd = 1
      call sdmn(m,n,c,Cv,kd,df)
      if ( Kf/=2 ) call rmn1(m,n,c,x,df,kd,R1f,R1d)
      if ( Kf>1 ) then
         call rmn2l(m,n,c,x,df,kd,R2f,R2d,id)
         if ( id>-8 ) call rmn2sp(m,n,c,x,Cv,df,kd,R2f,R2d)
      endif

   end subroutine rswfp