GNU Octave  4.2.1
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
cbunk.f
Go to the documentation of this file.
1  SUBROUTINE cbunk(Z, FNU, KODE, MR, N, Y, NZ, TOL, ELIM, ALIM)
2 C***BEGIN PROLOGUE CBUNK
3 C***REFER TO CBESK,CBESH
4 C
5 C CBUNK COMPUTES THE K BESSEL FUNCTION FOR FNU.GT.FNUL.
6 C ACCORDING TO THE UNIFORM ASYMPTOTIC EXPANSION FOR K(FNU,Z)
7 C IN CUNK1 AND THE EXPANSION FOR H(2,FNU,Z) IN CUNK2
8 C
9 C***ROUTINES CALLED CUNK1,CUNK2
10 C***END PROLOGUE CBUNK
11  COMPLEX Y, Z
12  REAL ALIM, AX, AY, ELIM, FNU, TOL, XX, YY
13  INTEGER KODE, MR, N, NZ
14  dimension y(n)
15  nz = 0
16  xx = REAL(z)
17  yy = aimag(z)
18  ax = abs(xx)*1.7321e0
19  ay = abs(yy)
20  IF (ay.GT.ax) go to 10
21 C-----------------------------------------------------------------------
22 C ASYMPTOTIC EXPANSION FOR K(FNU,Z) FOR LARGE FNU APPLIED IN
23 C -PI/3.LE.ARG(Z).LE.PI/3
24 C-----------------------------------------------------------------------
25  CALL cunk1(z, fnu, kode, mr, n, y, nz, tol, elim, alim)
26  go to 20
27  10 CONTINUE
28 C-----------------------------------------------------------------------
29 C ASYMPTOTIC EXPANSION FOR H(2,FNU,Z*EXP(M*HPI)) FOR LARGE FNU
30 C APPLIED IN PI/3.LT.ABS(ARG(Z)).LE.PI/2 WHERE M=+I OR -I
31 C AND HPI=PI/2
32 C-----------------------------------------------------------------------
33  CALL cunk2(z, fnu, kode, mr, n, y, nz, tol, elim, alim)
34  20 CONTINUE
35  RETURN
36  END
subroutine cunk1(Z, FNU, KODE, MR, N, Y, NZ, TOL, ELIM, ALIM)
Definition: cunk1.f:1
subroutine cbunk(Z, FNU, KODE, MR, N, Y, NZ, TOL, ELIM, ALIM)
Definition: cbunk.f:1
subroutine cunk2(Z, FNU, KODE, MR, N, Y, NZ, TOL, ELIM, ALIM)
Definition: cunk2.f:1
may be zero for pure relative error test tem the relative tolerance must be greater than or equal to
Definition: Quad-opts.cc:233
OCTAVE_EXPORT octave_value_list etc The functions then dimension(columns)
OCTAVE_EXPORT octave_value_list return the value of the option it must match the dimension of the state and the relative tolerance must also be a vector of the same length tem it must match the dimension of the state and the absolute tolerance must also be a vector of the same length The local error test applied at each integration step is xample roup abs(local error in x(i))<