Scrambler  1
specfun.f90
Go to the documentation of this file.
00001 !#########################################################################
00002 !               
00003 !    Copyright (C) 2003-2012 Department of Physics and Astronomy,
00004 !                            University of Rochester,
00005 !                            Rochester, NY
00006 !
00007 !    specfun.f90 is part of AstroBEAR.
00008 !
00009 !    AstroBEAR is free software: you can redistribute it and/or modify    
00010 !    it under the terms of the GNU General Public License as published by 
00011 !    the Free Software Foundation, either version 3 of the License, or    
00012 !    (at your option) any later version.
00013 !
00014 !    AstroBEAR is distributed in the hope that it will be useful, 
00015 !    but WITHOUT ANY WARRANTY; without even the implied warranty of
00016 !    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00017 !    GNU General Public License for more details.
00018 !
00019 !    You should have received a copy of the GNU General Public License
00020 !    along with AstroBEAR.  If not, see <http://www.gnu.org/licenses/>.
00021 !
00022 !#########################################################################
00023 ! Modified specfun library to be a module to be built along with AstroBEAR
00024 ! Can replace netlib.f
00025 MODULE SpecFun
00026 
00027    PUBLIC
00028 
00029 CONTAINS
00030 
00031    function besei0 ( x )
00032 
00033       !*****************************************************************************80
00034       !
00035       !! BESEI0 evaluates the exponentially scaled Bessel I0(X) function.
00036       !
00037       !  Discussion:
00038       !
00039       !    This routine computes approximate values for the modified Bessel
00040       !    function of the first kind of order zero multiplied by EXP(-ABS(X)).
00041       !
00042       !  Licensing:
00043       !
00044       !    This code is distributed under the GNU LGPL license.
00045       !
00046       !  Modified:
00047       !
00048       !    03 April 2007
00049       !
00050       !  Author:
00051       !
00052       !    Original FORTRAN77 version by William Cody.
00053       !    FORTRAN90 version by John Burkardt.
00054       !
00055       !  Parameters:
00056       !
00057       !    Input, real ( kind = 8 ) X, the argument of the function.
00058       !
00059       !    Output, real ( kind = 8 ) BESEI0, the value of the function.
00060       !
00061       implicit none
00062 
00063       real ( kind = 8 ) besei0
00064       integer ( kind = 4 ) jint
00065       real ( kind = 8 ) result
00066       real ( kind = 8 ) x
00067 
00068       jint = 2
00069       call calci0 ( x, result, jint )
00070       besei0 = result
00071 
00072       return
00073    end function besei0
00074    function besei1 ( x )
00075 
00076       !*****************************************************************************80
00077       !
00078       !! BESEI1 evaluates the exponentially scaled Bessel I1(X) function.
00079       !
00080       !  Discussion:
00081       !
00082       !    This routine computes approximate values for the
00083       !    modified Bessel function of the first kind of order one
00084       !    multiplied by EXP(-ABS(X)).
00085       !
00086       !  Licensing:
00087       !
00088       !    This code is distributed under the GNU LGPL license.
00089       !
00090       !  Modified:
00091       !
00092       !    03 April 2007
00093       !
00094       !  Author:
00095       !
00096       !    Original FORTRAN77 version by William Cody.
00097       !    FORTRAN90 version by John Burkardt.
00098       !
00099       !  Parameters:
00100       !
00101       !    Input, real ( kind = 8 ) X, the argument of the function.
00102       !
00103       !    Output, real ( kind = 8 ) BESEI1, the value of the function.
00104       !
00105       implicit none
00106 
00107       real ( kind = 8 ) besei1
00108       integer ( kind = 4 ) jint
00109       real ( kind = 8 ) result
00110       real ( kind = 8 ) x
00111 
00112       jint = 2
00113       call calci1 ( x, result, jint )
00114       besei1 = result
00115 
00116       return
00117    end function besei1
00118    function besek0 ( x )
00119 
00120       !*****************************************************************************80
00121       !
00122       !! BESEK0 evaluates the exponentially scaled Bessel K0(X) function.
00123       !
00124       !  Discussion:
00125       !
00126       !    This routine computes approximate values for the
00127       !    modified Bessel function of the second kind of order zero
00128       !    multiplied by the exponential function.
00129       !
00130       !  Licensing:
00131       !
00132       !    This code is distributed under the GNU LGPL license.
00133       !
00134       !  Modified:
00135       !
00136       !    03 April 2007
00137       !
00138       !  Author:
00139       !
00140       !    Original FORTRAN77 version by William Cody, Laura Stoltz.
00141       !    FORTRAN90 version by John Burkardt.
00142       !
00143       !  Parameters:
00144       !
00145       !    Input, real ( kind = 8 ) X, the argument of the function.
00146       !    0 < X.
00147       !
00148       !    Output, real ( kind = 8 ) BESK0, the value of the function.
00149       !
00150       implicit none
00151 
00152       real ( kind = 8 ) besek0
00153       integer ( kind = 4 ) jint
00154       real ( kind = 8 ) result
00155       real ( kind = 8 ) x
00156 
00157       jint = 2
00158       call calck0 ( x, result, jint )
00159       besek0 = result
00160 
00161       return
00162    end function besek0
00163    function besek1 ( x )
00164 
00165       !*****************************************************************************80
00166       !
00167       !! BESEK1 evaluates the exponentially scaled Bessel K1(X) function.
00168       !
00169       !  Discussion:
00170       !
00171       !    This routine computes approximate values for the
00172       !    modified Bessel function of the second kind of order one
00173       !    multiplied by the exponential function, for arguments
00174       !    XLEAST <= ARG <= XMAX.
00175       !
00176       !  Licensing:
00177       !
00178       !    This code is distributed under the GNU LGPL license.
00179       !
00180       !  Modified:
00181       !
00182       !    03 April 2007
00183       !
00184       !  Author:
00185       !
00186       !    Original FORTRAN77 version by William Cody.
00187       !    FORTRAN90 version by John Burkardt.
00188       !
00189       !  Parameters:
00190       !
00191       !    Input, real ( kind = 8 ) X, the argument of the function.
00192       !
00193       !    Output, real ( kind = 8 ) BESEK1, the value of the function.
00194       !
00195       implicit none
00196 
00197       real ( kind = 8 ) besek1
00198       integer ( kind = 4 ) jint
00199       real ( kind = 8 ) result
00200       real ( kind = 8 ) x
00201 
00202       jint = 2
00203       call calck1 ( x, result, jint )
00204       besek1 = result
00205 
00206       return
00207    end function besek1
00208    function besi0 ( x )
00209 
00210       !*****************************************************************************80
00211       !
00212       !! BESI0 evaluates the Bessel I0(X) function.
00213       !
00214       !  Discussion:
00215       !
00216       !    This routine computes approximate values for
00217       !    modified Bessel functions of the first kind of order zero for
00218       !    arguments ABS(ARG) <= XMAX.
00219       !
00220       !    See comments heading CALCI0.
00221       !
00222       !  Licensing:
00223       !
00224       !    This code is distributed under the GNU LGPL license.
00225       !
00226       !  Modified:
00227       !
00228       !    03 April 2007
00229       !
00230       !  Author:
00231       !
00232       !    Original FORTRAN77 version by William Cody.
00233       !    FORTRAN90 version by John Burkardt.
00234       !
00235       !  Parameters:
00236       !
00237       !    Input, real ( kind = 8 ) X, the argument of the function.
00238       !
00239       !    Output, real ( kind = 8 ) BESI0, the value of the function.
00240       !
00241       implicit none
00242 
00243       real ( kind = 8 ) besi0
00244       integer ( kind = 4 ) jint
00245       real ( kind = 8 ) result
00246       real ( kind = 8 ) x
00247 
00248       jint = 1
00249       call calci0 ( x, result, jint )
00250       besi0 = result
00251 
00252       return
00253    end function besi0
00254    function besi1 ( x )
00255 
00256       !*****************************************************************************80
00257       !
00258       !! BESI1 evaluates the Bessel I1(X) function.
00259       !
00260       !  Discussion:
00261       !
00262       !    This routine computes approximate values for
00263       !    modified Bessel functions of the first kind of order one for
00264       !    arguments ABS(ARG) <= XMAX.
00265       !
00266       !    See comments heading CALCI1.
00267       !
00268       !  Licensing:
00269       !
00270       !    This code is distributed under the GNU LGPL license.
00271       !
00272       !  Modified:
00273       !
00274       !    03 April 2007
00275       !
00276       !  Author:
00277       !
00278       !    Original FORTRAN77 version by William Cody.
00279       !    FORTRAN90 version by John Burkardt.
00280       !
00281       !  Parameters:
00282       !
00283       !    Input, real ( kind = 8 ) X, the argument of the function.
00284       !
00285       !    Output, real ( kind = 8 ) BESI1, the value of the function.
00286       !
00287       implicit none
00288 
00289       real ( kind = 8 ) besi1
00290       integer ( kind = 4 ) jint
00291       real ( kind = 8 ) result
00292       real ( kind = 8 ) x
00293 
00294       jint = 1
00295       call calci1 ( x, result, jint )
00296       besi1 = result
00297 
00298       return
00299    end function besi1
00300    function besj0 ( x )
00301 
00302       !*****************************************************************************80
00303       !
00304       !! BESJ0 evaluates the Bessel J0(X) function.
00305       !
00306       !  Discussion:
00307       !
00308       !    This routine computes approximate values for Bessel functions
00309       !    of the first kind of order zero for arguments  |X| <= XMAX
00310       !
00311       !    See comments heading CALJY0.
00312       !
00313       !  Licensing:
00314       !
00315       !    This code is distributed under the GNU LGPL license.
00316       !
00317       !  Modified:
00318       !
00319       !    03 April 2007
00320       !
00321       !  Author:
00322       !
00323       !    Original FORTRAN77 version by William Cody.
00324       !    FORTRAN90 version by John Burkardt.
00325       !
00326       !  Parameters:
00327       !
00328       !    Input, real ( kind = 8 ) X, the argument of the function.
00329       !
00330       !    Output, real ( kind = 8 ) BESJ0, the value of the function.
00331       !
00332       implicit none
00333 
00334       real ( kind = 8 ) besj0
00335       integer ( kind = 4 ) jint
00336       real ( kind = 8 ) result
00337       real ( kind = 8 ) x
00338 
00339       jint = 0
00340       call caljy0 ( x, result, jint )
00341       besj0 = result
00342 
00343       return
00344    end function besj0
00345    function besj1 ( x )
00346 
00347       !*****************************************************************************80
00348       !
00349       !! BESJ1 evaluates the Bessel J1(X) function.
00350       !
00351       !  Discussion:
00352       !
00353       !    This routine computes approximate values for Bessel functions
00354       !    of the first kind of order zero for arguments  |X| <= XMAX
00355       !
00356       !    See comments heading CALJY1.
00357       !
00358       !  Licensing:
00359       !
00360       !    This code is distributed under the GNU LGPL license.
00361       !
00362       !  Modified:
00363       !
00364       !    03 April 2007
00365       !
00366       !  Author:
00367       !
00368       !    Original FORTRAN77 version by William Cody.
00369       !    FORTRAN90 version by John Burkardt.
00370       !
00371       !  Parameters:
00372       !
00373       !    Input, real ( kind = 8 ) X, the argument of the function.
00374       !
00375       !    Output, real ( kind = 8 ) BESJ1, the value of the function.
00376       !
00377       implicit none
00378 
00379       real ( kind = 8 ) besj1
00380       integer ( kind = 4 ) jint
00381       real ( kind = 8 ) result
00382       real ( kind = 8 ) x
00383 
00384       jint = 0
00385       call caljy1 ( x, result, jint )
00386       besj1 = result
00387 
00388       return
00389    end function besj1
00390    function besk0 ( x )
00391 
00392       !*****************************************************************************80
00393       !
00394       !! BESK0 evaluates the Bessel K0(X) function.
00395       !
00396       !  Discussion:
00397       !
00398       !    This routine computes approximate values for the
00399       !    modified Bessel function of the second kind of order zero
00400       !    for arguments 0.0 < ARG <= XMAX.
00401       !
00402       !    See comments heading CALCK0.
00403       !
00404       !  Licensing:
00405       !
00406       !    This code is distributed under the GNU LGPL license.
00407       !
00408       !  Modified:
00409       !
00410       !    03 April 2007
00411       !
00412       !  Author:
00413       !
00414       !    Original FORTRAN77 version by William Cody, Laura Stoltz.
00415       !    FORTRAN90 version by John Burkardt.
00416       !
00417       !  Parameters:
00418       !
00419       !    Input, real ( kind = 8 ) X, the argument of the function.
00420       !
00421       !    Output, real ( kind = 8 ) BESK0, the value of the function.
00422       !
00423       implicit none
00424 
00425       real ( kind = 8 ) besk0
00426       integer ( kind = 4 ) jint
00427       real ( kind = 8 ) result
00428       real ( kind = 8 ) x
00429 
00430       jint = 1
00431       call calck0 ( x, result, jint )
00432       besk0 = result
00433 
00434       return
00435    end function besk0
00436    function besk1 ( x )
00437 
00438       !*****************************************************************************80
00439       !
00440       !! BESK1 evaluates the Bessel K1(X) function.
00441       !
00442       !  Discussion:
00443       !
00444       !    This routine computes approximate values for the
00445       !    modified Bessel function of the second kind of order one
00446       !    for arguments XLEAST <= ARG <= XMAX.
00447       !
00448       !  Licensing:
00449       !
00450       !    This code is distributed under the GNU LGPL license.
00451       !
00452       !  Modified:
00453       !
00454       !    03 April 2007
00455       !
00456       !  Author:
00457       !
00458       !    Original FORTRAN77 version by William Cody.
00459       !    FORTRAN90 version by John Burkardt.
00460       !
00461       !  Parameters:
00462       !
00463       !    Input, real ( kind = 8 ) X, the argument of the function.
00464       !
00465       !    Output, real ( kind = 8 ) BESK1, the value of the function.
00466       !
00467       implicit none
00468 
00469       real ( kind = 8 ) besk1
00470       integer ( kind = 4 ) jint
00471       real ( kind = 8 ) result
00472       real ( kind = 8 ) x
00473 
00474       jint = 1
00475       call calck1 ( x, result, jint )
00476       besk1 = result
00477 
00478       return
00479    end function besk1
00480    function besy0 ( x )
00481 
00482       !*****************************************************************************80
00483       !
00484       !! BESY0 evaluates the Bessel Y0(X) function.
00485       !
00486       !  Discussion:
00487       !
00488       !    This routine computes approximate values for Bessel functions
00489       !    of the second kind of order zero for arguments 0 < X <= XMAX.
00490       !
00491       !    See comments heading CALJY0.
00492       !
00493       !  Licensing:
00494       !
00495       !    This code is distributed under the GNU LGPL license.
00496       !
00497       !  Modified:
00498       !
00499       !    03 April 2007
00500       !
00501       !  Author:
00502       !
00503       !    Original FORTRAN77 version by William Cody.
00504       !    FORTRAN90 version by John Burkardt.
00505       !
00506       !  Parameters:
00507       !
00508       !    Input, real ( kind = 8 ) X, the argument of the function.
00509       !
00510       !    Output, real ( kind = 8 ) BESY0, the value of the function.
00511       !
00512       implicit none
00513 
00514       real ( kind = 8 ) besy0
00515       integer ( kind = 4 ) jint
00516       real ( kind = 8 ) result
00517       real ( kind = 8 ) x
00518 
00519       jint = 1
00520       call caljy0 ( x, result, jint )
00521       besy0 = result
00522 
00523       return
00524    end function besy0
00525    function besy1 ( x )
00526 
00527       !*****************************************************************************80
00528       !
00529       !! BESY1 evaluates the Bessel Y1(X) function.
00530       !
00531       !  Discussion:
00532       !
00533       !    This routine computes approximate values for Bessel functions
00534       !    of the second kind of order zero for arguments 0 < X <= XMAX.
00535       !
00536       !    See comments heading CALJY1.
00537       !
00538       !  Licensing:
00539       !
00540       !    This code is distributed under the GNU LGPL license.
00541       !
00542       !  Modified:
00543       !
00544       !    03 April 2007
00545       !
00546       !  Author:
00547       !
00548       !    Original FORTRAN77 version by William Cody.
00549       !    FORTRAN90 version by John Burkardt.
00550       !
00551       !  Parameters:
00552       !
00553       !    Input, real ( kind = 8 ) X, the argument of the function.
00554       !
00555       !    Output, real ( kind = 8 ) BESY1, the value of the function.
00556       !
00557       implicit none
00558 
00559       real ( kind = 8 ) besy1
00560       integer ( kind = 4 ) jint
00561       real ( kind = 8 ) result
00562       real ( kind = 8 ) x
00563 
00564       jint = 1
00565       call caljy1 ( x, result, jint )
00566       besy1 = result
00567 
00568       return
00569    end function besy1
00570    subroutine calcei ( arg, result, jint )
00571 
00572       !*****************************************************************************80
00573       !
00574       !! CALCEI computes various exponential integrals.
00575       !
00576       !  Discussion:
00577       !
00578       !    This routine computes the exponential integrals Ei(x),
00579       !    E1(x), and  exp(-x)*Ei(x) for real arguments x where
00580       !
00581       !           integral (from t=-oo to t=x) (exp(t)/t),  x > 0,
00582       !    Ei(x) =
00583       !          -integral (from t=-x to t=+oo) (exp(t)/t),  x < 0,
00584       !
00585       !    and where the first integral is a principal value integral.
00586       !
00587       !  Licensing:
00588       !
00589       !    This code is distributed under the GNU LGPL license.
00590       !
00591       !  Modified:
00592       !
00593       !    03 April 2007
00594       !
00595       !  Author:
00596       !
00597       !    Original FORTRAN77 version by William Cody.
00598       !    FORTRAN90 version by John Burkardt.
00599       !
00600       !  Reference:
00601       !
00602       !    William Cody, Henry Thacher,
00603       !    Rational Chebyshev Approximations for the Exponential
00604       !    Integral E1(x),
00605       !    Mathematics of Computation,
00606       !    Volume 22, Number 103, July 1968, pages 641-649.
00607       !
00608       !    William Cody, Henry Thacher,
00609       !    Chebyshev Approximations for the Exponential
00610       !    Integral Ei(x),
00611       !    Mathematics of Computation,
00612       !    Volume 23, Number 106, April 1969, pages 289-303.
00613       !
00614       !  Parameters:
00615       !
00616       !    Input, real ( kind = 8 ) ARG, the argument.  The argument must not
00617       !    be zero.  If JINT = 2, then the argument must be strictly positive.
00618       !
00619       !    Output, real ( kind = 8 ) RESULT, the value of the function,
00620       !    which depends on the input value of JINT:
00621       !    1, RESULT = EI ( ARG );
00622       !    2, RESULT = EONE ( ARG );
00623       !    3, RESULT = EXPEI ( ARG ).
00624       !
00625       !    Input, integer ( kind = 4 ) JINT, chooses the function to be computed.
00626       !    1, Ei(x);
00627       !    2, -Ei(-x);
00628       !    3, exp(-x)*Ei(x).
00629       !
00630       implicit none
00631 
00632       real ( kind = 8 ) a(7)
00633       real ( kind = 8 ) arg
00634       real ( kind = 8 ) b(6)
00635       real ( kind = 8 ) c(9)
00636       real ( kind = 8 ) d(9)
00637       real ( kind = 8 ) e(10)
00638       real ( kind = 8 ) ei
00639       real ( kind = 8 ) exp40
00640       real ( kind = 8 ) f(10)
00641       real ( kind = 8 ) four
00642       real ( kind = 8 ) fourty
00643       real ( kind = 8 ) frac
00644       real ( kind = 8 ) half
00645       integer ( kind = 4 ) i
00646       integer ( kind = 4 ) jint
00647       real ( kind = 8 ) one
00648       real ( kind = 8 ) p(10)
00649       real ( kind = 8 ) plg(4)
00650       real ( kind = 8 ) px(10)
00651       real ( kind = 8 ) p037
00652       real ( kind = 8 ) p1(10)
00653       real ( kind = 8 ) p2(10)
00654       real ( kind = 8 ) q(10)
00655       real ( kind = 8 ) qlg(4)
00656       real ( kind = 8 ) qx(10)
00657       real ( kind = 8 ) q1(9)
00658       real ( kind = 8 ) q2(9)
00659       real ( kind = 8 ) r(10)
00660       real ( kind = 8 ) result
00661       real ( kind = 8 ) s(9)
00662       real ( kind = 8 ) six
00663       real ( kind = 8 ) sump
00664       real ( kind = 8 ) sumq
00665       real ( kind = 8 ) t
00666       real ( kind = 8 ) three
00667       real ( kind = 8 ) twelve
00668       real ( kind = 8 ) two
00669       real ( kind = 8 ) two4
00670       real ( kind = 8 ) w
00671       real ( kind = 8 ) x
00672       real ( kind = 8 ) x0
00673       real ( kind = 8 ) x01
00674       real ( kind = 8 ) x02
00675       real ( kind = 8 ) x11
00676       real ( kind = 8 ) xbig
00677       real ( kind = 8 ) xinf
00678       real ( kind = 8 ) xmax
00679       real ( kind = 8 ) xmx0
00680       real ( kind = 8 ) y
00681       real ( kind = 8 ) ysq
00682       real ( kind = 8 ) zero
00683       !
00684       !  Mathematical constants
00685       !  EXP40 = exp(40)
00686       !  X0 = zero of Ei
00687       !  X01/X11 + X02 = zero of Ei to extra precision
00688       !
00689       data zero /0.0d0 /
00690       data p037 /0.037d0/
00691       data half / 0.5d0 /
00692       data one / 1.0d0 /
00693       data two / 2.0d0/
00694       data three /3.0d0 /
00695       data four / 4.0d0 /
00696       data six / 6.0d0 /
00697       data twelve / 12.d0 /
00698       data two4 / 24.0d0/
00699       data fourty / 40.0d0 /
00700       data exp40 / 2.3538526683701998541d17/
00701       data x01 /381.5d0 /
00702       data x11 / 1024.0d0 /
00703       data x02 /-5.1182968633365538008d-5/
00704       data x0 / 3.7250741078136663466d-1 /
00705       !
00706       !  Machine-dependent constants
00707       !
00708       data xinf /1.79d+308/
00709       data xmax /716.351d0/
00710       data xbig /701.84d0/
00711       !
00712       !  Coefficients  for -1.0 <= X < 0.0
00713       !
00714       data a /1.1669552669734461083368d2, 2.1500672908092918123209d3, &
00715            1.5924175980637303639884d4, 8.9904972007457256553251d4, &
00716            1.5026059476436982420737d5,-1.4815102102575750838086d5, &
00717            5.0196785185439843791020d0/
00718 
00719       data b /4.0205465640027706061433d1, 7.5043163907103936624165d2, &
00720            8.1258035174768735759855d3, 5.2440529172056355429883d4, &
00721            1.8434070063353677359298d5, 2.5666493484897117319268d5/
00722       !
00723       !  Coefficients for -4.0 <= X < -1.0
00724       !
00725       data c /3.828573121022477169108d-1, 1.107326627786831743809d+1, &
00726            7.246689782858597021199d+1, 1.700632978311516129328d+2, &
00727            1.698106763764238382705d+2, 7.633628843705946890896d+1, &
00728            1.487967702840464066613d+1, 9.999989642347613068437d-1, &
00729            1.737331760720576030932d-8/
00730       data d /8.258160008564488034698d-2, 4.344836335509282083360d+0, &
00731            4.662179610356861756812d+1, 1.775728186717289799677d+2, &
00732            2.953136335677908517423d+2, 2.342573504717625153053d+2, &
00733            9.021658450529372642314d+1, 1.587964570758947927903d+1, &
00734            1.000000000000000000000d+0/
00735       !
00736       !  Coefficients for X < -4.0
00737       !
00738       data e /1.3276881505637444622987d+2,3.5846198743996904308695d+4, &
00739            1.7283375773777593926828d+5,2.6181454937205639647381d+5, &
00740            1.7503273087497081314708d+5,5.9346841538837119172356d+4, &
00741            1.0816852399095915622498d+4,1.0611777263550331766871d03, &
00742            5.2199632588522572481039d+1,9.9999999999999999087819d-1/
00743       data f /3.9147856245556345627078d+4,2.5989762083608489777411d+5, &
00744            5.5903756210022864003380d+5,5.4616842050691155735758d+5, &
00745            2.7858134710520842139357d+5,7.9231787945279043698718d+4, &
00746            1.2842808586627297365998d+4,1.1635769915320848035459d+3, &
00747            5.4199632588522559414924d+1,1.0d0/
00748       !
00749       !  Coefficients for rational approximation to ln(x/a), |1-x/a| < .1
00750       !
00751       data plg /-2.4562334077563243311d+01,2.3642701335621505212d+02, &
00752            -5.4989956895857911039d+02,3.5687548468071500413d+02/
00753       data qlg /-3.5553900764052419184d+01,1.9400230218539473193d+02, &
00754            -3.3442903192607538956d+02,1.7843774234035750207d+02/
00755       !
00756       !  Coefficients for  0.0 < X < 6.0,
00757       !  ratio of Chebyshev polynomials
00758       !
00759       data p /-1.2963702602474830028590d01,-1.2831220659262000678155d03, &
00760            -1.4287072500197005777376d04,-1.4299841572091610380064d06, &
00761            -3.1398660864247265862050d05,-3.5377809694431133484800d08, &
00762            3.1984354235237738511048d08,-2.5301823984599019348858d10, &
00763            1.2177698136199594677580d10,-2.0829040666802497120940d11/
00764       data q / 7.6886718750000000000000d01,-5.5648470543369082846819d03, &
00765            1.9418469440759880361415d05,-4.2648434812177161405483d06, &
00766            6.4698830956576428587653d07,-7.0108568774215954065376d08, &
00767            5.4229617984472955011862d09,-2.8986272696554495342658d10, &
00768            9.8900934262481749439886d10,-8.9673749185755048616855d10/
00769       !
00770       !  J-fraction coefficients for 6.0 <= X < 12.0
00771       !
00772       data r/-2.645677793077147237806d00,-2.378372882815725244124d00, &
00773            -2.421106956980653511550d01, 1.052976392459015155422d01, &
00774            1.945603779539281810439d01,-3.015761863840593359165d01, &
00775            1.120011024227297451523d01,-3.988850730390541057912d00, &
00776            9.565134591978630774217d00, 9.981193787537396413219d-1/
00777       data s/ 1.598517957704779356479d-4, 4.644185932583286942650d00, &
00778            3.697412299772985940785d02,-8.791401054875438925029d00, &
00779            7.608194509086645763123d02, 2.852397548119248700147d01, &
00780            4.731097187816050252967d02,-2.369210235636181001661d02, &
00781            1.249884822712447891440d00/
00782       !
00783       !  J-fraction coefficients for 12.0 <= X < 24.0
00784       !
00785       data p1/-1.647721172463463140042d00,-1.860092121726437582253d01, &
00786            -1.000641913989284829961d01,-2.105740799548040450394d01, &
00787            -9.134835699998742552432d-1,-3.323612579343962284333d01, &
00788            2.495487730402059440626d01, 2.652575818452799819855d01, &
00789            -1.845086232391278674524d00, 9.999933106160568739091d-1/
00790       data q1/ 9.792403599217290296840d01, 6.403800405352415551324d01, &
00791            5.994932325667407355255d01, 2.538819315630708031713d02, &
00792            4.429413178337928401161d01, 1.192832423968601006985d03, &
00793            1.991004470817742470726d02,-1.093556195391091143924d01, &
00794            1.001533852045342697818d00/
00795       !
00796       !  J-fraction coefficients for  24 <= X.
00797       !
00798       data p2/ 1.75338801265465972390d02,-2.23127670777632409550d02, &
00799            -1.81949664929868906455d01,-2.79798528624305389340d01, &
00800            -7.63147701620253630855d00,-1.52856623636929636839d01, &
00801            -7.06810977895029358836d00,-5.00006640413131002475d00, &
00802            -3.00000000320981265753d00, 1.00000000000000485503d00/
00803       data q2/ 3.97845977167414720840d04, 3.97277109100414518365d00, &
00804            1.37790390235747998793d02, 1.17179220502086455287d02, &
00805            7.04831847180424675988d01,-1.20187763547154743238d01, &
00806            -7.99243595776339741065d00,-2.99999894040324959612d00, &
00807            1.99999999999048104167d00/
00808 
00809       x = arg
00810 
00811       if ( x == zero ) then
00812 
00813          ei = - xinf
00814 
00815          if ( jint == 2 ) then
00816             ei = -ei
00817          end if
00818          !
00819          !  Calculate EI for negative argument or for E1.
00820          !
00821       else if ( x < zero .or. jint == 2 ) then
00822 
00823          y = abs ( x )
00824 
00825          if ( y <= one ) then
00826 
00827             sump = a(7) * y + a(1)
00828             sumq = y + b(1)
00829             do i = 2, 6
00830                sump = sump * y + a(i)
00831                sumq = sumq * y + b(i)
00832             end do
00833             ei = log ( y ) - sump / sumq
00834 
00835             if ( jint == 3 ) then
00836                ei = ei * exp ( y )
00837             end if
00838 
00839          else if ( y <= four ) then
00840 
00841             w = one / y
00842             sump = c(1)
00843             sumq = d(1)
00844             do i = 2, 9
00845                sump = sump * w + c(i)
00846                sumq = sumq * w + d(i)
00847             end do
00848 
00849             ei = - sump / sumq
00850 
00851             if ( jint /= 3 ) then
00852                ei = ei * exp ( -y )
00853             end if
00854 
00855          else
00856 
00857             if ( xbig < y .and. jint < 3 ) then
00858 
00859                ei = zero
00860 
00861             else
00862 
00863                w = one / y
00864                sump = e(1)
00865                sumq = f(1)
00866                do i = 2, 10
00867                   sump = sump * w + e(i)
00868                   sumq = sumq * w + f(i)
00869                end do
00870 
00871                ei = -w * ( one - w * sump / sumq )
00872 
00873                if ( jint /= 3 ) then
00874                   ei = ei * exp ( - y )
00875                end if
00876 
00877             end if
00878 
00879          end if
00880 
00881          if ( jint == 2 ) then
00882             ei = -ei
00883          end if
00884          !
00885          !  To improve conditioning, rational approximations are expressed
00886          !  in terms of Chebyshev polynomials for 0 <= X < 6, and in
00887          !  continued fraction form for larger X.
00888          !
00889       else if ( x < six ) then
00890 
00891          t = x + x
00892          t = t / three - two
00893          px(1) = zero
00894          qx(1) = zero
00895          px(2) = p(1)
00896          qx(2) = q(1)
00897          do i = 2, 9
00898             px(i+1) = t * px(i) - px(i-1) + p(i)
00899             qx(i+1) = t * qx(i) - qx(i-1) + q(i)
00900          end do
00901          sump = half * t * px(10) - px(9) + p(10)
00902          sumq = half * t * qx(10) - qx(9) + q(10)
00903          frac = sump / sumq
00904          xmx0 = ( x - x01 / x11 ) - x02
00905 
00906          if ( p037 <= abs ( xmx0 ) ) then
00907 
00908             ei = log ( x / x0 ) + xmx0 * frac
00909 
00910             if ( jint == 3 ) then
00911                ei = exp ( - x ) * ei
00912             end if
00913             !
00914             !  Special approximation to ln(X/X0) for X close to X0.
00915             !
00916          else
00917 
00918             y = xmx0 / ( x + x0 )
00919             ysq = y * y
00920             sump = plg(1)
00921             sumq = ysq + qlg(1)
00922             do i = 2, 4
00923                sump = sump * ysq + plg(i)
00924                sumq = sumq * ysq + qlg(i)
00925             end do
00926             ei = ( sump / ( sumq * ( x + x0 ) ) + frac ) * xmx0
00927 
00928             if ( jint == 3 ) then
00929                ei = exp ( - x ) * ei
00930             end if
00931 
00932          end if
00933 
00934       else if ( x < twelve ) then
00935 
00936          frac = zero
00937          do i = 1, 9
00938             frac = s(i) / ( r(i) + x + frac )
00939          end do
00940 
00941          ei = ( r(10) + frac ) / x
00942          if ( jint /= 3 ) then
00943             ei = ei * exp ( x )
00944          end if
00945 
00946       else if ( x <= two4 ) then
00947 
00948          frac = zero
00949          do i = 1, 9
00950             frac = q1(i) / ( p1(i) + x + frac )
00951          end do
00952 
00953          ei = ( p1(10) + frac ) / x
00954 
00955          if ( jint /= 3 ) then
00956             ei = ei * exp ( x )
00957          end if
00958 
00959       else
00960 
00961          if ( xmax <= x .and. jint < 3 ) then
00962 
00963             ei = xinf
00964 
00965          else
00966 
00967             y = one / x
00968             frac = zero
00969             do i = 1, 9
00970                frac = q2(i) / ( p2(i) + x + frac )
00971             end do
00972             frac = p2(10) + frac
00973             ei = y + y * y * frac
00974 
00975             if ( jint /= 3 ) then
00976 
00977                if ( x <= xmax - two4 ) then
00978                   ei = ei * exp ( x )
00979                   !
00980                   !  Calculation reformulated to avoid premature overflow.
00981                   !
00982                else
00983                   ei = ( ei * exp ( x - fourty ) ) * exp40
00984                end if
00985 
00986             end if
00987          end if
00988       end if
00989 
00990       result = ei
00991 
00992       return
00993    end subroutine calcei
00994    subroutine calci0 ( arg, result, jint )
00995 
00996       !*****************************************************************************80
00997       !
00998       !! CALCI0 computes various I0 Bessel functions.
00999       !
01000       !  Discussion:
01001       !
01002       !    This routine computes modified Bessel functions of the first kind
01003       !    and order zero, I0(X) and EXP(-ABS(X))*I0(X), for real
01004       !    arguments X.
01005       !
01006       !    The main computation evaluates slightly modified forms of
01007       !    minimax approximations generated by Blair and Edwards, Chalk
01008       !    River (Atomic Energy of Canada Limited) Report AECL-4928,
01009       !    October, 1974.
01010       !
01011       !  Licensing:
01012       !
01013       !    This code is distributed under the GNU LGPL license.
01014       !
01015       !  Modified:
01016       !
01017       !    03 April 2007
01018       !
01019       !  Author:
01020       !
01021       !    Original FORTRAN77 version by William Cody, Laura Stoltz.
01022       !    FORTRAN90 version by John Burkardt.
01023       !
01024       !  Parameters:
01025       !
01026       !    Input, real ( kind = 8 ) ARG, the argument.  If JINT = 1, then
01027       !    the argument must be less than XMAX.
01028       !
01029       !    Output, real ( kind = 8 ) RESULT, the value of the function,
01030       !    which depends on the input value of JINT:
01031       !    1, RESULT = I0(x);
01032       !    2, RESULT = exp(-x) * I0(x);
01033       !
01034       !    Input, integer ( kind = 4 ) JINT, chooses the function to be computed.
01035       !    1, I0(x);
01036       !    2, exp(-x) * I0(x);
01037       !
01038       implicit none
01039 
01040       real ( kind = 8 ) a
01041       real ( kind = 8 ) arg
01042       real ( kind = 8 ) b
01043       real ( kind = 8 ) exp40
01044       real ( kind = 8 ) forty
01045       integer ( kind = 4 ) i
01046       integer ( kind = 4 ) jint
01047       real ( kind = 8 ) one
01048       real ( kind = 8 ) one5
01049       real ( kind = 8 ) p(15)
01050       real ( kind = 8 ) pp(8)
01051       real ( kind = 8 ) q(5)
01052       real ( kind = 8 ) qq(7)
01053       real ( kind = 8 ) result
01054       real ( kind = 8 ) rec15
01055       real ( kind = 8 ) sump
01056       real ( kind = 8 ) sumq
01057       real ( kind = 8 ) two25
01058       real ( kind = 8 ) x
01059       real ( kind = 8 ) xinf
01060       real ( kind = 8 ) xmax
01061       real ( kind = 8 ) xsmall
01062       real ( kind = 8 ) xx
01063       !
01064       !  Mathematical constants
01065       !
01066       data one /1.0d0/
01067       data one5 /15.0d0/
01068       data exp40 /2.353852668370199854d17/
01069       data forty /40.0d0/
01070       data rec15 /6.6666666666666666666d-2/
01071       data two25 /225.0d0/
01072       !
01073       !  Machine-dependent constants
01074       !
01075       data xsmall /5.55d-17/
01076       data xinf /1.79d308/
01077       data xmax /713.986d0/
01078       !
01079       !  Coefficients for XSMALL <= ABS(ARG) < 15.0
01080       !
01081       data  p/-5.2487866627945699800d-18,-1.5982226675653184646d-14, &
01082            -2.6843448573468483278d-11,-3.0517226450451067446d-08, &
01083            -2.5172644670688975051d-05,-1.5453977791786851041d-02, &
01084            -7.0935347449210549190d+00,-2.4125195876041896775d+03, &
01085            -5.9545626019847898221d+05,-1.0313066708737980747d+08, &
01086            -1.1912746104985237192d+10,-8.4925101247114157499d+11, &
01087            -3.2940087627407749166d+13,-5.5050369673018427753d+14, &
01088            -2.2335582639474375249d+15/
01089       data  q/-3.7277560179962773046d+03, 6.5158506418655165707d+06, &
01090            -6.5626560740833869295d+09, 3.7604188704092954661d+12, &
01091            -9.7087946179594019126d+14/
01092       !
01093       !  Coefficients for 15.0 <= ABS(ARG)
01094       !
01095       data pp/-3.9843750000000000000d-01, 2.9205384596336793945d+00, &
01096            -2.4708469169133954315d+00, 4.7914889422856814203d-01, &
01097            -3.7384991926068969150d-03,-2.6801520353328635310d-03, &
01098            9.9168777670983678974d-05,-2.1877128189032726730d-06/
01099       data qq/-3.1446690275135491500d+01, 8.5539563258012929600d+01, &
01100            -6.0228002066743340583d+01, 1.3982595353892851542d+01, &
01101            -1.1151759188741312645d+00, 3.2547697594819615062d-02, &
01102            -5.5194330231005480228d-04/
01103 
01104       x = abs ( arg )
01105 
01106       if ( x < xsmall ) then
01107 
01108          result = one
01109          !
01110          !  XSMALL <= ABS(ARG) < 15.0.
01111          !
01112       else if ( x < one5 ) then
01113 
01114          xx = x * x
01115          sump = p(1)
01116          do i = 2, 15
01117             sump = sump * xx + p(i)
01118          end do
01119          xx = xx - two25
01120 
01121          sumq = (((( &
01122               xx + q(1) ) &
01123               * xx + q(2) ) &
01124               * xx + q(3) ) &
01125               * xx + q(4) ) &
01126               * xx + q(5)
01127 
01128          result = sump / sumq
01129 
01130          if ( jint == 2 ) then
01131             result = result * exp ( - x )
01132          end if
01133 
01134       else if ( one5 <= x ) then
01135 
01136          if ( jint == 1 .and. xmax < x ) then
01137             result = xinf
01138          else
01139             !
01140             !  15.0 <= ABS(ARG).
01141             !
01142             xx = one / x - rec15
01143 
01144             sump = (((((( &
01145                  pp(1) &
01146                  * xx + pp(2) ) &
01147                  * xx + pp(3) ) &
01148                  * xx + pp(4) ) &
01149                  * xx + pp(5) ) &
01150                  * xx + pp(6) ) &
01151                  * xx + pp(7) ) &
01152                  * xx + pp(8)
01153 
01154             sumq = (((((( &
01155                  xx + qq(1) ) &
01156                  * xx + qq(2) ) &
01157                  * xx + qq(3) ) &
01158                  * xx + qq(4) ) &
01159                  * xx + qq(5) ) &
01160                  * xx + qq(6) ) &
01161                  * xx + qq(7)
01162 
01163             result = sump / sumq
01164 
01165             if ( jint == 2 ) then
01166                result = ( result - pp(1) ) / sqrt ( x )
01167             else
01168                !
01169                !  Calculation reformulated to avoid premature overflow.
01170                !
01171                if ( x .le.( xmax - one5 ) ) then
01172                   a = exp ( x )
01173                   b = one
01174                else
01175                   a = exp ( x - forty )
01176                   b = exp40
01177                end if
01178 
01179                result = ( ( result * a - pp(1) * a ) / sqrt ( x ) ) * b
01180 
01181             end if
01182 
01183          end if
01184 
01185       end if
01186 
01187       return
01188    end subroutine calci0
01189    subroutine calci1 ( arg, result, jint )
01190 
01191       !*****************************************************************************80
01192       !
01193       !! CALCI1 computes various I1 Bessel functions.
01194       !
01195       !  Discussion:
01196       !
01197       !    This routine computes modified Bessel functioons of the first kind
01198       !    and order one, I1(X) and EXP(-ABS(X))*I1(X), for real
01199       !    arguments X.
01200       !
01201       !    The main computation evaluates slightly modified forms of
01202       !    minimax approximations generated by Blair and Edwards, Chalk
01203       !    River (Atomic Energy of Canada Limited) Report AECL-4928,
01204       !    October, 1974.
01205       !
01206       !  Licensing:
01207       !
01208       !    This code is distributed under the GNU LGPL license.
01209       !
01210       !  Modified:
01211       !
01212       !    03 April 2007
01213       !
01214       !  Author:
01215       !
01216       !    Original FORTRAN77 version by William Cody, Laura Stoltz.
01217       !    FORTRAN90 version by John Burkardt.
01218       !
01219       !  Parameters:
01220       !
01221       !    Input, real ( kind = 8 ) ARG, the argument.  If JINT = 1, then
01222       !    the argument must be less than XMAX.
01223       !
01224       !    Output, real ( kind = 8 ) RESULT, the value of the function,
01225       !    which depends on the input value of JINT:
01226       !    1, RESULT = I1(x);
01227       !    2, RESULT = exp(-x) * I1(x);
01228       !
01229       !    Input, integer ( kind = 4 ) JINT, chooses the function to be computed.
01230       !    1, I1(x);
01231       !    2, exp(-x) * I1(x);
01232       !
01233       implicit none
01234 
01235       real ( kind = 8 ) a
01236       real ( kind = 8 ) arg
01237       real ( kind = 8 ) b
01238       real ( kind = 8 ) exp40
01239       real ( kind = 8 ) forty
01240       real ( kind = 8 ) half
01241       integer ( kind = 4 ) j
01242       integer ( kind = 4 ) jint
01243       real ( kind = 8 ) one
01244       real ( kind = 8 ) one5
01245       real ( kind = 8 ) p(15)
01246       real ( kind = 8 ) pbar
01247       real ( kind = 8 ) pp(8)
01248       real ( kind = 8 ) q(5)
01249       real ( kind = 8 ) qq(6)
01250       real ( kind = 8 ) rec15
01251       real ( kind = 8 ) result
01252       real ( kind = 8 ) sump
01253       real ( kind = 8 ) sumq
01254       real ( kind = 8 ) two25
01255       real ( kind = 8 ) x
01256       real ( kind = 8 ) xinf
01257       real ( kind = 8 ) xmax
01258       real ( kind = 8 ) xsmall
01259       real ( kind = 8 ) xx
01260       real ( kind = 8 ) zero
01261       !
01262       !  Mathematical constants
01263       !
01264       data one /1.0d0/
01265       data one5 /15.0d0/
01266       data exp40 /2.353852668370199854d17/
01267       data forty /40.0d0/
01268       data rec15 /6.6666666666666666666d-2/
01269       data two25 /225.0d0/
01270       data half /0.5d0/
01271       data zero /0.0d0/
01272       !
01273       !  Machine-dependent constants
01274       !
01275       data xsmall /5.55d-17/
01276       data xinf /1.79d308/
01277       data xmax /713.987d0/
01278       !
01279       !  Coefficients for XSMALL <= ABS(ARG) < 15.0
01280       !
01281       data p/-1.9705291802535139930d-19,-6.5245515583151902910d-16, &
01282            -1.1928788903603238754d-12,-1.4831904935994647675d-09, &
01283            -1.3466829827635152875d-06,-9.1746443287817501309d-04, &
01284            -4.7207090827310162436d-01,-1.8225946631657315931d+02, &
01285            -5.1894091982308017540d+04,-1.0588550724769347106d+07, &
01286            -1.4828267606612366099d+09,-1.3357437682275493024d+11, &
01287            -6.9876779648010090070d+12,-1.7732037840791591320d+14, &
01288            -1.4577180278143463643d+15/
01289       data q/-4.0076864679904189921d+03, 7.4810580356655069138d+06, &
01290            -8.0059518998619764991d+09, 4.8544714258273622913d+12, &
01291            -1.3218168307321442305d+15/
01292       !
01293       !  Coefficients for 15.0 <= ABS(ARG)
01294       !
01295       data pp/-6.0437159056137600000d-02, 4.5748122901933459000d-01, &
01296            -4.2843766903304806403d-01, 9.7356000150886612134d-02, &
01297            -3.2457723974465568321d-03,-3.6395264712121795296d-04, &
01298            1.6258661867440836395d-05,-3.6347578404608223492d-07/
01299       data qq/-3.8806586721556593450d+00, 3.2593714889036996297d+00, &
01300            -8.5017476463217924408d-01, 7.4212010813186530069d-02, &
01301            -2.2835624489492512649d-03, 3.7510433111922824643d-05/
01302       data pbar/3.98437500d-01/
01303 
01304       x = abs ( arg )
01305       !
01306       !  Return for ABS(ARG) < XSMALL.
01307       !
01308       if ( x < xsmall ) then
01309 
01310          result = half * x
01311          !
01312          !  XSMALL <= ABS(ARG) < 15.0.
01313          !
01314       else if ( x < one5 ) then
01315 
01316          xx = x * x
01317          sump = p(1)
01318          do j = 2, 15
01319             sump = sump * xx + p(j)
01320          end do
01321          xx = xx - two25
01322 
01323          sumq = (((( &
01324               xx + q(1) ) &
01325               * xx + q(2) ) &
01326               * xx + q(3) ) &
01327               * xx + q(4) ) &
01328               * xx + q(5)
01329 
01330          result = ( sump / sumq ) * x
01331 
01332          if ( jint == 2 ) then
01333             result = result * exp ( -x )
01334          end if
01335 
01336       else if ( jint == 1 .and. xmax < x ) then
01337 
01338          result = xinf
01339 
01340       else
01341          !
01342          !  15.0 <= ABS(ARG).
01343          !
01344          xx = one / x - rec15
01345 
01346          sump = (((((( &
01347               pp(1) &
01348               * xx + pp(2) ) &
01349               * xx + pp(3) ) &
01350               * xx + pp(4) ) &
01351               * xx + pp(5) ) &
01352               * xx + pp(6) ) &
01353               * xx + pp(7) ) &
01354               * xx + pp(8)
01355 
01356          sumq = ((((( &
01357               xx + qq(1) ) &
01358               * xx + qq(2) ) &
01359               * xx + qq(3) ) &
01360               * xx + qq(4) ) &
01361               * xx + qq(5) ) &
01362               * xx + qq(6)
01363 
01364          result = sump / sumq
01365 
01366          if ( jint /= 1 ) then
01367             result = ( result + pbar ) / sqrt ( x )
01368          else
01369             !
01370             !  Calculation reformulated to avoid premature overflow.
01371             !
01372             if ( xmax - one5 < x ) then
01373                a = exp ( x - forty )
01374                b = exp40
01375             else
01376                a = exp ( x )
01377                b = one
01378             end if
01379 
01380             result = ( ( result * a + pbar * a ) / sqrt ( x ) ) * b
01381 
01382          end if
01383       end if
01384 
01385       if ( arg < zero ) then
01386          result = -result
01387       end if
01388 
01389       return
01390    end subroutine calci1
01391    subroutine calck0 ( arg, result, jint )
01392 
01393       !*****************************************************************************80
01394       !
01395       !! CALCK0 computes various K0 Bessel functions.
01396       !
01397       !  Discussion:
01398       !
01399       !    This routine computes modified Bessel functions of the second kind
01400       !    and order zero, K0(X) and EXP(X)*K0(X), for real
01401       !    arguments X.
01402       !
01403       !    The main computation evaluates slightly modified forms of near
01404       !    minimax rational approximations generated by Russon and Blair,
01405       !    Chalk River (Atomic Energy of Canada Limited) Report AECL-3461,
01406       !    1969.
01407       !
01408       !  Licensing:
01409       !
01410       !    This code is distributed under the GNU LGPL license.
01411       !
01412       !  Modified:
01413       !
01414       !    03 April 2007
01415       !
01416       !  Author:
01417       !
01418       !    Original FORTRAN77 version by William Cody, Laura Stoltz.
01419       !    FORTRAN90 version by John Burkardt.
01420       !
01421       !  Parameters:
01422       !
01423       !    Input, real ( kind = 8 ) ARG, the argument.  0 < ARG is
01424       !    always required.  If JINT = 1, then the argument must also be
01425       !    less than XMAX.
01426       !
01427       !    Output, real ( kind = 8 ) RESULT, the value of the function,
01428       !    which depends on the input value of JINT:
01429       !    1, RESULT = K0(x);
01430       !    2, RESULT = exp(x) * K0(x);
01431       !
01432       !    Input, integer ( kind = 4 ) JINT, chooses the function to be computed.
01433       !    1, K0(x);
01434       !    2, exp(x) * K0(x);
01435       !
01436       implicit none
01437 
01438       integer ( kind = 4 ) i
01439       integer ( kind = 4 ) jint
01440       real ( kind = 8 ) arg
01441       real ( kind = 8 ) f(4)
01442       real ( kind = 8 ) g(3)
01443       real ( kind = 8 ) one
01444       real ( kind = 8 ) p(6)
01445       real ( kind = 8 ) pp(10)
01446       real ( kind = 8 ) q(2)
01447       real ( kind = 8 ) qq(10)
01448       real ( kind = 8 ) result
01449       real ( kind = 8 ) sumf
01450       real ( kind = 8 ) sumg
01451       real ( kind = 8 ) sump
01452       real ( kind = 8 ) sumq
01453       real ( kind = 8 ) temp
01454       real ( kind = 8 ) x
01455       real ( kind = 8 ) xinf
01456       real ( kind = 8 ) xmax
01457       real ( kind = 8 ) xsmall
01458       real ( kind = 8 ) xx
01459       real ( kind = 8 ) zero
01460       !
01461       !  Mathematical constants
01462       !
01463       data one /1.0d0/
01464       data zero /0.0d0/
01465       !
01466       !  Machine-dependent constants
01467       !
01468       data xsmall /1.11d-16/
01469       data xinf /1.79d+308/
01470       data xmax /705.342d0/
01471       !
01472       !  Coefficients for XSMALL <= ARG <= 1.0
01473       !
01474       data   p/ 5.8599221412826100000d-04, 1.3166052564989571850d-01, &
01475            1.1999463724910714109d+01, 4.6850901201934832188d+02, &
01476            5.9169059852270512312d+03, 2.4708152720399552679d+03/
01477       data   q/-2.4994418972832303646d+02, 2.1312714303849120380d+04/
01478       data   f/-1.6414452837299064100d+00,-2.9601657892958843866d+02, &
01479            -1.7733784684952985886d+04,-4.0320340761145482298d+05/
01480       data   g/-2.5064972445877992730d+02, 2.9865713163054025489d+04, &
01481            -1.6128136304458193998d+06/
01482       !
01483       !  Coefficients for  1.0 < ARG
01484       !
01485       data  pp/ 1.1394980557384778174d+02, 3.6832589957340267940d+03, &
01486            3.1075408980684392399d+04, 1.0577068948034021957d+05, &
01487            1.7398867902565686251d+05, 1.5097646353289914539d+05, &
01488            7.1557062783764037541d+04, 1.8321525870183537725d+04, &
01489            2.3444738764199315021d+03, 1.1600249425076035558d+02/
01490       data  qq/ 2.0013443064949242491d+02, 4.4329628889746408858d+03, &
01491            3.1474655750295278825d+04, 9.7418829762268075784d+04, &
01492            1.5144644673520157801d+05, 1.2689839587977598727d+05, &
01493            5.8824616785857027752d+04, 1.4847228371802360957d+04, &
01494            1.8821890840982713696d+03, 9.2556599177304839811d+01/
01495 
01496       x = arg
01497       !
01498       !  0.0 < ARG <= 1.0.
01499       !
01500       if ( zero < x ) then
01501 
01502          if ( x <= one ) then
01503 
01504             temp = log ( x )
01505 
01506             if ( x < xsmall ) then
01507                !
01508                !  Return for small ARG.
01509                !
01510                result = p(6) / q(2) - temp
01511 
01512             else
01513 
01514                xx = x * x
01515 
01516                sump = (((( &
01517                     p(1) &
01518                     * xx + p(2) ) &
01519                     * xx + p(3) ) &
01520                     * xx + p(4) ) &
01521                     * xx + p(5) ) &
01522                     * xx + p(6)
01523 
01524                sumq = ( xx + q(1) ) * xx + q(2)
01525                sumf = ( ( &
01526                     f(1) &
01527                     * xx + f(2) ) &
01528                     * xx + f(3) ) &
01529                     * xx + f(4)
01530 
01531                sumg = ( ( xx + g(1) ) * xx + g(2) ) * xx + g(3)
01532 
01533                result = sump / sumq - xx * sumf * temp / sumg - temp
01534 
01535                if ( jint == 2 ) then
01536                   result = result * exp ( x )
01537                end if
01538 
01539             end if
01540 
01541          else if ( jint == 1 .and. xmax < x ) then
01542             !
01543             !  Error return for XMAX < ARG.
01544             !
01545             result = zero
01546 
01547          else
01548             !
01549             !  1.0 < ARG.
01550             !
01551             xx = one / x
01552             sump = pp(1)
01553             do i = 2, 10
01554                sump = sump * xx + pp(i)
01555             end do
01556 
01557             sumq = xx
01558             do i = 1, 9
01559                sumq = ( sumq + qq(i) ) * xx
01560             end do
01561             sumq = sumq + qq(10)
01562             result = sump / sumq / sqrt ( x )
01563 
01564             if ( jint == 1 ) then
01565                result = result * exp ( -x )
01566             end if
01567 
01568          end if
01569 
01570       else
01571          !
01572          !  Error return for ARG <= 0.0.
01573          !
01574          result = xinf
01575 
01576       end if
01577 
01578       return
01579    end subroutine calck0
01580    subroutine calck1 ( arg, result, jint )
01581 
01582       !*****************************************************************************80
01583       !
01584       !! CALCK1 computes various K1 Bessel functions.
01585       !
01586       !  Discussion:
01587       !
01588       !    This routine computes modified Bessel functions of the second kind
01589       !    and order one, K1(X) and EXP(X)*K1(X), for real arguments X.
01590       !
01591       !    The main computation evaluates slightly modified forms of near
01592       !    minimax rational approximations generated by Russon and Blair,
01593       !    Chalk River (Atomic Energy of Canada Limited) Report AECL-3461,
01594       !    1969.
01595       !
01596       !  Licensing:
01597       !
01598       !    This code is distributed under the GNU LGPL license.
01599       !
01600       !  Modified:
01601       !
01602       !    03 April 2007
01603       !
01604       !  Author:
01605       !
01606       !    Original FORTRAN77 version by William Cody, Laura Stoltz.
01607       !    FORTRAN90 version by John Burkardt.
01608       !
01609       !  Parameters:
01610       !
01611       !    Input, real ( kind = 8 ) ARG, the argument.  XLEAST < ARG is
01612       !    always required.  If JINT = 1, then the argument must also be
01613       !    less than XMAX.
01614       !
01615       !    Output, real ( kind = 8 ) RESULT, the value of the function,
01616       !    which depends on the input value of JINT:
01617       !    1, RESULT = K1(x);
01618       !    2, RESULT = exp(x) * K1(x);
01619       !
01620       !    Input, integer ( kind = 4 ) JINT, chooses the function to be computed.
01621       !    1, K1(x);
01622       !    2, exp(x) * K1(x);
01623       !
01624       implicit none
01625 
01626       real ( kind = 8 ) arg
01627       real ( kind = 8 ) f(5)
01628       real ( kind = 8 ) g(3)
01629       integer ( kind = 4 ) i
01630       integer ( kind = 4 ) jint
01631       real ( kind = 8 ) one
01632       real ( kind = 8 ) p(5)
01633       real ( kind = 8 ) pp(11)
01634       real ( kind = 8 ) q(3)
01635       real ( kind = 8 ) qq(9)
01636       real ( kind = 8 ) result
01637       real ( kind = 8 ) sumf
01638       real ( kind = 8 ) sumg
01639       real ( kind = 8 ) sump
01640       real ( kind = 8 ) sumq
01641       real ( kind = 8 ) x
01642       real ( kind = 8 ) xinf
01643       real ( kind = 8 ) xmax
01644       real ( kind = 8 ) xleast
01645       real ( kind = 8 ) xsmall
01646       real ( kind = 8 ) xx
01647       real ( kind = 8 ) zero
01648       !
01649       !  Mathematical constants
01650       !
01651       data one /1.0d0/
01652       data zero /0.0d0/
01653       !
01654       !  Machine-dependent constants
01655       !
01656       data xleast /2.23d-308/
01657       data xsmall /1.11d-16/
01658       data xinf /1.79d+308/
01659       data xmax /705.343d+0/
01660       !
01661       !  Coefficients for  XLEAST <=  ARG  <= 1.0
01662       !
01663       data   p/ 4.8127070456878442310d-1, 9.9991373567429309922d+1, &
01664            7.1885382604084798576d+3, 1.7733324035147015630d+5, &
01665            7.1938920065420586101d+5/
01666       data   q/-2.8143915754538725829d+2, 3.7264298672067697862d+4, &
01667            -2.2149374878243304548d+6/
01668       data   f/-2.2795590826955002390d-1,-5.3103913335180275253d+1, &
01669            -4.5051623763436087023d+3,-1.4758069205414222471d+5, &
01670            -1.3531161492785421328d+6/
01671       data   g/-3.0507151578787595807d+2, 4.3117653211351080007d+4, &
01672            -2.7062322985570842656d+6/
01673       !
01674       !  Coefficients for  1.0 < ARG
01675       !
01676       data  pp/ 6.4257745859173138767d-2, 7.5584584631176030810d+0, &
01677            1.3182609918569941308d+2, 8.1094256146537402173d+2, &
01678            2.3123742209168871550d+3, 3.4540675585544584407d+3, &
01679            2.8590657697910288226d+3, 1.3319486433183221990d+3, &
01680            3.4122953486801312910d+2, 4.4137176114230414036d+1, &
01681            2.2196792496874548962d+0/
01682       data  qq/ 3.6001069306861518855d+1, 3.3031020088765390854d+2, &
01683            1.2082692316002348638d+3, 2.1181000487171943810d+3, &
01684            1.9448440788918006154d+3, 9.6929165726802648634d+2, &
01685            2.5951223655579051357d+2, 3.4552228452758912848d+1, &
01686            1.7710478032601086579d+0/
01687 
01688       x = arg
01689       !
01690       !  Error return for ARG < XLEAST.
01691       !
01692       if ( x < xleast ) then
01693 
01694          result = xinf
01695          !
01696          !  XLEAST <= ARG <= 1.0.
01697          !
01698       else if ( x <= one ) then
01699 
01700          if ( x < xsmall ) then
01701             !
01702             !  Return for small ARG.
01703             !
01704             result = one / x
01705 
01706          else
01707 
01708             xx = x * x
01709 
01710             sump = (((( &
01711                  p(1) &
01712                  * xx + p(2) ) &
01713                  * xx + p(3) ) &
01714                  * xx + p(4) ) &
01715                  * xx + p(5) ) &
01716                  * xx + q(3)
01717 
01718             sumq = (( &
01719                  xx + q(1) ) &
01720                  * xx + q(2) ) &
01721                  * xx + q(3)
01722 
01723             sumf = ((( &
01724                  f(1) &
01725                  * xx + f(2) ) &
01726                  * xx + f(3) ) &
01727                  * xx + f(4) ) &
01728                  * xx + f(5)
01729 
01730             sumg = (( &
01731                  xx + g(1) ) &
01732                  * xx + g(2) ) &
01733                  * xx + g(3)
01734 
01735             result = ( xx * log ( x ) * sumf / sumg + sump / sumq ) / x
01736 
01737             if ( jint == 2 ) then
01738                result = result * exp ( x )
01739             end if
01740 
01741          end if
01742 
01743       else if ( jint == 1 .and. xmax < x ) then
01744          !
01745          !  Error return for XMAX < ARG.
01746          !
01747          result = zero
01748 
01749       else
01750          !
01751          !  1.0 < ARG.
01752          !
01753          xx = one / x
01754 
01755          sump = pp(1)
01756          do i = 2, 11
01757             sump = sump * xx + pp(i)
01758          end do
01759 
01760          sumq = xx
01761          do i = 1, 8
01762             sumq = ( sumq + qq(i) ) * xx
01763          end do
01764          sumq = sumq + qq(9)
01765 
01766          result = sump / sumq / sqrt ( x )
01767 
01768          if ( jint == 1 ) then
01769             result = result * exp ( -x )
01770          end if
01771 
01772       end if
01773 
01774       return
01775    end subroutine calck1
01776    subroutine calerf ( arg, result, jint )
01777 
01778       !*****************************************************************************80
01779       !
01780       !! CALERF computes various forms of the error function.
01781       !
01782       !  Discussion:
01783       !
01784       !    This routine evaluates erf(x), erfc(x), and exp(x*x)*erfc(x)
01785       !    for a real argument x.
01786       !
01787       !  Licensing:
01788       !
01789       !    This code is distributed under the GNU LGPL license.
01790       !
01791       !  Modified:
01792       !
01793       !    03 April 2007
01794       !
01795       !  Author:
01796       !
01797       !    Original FORTRAN77 version by William Cody.
01798       !    FORTRAN90 version by John Burkardt.
01799       !
01800       !  Reference:
01801       !
01802       !    William Cody,
01803       !    Rational Chebyshev Approximations for the Error Function,
01804       !    Mathematics of Computation,
01805       !    Volume 23, Number 107, July 1969, pages 631-638.
01806       !
01807       !  Parameters:
01808       !
01809       !    Input, real ( kind = 8 ) ARG, the argument.  If JINT is 1, the
01810       !    argument must be less than XBIG.  If JINT is 2, the argument
01811       !    must lie between XNEG and XMAX.
01812       !
01813       !    Output, real ( kind = 8 ) RESULT, the value of the function,
01814       !    which depends on the input value of JINT:
01815       !    0, RESULT = erf(x);
01816       !    1, RESULT = erfc(x) = 1 - erf(x);
01817       !    2, RESULT = exp(x*x)*erfc(x) = exp(x*x) - erf(x*x)*erf(x).
01818       !
01819       !    Input, integer ( kind = 4 ) JINT, chooses the function to be computed.
01820       !    0, erf(x);
01821       !    1, erfc(x);
01822       !    2, exp(x*x)*erfc(x).
01823       !
01824       implicit none
01825 
01826       real ( kind = 8 ) a(5)
01827       real ( kind = 8 ) arg
01828       real ( kind = 8 ) b(4)
01829       real ( kind = 8 ) c(9)
01830       real ( kind = 8 ) d(8)
01831       real ( kind = 8 ) del
01832       real ( kind = 8 ) four
01833       real ( kind = 8 ) half
01834       integer ( kind = 4 ) i
01835       integer ( kind = 4 ) jint
01836       real ( kind = 8 ) one
01837       real ( kind = 8 ) p(6)
01838       real ( kind = 8 ) q(5)
01839       real ( kind = 8 ) result
01840       real ( kind = 8 ) sixten
01841       real ( kind = 8 ) sqrpi
01842       real ( kind = 8 ) two
01843       real ( kind = 8 ) thresh
01844       real ( kind = 8 ) x
01845       real ( kind = 8 ) xbig
01846       real ( kind = 8 ) xden
01847       real ( kind = 8 ) xhuge
01848       real ( kind = 8 ) xinf
01849       real ( kind = 8 ) xmax
01850       real ( kind = 8 ) xneg
01851       real ( kind = 8 ) xnum
01852       real ( kind = 8 ) xsmall
01853       real ( kind = 8 ) y
01854       real ( kind = 8 ) ysq
01855       real ( kind = 8 ) zero
01856       !
01857       !  Mathematical constants
01858       !
01859       data four /4.0d0/
01860       data one /1.0d0/
01861       data half /0.5d0/
01862       data two /2.0d0/
01863       data zero /0.0d0/
01864       data sqrpi /5.6418958354775628695d-1/
01865       data thresh /0.46875d0/
01866       data sixten /16.0d0/
01867       !
01868       !  Machine-dependent constants
01869       !
01870       data xinf /1.79d308 /
01871       data xneg / -26.628d0 /
01872       data xsmall /1.11d-16/
01873       data xbig /26.543d0 /
01874       data xhuge /6.71d7/
01875       data xmax /2.53d307/
01876       !
01877       !  Coefficients for approximation to  erf  in first interval
01878       !
01879       data a/3.16112374387056560d00,1.13864154151050156d02, &
01880            3.77485237685302021d02,3.20937758913846947d03, &
01881            1.85777706184603153d-1/
01882       data b/2.36012909523441209d01,2.44024637934444173d02, &
01883            1.28261652607737228d03,2.84423683343917062d03/
01884       !
01885       !  Coefficients for approximation to  erfc  in second interval
01886       !
01887       data c/5.64188496988670089d-1,8.88314979438837594d0, &
01888            6.61191906371416295d01,2.98635138197400131d02, &
01889            8.81952221241769090d02,1.71204761263407058d03, &
01890            2.05107837782607147d03,1.23033935479799725d03, &
01891            2.15311535474403846d-8/
01892       data d/1.57449261107098347d01,1.17693950891312499d02, &
01893            5.37181101862009858d02,1.62138957456669019d03, &
01894            3.29079923573345963d03,4.36261909014324716d03, &
01895            3.43936767414372164d03,1.23033935480374942d03/
01896       !
01897       !  Coefficients for approximation to  erfc  in third interval
01898       !
01899       data p/3.05326634961232344d-1,3.60344899949804439d-1, &
01900            1.25781726111229246d-1,1.60837851487422766d-2, &
01901            6.58749161529837803d-4,1.63153871373020978d-2/
01902       data q/2.56852019228982242d00,1.87295284992346047d00, &
01903            5.27905102951428412d-1,6.05183413124413191d-2, &
01904            2.33520497626869185d-3/
01905 
01906       x = arg
01907       y = abs ( x )
01908       !
01909       !  Evaluate erf for |X| <= 0.46875.
01910       !
01911       if ( y <= thresh ) then
01912 
01913          ysq = zero
01914          if ( xsmall < y ) then
01915             ysq = y * y
01916          end if
01917 
01918          xnum = a(5) * ysq
01919          xden = ysq
01920 
01921          do i = 1, 3
01922             xnum = ( xnum + a(i) ) * ysq
01923             xden = ( xden + b(i) ) * ysq
01924          end do
01925 
01926          result = x * ( xnum + a(4) ) / ( xden + b(4) )
01927 
01928          if ( jint /= 0 ) then
01929             result = one - result
01930          end if
01931 
01932          if ( jint == 2 ) then
01933             result = exp ( ysq ) * result
01934          end if
01935 
01936          return
01937          !
01938          !  Evaluate erfc for 0.46875 <= |X| <= 4.0.
01939          !
01940       else if ( y <= four ) then
01941 
01942          xnum = c(9) * y
01943          xden = y
01944 
01945          do i = 1, 7
01946             xnum = ( xnum + c(i) ) * y
01947             xden = ( xden + d(i) ) * y
01948          end do
01949 
01950          result = ( xnum + c(8) ) / ( xden + d(8) )
01951 
01952          if ( jint /= 2 ) then
01953             ysq = aint ( y * sixten ) / sixten
01954             del = ( y - ysq ) * ( y + ysq )
01955             result = exp ( -ysq * ysq ) * exp ( -del ) * result
01956          end if
01957          !
01958          !  Evaluate erfc for 4.0 < |X|.
01959          !
01960       else
01961 
01962          result = zero
01963 
01964          if ( xbig <= y ) then
01965 
01966             if ( jint /= 2 .or. xmax <= y ) then
01967                go to 300
01968             end if
01969 
01970             if ( xhuge <= y ) then
01971                result = sqrpi / y
01972                go to 300
01973             end if
01974 
01975          end if
01976 
01977          ysq = one / ( y * y )
01978          xnum = p(6) * ysq
01979          xden = ysq
01980          do i = 1, 4
01981             xnum = ( xnum + p(i) ) * ysq
01982             xden = ( xden + q(i) ) * ysq
01983          end do
01984 
01985          result = ysq * ( xnum + p(5) ) / ( xden + q(5) )
01986          result = ( sqrpi -  result ) / y
01987 
01988          if ( jint /= 2 ) then
01989             ysq = aint ( y * sixten ) / sixten
01990             del = ( y - ysq ) * ( y + ysq )
01991             result = exp ( -ysq * ysq ) * exp ( -del ) * result
01992          end if
01993 
01994       end if
01995       !
01996       !  Fix up for negative argument, erf, etc.
01997       !
01998 300   continue
01999 
02000       if ( jint == 0 ) then
02001 
02002          result = ( half - result ) + half
02003          if ( x < zero ) then
02004             result = -result
02005          end if
02006 
02007       else if ( jint == 1 ) then
02008 
02009          if ( x < zero ) then
02010             result = two - result
02011          end if
02012 
02013       else
02014 
02015          if ( x < zero ) then
02016 
02017             if ( x < xneg ) then
02018                result = xinf
02019             else
02020                ysq = aint ( x * sixten ) / sixten
02021                del = ( x - ysq ) * ( x + ysq )
02022                y = exp ( ysq * ysq ) * exp ( del )
02023                result = ( y + y ) - result
02024             end if
02025 
02026          end if
02027 
02028       end if
02029 
02030       return
02031    end subroutine calerf
02032    subroutine caljy0 ( arg, result, jint )
02033 
02034       !*****************************************************************************80
02035       !
02036       !! CALJY0 computes various J0 and Y0 Bessel functions.
02037       !
02038       !  Discussion:
02039       !
02040       !    This routine computes zero-order Bessel functions of the first and
02041       !    second kind (J0 and Y0), for real arguments X, where 0 < X <= XMAX
02042       !    for Y0, and |X| <= XMAX for J0.
02043       !
02044       !  Licensing:
02045       !
02046       !    This code is distributed under the GNU LGPL license.
02047       !
02048       !  Modified:
02049       !
02050       !    03 April 2007
02051       !
02052       !  Author:
02053       !
02054       !    Original FORTRAN77 version by William Cody.
02055       !    FORTRAN90 version by John Burkardt.
02056       !
02057       !  Reference:
02058       !
02059       !    John Hart, Ward Cheney, Charles Lawson, Hans Maehly,
02060       !    Charles Mesztenyi, John Rice, Henry Thatcher,
02061       !    Christoph Witzgall,
02062       !    Computer Approximations,
02063       !    Wiley, 1968,
02064       !    LC: QA297.C64.
02065       !
02066       !  Parameters:
02067       !
02068       !    Input, real ( kind = 8 ) ARG, the argument.  If JINT = 0, ARG
02069       !    must satisfy
02070       !     -XMAX < ARG < XMAX;
02071       !    If JINT = 1, then ARG must satisfy
02072       !      0 < ARG < XMAX.
02073       !
02074       !    Output, real ( kind = 8 ) RESULT, the value of the function,
02075       !    which depends on the input value of JINT:
02076       !    0, RESULT = J0(x);
02077       !    1, RESULT = Y0(x);
02078       !
02079       !    Input, integer ( kind = 4 ) JINT, chooses the function to be computed.
02080       !    0, J0(x);
02081       !    1, Y0(x);
02082       !
02083       implicit none
02084 
02085       integer ( kind = 4 ) i
02086       integer ( kind = 4 ) jint
02087       real ( kind = 8 ) arg
02088       real ( kind = 8 ) ax
02089       real ( kind = 8 ) cons
02090       real ( kind = 8 ) down
02091       real ( kind = 8 ) eight
02092       real ( kind = 8 ) five5
02093       real ( kind = 8 ) four
02094       real ( kind = 8 ) one
02095       real ( kind = 8 ) oneov8
02096       real ( kind = 8 ) pi2
02097       real ( kind = 8 ) pj0(7)
02098       real ( kind = 8 ) pj1(8)
02099       real ( kind = 8 ) plg(4)
02100       real ( kind = 8 ) prod
02101       real ( kind = 8 ) py0(6)
02102       real ( kind = 8 ) py1(7)
02103       real ( kind = 8 ) py2(8)
02104       real ( kind = 8 ) p0(6)
02105       real ( kind = 8 ) p1(6)
02106       real ( kind = 8 ) p17
02107       real ( kind = 8 ) qj0(5)
02108       real ( kind = 8 ) qj1(7)
02109       real ( kind = 8 ) qlg(4)
02110       real ( kind = 8 ) qy0(5)
02111       real ( kind = 8 ) qy1(6)
02112       real ( kind = 8 ) qy2(7)
02113       real ( kind = 8 ) q0(5)
02114       real ( kind = 8 ) q1(5)
02115       real ( kind = 8 ) resj
02116       real ( kind = 8 ) result
02117       real ( kind = 8 ) r0
02118       real ( kind = 8 ) r1
02119       real ( kind = 8 ) sixty4
02120       real ( kind = 8 ) three
02121       real ( kind = 8 ) twopi
02122       real ( kind = 8 ) twopi1
02123       real ( kind = 8 ) twopi2
02124       real ( kind = 8 ) two56
02125       real ( kind = 8 ) up
02126       real ( kind = 8 ) w
02127       real ( kind = 8 ) wsq
02128       real ( kind = 8 ) xden
02129       real ( kind = 8 ) xinf
02130       real ( kind = 8 ) xmax
02131       real ( kind = 8 ) xnum
02132       real ( kind = 8 ) xsmall
02133       real ( kind = 8 ) xj0
02134       real ( kind = 8 ) xj1
02135       real ( kind = 8 ) xj01
02136       real ( kind = 8 ) xj02
02137       real ( kind = 8 ) xj11
02138       real ( kind = 8 ) xj12
02139       real ( kind = 8 ) xy
02140       real ( kind = 8 ) xy0
02141       real ( kind = 8 ) xy01
02142       real ( kind = 8 ) xy02
02143       real ( kind = 8 ) xy1
02144       real ( kind = 8 ) xy11
02145       real ( kind = 8 ) xy12
02146       real ( kind = 8 ) xy2
02147       real ( kind = 8 ) xy21
02148       real ( kind = 8 ) xy22
02149       real ( kind = 8 ) z
02150       real ( kind = 8 ) zero
02151       real ( kind = 8 ) zsq
02152       !
02153       !  Mathematical constants
02154       !  CONS = ln(.5) + Euler's gamma
02155       !
02156       data zero / 0.0d0 /
02157       data one /1.0d0 /
02158       data three /3.0d0 /
02159       data four /4.0d0 /
02160       data eight /8.0d0/
02161       data five5 / 5.5d0 /
02162       data sixty4 /64.0d0 /
02163       data oneov8 /0.125d0 /
02164       data p17 /1.716d-1/
02165       data two56 /256.0d0/
02166       data cons / -1.1593151565841244881d-1/
02167       data pi2 /6.3661977236758134308d-1/
02168       data twopi /6.2831853071795864769d0/
02169       data twopi1 /6.28125d0 /
02170       data twopi2 / 1.9353071795864769253d-3/
02171       !
02172       !  Machine-dependent constants
02173       !
02174       data xmax /1.07d+09/
02175       data xsmall /9.31d-10/
02176       data xinf /1.7d+38/
02177       !
02178       !  Zeroes of Bessel functions
02179       !
02180       data xj0 /2.4048255576957727686d+0/
02181       data xj1 /5.5200781102863106496d+0/
02182       data xy0 /8.9357696627916752158d-1/
02183       data xy1 /3.9576784193148578684d+0/
02184       data xy2 /7.0860510603017726976d+0/
02185       data xj01 / 616.0d+0/
02186       data xj02 /-1.4244423042272313784d-03/
02187       data xj11 /1413.0d+0/
02188       data xj12 / 5.4686028631064959660d-04/
02189       data xy01 / 228.0d+0/
02190       data xy02 / 2.9519662791675215849d-03/
02191       data xy11 /1013.0d+0/
02192       data xy12 / 6.4716931485786837568d-04/
02193       data xy21 /1814.0d+0/
02194       data xy22 / 1.1356030177269762362d-04/
02195       !
02196       !  Coefficients for rational approximation to ln(x/a)
02197       !
02198       data plg/-2.4562334077563243311d+01,2.3642701335621505212d+02, &
02199            -5.4989956895857911039d+02,3.5687548468071500413d+02/
02200       data qlg/-3.5553900764052419184d+01,1.9400230218539473193d+02, &
02201            -3.3442903192607538956d+02,1.7843774234035750207d+02/
02202       !
02203       !  Coefficients for rational approximation of
02204       !  J0(X) / (X**2 - XJ0**2),  XSMALL < |X| <= 4.0
02205       !
02206       data pj0/6.6302997904833794242d+06,-6.2140700423540120665d+08, &
02207            2.7282507878605942706d+10,-4.1298668500990866786d+11, &
02208            -1.2117036164593528341d-01, 1.0344222815443188943d+02, &
02209            -3.6629814655107086448d+04/
02210       data qj0/4.5612696224219938200d+05, 1.3985097372263433271d+08, &
02211            2.6328198300859648632d+10, 2.3883787996332290397d+12, &
02212            9.3614022392337710626d+02/
02213       !
02214       !  Coefficients for rational approximation of
02215       !  J0(X) / (X**2 - XJ1**2), 4.0 < |X| <= 8.0
02216       !
02217       data pj1/4.4176707025325087628d+03, 1.1725046279757103576d+04, &
02218            1.0341910641583726701d+04,-7.2879702464464618998d+03, &
02219            -1.2254078161378989535d+04,-1.8319397969392084011d+03, &
02220            4.8591703355916499363d+01, 7.4321196680624245801d+02/
02221       data qj1/3.3307310774649071172d+02,-2.9458766545509337327d+03, &
02222            1.8680990008359188352d+04,-8.4055062591169562211d+04, &
02223            2.4599102262586308984d+05,-3.5783478026152301072d+05, &
02224            -2.5258076240801555057d+01/
02225       !
02226       !  Coefficients for rational approximation of
02227       !  (Y0(X) - 2 LN(X/XY0) J0(X)) / (X**2 - XY0**2),
02228       !  XSMALL < |X| <= 3.0
02229       !
02230       data py0/1.0102532948020907590d+04,-2.1287548474401797963d+06, &
02231            2.0422274357376619816d+08,-8.3716255451260504098d+09, &
02232            1.0723538782003176831d+11,-1.8402381979244993524d+01/
02233       data qy0/6.6475986689240190091d+02, 2.3889393209447253406d+05, &
02234            5.5662956624278251596d+07, 8.1617187777290363573d+09, &
02235            5.8873865738997033405d+11/
02236       !
02237       !  Coefficients for rational approximation of
02238       !  (Y0(X) - 2 LN(X/XY1) J0(X)) / (X**2 - XY1**2),
02239       !  3.0 < |X| <= 5.5
02240       !
02241       data py1/-1.4566865832663635920d+04, 4.6905288611678631510d+06, &
02242            -6.9590439394619619534d+08, 4.3600098638603061642d+10, &
02243            -5.5107435206722644429d+11,-2.2213976967566192242d+13, &
02244            1.7427031242901594547d+01/
02245       data qy1/ 8.3030857612070288823d+02, 4.0669982352539552018d+05, &
02246            1.3960202770986831075d+08, 3.4015103849971240096d+10, &
02247            5.4266824419412347550d+12, 4.3386146580707264428d+14/
02248       !
02249       !  Coefficients for rational approximation of
02250       !  (Y0(X) - 2 LN(X/XY2) J0(X)) / (X**2 - XY2**2),
02251       !  5.5 < |X| <= 8.0
02252       !
02253       data py2/ 2.1363534169313901632d+04,-1.0085539923498211426d+07, &
02254            2.1958827170518100757d+09,-1.9363051266772083678d+11, &
02255            -1.2829912364088687306d+11, 6.7016641869173237784d+14, &
02256            -8.0728726905150210443d+15,-1.7439661319197499338d+01/
02257       data qy2/ 8.7903362168128450017d+02, 5.3924739209768057030d+05, &
02258            2.4727219475672302327d+08, 8.6926121104209825246d+10, &
02259            2.2598377924042897629d+13, 3.9272425569640309819d+15, &
02260            3.4563724628846457519d+17/
02261       !
02262       !  Coefficients for Hart,s approximation, 8.0 < |X|.
02263       !
02264       data p0/3.4806486443249270347d+03, 2.1170523380864944322d+04, &
02265            4.1345386639580765797d+04, 2.2779090197304684302d+04, &
02266            8.8961548424210455236d-01, 1.5376201909008354296d+02/
02267       data q0/3.5028735138235608207d+03, 2.1215350561880115730d+04, &
02268            4.1370412495510416640d+04, 2.2779090197304684318d+04, &
02269            1.5711159858080893649d+02/
02270       data p1/-2.2300261666214198472d+01,-1.1183429920482737611d+02, &
02271            -1.8591953644342993800d+02,-8.9226600200800094098d+01, &
02272            -8.8033303048680751817d-03,-1.2441026745835638459d+00/
02273       data q1/1.4887231232283756582d+03, 7.2642780169211018836d+03, &
02274            1.1951131543434613647d+04, 5.7105024128512061905d+03, &
02275            9.0593769594993125859d+01/
02276       !
02277       !  Check for error conditions.
02278       !
02279       ax = abs ( arg )
02280 
02281       if ( jint == 1 .and. arg <= zero ) then
02282          result = -xinf
02283          return
02284       else if ( xmax < ax ) then
02285          result = zero
02286          return
02287       end if
02288 
02289       if ( eight < ax ) then
02290          go to 800
02291       end if
02292 
02293       if ( ax <= xsmall ) then
02294          if ( jint == 0 ) then
02295             result = one
02296          else
02297             result = pi2 * ( log ( ax ) + cons )
02298          end if
02299          return
02300       end if
02301       !
02302       !  Calculate J0 for appropriate interval, preserving
02303       !  accuracy near the zero of J0.
02304       !
02305       zsq = ax * ax
02306 
02307       if ( ax <= four ) then
02308          xnum = ( pj0(5) * zsq + pj0(6) ) * zsq + pj0(7)
02309          xden = zsq + qj0(5)
02310          do i = 1, 4
02311             xnum = xnum * zsq + pj0(i)
02312             xden = xden * zsq + qj0(i)
02313          end do
02314          prod = ( ( ax - xj01 / two56 ) - xj02 ) * ( ax + xj0 )
02315       else
02316          wsq = one - zsq / sixty4
02317          xnum = pj1(7) * wsq + pj1(8)
02318          xden = wsq + qj1(7)
02319          do i = 1, 6
02320             xnum = xnum * wsq + pj1(i)
02321             xden = xden * wsq + qj1(i)
02322          end do
02323          prod = ( ax + xj1 ) * ( ( ax - xj11 / two56 ) - xj12 )
02324       end if
02325 
02326       result = prod * xnum / xden
02327 
02328       if ( jint == 0 ) then
02329          return
02330       end if
02331       !
02332       !  Calculate Y0.  First find  RESJ = pi/2 ln(x/xn) J0(x),
02333       !  where xn is a zero of Y0.
02334       !
02335       if ( ax <= three ) then
02336          up = ( ax - xy01 / two56 ) - xy02
02337          xy = xy0
02338       else if ( ax <= five5 ) then
02339          up = ( ax - xy11 / two56 ) - xy12
02340          xy = xy1
02341       else
02342          up = ( ax - xy21 / two56 ) - xy22
02343          xy = xy2
02344       end if
02345 
02346       down = ax + xy
02347 
02348       if ( abs ( up ) < p17 * down ) then
02349          w = up / down
02350          wsq = w * w
02351          xnum = plg(1)
02352          xden = wsq + qlg(1)
02353          do i = 2, 4
02354             xnum = xnum * wsq + plg(i)
02355             xden = xden * wsq + qlg(i)
02356          end do
02357          resj = pi2 * result * w * xnum / xden
02358       else
02359          resj = pi2 * result * log ( ax / xy )
02360       end if
02361       !
02362       !  Now calculate Y0 for appropriate interval, preserving
02363       !  accuracy near the zero of Y0.
02364       !
02365       if ( ax <= three ) then
02366          xnum = py0(6) * zsq + py0(1)
02367          xden = zsq + qy0(1)
02368          do i = 2, 5
02369             xnum = xnum * zsq + py0(i)
02370             xden = xden * zsq + qy0(i)
02371          end do
02372       else if ( ax <= five5 ) then
02373          xnum = py1(7) * zsq + py1(1)
02374          xden = zsq + qy1(1)
02375          do i = 2, 6
02376             xnum = xnum * zsq + py1(i)
02377             xden = xden * zsq + qy1(i)
02378          end do
02379       else
02380          xnum = py2(8) * zsq + py2(1)
02381          xden = zsq + qy2(1)
02382          do i = 2, 7
02383             xnum = xnum * zsq + py2(i)
02384             xden = xden * zsq + qy2(i)
02385          end do
02386       end if
02387 
02388       result = resj + up * down * xnum / xden
02389 
02390       return
02391       !
02392       !  Calculate J0 or Y0 for 8.0 < |ARG|.
02393       !
02394 800   continue
02395 
02396       z = eight / ax
02397       w = ax / twopi
02398       w = aint ( w ) + oneov8
02399       w = ( ax - w * twopi1 ) - w * twopi2
02400       zsq = z * z
02401       xnum = p0(5) * zsq + p0(6)
02402       xden = zsq + q0(5)
02403       up = p1(5) * zsq + p1(6)
02404       down = zsq + q1(5)
02405 
02406       do i = 1, 4
02407          xnum = xnum * zsq + p0(i)
02408          xden = xden * zsq + q0(i)
02409          up = up * zsq + p1(i)
02410          down = down * zsq + q1(i)
02411       end do
02412 
02413       r0 = xnum / xden
02414       r1 = up / down
02415 
02416       if ( jint == 0 ) then
02417          result = sqrt ( pi2 / ax ) &
02418               * ( r0 * cos ( w ) - z * r1 * sin ( w ) )
02419       else
02420          result = sqrt ( pi2 / ax ) &
02421               * ( r0 * sin ( w ) + z * r1 * cos ( w ) )
02422       end if
02423 
02424       return
02425    end subroutine caljy0
02426    subroutine caljy1 ( arg, result, jint )
02427 
02428       !*****************************************************************************80
02429       !
02430       !! CALJY1 computes various J1 and Y1 Bessel functions.
02431       !
02432       !  Discussion:
02433       !
02434       !    This routine computes first-order Bessel functions of the first and
02435       !    second kind (J1 and Y1), for real arguments X, where 0 < X <= XMAX
02436       !    for Y1, and |X| <= XMAX for J1.
02437       !
02438       !  Licensing:
02439       !
02440       !    This code is distributed under the GNU LGPL license.
02441       !
02442       !  Modified:
02443       !
02444       !    03 April 2007
02445       !
02446       !  Author:
02447       !
02448       !    Original FORTRAN77 version by William Cody.
02449       !    FORTRAN90 version by John Burkardt.
02450       !
02451       !  Reference:
02452       !
02453       !    John Hart, Ward Cheney, Charles Lawson, Hans Maehly,
02454       !    Charles Mesztenyi, John Rice, Henry Thatcher,
02455       !    Christoph Witzgall,
02456       !    Computer Approximations,
02457       !    Wiley, 1968,
02458       !    LC: QA297.C64.
02459       !
02460       !  Parameters:
02461       !
02462       !    Input, real ( kind = 8 ) ARG, the argument.  If JINT = 0, ARG
02463       !    must satisfy
02464       !     -XMAX < ARG < XMAX;
02465       !    If JINT = 1, then ARG must satisfy
02466       !      0 < ARG < XMAX.
02467       !
02468       !    Output, real ( kind = 8 ) RESULT, the value of the function,
02469       !    which depends on the input value of JINT:
02470       !    0, RESULT = J1(x);
02471       !    1, RESULT = Y1(x);
02472       !
02473       !    Input, integer ( kind = 4 ) JINT, chooses the function to be computed.
02474       !    0, J1(x);
02475       !    1, Y1(x);
02476       !
02477       implicit none
02478 
02479       real ( kind = 8 ) arg
02480       real ( kind = 8 ) ax
02481       real ( kind = 8 ) down
02482       real ( kind = 8 ) eight
02483       real ( kind = 8 ) four
02484       real ( kind = 8 ) half
02485       integer ( kind = 4 ) i
02486       integer ( kind = 4 ) jint
02487       real ( kind = 8 ) pi2
02488       real ( kind = 8 ) pj0(7)
02489       real ( kind = 8 ) pj1(8)
02490       real ( kind = 8 ) plg(4)
02491       real ( kind = 8 ) prod
02492       real ( kind = 8 ) py0(7)
02493       real ( kind = 8 ) py1(9)
02494       real ( kind = 8 ) p0(6)
02495       real ( kind = 8 ) p1(6)
02496       real ( kind = 8 ) p17
02497       real ( kind = 8 ) qj0(5)
02498       real ( kind = 8 ) qj1(7)
02499       real ( kind = 8 ) qlg(4)
02500       real ( kind = 8 ) qy0(6)
02501       real ( kind = 8 ) qy1(8)
02502       real ( kind = 8 ) q0(6)
02503       real ( kind = 8 ) q1(6)
02504       real ( kind = 8 ) resj
02505       real ( kind = 8 ) result
02506       real ( kind = 8 ) rtpi2
02507       real ( kind = 8 ) r0
02508       real ( kind = 8 ) r1
02509       real ( kind = 8 ) throv8
02510       real ( kind = 8 ) twopi
02511       real ( kind = 8 ) twopi1
02512       real ( kind = 8 ) twopi2
02513       real ( kind = 8 ) two56
02514       real ( kind = 8 ) up
02515       real ( kind = 8 ) w
02516       real ( kind = 8 ) wsq
02517       real ( kind = 8 ) xden
02518       real ( kind = 8 ) xinf
02519       real ( kind = 8 ) xmax
02520       real ( kind = 8 ) xnum
02521       real ( kind = 8 ) xsmall
02522       real ( kind = 8 ) xj0
02523       real ( kind = 8 ) xj1
02524       real ( kind = 8 ) xj01
02525       real ( kind = 8 ) xj02
02526       real ( kind = 8 ) xj11
02527       real ( kind = 8 ) xj12
02528       real ( kind = 8 ) xy
02529       real ( kind = 8 ) xy0
02530       real ( kind = 8 ) xy01
02531       real ( kind = 8 ) xy02
02532       real ( kind = 8 ) xy1
02533       real ( kind = 8 ) xy11
02534       real ( kind = 8 ) xy12
02535       real ( kind = 8 ) z
02536       real ( kind = 8 ) zero
02537       real ( kind = 8 ) zsq
02538       !
02539       !  Mathematical constants
02540       !
02541       data eight /8.0d0/
02542       data four /4.0d0/
02543       data half /0.5d0/
02544       data throv8 /0.375d0/
02545       data pi2 /6.3661977236758134308d-1/
02546       data p17 /1.716d-1/
02547       data twopi /6.2831853071795864769d+0/
02548       data zero /0.0d0/
02549       data twopi1 /6.28125d0/
02550       data twopi2 /1.9353071795864769253d-03/
02551       data two56 /256.0d+0/
02552       data rtpi2 /7.9788456080286535588d-1/
02553       !
02554       !  Machine-dependent constants
02555       !
02556       data xmax /1.07d+09/
02557       data xsmall /9.31d-10/
02558       data xinf /1.7d+38/
02559       !
02560       !  Zeroes of Bessel functions
02561       !
02562       data xj0 /3.8317059702075123156d+0/
02563       data xj1 /7.0155866698156187535d+0/
02564       data xy0 /2.1971413260310170351d+0/
02565       data xy1 /5.4296810407941351328d+0/
02566       data xj01 / 981.0d+0/
02567       data xj02 /-3.2527979248768438556d-04/
02568       data xj11 /1796.0d+0/
02569       data xj12 /-3.8330184381246462950d-05/
02570       data xy01 / 562.0d+0/
02571       data xy02 / 1.8288260310170351490d-03/
02572       data xy11 /1390.0d+0/
02573       data xy12 /-6.4592058648672279948d-06/
02574       !
02575       !  Coefficients for rational approximation to ln(x/a)
02576       !
02577       data plg/-2.4562334077563243311d+01,2.3642701335621505212d+02, &
02578            -5.4989956895857911039d+02,3.5687548468071500413d+02/
02579       data qlg/-3.5553900764052419184d+01,1.9400230218539473193d+02, &
02580            -3.3442903192607538956d+02,1.7843774234035750207d+02/
02581       !
02582       !  Coefficients for rational approximation of
02583       !  J1(X) / (X * (X**2 - XJ0**2)), XSMALL < |X| <=  4.0
02584       !
02585       data pj0/9.8062904098958257677d+05,-1.1548696764841276794d+08, &
02586            6.6781041261492395835d+09,-1.4258509801366645672d+11, &
02587            -4.4615792982775076130d+03, 1.0650724020080236441d+01, &
02588            -1.0767857011487300348d-02/
02589       data qj0/5.9117614494174794095d+05, 2.0228375140097033958d+08, &
02590            4.2091902282580133541d+10, 4.1868604460820175290d+12, &
02591            1.0742272239517380498d+03/
02592       !
02593       !  Coefficients for rational approximation of
02594       !  J1(X) / (X * (X**2 - XJ1**2)), 4.0 < |X| <= 8.0
02595       !
02596       data pj1/4.6179191852758252280d+00,-7.1329006872560947377d+03, &
02597            4.5039658105749078904d+06,-1.4437717718363239107d+09, &
02598            2.3569285397217157313d+11,-1.6324168293282543629d+13, &
02599            1.1357022719979468624d+14, 1.0051899717115285432d+15/
02600       data qj1/1.1267125065029138050d+06, 6.4872502899596389593d+08, &
02601            2.7622777286244082666d+11, 8.4899346165481429307d+13, &
02602            1.7128800897135812012d+16, 1.7253905888447681194d+18, &
02603            1.3886978985861357615d+03/
02604       !
02605       !  Coefficients for rational approximation of
02606       !  (Y1(X) - 2 LN(X/XY0) J1(X)) / (X**2 - XY0**2),
02607       !  XSMALL < |X| <=  4.0
02608       !
02609       data py0/2.2157953222280260820d+05,-5.9157479997408395984d+07, &
02610            7.2144548214502560419d+09,-3.7595974497819597599d+11, &
02611            5.4708611716525426053d+12, 4.0535726612579544093d+13, &
02612            -3.1714424660046133456d+02/
02613       data qy0/8.2079908168393867438d+02, 3.8136470753052572164d+05, &
02614            1.2250435122182963220d+08, 2.7800352738690585613d+10, &
02615            4.1272286200406461981d+12, 3.0737873921079286084d+14/
02616       !
02617       !  Coefficients for rational approximation of
02618       !  (Y1(X) - 2 LN(X/XY1) J1(X)) / (X**2 - XY1**2),
02619       !  4.0 < |X| <= 8.0
02620       !
02621       data py1/ 1.9153806858264202986d+06,-1.1957961912070617006d+09, &
02622            3.7453673962438488783d+11,-5.9530713129741981618d+13, &
02623            4.0686275289804744814d+15,-2.3638408497043134724d+16, &
02624            -5.6808094574724204577d+18, 1.1514276357909013326d+19, &
02625            -1.2337180442012953128d+03/
02626       data qy1/ 1.2855164849321609336d+03, 1.0453748201934079734d+06, &
02627            6.3550318087088919566d+08, 3.0221766852960403645d+11, &
02628            1.1187010065856971027d+14, 3.0837179548112881950d+16, &
02629            5.6968198822857178911d+18, 5.3321844313316185697d+20/
02630       !
02631       !  Coefficients for Hart's approximation, 8.0 < |X|.
02632       !
02633       data p0/-1.0982405543459346727d+05,-1.5235293511811373833d+06, &
02634            -6.6033732483649391093d+06,-9.9422465050776411957d+06, &
02635            -4.4357578167941278571d+06,-1.6116166443246101165d+03/
02636       data q0/-1.0726385991103820119d+05,-1.5118095066341608816d+06, &
02637            -6.5853394797230870728d+06,-9.9341243899345856590d+06, &
02638            -4.4357578167941278568d+06,-1.4550094401904961825d+03/
02639       data p1/ 1.7063754290207680021d+03, 1.8494262873223866797d+04, &
02640            6.6178836581270835179d+04, 8.5145160675335701966d+04, &
02641            3.3220913409857223519d+04, 3.5265133846636032186d+01/
02642       data q1/ 3.7890229745772202641d+04, 4.0029443582266975117d+05, &
02643            1.4194606696037208929d+06, 1.8194580422439972989d+06, &
02644            7.0871281941028743574d+05, 8.6383677696049909675d+02/
02645       !
02646       !  Check for error conditions.
02647       !
02648       ax = abs ( arg )
02649 
02650       if ( jint == 1 .and. ( arg <= zero .or. &
02651            ( arg < half .and. ax * xinf < pi2 ) ) ) then
02652          result = -xinf
02653          return
02654       else if ( xmax < ax ) then
02655          result = zero
02656          return
02657       end if
02658 
02659       if ( eight < ax ) then
02660          go to 800
02661       else if ( ax <= xsmall ) then
02662          if ( jint == 0 ) then
02663             result = arg * half
02664          else
02665             result = -pi2 / ax
02666          end if
02667          return
02668       end if
02669       !
02670       !  Calculate J1 for appropriate interval, preserving
02671       !  accuracy near the zero of J1.
02672       !
02673       zsq = ax * ax
02674 
02675       if ( ax <= four ) then
02676          xnum = ( pj0(7) * zsq + pj0(6) ) * zsq + pj0(5)
02677          xden = zsq + qj0(5)
02678          do i = 1, 4
02679             xnum = xnum * zsq + pj0(i)
02680             xden = xden * zsq + qj0(i)
02681          end do
02682          prod = arg * ( ( ax - xj01 / two56 ) - xj02 ) * ( ax + xj0 )
02683       else
02684          xnum = pj1(1)
02685          xden = ( zsq + qj1(7) ) * zsq + qj1(1)
02686          do i = 2, 6
02687             xnum = xnum * zsq + pj1(i)
02688             xden = xden * zsq + qj1(i)
02689          end do
02690          xnum = xnum * ( ax - eight ) * ( ax + eight ) + pj1(7)
02691          xnum = xnum * ( ax - four ) * ( ax + four ) + pj1(8)
02692          prod = arg * ( ( ax - xj11 / two56 ) - xj12 ) * ( ax + xj1 )
02693       end if
02694 
02695       result = prod * ( xnum / xden )
02696 
02697       if ( jint == 0 ) then
02698          return
02699       end if
02700       !
02701       !  Calculate Y1.  First find RESJ = pi/2 ln(x/xn) J1(x),
02702       !  where xn is a zero of Y1.
02703       !
02704       if ( ax <= four ) then
02705          up = ( ax - xy01 / two56 ) - xy02
02706          xy = xy0
02707       else
02708          up = ( ax - xy11 / two56 ) - xy12
02709          xy = xy1
02710       end if
02711 
02712       down = ax + xy
02713 
02714       if ( abs ( up ) < p17 * down ) then
02715          w = up / down
02716          wsq = w * w
02717          xnum = plg(1)
02718          xden = wsq + qlg(1)
02719          do i = 2, 4
02720             xnum = xnum * wsq + plg(i)
02721             xden = xden * wsq + qlg(i)
02722          end do
02723          resj = pi2 * result * w * xnum / xden
02724       else
02725          resj = pi2 * result * log ( ax / xy )
02726       end if
02727       !
02728       !  Now calculate Y1 for appropriate interval, preserving
02729       !  accuracy near the zero of Y1.
02730       !
02731       if ( ax <= four ) then
02732          xnum = py0(7) * zsq + py0(1)
02733          xden = zsq + qy0(1)
02734          do i = 2, 6
02735             xnum = xnum * zsq + py0(i)
02736             xden = xden * zsq + qy0(i)
02737          end do
02738       else
02739          xnum = py1(9) * zsq + py1(1)
02740          xden = zsq + qy1(1)
02741          do i = 2, 8
02742             xnum = xnum * zsq + py1(i)
02743             xden = xden * zsq + qy1(i)
02744          end do
02745       end if
02746 
02747       result = resj + ( up * down / ax ) * xnum / xden
02748       return
02749       !
02750       !  Calculate J1 or Y1 for 8.0 < |ARG|.
02751       !
02752 800   continue
02753 
02754       z = eight / ax
02755       w = aint ( ax / twopi ) + throv8
02756       w = ( ax - w * twopi1 ) - w * twopi2
02757       zsq = z * z
02758       xnum = p0(6)
02759       xden = zsq + q0(6)
02760       up = p1(6)
02761       down = zsq + q1(6)
02762 
02763       do i = 1, 5
02764          xnum = xnum * zsq + p0(i)
02765          xden = xden * zsq + q0(i)
02766          up = up * zsq + p1(i)
02767          down = down * zsq + q1(i)
02768       end do
02769 
02770       r0 = xnum / xden
02771       r1 = up / down
02772 
02773       if ( jint == 0 ) then
02774          result = ( rtpi2 / sqrt ( ax ) ) &
02775               * ( r0 * cos ( w ) - z * r1 * sin ( w ) )
02776       else
02777          result = ( rtpi2 / sqrt ( ax ) ) &
02778               * ( r0 * sin ( w ) + z * r1 * cos ( w ) )
02779       end if
02780 
02781       if ( jint == 0 .and. arg < zero ) then
02782          result = -result
02783       end if
02784 
02785       return
02786    end subroutine caljy1
02787    function daw ( xx )
02788 
02789       !*****************************************************************************80
02790       !
02791       !! DAW evaluates Dawson's integral function.
02792       !
02793       !  Discussion:
02794       !
02795       !    This routine evaluates Dawson's integral,
02796       !
02797       !      F(x) = exp ( - x * x ) * Integral ( 0 <= t <= x ) exp ( t * t ) dt
02798       !
02799       !    for a real argument x.
02800       !
02801       !  Licensing:
02802       !
02803       !    This code is distributed under the GNU LGPL license.
02804       !
02805       !  Modified:
02806       !
02807       !    03 April 2007
02808       !
02809       !  Author:
02810       !
02811       !    Original FORTRAN77 version by William Cody.
02812       !    FORTRAN90 version by John Burkardt.
02813       !
02814       !  Reference:
02815       !
02816       !    William Cody, Kathleen Paciorek, Henry Thacher,
02817       !    Chebyshev Approximations for Dawson's Integral,
02818       !    Mathematics of Computation,
02819       !    Volume 24, Number 109, January 1970, pages 171-178.
02820       !
02821       !  Parameters:
02822       !
02823       !    Input, real ( kind = 8 ) XX, the argument of the function.
02824       !
02825       !    Output, real ( kind = 8 ) DAW, the value of the function.
02826       !
02827       implicit none
02828 
02829       real ( kind = 8 ) daw
02830       real ( kind = 8 ) frac
02831       real ( kind = 8 ) half
02832       integer ( kind = 4 ) i
02833       real ( kind = 8 ) one
02834       real ( kind = 8 ) one225
02835       real ( kind = 8 ) p1(10)
02836       real ( kind = 8 ) p2(10)
02837       real ( kind = 8 ) p3(10)
02838       real ( kind = 8 ) p4(10)
02839       real ( kind = 8 ) q1(10)
02840       real ( kind = 8 ) q2(9)
02841       real ( kind = 8 ) q3(9)
02842       real ( kind = 8 ) q4(9)
02843       real ( kind = 8 ) six25
02844       real ( kind = 8 ) sump
02845       real ( kind = 8 ) sumq
02846       real ( kind = 8 ) two5
02847       real ( kind = 8 ) w2
02848       real ( kind = 8 ) x
02849       real ( kind = 8 ) xlarge
02850       real ( kind = 8 ) xmax
02851       real ( kind = 8 ) xsmall
02852       real ( kind = 8 ) xx
02853       real ( kind = 8 ) y
02854       real ( kind = 8 ) zero
02855       !
02856       !  Mathematical constants.
02857       !
02858       data zero / 0.0D+00 /
02859       data half / 0.5D+00 /
02860       data one / 1.0D+00 /
02861       data six25 /6.25D+00 /
02862       data one225 /12.25d0 /
02863       data two5 /25.0d0/
02864       !
02865       !  Machine-dependent constants
02866       !
02867       data xsmall /1.05d-08/
02868       data xlarge /9.49d+07/
02869       data xmax /2.24d+307/
02870       !
02871       !  Coefficients for R(9,9) approximation for  |x| < 2.5
02872       !
02873       data p1/-2.69020398788704782410d-12, 4.18572065374337710778d-10, &
02874            -1.34848304455939419963d-08, 9.28264872583444852976d-07, &
02875            -1.23877783329049120592d-05, 4.07205792429155826266d-04, &
02876            -2.84388121441008500446d-03, 4.70139022887204722217d-02, &
02877            -1.38868086253931995101d-01, 1.00000000000000000004d+00/
02878       data q1/ 1.71257170854690554214d-10, 1.19266846372297253797d-08, &
02879            4.32287827678631772231d-07, 1.03867633767414421898d-05, &
02880            1.78910965284246249340d-04, 2.26061077235076703171d-03, &
02881            2.07422774641447644725d-02, 1.32212955897210128811d-01, &
02882            5.27798580412734677256d-01, 1.00000000000000000000d+00/
02883       !
02884       !  Coefficients for R(9,9) approximation in J-fraction form
02885       !  for  x in [2.5, 3.5)
02886       !
02887       data p2/-1.70953804700855494930d+00,-3.79258977271042880786d+01, &
02888            2.61935631268825992835d+01, 1.25808703738951251885d+01, &
02889            -2.27571829525075891337d+01, 4.56604250725163310122d+00, &
02890            -7.33080089896402870750d+00, 4.65842087940015295573d+01, &
02891            -1.73717177843672791149d+01, 5.00260183622027967838d-01/
02892       data q2/ 1.82180093313514478378d+00, 1.10067081034515532891d+03, &
02893            -7.08465686676573000364d+00, 4.53642111102577727153d+02, &
02894            4.06209742218935689922d+01, 3.02890110610122663923d+02, &
02895            1.70641269745236227356d+02, 9.51190923960381458747d+02, &
02896            2.06522691539642105009d-01/
02897       !
02898       !  Coefficients for R(9,9) approximation in J-fraction form
02899       !  for  x in [3.5, 5.0]
02900       !
02901       data p3/-4.55169503255094815112d+00,-1.86647123338493852582d+01, &
02902            -7.36315669126830526754d+00,-6.68407240337696756838d+01, &
02903            4.84507265081491452130d+01, 2.69790586735467649969d+01, &
02904            -3.35044149820592449072d+01, 7.50964459838919612289d+00, &
02905            -1.48432341823343965307d+00, 4.99999810924858824981d-01/
02906       data q3/ 4.47820908025971749852d+01, 9.98607198039452081913d+01, &
02907            1.40238373126149385228d+01, 3.48817758822286353588d+03, &
02908            -9.18871385293215873406d+00, 1.24018500009917163023d+03, &
02909            -6.88024952504512254535d+01,-2.31251575385145143070d+00, &
02910            2.50041492369922381761d-01/
02911       !
02912       !  Coefficients for R(9,9) approximation in J-fraction form
02913       !  for 5.0 < |x|.
02914       !
02915       data p4/-8.11753647558432685797d+00,-3.84043882477454453430d+01, &
02916            -2.23787669028751886675d+01,-2.88301992467056105854d+01, &
02917            -5.99085540418222002197d+00,-1.13867365736066102577d+01, &
02918            -6.52828727526980741590d+00,-4.50002293000355585708d+00, &
02919            -2.50000000088955834952d+00, 5.00000000000000488400d-01/
02920       data q4/ 2.69382300417238816428d+02, 5.04198958742465752861d+01, &
02921            6.11539671480115846173d+01, 2.08210246935564547889d+02, &
02922            1.97325365692316183531d+01,-1.22097010558934838708d+01, &
02923            -6.99732735041547247161d+00,-2.49999970104184464568d+00, &
02924            7.49999999999027092188d-01/
02925 
02926       x = xx
02927 
02928       if ( xlarge < abs ( x ) ) then
02929 
02930          if ( abs ( x ) <= xmax ) then
02931             daw = half / x
02932          else
02933             daw = zero
02934          end if
02935 
02936       else if ( abs ( x ) < xsmall ) then
02937 
02938          daw = x
02939 
02940       else
02941 
02942          y = x * x
02943          !
02944          !  ABS(X) < 2.5.
02945          !
02946          if ( y < six25 ) then
02947 
02948             sump = p1(1)
02949             sumq = q1(1)
02950             do i = 2, 10
02951                sump = sump * y + p1(i)
02952                sumq = sumq * y + q1(i)
02953             end do
02954 
02955             daw = x * sump / sumq
02956             !
02957             !  2.5 <= ABS(X) < 3.5.
02958             !
02959          else if ( y < one225 ) then
02960 
02961             frac = zero
02962             do i = 1, 9
02963                frac = q2(i) / ( p2(i) + y + frac )
02964             end do
02965 
02966             daw = ( p2(10) + frac ) / x
02967             !
02968             !  3.5 <= ABS(X) < 5.0.
02969             !
02970          else if ( y < two5 ) then
02971 
02972             frac = zero
02973             do i = 1, 9
02974                frac = q3(i) / ( p3(i) + y + frac )
02975             end do
02976 
02977             daw = ( p3(10) + frac ) / x
02978 
02979          else
02980             !
02981             !  5.0 <= ABS(X) <= XLARGE.
02982             !
02983             w2 = one / x / x
02984 
02985             frac = zero
02986             do i = 1, 9
02987                frac = q4(i) / ( p4(i) + y + frac )
02988             end do
02989             frac = p4(10) + frac
02990 
02991             daw = ( half + half * w2 * frac ) / x
02992 
02993          end if
02994 
02995       end if
02996 
02997       return
02998    end function daw
02999    function dlgama ( x )
03000 
03001       !*****************************************************************************80
03002       !
03003       !! DLGAMA evaluates log ( Gamma ( X ) ) for a real argument.
03004       !
03005       !  Discussion:
03006       !
03007       !    This routine calculates the LOG(GAMMA) function for a positive real
03008       !    argument X.  Computation is based on an algorithm outlined in
03009       !    references 1 and 2.  The program uses rational functions that
03010       !    theoretically approximate LOG(GAMMA) to at least 18 significant
03011       !    decimal digits.  The approximation for X > 12 is from reference
03012       !    3, while approximations for X < 12.0 are similar to those in
03013       !    reference 1, but are unpublished.
03014       !
03015       !  Licensing:
03016       !
03017       !    This code is distributed under the GNU LGPL license.
03018       !
03019       !  Modified:
03020       !
03021       !    03 April 2007
03022       !
03023       !  Author:
03024       !
03025       !    Original FORTRAN77 version by William Cody, Laura Stoltz.
03026       !    FORTRAN90 version by John Burkardt.
03027       !
03028       !  Reference:
03029       !
03030       !    William Cody, Kenneth Hillstrom,
03031       !    Chebyshev Approximations for the Natural Logarithm of the
03032       !    Gamma Function,
03033       !    Mathematics of Computation,
03034       !    Volume 21, Number 98, April 1967, pages 198-203.
03035       !
03036       !    Kenneth Hillstrom,
03037       !    ANL/AMD Program ANLC366S, DGAMMA/DLGAMA,
03038       !    May 1969.
03039       !
03040       !    John Hart, Ward Cheney, Charles Lawson, Hans Maehly,
03041       !    Charles Mesztenyi, John Rice, Henry Thatcher,
03042       !    Christoph Witzgall,
03043       !    Computer Approximations,
03044       !    Wiley, 1968,
03045       !    LC: QA297.C64.
03046       !
03047       !  Parameters:
03048       !
03049       !    Input, real ( kind = 8 ) X, the argument of the function.
03050       !
03051       !    Output, real ( kind = 8 ) DLGAMA, the value of the function.
03052       !
03053       implicit none
03054 
03055       real ( kind = 8 ) c(7)
03056       real ( kind = 8 ) corr
03057       real ( kind = 8 ) d1
03058       real ( kind = 8 ) d2
03059       real ( kind = 8 ) d4
03060       real ( kind = 8 ) dlgama
03061       real ( kind = 8 ) eps
03062       real ( kind = 8 ) frtbig
03063       real ( kind = 8 ) four
03064       real ( kind = 8 ) half
03065       integer ( kind = 4 ) i
03066       real ( kind = 8 ) one
03067       real ( kind = 8 ) pnt68
03068       real ( kind = 8 ) p1(8)
03069       real ( kind = 8 ) p2(8)
03070       real ( kind = 8 ) p4(8)
03071       real ( kind = 8 ) q1(8)
03072       real ( kind = 8 ) q2(8)
03073       real ( kind = 8 ) q4(8)
03074       real ( kind = 8 ) res
03075       real ( kind = 8 ) sqrtpi
03076       real ( kind = 8 ) thrhal
03077       real ( kind = 8 ) twelve
03078       real ( kind = 8 ) two
03079       real ( kind = 8 ) x
03080       real ( kind = 8 ) xbig
03081       real ( kind = 8 ) xden
03082       real ( kind = 8 ) xinf
03083       real ( kind = 8 ) xm1
03084       real ( kind = 8 ) xm2
03085       real ( kind = 8 ) xm4
03086       real ( kind = 8 ) xnum
03087       real ( kind = 8 ) y
03088       real ( kind = 8 ) ysq
03089       real ( kind = 8 ) zero
03090       !
03091       !  Mathematical constants
03092       !
03093       data one /1.0d0/
03094       data half /0.5d0/
03095       data twelve /12.0d0/
03096       data zero /0.0d0/
03097       data four /4.0d0/
03098       data thrhal /1.5d0/
03099       data two /2.0d0/
03100       data pnt68 /0.6796875d0/
03101       data sqrtpi /0.9189385332046727417803297d0/
03102       !
03103       !  Machine dependent parameters
03104       !
03105       data xbig /2.55d305/
03106       data xinf /1.79d308/
03107       data eps /2.22d-16/
03108       data frtbig /2.25d76/
03109       !
03110       !  Numerator and denominator coefficients for rational minimax
03111       !  approximation over (0.5,1.5).
03112       !
03113       data d1/-5.772156649015328605195174d-1/
03114       data p1/4.945235359296727046734888d0,2.018112620856775083915565d2, &
03115            2.290838373831346393026739d3,1.131967205903380828685045d4, &
03116            2.855724635671635335736389d4,3.848496228443793359990269d4, &
03117            2.637748787624195437963534d4,7.225813979700288197698961d3/
03118       data q1/6.748212550303777196073036d1,1.113332393857199323513008d3, &
03119            7.738757056935398733233834d3,2.763987074403340708898585d4, &
03120            5.499310206226157329794414d4,6.161122180066002127833352d4, &
03121            3.635127591501940507276287d4,8.785536302431013170870835d3/
03122       !
03123       !  Numerator and denominator coefficients for rational minimax
03124       !  Approximation over (1.5,4.0).
03125       !
03126       data d2/4.227843350984671393993777d-1/
03127       data p2/4.974607845568932035012064d0,5.424138599891070494101986d2, &
03128            1.550693864978364947665077d4,1.847932904445632425417223d5, &
03129            1.088204769468828767498470d6,3.338152967987029735917223d6, &
03130            5.106661678927352456275255d6,3.074109054850539556250927d6/
03131       data q2/1.830328399370592604055942d2,7.765049321445005871323047d3, &
03132            1.331903827966074194402448d5,1.136705821321969608938755d6, &
03133            5.267964117437946917577538d6,1.346701454311101692290052d7, &
03134            1.782736530353274213975932d7,9.533095591844353613395747d6/
03135       !
03136       !  Numerator and denominator coefficients for rational minimax
03137       !  Approximation over (4.0,12.0).
03138       !
03139       data d4/1.791759469228055000094023d0/
03140       data p4/1.474502166059939948905062d4,2.426813369486704502836312d6, &
03141            1.214755574045093227939592d8,2.663432449630976949898078d9, &
03142            2.940378956634553899906876d10,1.702665737765398868392998d11, &
03143            4.926125793377430887588120d11,5.606251856223951465078242d11/
03144       data q4/2.690530175870899333379843d3,6.393885654300092398984238d5, &
03145            4.135599930241388052042842d7,1.120872109616147941376570d9, &
03146            1.488613728678813811542398d10,1.016803586272438228077304d11, &
03147            3.417476345507377132798597d11,4.463158187419713286462081d11/
03148       !
03149       !  Coefficients for minimax approximation over (12, INF).
03150       !
03151       data c/-1.910444077728d-03,8.4171387781295d-04, &
03152            -5.952379913043012d-04,7.93650793500350248d-04, &
03153            -2.777777777777681622553d-03,8.333333333333333331554247d-02, &
03154            5.7083835261d-03/
03155 
03156       y = x
03157 
03158       if ( zero < y .and. y <= xbig ) then
03159 
03160          if ( y <= eps ) then
03161 
03162             res = -log ( y )
03163             !
03164             !  EPS < X <= 1.5.
03165             !
03166          else if ( y <= thrhal ) then
03167 
03168             if ( y < pnt68 ) then
03169                corr = -log ( y )
03170                xm1 = y
03171             else
03172                corr = zero
03173                xm1 = ( y - half ) - half
03174             end if
03175 
03176             if ( y <= half .or. pnt68 <= y ) then
03177 
03178                xden = one
03179                xnum = zero
03180                do i = 1, 8
03181                   xnum = xnum * xm1 + p1(i)
03182                   xden = xden * xm1 + q1(i)
03183                end do
03184 
03185                res = corr + ( xm1 * ( d1 + xm1 * ( xnum / xden ) ) )
03186 
03187             else
03188 
03189                xm2 = ( y - half ) - half
03190                xden = one
03191                xnum = zero
03192                do i = 1, 8
03193                   xnum = xnum * xm2 + p2(i)
03194                   xden = xden * xm2 + q2(i)
03195                end do
03196 
03197                res = corr + xm2 * (d2 + xm2 * ( xnum / xden ) )
03198 
03199             end if
03200             !
03201             !  1.5 < X <= 4.0.
03202             !
03203          else if ( y <= four ) then
03204 
03205             xm2 = y - two
03206             xden = one
03207             xnum = zero
03208             do i = 1, 8
03209                xnum = xnum * xm2 + p2(i)
03210                xden = xden * xm2 + q2(i)
03211             end do
03212 
03213             res = xm2 * ( d2 + xm2 * ( xnum / xden ) )
03214             !
03215             !  4.0 < X <= 12.0.
03216             !
03217          else if ( y <= twelve ) then
03218 
03219             xm4 = y - four
03220             xden = -one
03221             xnum = zero
03222             do i = 1, 8
03223                xnum = xnum * xm4 + p4(i)
03224                xden = xden * xm4 + q4(i)
03225             end do
03226 
03227             res = d4 + xm4 * ( xnum / xden )
03228             !
03229             !  Evaluate for 12 <= argument.
03230             !
03231          else
03232 
03233             res = zero
03234 
03235             if ( y <= frtbig ) then
03236 
03237                res = c(7)
03238                ysq = y * y
03239 
03240                do i = 1, 6
03241                   res = res / ysq + c(i)
03242                end do
03243 
03244             end if
03245 
03246             res = res / y
03247             corr = log ( y )
03248             res = res + sqrtpi - half * corr
03249             res = res + y * ( corr - one )
03250 
03251          end if
03252          !
03253          !  Return for bad arguments.
03254          !
03255       else
03256 
03257          res = xinf
03258 
03259       end if
03260       !
03261       !  Final adjustments and return.
03262       !
03263       dlgama = res
03264 
03265       return
03266    end function dlgama
03267    subroutine dsubn ( x, nmax, xmax, d )
03268 
03269       !*****************************************************************************80
03270       !
03271       !! DSUBN evaluates derivatives of Ei(X).
03272       !
03273       !  Discussion:
03274       !
03275       !    Translation of Gautschi CACM Algorithm 282 for derivatives of Ei(x).
03276       !
03277       !  Licensing:
03278       !
03279       !    This code is distributed under the GNU LGPL license.
03280       !
03281       !  Modified:
03282       !
03283       !    03 April 2007
03284       !
03285       !  Author:
03286       !
03287       !    Original FORTRAN77 version by William Cody.
03288       !    FORTRAN90 version by John Burkardt.
03289       !
03290       !  Reference:
03291       !
03292       !    Walter Gautschi,
03293       !    Algorithm 282:
03294       !    Derivatives of EXP(X)/X, COS(X)/X, and SIN(X)/X,
03295       !    Communications of the ACM,
03296       !    Volume 9, April 1966, page 272.
03297       !
03298       !  Parameters:
03299       !
03300       !    Input, real ( kind = 8 ) X, the argument.
03301       !
03302       !    Input, integer ( kind = 4 ) NMAX, the number of derivatives to compute.
03303       !
03304       !    Input, real ( kind = 8 ) XMAX, the largest finite floating-point number,
03305       !    as computed, for instance, by MACHAR.
03306       !
03307       !    Output, real ( kind = 8 ) D(0:NMAX), the value and first NMAX derivatives
03308       !    of Ei(X).
03309       !
03310       implicit none
03311 
03312       integer ( kind = 4 ) nmax
03313 
03314       real ( kind = 8 ) b0
03315       real ( kind = 8 ) b1
03316       real ( kind = 8 ) b2
03317       real ( kind = 8 ) b3
03318       real ( kind = 8 ) b4
03319       real ( kind = 8 ) b5
03320       real ( kind = 8 ) b6
03321       logical bool1
03322       logical bool2
03323       real ( kind = 8 ) c0
03324       real ( kind = 8 ) c1
03325       real ( kind = 8 ) d(0:nmax)
03326       real ( kind = 8 ) e
03327       real ( kind = 8 ) en
03328       integer ( kind = 4 ) j
03329       integer ( kind = 4 ) lim
03330       integer ( kind = 4 ) mini
03331       integer ( kind = 4 ) n
03332       integer ( kind = 4 ) n0
03333       integer ( kind = 4 ) n1
03334       real ( kind = 8 ) one
03335       real ( kind = 8 ) p
03336       real ( kind = 8 ) q
03337       real ( kind = 8 ) t
03338       real ( kind = 8 ) ten
03339       real ( kind = 8 ) two
03340       real ( kind = 8 ) x
03341       real ( kind = 8 ) xmax
03342       real ( kind = 8 ) x1
03343       real ( kind = 8 ) z
03344       real ( kind = 8 ) zero
03345 
03346       data zero /0.0d0/
03347       data one /1.0d0/
03348       data two /2.0d0/
03349       data ten /10.0d0/
03350       data c0 /2.7183d0/
03351       data c1 /4.67452d1/
03352       data b0 /5.7941d-5/
03353       data b1 /-1.76148d-3/
03354       data b2 /2.08645d-2/
03355       data b3 /-1.29013d-1/
03356       data b4 /8.5777d-1/
03357       data b5 /1.0125d0/
03358       data b6 /7.75d-1/
03359 
03360       x1 = abs ( x )
03361       n0 = int ( x1 )
03362       e = exp ( x )
03363       d(0) = e / x
03364       bool1 = ( x < zero ) .or. ( x1 <= two )
03365       bool2 = n0 < nmax
03366       mini = min ( n0, nmax )
03367 
03368       if ( bool1 ) then
03369          lim = nmax
03370       else
03371          lim = mini
03372       end if
03373 
03374       n = 1
03375       en = one
03376 
03377 50    continue
03378 
03379       d(n) = ( e - en * d(n-1) ) / x
03380       n = n + 1
03381       en = en + one
03382 
03383       if ( x1 < one ) then
03384 
03385          if ( abs ( d(n-1) ) < abs ( xmax * x / en ) .and. &
03386               n <= lim ) then
03387             go to 50
03388          end if
03389 
03390       else
03391 
03392          if ( abs ( d(n-1) / x ) < xmax / en .and. n <= lim ) then
03393             go to 50
03394          end if
03395 
03396       end if
03397 
03398       do j = n, lim
03399          d(n) = zero
03400       end do
03401 
03402       if ( .not. bool1 .and. bool2 ) then
03403 
03404          t = ( x1 + c1 ) / ( c0 * x1 )
03405 
03406          if ( t < ten ) then
03407             t = (((( &
03408                  b0 &
03409                  * t + b1 ) &
03410                  * t + b2 ) &
03411                  * t + b3 ) &
03412                  * t + b4 ) &
03413                  * t + b5
03414          else
03415             z = log ( t ) - b6
03416             p = ( b6 - log ( z ) ) / ( one + z )
03417             p = one / ( one + p )
03418             t = t * p / z
03419          end if
03420 
03421          n1 = c0 * x1 * t - one
03422          if ( n1 < nmax ) then
03423             n1 = nmax
03424          end if
03425 
03426          q = one / x
03427          en = one
03428          do n = 1, n1+1
03429             q = -en * q / x
03430             en = en + one
03431          end do
03432 
03433          do n = n1, n0+1, -1
03434 
03435             en = en - one
03436             q = ( e - x * q ) / en
03437 
03438             if ( n <= nmax ) then
03439                d(n) = q
03440             end if
03441 
03442          end do
03443 
03444       end if
03445 
03446       return
03447    end subroutine dsubn
03448    function ei ( x )
03449 
03450       !*****************************************************************************80
03451       !
03452       !! EI evaluates the exponential integral Ei(X).
03453       !
03454       !  Discussion:
03455       !
03456       !    This routine computes approximate values for the
03457       !    exponential integral Ei(x), where x is real.
03458       !
03459       !  Licensing:
03460       !
03461       !    This code is distributed under the GNU LGPL license.
03462       !
03463       !  Modified:
03464       !
03465       !    03 April 2007
03466       !
03467       !  Author:
03468       !
03469       !    Original FORTRAN77 version by William Cody.
03470       !    FORTRAN90 version by John Burkardt.
03471       !
03472       !  Parameters:
03473       !
03474       !    Input, real ( kind = 8 ) X, the argument of the function.
03475       !
03476       !    Output, real ( kind = 8 ) EI, the value of the function.
03477       !
03478       implicit none
03479 
03480       real ( kind = 8 ) ei
03481       integer ( kind = 4 ) jint
03482       real ( kind = 8 ) result
03483       real ( kind = 8 ) x
03484 
03485       jint = 1
03486       call calcei ( x, result, jint )
03487       ei = result
03488 
03489       return
03490    end function ei
03491    function eone ( x )
03492 
03493       !*****************************************************************************80
03494       !
03495       !! EONE evaluates the exponential integral E1(X).
03496       !
03497       !  Discussion:
03498       !
03499       !    This routine computes approximate values for the
03500       !    exponential integral E1(x), where x is real.
03501       !
03502       !  Licensing:
03503       !
03504       !    This code is distributed under the GNU LGPL license.
03505       !
03506       !  Modified:
03507       !
03508       !    03 April 2007
03509       !
03510       !  Author:
03511       !
03512       !    Original FORTRAN77 version by William Cody.
03513       !    FORTRAN90 version by John Burkardt.
03514       !
03515       !  Parameters:
03516       !
03517       !    Input, real ( kind = 8 ) X, the argument of the function.
03518       !
03519       !    Output, real ( kind = 8 ) EONE, the value of the function.
03520       !
03521       implicit none
03522 
03523       real ( kind = 8 ) eone
03524       integer ( kind = 4 ) jint
03525       real ( kind = 8 ) result
03526       real ( kind = 8 ) x
03527 
03528       jint = 2
03529       call calcei ( x, result, jint )
03530       eone = result
03531 
03532       return
03533    end function eone
03534    function expei ( x )
03535 
03536       !*****************************************************************************80
03537       !
03538       !! EXPEI evaluates the scaled exponential integral exp(-X) * Ei(X).
03539       !
03540       !  Discussion:
03541       !
03542       !    This routine computes approximate values for the
03543       !    function  exp(-x) * Ei(x), where  Ei(x)  is the exponential
03544       !    integral, and  x  is real.
03545       !
03546       !  Licensing:
03547       !
03548       !    This code is distributed under the GNU LGPL license.
03549       !
03550       !  Modified:
03551       !
03552       !    03 April 2007
03553       !
03554       !  Author:
03555       !
03556       !    Original FORTRAN77 version by William Cody.
03557       !    FORTRAN90 version by John Burkardt.
03558       !
03559       !  Parameters:
03560       !
03561       !    Input, real ( kind = 8 ) X, the argument of the function.
03562       !
03563       !    Output, real ( kind = 8 ) EXPEI, the value of the function.
03564       !
03565       implicit none
03566 
03567       real ( kind = 8 ) expei
03568       integer ( kind = 4 ) jint
03569       real ( kind = 8 ) result
03570       real ( kind = 8 ) x
03571 
03572       jint = 3
03573       call calcei ( x, result, jint )
03574       expei = result
03575 
03576       return
03577    end function expei
03578    subroutine machar ( ibeta, it, irnd, ngrd, machep, negep, iexp, &
03579         minexp, maxexp, eps, epsneg, xmin, xmax )
03580 
03581       !*****************************************************************************80
03582       !
03583       !! MACHAR determines various machine arithmetic parameters.
03584       !
03585       !  Discussion:
03586       !
03587       !    This routine is intended to determine the parameters
03588       !    of the floating-point arithmetic system specified below.
03589       !
03590       !  Licensing:
03591       !
03592       !    This code is distributed under the GNU LGPL license.
03593       !
03594       !  Modified:
03595       !
03596       !    03 April 2007
03597       !
03598       !  Author:
03599       !
03600       !    Original FORTRAN77 version by William Cody.
03601       !    FORTRAN90 version by John Burkardt.
03602       !
03603       !  Reference:
03604       !
03605       !    William Cody, William Waite,
03606       !    Software Manual for the Elementary Functions,
03607       !    Prentice Hall, 1980,
03608       !    ISBN: 0138220646,
03609       !    LC: QA331.C635.
03610       !
03611       !    Michael Malcolm,
03612       !    Algorithms to Reveal Properties of Floating Point Arithmetic,
03613       !    Communications of the ACM,
03614       !    Volume 15, Number 11, November 1972, pages 949-951.
03615       !
03616       !  Parameters:
03617       !
03618       !    Output, integer ( kind = 4 ) IBETA, the radix for the floating-point
03619       !    representation.
03620       !
03621       !    Output, integer ( kind = 4 ) IT, the number of base IBETA digits in
03622       !    the floating-point significand.
03623       !
03624       !    Output, integer ( kind = 4 ) IRND.
03625       !    * 0, if floating-point addition chops;
03626       !    * 1, if floating-point addition rounds, but not in the IEEE style;
03627       !    * 2, if floating-point addition rounds in the IEEE style;
03628       !    * 3, if floating-point addition chops, and there is partial underflow;
03629       !    * 4, if floating-point addition rounds, but not in the IEEE style,
03630       !      and there is partial underflow
03631       !    * 5, if floating-point addition rounds in the IEEE style, and there
03632       !      is partial underflow
03633       !
03634       !    Output, integer ( kind = 4 ) NGRD, the number of guard digits for
03635       !    multiplication with truncating arithmetic.  It is
03636       !    * 0, if floating-point arithmetic rounds, or if it truncates and only
03637       !      IT base IBETA digits participate in the post-normalization shift of the
03638       !      floating-point significand in multiplication;
03639       !    * 1, if floating-point arithmetic truncates and more than IT base IBETA
03640       !      digits participate in the post-normalization shift of the floating-point
03641       !      significand in multiplication.
03642       !
03643       !    Output, integer ( kind = 4 ) MACHEP. the largest negative integer such that
03644       !    1.0+FLOAT(IBETA)**MACHEP .NE. 1.0, except that MACHEP is bounded below
03645       !    by  -(IT+3).
03646       !
03647       !    Output, integer ( kind = 4 ) NEGEPS. the largest negative integer such that
03648       !    1.0-FLOAT(IBETA)**NEGEPS .NE. 1.0, except that NEGEPS is bounded below
03649       !    by  -(IT+3).
03650       !
03651       !    Output, integer ( kind = 4 ) IEXP, the number of bits (decimal places if
03652       !    IBETA = 10) reserved for the representation of the exponent, including the
03653       !    bias or sign, of a floating-point number.
03654       !
03655       !    Output, integer ( kind = 4 ) MINEXP, largest magnitude negative integer
03656       !    such that FLOAT(IBETA)**MINEXP is positive and normalized.
03657       !
03658       !    Output, integer ( kind = 4 ) MAXEXP, the smallest positive power of BETA
03659       !    that overflows.
03660       !
03661       !    Output, real ( kind = 8 ) EPS, the value FLOAT(IBETA)**MACHEP.
03662       !
03663       !    Output, real ( kind = 8 ) EPSNEG, the value FLOAT(IBETA)**NEGEPS.
03664       !
03665       !    Output, real ( kind = 8 ) XMIN, the smallest non-vanishing normalized
03666       !    floating-point power of the radix, that is, XMIN = FLOAT(IBETA)**MINEXP.
03667       !
03668       !    Output, real ( kind = 8 ) XMAX, the largest finite floating-point number.
03669       !    In particular XMAX = (1.0-EPSNEG)*FLOAT(IBETA)**MAXEXP.
03670       !    On some machines  XMAX  will be only the second, or perhaps third,
03671       !    largest number, being too small by 1 or 2 units in the last digit of
03672       !    the significand.
03673       !
03674       implicit none
03675 
03676       real ( kind = 8 ) a
03677       real ( kind = 8 ) b
03678       real ( kind = 8 ) beta
03679       real ( kind = 8 ) betain
03680       real ( kind = 8 ) betah
03681       real ( kind = 8 ) eps
03682       real ( kind = 8 ) epsneg
03683       integer ( kind = 4 ) i
03684       integer ( kind = 4 ) ibeta
03685       integer ( kind = 4 ) iexp
03686       integer ( kind = 4 ) irnd
03687       integer ( kind = 4 ) it
03688       integer ( kind = 4 ) itemp
03689       integer ( kind = 4 ) iz
03690       integer ( kind = 4 ) j
03691       integer ( kind = 4 ) k
03692       integer ( kind = 4 ) machep
03693       integer ( kind = 4 ) maxexp
03694       integer ( kind = 4 ) minexp
03695       integer ( kind = 4 ) mx
03696       integer ( kind = 4 ) negep
03697       integer ( kind = 4 ) ngrd
03698       integer ( kind = 4 ) nxres
03699       real ( kind = 8 ) one
03700       real ( kind = 8 ) t
03701       real ( kind = 8 ) temp
03702       real ( kind = 8 ) tempa
03703       real ( kind = 8 ) temp1
03704       real ( kind = 8 ) two
03705       real ( kind = 8 ) xmax
03706       real ( kind = 8 ) xmin
03707       real ( kind = 8 ) y
03708       real ( kind = 8 ) z
03709       real ( kind = 8 ) zero
03710 
03711       one = real ( 1, kind = 8 )
03712       two = one + one
03713       zero = one - one
03714       !
03715       !  Determine IBETA, BETA ala Malcolm.
03716       !
03717       a = one
03718 10    continue
03719       a = a + a
03720       temp = a + one
03721       temp1 = temp - a
03722       if ( temp1 - one == zero ) then
03723          go to 10
03724       end if
03725 
03726       b = one
03727 20    continue
03728       b = b + b
03729       temp = a + b
03730       itemp = int ( temp - a )
03731       if ( itemp == 0 ) then
03732          go to 20
03733       end if
03734 
03735       ibeta = itemp
03736       beta = real ( ibeta, kind = 8 )
03737       !
03738       !  Determine IT, IRND.
03739       !
03740       it = 0
03741       b = one
03742 100   continue
03743       it = it + 1
03744       b = b * beta
03745       temp = b + one
03746       temp1 = temp - b
03747       if ( temp1 - one == zero ) then
03748          go to 100
03749       end if
03750 
03751       irnd = 0
03752       betah = beta / two
03753       temp = a + betah
03754 
03755       if ( temp - a /= zero ) then
03756          irnd = 1
03757       end if
03758 
03759       tempa = a + beta
03760       temp = tempa + betah
03761 
03762       if ( irnd == 0 .and. temp - tempa /= zero ) then
03763          irnd = 2
03764       end if
03765       !
03766       !  Determine NEGEP, EPSNEG.
03767       !
03768       negep = it + 3
03769       betain = one / beta
03770       a = one
03771       do i = 1, negep
03772          a = a * betain
03773       end do
03774 
03775       b = a
03776 210   continue
03777       temp = one - a
03778       if ( temp - one == zero ) then
03779          a = a * beta
03780          negep = negep - 1
03781          go to 210
03782       end if
03783 
03784       negep = -negep
03785       epsneg = a
03786       !
03787       !  Determine MACHEP, EPS.
03788       !
03789       machep = -it - 3
03790       a = b
03791 300   continue
03792       temp = one + a
03793       if ( temp - one == zero ) then
03794          a = a * beta
03795          machep = machep + 1
03796          go to 300
03797       end if
03798 
03799       eps = a
03800       !
03801       !  Determine NGRD.
03802       !
03803       ngrd = 0
03804       temp = one + eps
03805 
03806       if ( irnd == 0 .and. temp * one - one /= zero ) then
03807          ngrd = 1
03808       end if
03809       !
03810       !  Determine IEXP, MINEXP, XMIN.
03811       !
03812       !  Loop to determine largest I and K = 2**I such that
03813       !    (1/BETA) ** (2**(I))
03814       !  does not underflow.
03815       !  Exit from loop is signaled by an underflow.
03816       !
03817       i = 0
03818       k = 1
03819       z = betain
03820       t = one + eps
03821       nxres = 0
03822 400   continue
03823       y = z
03824       z = y * y
03825       !
03826       !  Check for underflow here.
03827       !
03828       a = z * one
03829       temp = z * t
03830 
03831       if ( a + a == zero .or. y <= abs ( z ) ) then
03832          go to 410
03833       end if
03834 
03835       temp1 = temp * betain
03836       if ( temp1 * beta == z ) then
03837          go to 410
03838       end if
03839 
03840       i = i + 1
03841       k = k + k
03842       go to 400
03843 
03844 410   continue
03845       if ( ibeta == 10 ) then
03846          go to 420
03847       end if
03848 
03849       iexp = i + 1
03850       mx = k + k
03851       go to 450
03852       !
03853       !  This segment is for decimal machines only.
03854       !
03855 420   continue
03856 
03857       iexp = 2
03858 
03859       iz = ibeta
03860 
03861 430   continue
03862 
03863       if ( iz <= k ) then
03864          iz = iz * ibeta
03865          iexp = iexp + 1
03866          go to 430
03867       end if
03868 
03869       mx = iz + iz - 1
03870       !
03871       !  Loop to determine MINEXP, XMIN.
03872       !  Exit from loop is signaled by an underflow.
03873       !
03874 450   continue
03875 
03876       xmin = y
03877       y = y * betain
03878       !
03879       !  Check for underflow here.
03880       !
03881       a = y * one
03882       temp = y * t
03883       if ( a + a == zero .or. xmin <= abs ( y ) ) then
03884          go to 460
03885       end if
03886 
03887       k = k + 1
03888       temp1 = temp * betain
03889       if ( temp1 * beta /= y .or. temp == y ) then
03890          go to 450
03891       else
03892          nxres = 3
03893          xmin = y
03894       end if
03895 
03896 460   continue
03897 
03898       minexp = -k
03899       !
03900       !  Determine MAXEXP, XMAX.
03901       !
03902       if ( k + k - 3 < mx .or. ibeta == 10 ) then
03903          go to 500
03904       end if
03905 
03906       mx = mx + mx
03907       iexp = iexp + 1
03908 
03909 500   continue
03910 
03911       maxexp = mx + minexp
03912       !
03913       !  Adjust IRND to reflect partial underflow.
03914       !
03915       irnd = irnd + nxres
03916       !
03917       !  Adjust for IEEE-style machines.
03918       !
03919       if ( 2 <= irnd ) then
03920          maxexp = maxexp - 2
03921       end if
03922       !
03923       !  Adjust for machines with implicit leading bit in binary
03924       !  significand, and machines with radix point at extreme
03925       !  right of significand.
03926       !
03927       i = maxexp + minexp
03928 
03929       if ( ibeta == 2 .and. i == 0 ) then
03930          maxexp = maxexp - 1
03931       end if
03932 
03933       if ( 20 < i ) then
03934          maxexp = maxexp - 1
03935       end if
03936 
03937       if ( a /= y ) then
03938          maxexp = maxexp - 2
03939       end if
03940 
03941       xmax = one - epsneg
03942 
03943       if ( xmax * one /= xmax ) then
03944          xmax = one - beta * epsneg
03945       end if
03946 
03947       xmax = xmax / ( beta * beta * beta * xmin )
03948       i = maxexp + minexp + 3
03949 
03950       do j = 1, i
03951          if ( ibeta == 2 ) then
03952             xmax = xmax + xmax
03953          else
03954             xmax = xmax * beta
03955          end if
03956       end do
03957 
03958       return
03959    end subroutine machar
03960    function r8_erf ( x )
03961 
03962       !*****************************************************************************80
03963       !
03964       !! R8_ERF evaluates the error function.
03965       !
03966       !  Discussion:
03967       !
03968       !    This routine computes approximate values for erf(x).
03969       !
03970       !    This routine was renamed from "DERF" to "R8_ERF" to avoid being
03971       !    overshadowed by corresponding functions supplied by some compilers.
03972       !
03973       !    See comments heading CALERF.
03974       !
03975       !  Licensing:
03976       !
03977       !    This code is distributed under the GNU LGPL license.
03978       !
03979       !  Modified:
03980       !
03981       !    23 January 2008
03982       !
03983       !  Author:
03984       !
03985       !    Original FORTRAN77 version by William Cody.
03986       !    FORTRAN90 version by John Burkardt.
03987       !
03988       !  Parameters:
03989       !
03990       !    Input, real ( kind = 8 ) X, the argument of the function.
03991       !
03992       !    Output, real ( kind = 8 ) R8_ERF, the value of the function.
03993       !
03994       implicit none
03995 
03996       integer ( kind = 4 ) jint
03997       real ( kind = 8 ) r8_erf
03998       real ( kind = 8 ) result
03999       real ( kind = 8 ) x
04000 
04001       jint = 0
04002       call calerf ( x, result, jint )
04003       r8_erf = result
04004 
04005       return
04006    end function r8_erf
04007    function r8_erfc ( x )
04008 
04009       !*****************************************************************************80
04010       !
04011       !! R8_ERFC evaluates the complementary error function.
04012       !
04013       !  Discussion:
04014       !
04015       !    This routine computes approximate values for erfc(x).
04016       !
04017       !    This routine was renamed from "DERFC" to "R8_ERFC" to avoid being
04018       !    overshadowed by corresponding functions supplied by some compilers.
04019       !
04020       !    See comments heading CALERF.
04021       !
04022       !  Licensing:
04023       !
04024       !    This code is distributed under the GNU LGPL license.
04025       !
04026       !  Modified:
04027       !
04028       !    23 January 2008
04029       !
04030       !  Author:
04031       !
04032       !    Original FORTRAN77 version by William Cody.
04033       !    FORTRAN90 version by John Burkardt.
04034       !
04035       !  Parameters:
04036       !
04037       !    Input, real ( kind = 8 ) X, the argument of the function.
04038       !
04039       !    Output, real ( kind = 8 ) R8_ERFC, the value of the function.
04040       !
04041       implicit none
04042 
04043       integer ( kind = 4 ) jint
04044       real ( kind = 8 ) r8_erfc
04045       real ( kind = 8 ) result
04046       real ( kind = 8 ) x
04047 
04048       jint = 1
04049       call calerf ( x, result, jint )
04050       r8_erfc = result
04051 
04052       return
04053    end function r8_erfc
04054    function r8_erfcx ( x )
04055 
04056       !*****************************************************************************80
04057       !
04058       !! R8_ERFCX evaluates the exponentially scaled complementary error function.
04059       !
04060       !  Discussion:
04061       !
04062       !    This routine computes approximate values for exp(x*x) * erfc(x).
04063       !
04064       !    This routine was renamed from "DERFCX" to "R8_ERFCX" to match the
04065       !    renamings of DERF and DERFC.
04066       !
04067       !    See comments heading CALERF.
04068       !
04069       !  Licensing:
04070       !
04071       !    This code is distributed under the GNU LGPL license.
04072       !
04073       !  Modified:
04074       !
04075       !    23 January 2008
04076       !
04077       !  Author:
04078       !
04079       !    Original FORTRAN77 version by William Cody.
04080       !    FORTRAN90 version by John Burkardt.
04081       !
04082       !  Parameters:
04083       !
04084       !    Input, real ( kind = 8 ) X, the argument of the function.
04085       !
04086       !    Output, real ( kind = 8 ) R8_ERFCX, the value of the function.
04087       !
04088       implicit none
04089 
04090       integer ( kind = 4 ) jint
04091       real ( kind = 8 ) r8_erfcx
04092       real ( kind = 8 ) result
04093       real ( kind = 8 ) x
04094 
04095       jint = 2
04096       call calerf ( x, result, jint )
04097       r8_erfcx = result
04098 
04099       return
04100    end function r8_erfcx
04101    function r8_gamma ( x )
04102 
04103       !*****************************************************************************80
04104       !
04105       !! R8_GAMMA evaluates Gamma(X) for a real argument.
04106       !
04107       !  Discussion:
04108       !
04109       !    This function was originally named DGAMMA.
04110       !
04111       !    However, a number of Fortran compilers now include a library
04112       !    function of this name.  To avoid conflicts, this function was
04113       !    renamed R8_GAMMA.
04114       !
04115       !    This routine calculates the GAMMA function for a real argument X.
04116       !    Computation is based on an algorithm outlined in reference 1.
04117       !    The program uses rational functions that approximate the GAMMA
04118       !    function to at least 20 significant decimal digits.  Coefficients
04119       !    for the approximation over the interval (1,2) are unpublished.
04120       !    Those for the approximation for 12 <= X are from reference 2.
04121       !
04122       !  Licensing:
04123       !
04124       !    This code is distributed under the GNU LGPL license.
04125       !
04126       !  Modified:
04127       !
04128       !    18 January 2008
04129       !
04130       !  Author:
04131       !
04132       !    Original FORTRAN77 version by William Cody, Laura Stoltz.
04133       !    FORTRAN90 version by John Burkardt.
04134       !
04135       !  Reference:
04136       !
04137       !    William Cody,
04138       !    An Overview of Software Development for Special Functions,
04139       !    in Numerical Analysis Dundee, 1975,
04140       !    edited by GA Watson,
04141       !    Lecture Notes in Mathematics 506,
04142       !    Springer, 1976.
04143       !
04144       !    John Hart, Ward Cheney, Charles Lawson, Hans Maehly,
04145       !    Charles Mesztenyi, John Rice, Henry Thatcher,
04146       !    Christoph Witzgall,
04147       !    Computer Approximations,
04148       !    Wiley, 1968,
04149       !    LC: QA297.C64.
04150       !
04151       !  Parameters:
04152       !
04153       !    Input, real ( kind = 8 ) X, the argument of the function.
04154       !
04155       !    Output, real ( kind = 8 ) R8_GAMMA, the value of the function.
04156       !
04157       implicit none
04158       !
04159       !  Coefficients for minimax approximation over (12, INF).
04160       !
04161       real ( kind = 8 ), dimension ( 7 ) :: c = (/ 
04162            -1.910444077728D-03, 
04163            8.4171387781295D-04, 
04164            -5.952379913043012D-04, 
04165            7.93650793500350248D-04, 
04166            -2.777777777777681622553D-03, 
04167            8.333333333333333331554247D-02, 
04168            5.7083835261D-03 /)
04169       real ( kind = 8 ) eps
04170       real ( kind = 8 ) fact
04171       real ( kind = 8 ) half
04172       integer ( kind = 4 ) i
04173       integer ( kind = 4 ) n
04174       real ( kind = 8 ) one
04175       real ( kind = 8 ) p(8)
04176       logical parity
04177       real ( kind = 8 ) pi
04178       real ( kind = 8 ) q(8)
04179       real ( kind = 8 ) r8_gamma
04180       real ( kind = 8 ) res
04181       real ( kind = 8 ) sqrtpi
04182       real ( kind = 8 ) sum
04183       real ( kind = 8 ) twelve
04184       real ( kind = 8 ) two
04185       real ( kind = 8 ) x
04186       real ( kind = 8 ) xbig
04187       real ( kind = 8 ) xden
04188       real ( kind = 8 ) xinf
04189       real ( kind = 8 ) xminin
04190       real ( kind = 8 ) xnum
04191       real ( kind = 8 ) y
04192       real ( kind = 8 ) y1
04193       real ( kind = 8 ) ysq
04194       real ( kind = 8 ) z
04195       real ( kind = 8 ) zero
04196       !
04197       !  Mathematical constants
04198       !
04199       data one /1.0D+00 /
04200       data half /0.5D+00/
04201       data twelve /12.0D+00/
04202       data two /2.0D+00 /
04203       data zero /0.0D+00/
04204       data sqrtpi /0.9189385332046727417803297D+00/
04205       data pi /3.1415926535897932384626434D+00/
04206       !
04207       !  Machine dependent parameters
04208       !
04209       data xbig / 171.624D+00 /
04210       data xminin / 2.23D-308 /
04211       data eps / 2.22D-16 /
04212       data xinf /1.79D+308/
04213       !
04214       !  Numerator and denominator coefficients for rational minimax
04215       !  approximation over (1,2).
04216       !
04217       data p / -1.71618513886549492533811D+0,   2.47656508055759199108314D+01, &
04218            -3.79804256470945635097577D+02,  6.29331155312818442661052D+02, &
04219            8.66966202790413211295064D+02, -3.14512729688483675254357D+04, &
04220            -3.61444134186911729807069D+04,  6.64561438202405440627855D+04 /
04221 
04222       data q / -3.08402300119738975254353D+01,  3.15350626979604161529144D+02, &
04223            -1.01515636749021914166146D+03, -3.10777167157231109440444D+03, &
04224            2.25381184209801510330112D+04,  4.75584627752788110767815D+03, &
04225            -1.34659959864969306392456D+05, -1.15132259675553483497211D+05 /
04226 
04227       parity = .false.
04228       fact = one
04229       n = 0
04230       y = x
04231       !
04232       !  Argument is negative.
04233       !
04234       if ( y <= zero ) then
04235 
04236          y = - x
04237          y1 = aint ( y )
04238          res = y - y1
04239 
04240          if ( res /= zero ) then
04241 
04242             if ( y1 /= aint ( y1 * half ) * two ) then
04243                parity = .true.
04244             end if
04245 
04246             fact = - pi / sin ( pi * res )
04247             y = y + one
04248 
04249          else
04250 
04251             res = xinf
04252             r8_gamma = res
04253             return
04254 
04255          end if
04256 
04257       end if
04258       !
04259       !  Argument is positive.
04260       !
04261       if ( y < eps ) then
04262          !
04263          !  Argument < EPS.
04264          !
04265          if ( xminin <= y ) then
04266             res = one / y
04267          else
04268             res = xinf
04269             r8_gamma = res
04270             return
04271          end if
04272 
04273       else if ( y < twelve ) then
04274 
04275          y1 = y
04276          !
04277          !  0.0 < argument < 1.0.
04278          !
04279          if ( y < one ) then
04280 
04281             z = y
04282             y = y + one
04283             !
04284             !  1.0 < argument < 12.0.
04285             !  Reduce argument if necessary.
04286             !
04287          else
04288 
04289             n = int ( y ) - 1
04290             y = y - real ( n, kind = 8 )
04291             z = y - one
04292 
04293          end if
04294          !
04295          !  Evaluate approximation for 1.0 < argument < 2.0.
04296          !
04297          xnum = zero
04298          xden = one
04299          do i = 1, 8
04300             xnum = ( xnum + p(i) ) * z
04301             xden = xden * z + q(i)
04302          end do
04303 
04304          res = xnum / xden + one
04305          !
04306          !  Adjust result for case  0.0 < argument < 1.0.
04307          !
04308          if ( y1 < y ) then
04309 
04310             res = res / y1
04311             !
04312             !  Adjust result for case 2.0 < argument < 12.0.
04313             !
04314          else if ( y < y1 ) then
04315 
04316             do i = 1, n
04317                res = res * y
04318                y = y + one
04319             end do
04320 
04321          end if
04322 
04323       else
04324          !
04325          !  Evaluate for 12.0 <= argument.
04326          !
04327          if ( y <= xbig ) then
04328 
04329             ysq = y * y
04330             sum = c(7)
04331             do i = 1, 6
04332                sum = sum / ysq + c(i)
04333             end do
04334             sum = sum / y - y + sqrtpi
04335             sum = sum + ( y - half ) * log ( y )
04336             res = exp ( sum )
04337 
04338          else
04339 
04340             res = huge ( res )
04341             r8_gamma = res
04342             return
04343 
04344          end if
04345 
04346       end if
04347       !
04348       !  Final adjustments and return.
04349       !
04350       if ( parity ) then
04351          res = - res
04352       end if
04353 
04354       if ( fact /= one ) then
04355          res = fact / res
04356       end if
04357 
04358       r8_gamma = res
04359 
04360       return
04361    end function r8_gamma
04362    function r8_psi ( xx )
04363 
04364       !*****************************************************************************80
04365       !
04366       !! R8_PSI evaluates the function Psi(X).
04367       !
04368       !  Discussion:
04369       !
04370       !    This routine evaluates the logarithmic derivative of the
04371       !    Gamma function,
04372       !
04373       !      PSI(X) = d/dX ( GAMMA(X) ) / GAMMA(X)
04374       !             = d/dX LN ( GAMMA(X) )
04375       !
04376       !    for real X, where either
04377       !
04378       !      - XMAX1 < X < - XMIN, and X is not a negative integer,
04379       !
04380       !    or
04381       !
04382       !      XMIN < X.
04383       !
04384       !  Licensing:
04385       !
04386       !    This code is distributed under the GNU LGPL license.
04387       !
04388       !  Modified:
04389       !
04390       !    23 January 2008
04391       !
04392       !  Author:
04393       !
04394       !    Original FORTRAN77 version by William Cody.
04395       !    FORTRAN90 version by John Burkardt.
04396       !
04397       !  Reference:
04398       !
04399       !    William Cody, Anthony Strecok, Henry Thacher,
04400       !    Chebyshev Approximations for the Psi Function,
04401       !    Mathematics of Computation,
04402       !    Volume 27, Number 121, January 1973, pages 123-127.
04403       !
04404       !  Parameters:
04405       !
04406       !    Input, real ( kind = 8 ) XX, the argument of the function.
04407       !
04408       !    Output, real ( kind = 8 ) R8_PSI, the value of the function.
04409       !
04410       implicit none
04411 
04412       real ( kind = 8 ) aug
04413       real ( kind = 8 ) den
04414       real ( kind = 8 ), parameter :: four = 4.0D+00
04415       real ( kind = 8 ), parameter :: fourth = 0.25D+00
04416       real ( kind = 8 ), parameter :: half = 0.5D+00
04417       integer ( kind = 4 ) i
04418       integer ( kind = 4 ) n
04419       integer ( kind = 4 ) nq
04420       real ( kind = 8 ), parameter :: one = 1.0D+00
04421       real ( kind = 8 ), dimension ( 9 ) :: p1 = (/ 
04422            4.5104681245762934160D-03, 
04423            5.4932855833000385356D+00, 
04424            3.7646693175929276856D+02, 
04425            7.9525490849151998065D+03, 
04426            7.1451595818951933210D+04, 
04427            3.0655976301987365674D+05, 
04428            6.3606997788964458797D+05, 
04429            5.8041312783537569993D+05, 
04430            1.6585695029761022321D+05 /)
04431       real ( kind = 8 ), dimension ( 7 ) :: p2 = (/ 
04432            -2.7103228277757834192D+00, 
04433            -1.5166271776896121383D+01, 
04434            -1.9784554148719218667D+01, 
04435            -8.8100958828312219821D+00, 
04436            -1.4479614616899842986D+00, 
04437            -7.3689600332394549911D-02, 
04438            -6.5135387732718171306D-21 /)
04439       real ( kind = 8 ), parameter :: piov4 = 0.78539816339744830962D+00
04440       real ( kind = 8 ), dimension ( 8 ) :: q1 = (/ 
04441            9.6141654774222358525D+01, 
04442            2.6287715790581193330D+03, 
04443            2.9862497022250277920D+04, 
04444            1.6206566091533671639D+05, 
04445            4.3487880712768329037D+05, 
04446            5.4256384537269993733D+05, 
04447            2.4242185002017985252D+05, 
04448            6.4155223783576225996D-08 /)
04449       real ( kind = 8 ), dimension ( 6 ) :: q2 = (/ 
04450            4.4992760373789365846D+01, 
04451            2.0240955312679931159D+02, 
04452            2.4736979003315290057D+02, 
04453            1.0742543875702278326D+02, 
04454            1.7463965060678569906D+01, 
04455            8.8427520398873480342D-01 /)
04456       real ( kind = 8 ) r8_psi
04457       real ( kind = 8 ) sgn
04458       real ( kind = 8 ), parameter :: three = 3.0D+00
04459       real ( kind = 8 ) upper
04460       real ( kind = 8 ) w
04461       real ( kind = 8 ) x
04462       real ( kind = 8 ), parameter :: x01 = 187.0D+00
04463       real ( kind = 8 ), parameter :: x01d = 128.0D+00
04464       real ( kind = 8 ), parameter :: x02 = 6.9464496836234126266D-04
04465       real ( kind = 8 ), parameter :: xinf = 1.70D+38
04466       real ( kind = 8 ), parameter :: xlarge = 2.04D+15
04467       real ( kind = 8 ), parameter :: xmax1 = 3.60D+16
04468       real ( kind = 8 ), parameter :: xmin1 = 5.89D-39
04469       real ( kind = 8 ), parameter :: xsmall = 2.05D-09
04470       real ( kind = 8 ) xx
04471       real ( kind = 8 ) z
04472       real ( kind = 8 ), parameter :: zero = 0.0D+00
04473 
04474       x = xx
04475       w = abs ( x )
04476       aug = zero
04477       !
04478       !  Check for valid arguments, then branch to appropriate algorithm.
04479       !
04480       if ( xmax1 <= - x .or. w < xmin1 ) then
04481 
04482          if ( zero < x ) then
04483             r8_psi = - huge ( x )
04484          else
04485             r8_psi = huge ( x )
04486          end if
04487 
04488          return
04489       end if
04490 
04491       if ( x < half ) then
04492          !
04493          !  X < 0.5, use reflection formula: psi(1-x) = psi(x) + pi * cot(pi*x)
04494          !  Use 1/X for PI*COTAN(PI*X)  when  XMIN1 < |X| <= XSMALL.
04495          !
04496          if ( w <= xsmall ) then
04497 
04498             aug = - one / x
04499             !
04500             !  Argument reduction for cotangent.
04501             !
04502          else
04503 
04504             if ( x < zero ) then
04505                sgn = piov4
04506             else
04507                sgn = - piov4
04508             end if
04509 
04510             w = w - real ( int ( w ), kind = 8 )
04511             nq = int ( w * four )
04512             w = four * ( w - real ( nq, kind = 8 ) * fourth )
04513             !
04514             !  W is now related to the fractional part of 4.0 * X.
04515             !  Adjust argument to correspond to values in the first
04516             !  quadrant and determine the sign.
04517             !
04518             n = nq / 2
04519 
04520             if ( n + n /= nq ) then
04521                w = one - w
04522             end if
04523 
04524             z = piov4 * w
04525 
04526             if ( mod ( n, 2 ) /= 0 ) then
04527                sgn = - sgn
04528             end if
04529             !
04530             !  Determine the final value for  -pi * cotan(pi*x).
04531             !
04532             n = ( nq + 1 ) / 2
04533             if ( mod ( n, 2 ) == 0 ) then
04534                !
04535                !  Check for singularity.
04536                !
04537                if ( z == zero ) then
04538 
04539                   if ( zero < x ) then
04540                      r8_psi = - huge ( x )
04541                   else
04542                      r8_psi = huge ( x )
04543                   end if
04544 
04545                   return
04546                end if
04547 
04548                aug = sgn * ( four / tan ( z ) )
04549 
04550             else
04551 
04552                aug = sgn * ( four * tan ( z ) )
04553 
04554             end if
04555 
04556          end if
04557 
04558          x = one - x
04559 
04560       end if
04561       !
04562       !  0.5 <= X <= 3.0.
04563       !
04564       if ( x <= three ) then
04565 
04566          den = x
04567          upper = p1(1) * x
04568          do i = 1, 7
04569             den = ( den + q1(i) ) * x
04570             upper = ( upper + p1(i+1) ) * x
04571          end do
04572          den = ( upper + p1(9) ) / ( den + q1(8) )
04573          x = ( x - x01 / x01d ) - x02
04574          r8_psi = den * x + aug
04575          return
04576 
04577       end if
04578       !
04579       !  3.0 < X.
04580       !
04581       if ( x < xlarge ) then
04582          w = one / ( x * x )
04583          den = w
04584          upper = p2(1) * w
04585          do i = 1, 5
04586             den = ( den + q2(i) ) * w
04587             upper = ( upper + p2(i+1) ) * w
04588          end do
04589          aug = ( upper + p2(7) ) / ( den + q2(6) ) - half / x + aug
04590       end if
04591 
04592       r8_psi = aug + log ( x )
04593 
04594       return
04595    end function r8_psi
04596    function ren ( k )
04597 
04598       !*****************************************************************************80
04599       !
04600       !! REN is a random number generator.
04601       !
04602       !  Discussion:
04603       !
04604       !    This routine is intended for use on computers with
04605       !    fixed point wordlength of at least 29 bits.  It is
04606       !    best if the floating-point significand has at most
04607       !    29 bits.
04608       !
04609       !  Licensing:
04610       !
04611       !    This code is distributed under the GNU LGPL license.
04612       !
04613       !  Modified:
04614       !
04615       !    03 April 2007
04616       !
04617       !  Author:
04618       !
04619       !    Original FORTRAN77 version by William Cody.
04620       !    FORTRAN90 version by John Burkardt.
04621       !
04622       !  Reference:
04623       !
04624       !    Malcolm Pike, David Hill,
04625       !    Algorithm 266:
04626       !    Pseudo-Random Numbers,
04627       !    Communications of the ACM,
04628       !    Volume 8, Number 10, October 1965, page 605.
04629       !
04630       !  Parameters:
04631       !
04632       !    Input, integer ( kind = 4 ) K, a seed for the random number generator.
04633       !    (Not actually used.)
04634       !
04635       !    Output, real ( kind = 8 ) REN, a pseudorandom value.
04636       !
04637       implicit none
04638 
04639       real ( kind = 8 ), parameter :: c1 = 2796203.0D+00
04640       real ( kind = 8 ), parameter :: c2 = 1.0D-06
04641       real ( kind = 8 ), parameter :: c3 = 1.0D-12
04642       integer ( kind = 4 ) iy
04643       integer ( kind = 4 ) k
04644       real ( kind = 8 ) ren
04645 
04646       save iy
04647 
04648       data iy / 100001 /
04649 
04650       iy = iy * 125
04651       iy = iy - ( iy / 2796203 ) * 2796203
04652       ren = real ( iy, kind = 8 ) / c1 * ( 1.0D+00 + c2 + c3 )
04653 
04654       return
04655    end function ren
04656    subroutine ribesl ( x, alpha, nb, ize, b, ncalc )
04657 
04658       !*****************************************************************************80
04659       !
04660       !! RIBESL calculates I Bessel function with non-integer orders.
04661       !
04662       !  Discussion:
04663       !
04664       !    This routine calculates Bessel functions I SUB(N+ALPHA) (X)
04665       !    for non-negative argument X, and non-negative order N+ALPHA,
04666       !    with or without exponential scaling.
04667       !
04668       !    This program is based on a program written by David
04669       !    Sookne that computes values of the Bessel functions J or
04670       !    I of real argument and integer order.  Modifications include
04671       !    the restriction of the computation to the I Bessel function
04672       !    of non-negative real argument, the extension of the computation
04673       !    to arbitrary positive order, the inclusion of optional
04674       !    exponential scaling, and the elimination of most underflow.
04675       !
04676       !    In case of an error,  NCALC .NE. NB, and not all I's are
04677       !    calculated to the desired accuracy.
04678       !
04679       !    NCALC < 0:  An argument is out of range. For example,
04680       !    NB <= 0, IZE is not 1 or 2, or IZE=1 and ABS(X) .GE. EXPARG.
04681       !    In this case, the B-vector is not calculated, and NCALC is
04682       !    set to MIN0(NB,0)-1 so that NCALC .NE. NB.
04683       !
04684       !    NB .GT. NCALC .GT. 0: Not all requested function values could
04685       !    be calculated accurately.  This usually occurs because NB is
04686       !    much larger than ABS(X).  In this case, B(N) is calculated
04687       !    to the desired accuracy for N <= NCALC, but precision
04688       !    is lost for NCALC < N <= NB.  If B(N) does not vanish
04689       !    for N .GT. NCALC (because it is too small to be represented),
04690       !    and B(N)/B(NCALC) = 10**(-K), then only the first NSIG-K
04691       !    significant figures of B(N) can be trusted.
04692       !
04693       !  Licensing:
04694       !
04695       !    This code is distributed under the GNU LGPL license.
04696       !
04697       !  Modified:
04698       !
04699       !    03 April 2007
04700       !
04701       !  Author:
04702       !
04703       !    Original FORTRAN77 version by William Cody.
04704       !    FORTRAN90 version by John Burkardt.
04705       !
04706       !  Reference:
04707       !
04708       !    Frank Olver, David Sookne,
04709       !    A Note on Backward Recurrence Algorithms,
04710       !    Mathematics of Computation,
04711       !    Volume 26, 1972, pages 941-947.
04712       !
04713       !    David Sookne,
04714       !    Bessel Functions of Real Argument and Integer Order,
04715       !    NBS Journal of Research B,
04716       !    Volume 77B, 1973, pages 125-132.
04717       !
04718       !    William Cody,
04719       !    Algorithm 597:
04720       !    Sequence of Modified Bessel Functions of the First Kind,
04721       !    ACM Transactions of Mathematical Software,
04722       !    Volume 9, Number 2, June 1983, pages 242-245.
04723       !
04724       !  Parameters:
04725       !
04726       !    Input, real ( kind = 8 ) X, the argument for which the functions
04727       !    are to be calculated.
04728       !
04729       !    Input, real ( kind = 8 ) ALPHA,the fractional part of the order
04730       !    for which the functions are to be calculated.
04731       !    0 <= ALPHA < 1.0.
04732       !
04733       !    Input, integer ( kind = 4 ) NB, the number of functions to be calculated.
04734       !    The first function calculated is of order ALPHA, and the
04735       !    last is of order (NB - 1 + ALPHA).  1 <= NB.
04736       !
04737       !    Input, integer IZE, scaling option.
04738       !    1, unscaled I's are to calculated,
04739       !    2, exponentially scaled I's are to be calculated.
04740       !
04741       !    Output, real ( kind = 8 ) B(NB), the values of the functions
04742       !    I(ALPHA,X) through I(NB-1+ALPHA,X), with scaling if requested.
04743       !
04744       !    Output, integer ( kind = 4 ) NCALC, error indicator.  If NCALC = NB, then all the
04745       !    requested values were calculated to the desired accuracy.
04746       !
04747       !  Local Parameeters:
04748       !
04749       !    BETA, radix for the floating-point system
04750       !
04751       !    MINEXP, the smallest representable power of beta
04752       !
04753       !    MAXEXP, the smallest power of BETA that overflows
04754       !
04755       !    IT, the number of bits in the mantissa of a working precision
04756       !    variable.
04757       !
04758       !    NSIG, the decimal significance desired.  Should be set to
04759       !    INT(LOG10(2)*IT+1).  Setting NSIG lower will result
04760       !    in decreased accuracy while setting NSIG higher will
04761       !    increase CPU time without increasing accuracy.  The
04762       !    truncation error is limited to a relative error of
04763       !    T=.5*10**(-NSIG).
04764       !
04765       !    ENTEN = 10.0**K, where K is the largest integer such that
04766       !    ENTEN is machine-representable in working precision
04767       !
04768       !    ENSIG = 10.0**NSIG
04769       !
04770       !    RTNSIG = 10.0**(-K) for the smallest integer K such that
04771       !    NSIG/4 <= K.
04772       !
04773       !    ENMTEN, the smallest ABS(X) such that X/4 does not underflow
04774       !
04775       !    XLARGE, the upper limit on the magnitude of X when IZE=2.  Bear
04776       !    in mind that if ABS(X)=N, then at least N iterations
04777       !    of the backward recursion will be executed.  The value
04778       !    of 10000.0 is used on every machine.
04779       !
04780       !    EXPARG, the largest working precision argument that the library
04781       !    EXP routine can handle and upper limit on the
04782       !    magnitude of X when IZE=1; approximately LOG(BETA**MAXEXP).
04783       !
04784       implicit none
04785 
04786       integer ( kind = 4 ) nb
04787 
04788       real ( kind = 8 ) alpha
04789       real ( kind = 8 ) b(nb)
04790       real ( kind = 8 ) const
04791 !      real ( kind = 8 ) r8_gamma
04792       real ( kind = 8 ) em
04793       real ( kind = 8 ) empal
04794       real ( kind = 8 ) emp2al
04795       real ( kind = 8 ) en
04796       real ( kind = 8 ) enmten
04797       real ( kind = 8 ) ensig
04798       real ( kind = 8 ) enten
04799       real ( kind = 8 ) exparg
04800       real ( kind = 8 ) half
04801       real ( kind = 8 ) halfx
04802       integer ( kind = 4 ) ize
04803       integer ( kind = 4 ) k
04804       integer ( kind = 4 ) l
04805       integer ( kind = 4 ) magx
04806       integer ( kind = 4 ) n
04807       integer ( kind = 4 ) nbmx
04808       integer ( kind = 4 ) ncalc
04809       integer ( kind = 4 ) nend
04810       integer ( kind = 4 ) nsig
04811       integer ( kind = 4 ) nstart
04812       real ( kind = 8 ) one
04813       real ( kind = 8 ) p
04814       real ( kind = 8 ) plast
04815       real ( kind = 8 ) pold
04816       real ( kind = 8 ) psave
04817       real ( kind = 8 ) psavel
04818       real ( kind = 8 ) rtnsig
04819       real ( kind = 8 ) sum
04820       real ( kind = 8 ) tempa
04821       real ( kind = 8 ) tempb
04822       real ( kind = 8 ) tempc
04823       real ( kind = 8 ) test
04824       real ( kind = 8 ) tover
04825       real ( kind = 8 ) two
04826       real ( kind = 8 ) x
04827       real ( kind = 8 ) xlarge
04828       real ( kind = 8 ) zero
04829 
04830       data const / 1.585D+00 /
04831       data enmten / 8.9D-308 /
04832       data ensig / 1.0D+16 /
04833       data enten / 1.0D+308 /
04834       data exparg / 709.0D+00 /
04835       data half / 0.5D+00 /
04836       data nsig / 16 /
04837       data one / 1.0D+00 /
04838       data rtnsig / 1.0D-04 /
04839       data two / 2.0D+00 /
04840       data xlarge / 1.0D+04 /
04841       data zero / 0.0D+00 /
04842       !
04843       !  Check the input.
04844       !
04845       if ( nb <= 0 ) then
04846          ncalc = -1
04847          write ( *, '(a)' ) ' '
04848          write ( *, '(a)' ) 'RIBESL - Fatal error!'
04849          write ( *, '(a)' ) '  NB <= 0.'
04850          return
04851       end if
04852 
04853       if ( x < 0.0D+00 ) then
04854          ncalc = -1
04855          write ( *, '(a)' ) ' '
04856          write ( *, '(a)' ) 'RIBESL - Fatal error!'
04857          write ( *, '(a)' ) '  X < 0.0.'
04858          return
04859       end if
04860 
04861       if ( alpha < 0.0D+00 ) then
04862          ncalc = -1
04863          write ( *, '(a)' ) ' '
04864          write ( *, '(a)' ) 'RIBESL - Fatal error!'
04865          write ( *, '(a)' ) '  ALPHA < 0.'
04866          return
04867       end if
04868 
04869       if ( 1.0D+00 < alpha ) then
04870          ncalc = -1
04871          write ( *, '(a)' ) ' '
04872          write ( *, '(a)' ) 'RIBESL - Fatal error!'
04873          write ( *, '(a)' ) '  1 < ALPHA.'
04874          return
04875       end if
04876       !
04877       !  Check for X, NB, OR IZE out of range.
04878       !
04879       if ( &
04880            ( ize == 1 .and. x <= exparg ) .or. &
04881            ( ize == 2 .and. x <= xlarge ) ) then
04882          !
04883          !  Use 2-term ascending series for small X.
04884          !
04885          ncalc = nb
04886          magx = int ( x )
04887          !
04888          !  Initialize the forward sweep, the P-sequence of Olver.
04889          !
04890          if ( rtnsig <= x ) then
04891 
04892             nbmx = nb - magx
04893             n = magx + 1
04894             en = real ( n + n, kind = 8 ) + ( alpha + alpha )
04895             plast = one
04896             p = en / x
04897             !
04898             !  Calculate general significance test.
04899             !
04900             test = ensig + ensig
04901 
04902             if ( 5 * nsig < 2 * magx ) then
04903                test = sqrt ( test * p )
04904             else
04905                test = test / const**magx
04906             end if
04907 
04908             if ( 3 <= nbmx ) then
04909                !
04910                !  Calculate P-sequence until N = NB-1.  Check for possible overflow.
04911                !
04912                tover = enten / ensig
04913                nstart = magx + 2
04914                nend = nb - 1
04915 
04916                do k = nstart, nend
04917 
04918                   n = k
04919                   en = en + two
04920                   pold = plast
04921                   plast = p
04922                   p = en * plast / x + pold
04923                   !
04924                   !  To avoid overflow, divide P-sequence by TOVER.  Calculate
04925                   !  P-sequence until 1 < ABS(P).
04926                   !
04927                   if ( tover < p ) then
04928 
04929                      tover = enten
04930                      p = p / tover
04931                      plast = plast / tover
04932                      psave = p
04933                      psavel = plast
04934                      nstart = n + 1
04935 
04936 60                   continue
04937 
04938                      n = n + 1
04939                      en = en + two
04940                      pold = plast
04941                      plast = p
04942                      p = en * plast / x + pold
04943 
04944                      if ( p <= one ) then
04945                         go to 60
04946                      end if
04947 
04948                      tempb = en / x
04949                      !
04950                      !  Calculate backward test, and find NCALC, the highest N
04951                      !  such that the test is passed.
04952                      !
04953                      test = pold * plast / ensig
04954                      test = test * ( half - half / ( tempb * tempb ) )
04955                      p = plast * tover
04956                      n = n - 1
04957                      en = en - two
04958                      nend = min0 ( nb, n )
04959 
04960                      do l = nstart, nend
04961 
04962                         ncalc = l
04963                         pold = psavel
04964                         psavel = psave
04965                         psave = en * psavel / x + pold
04966 
04967                         if ( test < psave * psavel ) then
04968                            go to 90
04969                         end if
04970 
04971                      end do
04972 
04973                      ncalc = nend + 1
04974 
04975 90                   continue
04976 
04977                      ncalc = ncalc - 1
04978                      go to 120
04979 
04980                   end if
04981 
04982                end do
04983 
04984                n = nend
04985                en = real ( n + n, kind = 8 ) + ( alpha + alpha )
04986                !
04987                !  Calculate special significance test for 2 < NBMX.
04988                !
04989                test = max ( test, sqrt ( plast * ensig ) * sqrt ( p + p ) )
04990 
04991             end if
04992             !
04993             !  Calculate P-sequence until significance test passed.
04994             !
04995 110         continue
04996 
04997             n = n + 1
04998             en = en + two
04999             pold = plast
05000             plast = p
05001             p = en * plast / x + pold
05002             if ( p < test ) then
05003                go to 110
05004             end if
05005             !
05006             !  Initialize the backward recursion and the normalization sum.
05007             !
05008 120         continue
05009 
05010             n = n + 1
05011             en = en + two
05012             tempb = zero
05013             tempa = one / p
05014             em = real ( n, kind = 8 ) - one
05015             empal = em + alpha
05016             emp2al = ( em - one ) + ( alpha + alpha )
05017             sum = tempa * empal * emp2al / em
05018             nend = n - nb
05019             !
05020             !  N < NB, so store B(N) and set higher orders to zero.
05021             !
05022             if ( nend < 0 ) then
05023 
05024                b(n) = tempa
05025                nend = -nend
05026                do l = 1, nend
05027                   b(n+l) = zero
05028                end do
05029 
05030             else
05031                !
05032                !  Recur backward via difference equation, calculating (but
05033                !  not storing) B(N), until N = NB.
05034                !
05035                if ( 0 < nend ) then
05036 
05037                   do l = 1, nend
05038 
05039                      n = n - 1
05040                      en = en - two
05041                      tempc = tempb
05042                      tempb = tempa
05043                      tempa = ( en * tempb ) / x + tempc
05044                      em = em - one
05045                      emp2al = emp2al - one
05046 
05047                      if ( n == 1 ) then
05048                         go to 150
05049                      end if
05050 
05051                      if ( n == 2 ) then
05052                         emp2al = one
05053                      end if
05054 
05055                      empal = empal - one
05056                      sum = ( sum + tempa * empal ) * emp2al / em
05057 
05058                   end do
05059 
05060                end if
05061                !
05062                !  Store B(NB).
05063                !
05064 150            continue
05065 
05066                b(n) = tempa
05067 
05068                if ( nb <= 1 ) then
05069                   sum = ( sum + sum ) + tempa
05070                   go to 230
05071                end if
05072                !
05073                !  Calculate and Store B(NB-1).
05074                !
05075                n = n - 1
05076                en = en - two
05077                b(n)  = ( en * tempa ) / x + tempb
05078 
05079                if ( n == 1 ) then
05080                   go to 220
05081                end if
05082 
05083                em = em - one
05084                emp2al = emp2al - one
05085 
05086                if ( n == 2 ) then
05087                   emp2al = one
05088                end if
05089 
05090                empal = empal - one
05091                sum = ( sum + b(n) * empal ) * emp2al / em
05092 
05093             end if
05094 
05095             nend = n - 2
05096             !
05097             !  Calculate via difference equation and store B(N), until N = 2.
05098             !
05099             if ( 0 < nend ) then
05100 
05101                do l = 1, nend
05102                   n = n - 1
05103                   en = en - two
05104                   b(n) = ( en * b(n+1) ) / x + b(n+2)
05105                   em = em - one
05106                   emp2al = emp2al - one
05107                   if ( n == 2 ) then
05108                      emp2al = one
05109                   end if
05110                   empal = empal - one
05111                   sum = ( sum + b(n) * empal ) * emp2al / em
05112                end do
05113 
05114             end if
05115             !
05116             !  Calculate B(1).
05117             !
05118             b(1) = two * empal * b(2) / x + b(3)
05119 
05120 220         continue
05121 
05122             sum = ( sum + sum ) + b(1)
05123             !
05124             !  Normalize.  Divide all B(N) by sum.
05125             !
05126 230         continue
05127 
05128             if ( alpha /= zero ) then
05129                sum = sum * r8_gamma ( one + alpha ) * ( x * half )**(-alpha)
05130             end if
05131 
05132             if ( ize == 1 ) then
05133                sum = sum * exp ( -x )
05134             end if
05135 
05136             tempa = enmten
05137 
05138             if ( one < sum ) then
05139                tempa = tempa * sum
05140             end if
05141 
05142             do n = 1, nb
05143                if ( b(n) < tempa ) then
05144                   b(n) = zero
05145                end if
05146                b(n) = b(n) / sum
05147             end do
05148 
05149             return
05150             !
05151             !  Two-term ascending series for small X.
05152             !
05153          else
05154 
05155             tempa = one
05156             empal = one + alpha
05157             halfx = zero
05158 
05159             if ( enmten < x ) then
05160                halfx = half * x
05161             end if
05162 
05163             if ( alpha /= zero ) then
05164                tempa = halfx**alpha / r8_gamma ( empal )
05165             end if
05166 
05167             if ( ize == 2 ) then
05168                tempa = tempa * exp ( -x )
05169             end if
05170 
05171             tempb = zero
05172 
05173             if ( one < x + one ) then
05174                tempb = halfx * halfx
05175             end if
05176 
05177             b(1) = tempa + tempa * tempb / empal
05178 
05179             if ( x /= zero .and. b(1) == zero ) then
05180                ncalc = 0
05181             end if
05182 
05183             if ( 1 < nb ) then
05184 
05185                if ( x == zero ) then
05186 
05187                   do n = 2, nb
05188                      b(n) = zero
05189                   end do
05190 
05191                else
05192                   !
05193                   !  Calculate higher-order functions.
05194                   !
05195                   tempc = halfx
05196                   tover = ( enmten + enmten ) / x
05197 
05198                   if ( tempb /= zero ) then
05199                      tover = enmten / tempb
05200                   end if
05201 
05202                   do n = 2, nb
05203 
05204                      tempa = tempa / empal
05205                      empal = empal + one
05206                      tempa = tempa * tempc
05207 
05208                      if ( tempa <= tover * empal ) then
05209                         tempa = zero
05210                      end if
05211 
05212                      b(n) = tempa + tempa * tempb / empal
05213 
05214                      if ( b(n) == zero .and. n < ncalc ) then
05215                         ncalc = n - 1
05216                      end if
05217 
05218                   end do
05219                end if
05220             end if
05221          end if
05222 
05223       else
05224 
05225          ncalc = min ( nb, 0 ) - 1
05226 
05227       end if
05228 
05229       return
05230    end subroutine ribesl
05231    subroutine rjbesl ( x, alpha, nb, b, ncalc )
05232 
05233       !*****************************************************************************80
05234       !
05235       !! RJBESL calculates J Bessel function with non-integer orders.
05236       !
05237       !  Discussion:
05238       !
05239       !    This routine calculates Bessel functions J sub(N+ALPHA) (X)
05240       !    for non-negative argument X, and non-negative order N+ALPHA.
05241       !
05242       !    This program is based on a program written by David Sookne
05243       !    that computes values of the Bessel functions J or I of real
05244       !    argument and integer order.  Modifications include the restriction
05245       !    of the computation to the J Bessel function of non-negative real
05246       !    argument, the extension of the computation to arbitrary positive
05247       !    order, and the elimination of most underflow.
05248       !
05249       !    In case of an error,  NCALC .NE. NB, and not all J's are
05250       !    calculated to the desired accuracy.
05251       !
05252       !    NCALC < 0:  An argument is out of range. For example,
05253       !    NBES <= 0, ALPHA < 0 or .GT. 1, or X is too large.
05254       !    In this case, B(1) is set to zero, the remainder of the
05255       !    B-vector is not calculated, and NCALC is set to
05256       !    MIN(NB,0)-1 so that NCALC .NE. NB.
05257       !
05258       !    NB .GT. NCALC .GT. 0: Not all requested function values could
05259       !    be calculated accurately.  This usually occurs because NB is
05260       !    much larger than ABS(X).  In this case, B(N) is calculated
05261       !    to the desired accuracy for N <= NCALC, but precision
05262       !    is lost for NCALC < N <= NB.  If B(N) does not vanish
05263       !    for N .GT. NCALC (because it is too small to be represented),
05264       !    and B(N)/B(NCALC) = 10**(-K), then only the first NSIG-K
05265       !    significant figures of B(N) can be trusted.
05266       !
05267       !  Licensing:
05268       !
05269       !    This code is distributed under the GNU LGPL license.
05270       !
05271       !  Modified:
05272       !
05273       !    03 April 2007
05274       !
05275       !  Author:
05276       !
05277       !    Original FORTRAN77 version by William Cody.
05278       !    FORTRAN90 version by John Burkardt.
05279       !
05280       !  Reference:
05281       !
05282       !    Frank Olver, David Sookne,
05283       !    A Note on Backward Recurrence Algorithms,
05284       !    Mathematics of Computation,
05285       !    Volume 26, 1972, pages 941-947.
05286       !
05287       !    David Sookne,
05288       !    Bessel Functions of Real Argument and Integer Order,
05289       !    NBS Journal of Res. B,
05290       !    Volume 77B, 1973, pages 125-132.
05291       !
05292       !  Parameters:
05293       !
05294       !    Input, real ( kind = 8 ) X, the argument for which the
05295       !    J's are to be calculated.
05296       !
05297       !    Input, real ( kind = 8 ) ALPHA, the fractional part of order for which
05298       !    the J's or exponentially scaled J's (J*exp(X)) are to be calculated.
05299       !    0 <= ALPHA < 1.0.
05300       !
05301       !    Input, integer ( kind = 4 ) NB, the number of functions to be calculated.
05302       !    0 < NB.  The first function calculated is of order ALPHA, and the
05303       !    last is of order (NB - 1 + ALPHA).
05304       !
05305       !    Output, real ( kind = 8 ) B(NB).  If RJBESL terminates normally, with
05306       !    NCALC = NB, then B contains the functions J/ALPHA/(X) through
05307       !    J/NB-1+ALPHA/(X), or the corresponding exponentially scaled functions.
05308       !
05309       !    Output, integer ( kind = 4 ) NCALC, error indicator.  If NCALC = NB, then all the
05310       !    requested values were calculated to the desired accuracy.
05311       !
05312       !  Local Parameters:
05313       !
05314       !    IT, the number of bits in the mantissa of a working precision
05315       !    variable.
05316       !
05317       !    NSIG, the decimal significance desired.  Should be set to
05318       !    INT(LOG10(2)*IT+1).  Setting NSIG lower will result
05319       !    in decreased accuracy while setting NSIG higher will
05320       !    increase CPU time without increasing accuracy.  The
05321       !    truncation error is limited to a relative error of
05322       !    T=.5*10**(-NSIG).
05323       !
05324       !    ENTEN = 10.0**K, where K is the largest integer such that
05325       !    ENTEN is machine-representable in working precision
05326       !
05327       !    ENSIG = 10.0**NSIG
05328       !
05329       !    RTNSIG = 10.0**(-K) for the smallest integer K such that
05330       !    K .GE. NSIG/4
05331       !
05332       !    ENMTEN, the smallest ABS(X) such that X/4 does not underflow
05333       !
05334       !    XLARGE, the upper limit on the magnitude of X.  If ABS(X)=N,
05335       !    then at least N iterations of the backward recursion
05336       !    will be executed.  The value of 10000.0 is used on
05337       !    every machine.
05338       !
05339       implicit none
05340 
05341       integer ( kind = 4 ) nb
05342 
05343       real ( kind = 8 ) alpha
05344       real ( kind = 8 ) alpem
05345       real ( kind = 8 ) alp2em
05346       real ( kind = 8 ) b(nb)
05347       real ( kind = 8 ) capp
05348       real ( kind = 8 ) capq
05349 !      real ( kind = 8 ) r8_gamma
05350       real ( kind = 8 ) eighth
05351       real ( kind = 8 ) em
05352       real ( kind = 8 ) en
05353       real ( kind = 8 ) enmten
05354       real ( kind = 8 ) ensig
05355       real ( kind = 8 ) enten
05356       real ( kind = 8 ) fact(25)
05357       real ( kind = 8 ) four
05358       real ( kind = 8 ) gnu
05359       real ( kind = 8 ) half
05360       real ( kind = 8 ) halfx
05361       integer ( kind = 4 ) i
05362       integer ( kind = 4 ) j
05363       integer ( kind = 4 ) k
05364       integer ( kind = 4 ) l
05365       integer ( kind = 4 ) m
05366       integer ( kind = 4 ) magx
05367       integer ( kind = 4 ) n
05368       integer ( kind = 4 ) nbmx
05369       integer ( kind = 4 ) ncalc
05370       integer ( kind = 4 ) nend
05371       integer ( kind = 4 ) nstart
05372       real ( kind = 8 ) one
05373       real ( kind = 8 ) one30
05374       real ( kind = 8 ) p
05375       real ( kind = 8 ) pi2
05376       real ( kind = 8 ) plast
05377       real ( kind = 8 ) pold
05378       real ( kind = 8 ) psave
05379       real ( kind = 8 ) psavel
05380       real ( kind = 8 ) rtnsig
05381       real ( kind = 8 ) s
05382       real ( kind = 8 ) sum
05383       real ( kind = 8 ) t
05384       real ( kind = 8 ) t1
05385       real ( kind = 8 ) tempa
05386       real ( kind = 8 ) tempb
05387       real ( kind = 8 ) tempc
05388       real ( kind = 8 ) test
05389       real ( kind = 8 ) three
05390       real ( kind = 8 ) three5
05391       real ( kind = 8 ) tover
05392       real ( kind = 8 ) two
05393       real ( kind = 8 ) twofiv
05394       real ( kind = 8 ) twopi1
05395       real ( kind = 8 ) twopi2
05396       real ( kind = 8 ) x
05397       real ( kind = 8 ) xc
05398       real ( kind = 8 ) xin
05399       real ( kind = 8 ) xk
05400       real ( kind = 8 ) xlarge
05401       real ( kind = 8 ) xm
05402       real ( kind = 8 ) vcos
05403       real ( kind = 8 ) vsin
05404       real ( kind = 8 ) z
05405       real ( kind = 8 ) zero
05406       !
05407       !  Mathematical constants
05408       !
05409       !   PI2    - 2 / PI
05410       !   TWOPI1 - first few significant digits of 2 * PI
05411       !   TWOPI2 - (2*PI - TWOPI) to working precision, i.e.,
05412       !            TWOPI1 + TWOPI2 = 2 * PI to extra precision.
05413       !
05414       data pi2 / 0.636619772367581343075535d0 /
05415       data twopi1 / 6.28125d0 /
05416       data twopi2 / 1.935307179586476925286767d-3 /
05417       data zero /0.0d0 /
05418       data eighth / 0.125d0 /
05419       data half / 0.5d0 /
05420       data one / 1.0d0/
05421       data two /2.0d0 /
05422       data three / 3.0d0 /
05423       data four / 4.0d0 /
05424       data twofiv /25.0d0/
05425       data one30 /130.0d0 /
05426       data three5 / 35.0d0/
05427       !
05428       !  Machine-dependent parameters
05429       !
05430       data enten /1.0d38 /
05431       data ensig / 1.0d17 /
05432       data rtnsig / 1.0d-4/
05433       data enmten /1.2d-37 /
05434       data xlarge / 1.0d4/
05435       !
05436       !  Factorial(N)
05437       !
05438       data fact / &
05439            1.0d0, &
05440            1.0d0, &
05441            2.0d0, &
05442            6.0d0, &
05443            24.0d0, &
05444            1.2d2, &
05445            7.2d2, &
05446            5.04d3, &
05447            4.032d4, &
05448            3.6288d5,3.6288d6,3.99168d7,4.790016d8,6.2270208d9, &
05449            8.71782912d10,1.307674368d12,2.0922789888d13,3.55687428096d14, &
05450            6.402373705728d15,1.21645100408832d17,2.43290200817664d18, &
05451            5.109094217170944d19,1.12400072777760768d21, &
05452            2.585201673888497664d22, &
05453            6.2044840173323943936d23/
05454       !
05455       !  Check for out of range arguments.
05456       !
05457       magx = int ( x )
05458 
05459       if ( &
05460            0 < nb .and. &
05461            zero <= x .and. &
05462            x <= xlarge .and. &
05463            zero <= alpha .and. &
05464            alpha < one ) then
05465          !
05466          !  Initialize result array to zero.
05467          !
05468          ncalc = nb
05469          do i = 1, nb
05470             b(i) = zero
05471          end do
05472          !
05473          !  Branch to use 2-term ascending series for small X and asymptotic
05474          !  form for large X when NB is not too large.
05475          !
05476          if ( x < rtnsig ) then
05477             !
05478             !  Two-term ascending series for small X.
05479             !
05480             tempa = one
05481             alpem = one + alpha
05482             halfx = zero
05483 
05484             if ( enmten < x ) then
05485                halfx = half * x
05486             end if
05487 
05488             if ( alpha /= zero ) then
05489                tempa = halfx**alpha / ( alpha * r8_gamma ( alpha ) )
05490             end if
05491 
05492             tempb = zero
05493 
05494             if ( one < x + one ) then
05495                tempb = -halfx * halfx
05496             end if
05497 
05498             b(1) = tempa + tempa * tempb / alpem
05499 
05500             if ( x /= zero .and. b(1) == zero ) then
05501                ncalc = 0
05502             end if
05503 
05504             if ( nb /= 1 ) then
05505 
05506                if ( x <= zero ) then
05507 
05508                   do n = 2, nb
05509                      b(n) = zero
05510                   end do
05511                   !
05512                   !  Calculate higher order functions.
05513                   !
05514                else
05515 
05516                   tempc = halfx
05517                   tover = ( enmten + enmten ) / x
05518 
05519                   if ( tempb /= zero ) then
05520                      tover = enmten / tempb
05521                   end if
05522 
05523                   do n = 2, nb
05524 
05525                      tempa = tempa / alpem
05526                      alpem = alpem + one
05527                      tempa = tempa * tempc
05528 
05529                      if ( tempa <= tover * alpem ) then
05530                         tempa = zero
05531                      end if
05532 
05533                      b(n) = tempa + tempa * tempb / alpem
05534 
05535                      if ( b(n) == zero .and. n < ncalc ) then
05536                         ncalc = n - 1
05537                      end if
05538 
05539                   end do
05540 
05541                end if
05542             end if
05543             !
05544             !  Asymptotic series for 21 < X.
05545             !
05546          else if ( twofiv < x .and. nb <= magx + 1 ) then
05547 
05548             xc = sqrt ( pi2 / x )
05549             xin = ( eighth / x )**2
05550             m = 11
05551 
05552             if ( three5 <= x ) then
05553                m = 8
05554             end if
05555 
05556             if ( one30 <= x ) then
05557                m = 4
05558             end if
05559 
05560             xm = four * real ( m, kind = 8 )
05561             !
05562             !  Argument reduction for SIN and COS routines.
05563             !
05564             t = aint ( x / ( twopi1 + twopi2 ) + half )
05565             z = ( ( x - t * twopi1 ) - t * twopi2 ) &
05566                  - ( alpha + half ) / pi2
05567             vsin = sin ( z )
05568             vcos = cos ( z )
05569             gnu = alpha + alpha
05570 
05571             do i = 1, 2
05572 
05573                s = ( ( xm - one ) - gnu ) * ( ( xm - one ) + gnu ) &
05574                     * xin * half
05575                t = ( gnu - ( xm - three ) ) * ( gnu + ( xm - three ) )
05576                capp = s * t / fact(2*m+1)
05577                t1 = ( gnu - ( xm + one ) ) * ( gnu + ( xm + one ) )
05578                capq = s * t1 / fact(2*m+2)
05579                xk = xm
05580                k = m + m
05581                t1 = t
05582 
05583                do j = 2, m
05584                   xk = xk - four
05585                   s = ( ( xk - one ) - gnu ) * ( ( xk - one ) + gnu )
05586                   t = ( gnu - ( xk - three ) ) * ( gnu + ( xk - three ) )
05587                   capp = ( capp + one / fact(k-1) ) * s * t * xin
05588                   capq = ( capq + one / fact(k) ) * s * t1 * xin
05589                   k = k - 2
05590                   t1 = t
05591                end do
05592 
05593                capp = capp + one
05594                capq = ( capq + one ) * ( gnu * gnu - one ) * ( eighth / x )
05595                b(i) = xc * ( capp * vcos - capq * vsin )
05596 
05597                if ( nb == 1 ) then
05598                   return
05599                end if
05600 
05601                t = vsin
05602                vsin = -vcos
05603                vcos = t
05604                gnu = gnu + two
05605 
05606             end do
05607             !
05608             !  If 2 < NB, compute J(X,ORDER+I)  I = 2, NB-1.
05609             !
05610             if ( 2 < nb ) then
05611                gnu = alpha + alpha + two
05612                do j = 3, nb
05613                   b(j) = gnu * b(j-1) / x - b(j-2)
05614                   gnu = gnu + two
05615                end do
05616             end if
05617             !
05618             !  Use recurrence to generate results.  First initialize the
05619             !  calculation of P's.
05620             !
05621          else
05622 
05623             nbmx = nb - magx
05624             n = magx + 1
05625             en = real ( n + n, kind = 8 ) + ( alpha + alpha )
05626             plast = one
05627             p = en / x
05628             !
05629             !  Calculate general significance test.
05630             !
05631             test = ensig + ensig
05632             !
05633             !  Calculate P's until N = NB-1.  Check for possible overflow.
05634             !
05635             if ( 3 <= nbmx ) then
05636 
05637                tover = enten / ensig
05638                nstart = magx + 2
05639                nend = nb - 1
05640                en = real ( nstart + nstart, kind = 8 ) - two + ( alpha + alpha )
05641 
05642                do k = nstart, nend
05643 
05644                   n = k
05645                   en = en + two
05646                   pold = plast
05647                   plast = p
05648                   p = en * plast / x - pold
05649                   !
05650                   !  To avoid overflow, divide P's by TOVER.  Calculate P's until
05651                   !  1 < ABS(P).
05652                   !
05653                   if ( tover < p ) then
05654 
05655                      tover = enten
05656                      p = p / tover
05657                      plast = plast / tover
05658                      psave = p
05659                      psavel = plast
05660                      nstart = n + 1
05661 
05662 100                  continue
05663 
05664                      n = n + 1
05665                      en = en + two
05666                      pold = plast
05667                      plast = p
05668                      p = en * plast / x - pold
05669 
05670                      if ( p <= one ) then
05671                         go to 100
05672                      end if
05673 
05674                      tempb = en / x
05675                      !
05676                      !  Calculate backward test and find NCALC, the highest N such that
05677                      !  the test is passed.
05678                      !
05679                      test = pold * plast &
05680                           * ( half - half / ( tempb * tempb ) )
05681                      test = test / ensig
05682                      p = plast * tover
05683                      n = n - 1
05684                      en = en - two
05685                      nend = min ( nb, n )
05686 
05687                      do l = nstart, nend
05688                         pold = psavel
05689                         psavel = psave
05690                         psave = en * psavel / x - pold
05691                         if ( test < psave * psavel ) then
05692                            ncalc = l - 1
05693                            go to 190
05694                         end if
05695                      end do
05696 
05697                      ncalc = nend
05698                      go to 190
05699 
05700                   end if
05701 
05702                end do
05703 
05704                n = nend
05705                en = real ( n + n, kind = 8 ) + ( alpha + alpha )
05706                !
05707                !  Calculate special significance test for 2 < NBMX.
05708                !
05709                test = max ( test, sqrt ( plast * ensig ) * sqrt ( p + p ) )
05710 
05711             end if
05712             !
05713             !  Calculate P's until significance test passes.
05714             !
05715 140         continue
05716 
05717             n = n + 1
05718             en = en + two
05719             pold = plast
05720             plast = p
05721             p = en * plast / x - pold
05722             if ( p < test ) then
05723                go to 140
05724             end if
05725             !
05726             !  Initialize the backward recursion and the normalization sum.
05727             !
05728 190         continue
05729 
05730             n = n + 1
05731             en = en + two
05732             tempb = zero
05733             tempa = one / p
05734             m = 2 * n - 4 * ( n / 2 )
05735             sum = zero
05736             em = real ( n / 2, kind = 8 )
05737             alpem = ( em - one ) + alpha
05738             alp2em = ( em + em ) + alpha
05739 
05740             if ( m /= 0 ) then
05741                sum = tempa * alpem * alp2em / em
05742             end if
05743 
05744             nend = n - nb
05745             !
05746             !  Recur backward via difference equation, calculating (but not
05747             !  storing) B(N), until N = NB.
05748             !
05749             if ( 0 < nend ) then
05750 
05751                do l = 1, nend
05752 
05753                   n = n - 1
05754                   en = en - two
05755                   tempc = tempb
05756                   tempb = tempa
05757                   tempa = ( en * tempb ) / x - tempc
05758                   m = 2 - m
05759 
05760                   if ( m /= 0 ) then
05761                      em = em - one
05762                      alp2em = ( em + em ) + alpha
05763                      if ( n == 1 ) then
05764                         go to 210
05765                      end if
05766                      alpem = ( em - one ) + alpha
05767                      if ( alpem == zero ) then
05768                         alpem = one
05769                      end if
05770                      sum = ( sum + tempa * alp2em ) * alpem / em
05771                   end if
05772 
05773                end do
05774 
05775             end if
05776             !
05777             !  Store B(NB).
05778             !
05779 210         continue
05780 
05781             b(n) = tempa
05782 
05783             if ( 0 <= nend ) then
05784 
05785                if ( nb <= 1 ) then
05786 
05787                   alp2em = alpha
05788                   if ( alpha + one == one ) then
05789                      alp2em = one
05790                   end if
05791                   sum = sum + b(1) * alp2em
05792                   go to 250
05793 
05794                else
05795                   !
05796                   !  Calculate and store B(NB-1).
05797                   !
05798                   n = n - 1
05799                   en = en - two
05800                   b(n) = ( en * tempa ) / x - tempb
05801 
05802                   if ( n == 1 ) then
05803                      go to 240
05804                   end if
05805 
05806                   m = 2 - m
05807 
05808                   if ( m /= 0 ) then
05809                      em = em - one
05810                      alp2em = ( em + em ) + alpha
05811                      alpem = ( em - one ) + alpha
05812                      if ( alpem == zero ) then
05813                         alpem = one
05814                      end if
05815                      sum = ( sum + b(n) * alp2em ) * alpem / em
05816                   end if
05817 
05818                end if
05819 
05820             end if
05821 
05822             nend = n - 2
05823             !
05824             !  Calculate via difference equation and store B(N), until N = 2.
05825             !
05826             if ( nend /= 0 ) then
05827 
05828                do l = 1, nend
05829                   n = n - 1
05830                   en = en - two
05831                   b(n) = ( en * b(n+1) ) / x - b(n+2)
05832                   m = 2 - m
05833                   if ( m /= 0 ) then
05834                      em = em - one
05835                      alp2em = ( em + em ) + alpha
05836                      alpem = ( em - one ) + alpha
05837                      if ( alpem == zero ) then
05838                         alpem = one
05839                      end if
05840                      sum = ( sum + b(n) * alp2em ) * alpem / em
05841                   end if
05842                end do
05843 
05844             end if
05845             !
05846             !  Calculate B(1).
05847             !
05848             b(1) = two * ( alpha + one ) * b(2) / x - b(3)
05849 
05850 240         continue
05851 
05852             em = em - one
05853             alp2em = ( em + em ) + alpha
05854 
05855             if ( alp2em == zero ) then
05856                alp2em = one
05857             end if
05858 
05859             sum = sum + b(1) * alp2em
05860             !
05861             !  Normalize.  Divide all B(N) by sum.
05862             !
05863 250         continue
05864 
05865             if ( alpha + one /= one ) then
05866                sum = sum * r8_gamma ( alpha ) * ( x * half )**( -alpha )
05867             end if
05868 
05869             tempa = enmten
05870             if ( one < sum ) then
05871                tempa = tempa * sum
05872             end if
05873 
05874             do n = 1, nb
05875                if ( abs ( b(n) ) < tempa ) then
05876                   b(n) = zero
05877                end if
05878                b(n) = b(n) / sum
05879             end do
05880 
05881          end if
05882          !
05883          !  Error return: X, NB, or ALPHA is out of range.
05884          !
05885       else
05886          b(1) = zero
05887          ncalc = min ( nb, 0 ) - 1
05888       end if
05889 
05890       return
05891    end subroutine rjbesl
05892    subroutine rkbesl ( x, alpha, nb, ize, bk, ncalc )
05893 
05894       !*****************************************************************************80
05895       !
05896       !! RKBESL calculates K Bessel function with non-integer orders.
05897       !
05898       !  Discussion:
05899       !
05900       !    This routine calculates modified Bessel functions of the second
05901       !    kind, K SUB(N+ALPHA) (X), for non-negative argument X, and
05902       !    non-negative order N+ALPHA, with or without exponential scaling.
05903       !
05904       !    This program is based on a program written by J. B. Campbell
05905       !    that computes values of the Bessel functions K of real
05906       !    argument and real order.  Modifications include the addition
05907       !    of non-scaled functions, parameterization of machine
05908       !    dependencies, and the use of more accurate approximations
05909       !    for SINH and SIN.
05910       !
05911       !    In case of an error, NCALC .NE. NB, and not all K's are
05912       !    calculated to the desired accuracy.
05913       !
05914       !    NCALC < -1:  An argument is out of range. For example,
05915       !    NB <= 0, IZE is not 1 or 2, or IZE=1 and ABS(X) .GE.
05916       !    XMAX.  In this case, the B-vector is not calculated,
05917       !    and NCALC is set to MIN0(NB,0)-2  so that NCALC .NE. NB.
05918       !
05919       !    NCALC = -1:  Either  K(ALPHA,X) .GE. XINF  or
05920       !    K(ALPHA+NB-1,X)/K(ALPHA+NB-2,X) .GE. XINF.  In this case,
05921       !    the B-vector is not calculated.  Note that again
05922       !    NCALC .NE. NB.
05923       !
05924       !    0 < NCALC < NB: Not all requested function values could
05925       !    be calculated accurately.  BK(I) contains correct function
05926       !    values for I <= NCALC, and contains the ratios
05927       !    K(ALPHA+I-1,X)/K(ALPHA+I-2,X) for the rest of the array.
05928       !
05929       !  Licensing:
05930       !
05931       !    This code is distributed under the GNU LGPL license.
05932       !
05933       !  Modified:
05934       !
05935       !    03 April 2007
05936       !
05937       !  Author:
05938       !
05939       !    Original FORTRAN77 version by William Cody, Laura Stoltz.
05940       !    FORTRAN90 version by John Burkardt.
05941       !
05942       !  Reference:
05943       !
05944       !    JB Campbell,
05945       !    On Temme's Algorithm for the Modified Bessel Functions of the
05946       !    Third Kind,
05947       !    ACM Transactions on Mathematical Software,
05948       !    Volume 6, Number 4, December 1980, pages 581-586.
05949       !
05950       !    JB Campbell,
05951       !    A FORTRAN IV Subroutine for the Modified Bessel Functions of
05952       !    the Third Kind of Real Order and Real Argument,
05953       !    Report NRC/ERB-925,
05954       !    National Research Council, Canada.
05955       !
05956       !  Parameters:
05957       !
05958       !    Input, real ( kind = 8 ) X, the non-negative argument for which
05959       !    K's or exponentially scaled K's (K*EXP(X))
05960       !    are to be calculated.  If K's are to be calculated,
05961       !    X must not be greater than XMAX.
05962       !
05963       !    Input, real ( kind = 8 ) ALPHA, the fractional part of order for which
05964       !    K's or exponentially scaled K's (K*EXP(X)) are to be calculated.
05965       !    0 <= ALPHA < 1.0.
05966       !
05967       !    Input, integer ( kind = 4 ) NB, the number of functions to be calculated, 
05968       !    NB .GT. 0.  The first function calculated is of order ALPHA, and the
05969       !    last is of order (NB - 1 + ALPHA).
05970       !
05971       !    Input, integer ( kind = 4 ) IZE, scaling option.
05972       !    1, unscaled functions are to calculated,
05973       !    2, exponentially scaled functions are to be calculated.
05974       !
05975       !    Output, real ( kind = 8 ) BK(NB), the results.  If the routine
05976       !    terminates normally, with NCALC = NB, the vector BK contains the
05977       !    functions K(ALPHA,X), ... , K(NB-1+ALPHA,X), or the corresponding
05978       !    exponentially scaled functions.
05979       !    If (0 < NCALC < NB), BK(I) contains correct function
05980       !    values for I <= NCALC, and contains the ratios
05981       !    K(ALPHA+I-1,X)/K(ALPHA+I-2,X) for the rest of the array.
05982       !
05983       !    Output, integer ( kind = 4 ) NCALC, error indicator.  If NCALC = NB, then 
05984       !    all the requested values were calculated to the desired accuracy.
05985       !
05986       implicit none
05987 
05988       real ( kind = 8 ) a
05989       real ( kind = 8 ) alpha
05990       real ( kind = 8 ) blpha
05991       real ( kind = 8 ) bk(1)
05992       real ( kind = 8 ) bk1
05993       real ( kind = 8 ) bk2
05994       real ( kind = 8 ) c
05995       real ( kind = 8 ) d
05996       real ( kind = 8 ) dm
05997       real ( kind = 8 ) d1
05998       real ( kind = 8 ) d2
05999       real ( kind = 8 ) d3
06000       real ( kind = 8 ) enu
06001       real ( kind = 8 ) eps
06002       real ( kind = 8 ) estf(7)
06003       real ( kind = 8 ) estm(6)
06004       real ( kind = 8 ) ex
06005       real ( kind = 8 ) four
06006       real ( kind = 8 ) f0
06007       real ( kind = 8 ) f1
06008       real ( kind = 8 ) f2
06009       real ( kind = 8 ) half
06010       integer ( kind = 4 ) i
06011       integer ( kind = 4 ) iend
06012       integer ( kind = 4 ) itemp
06013       integer ( kind = 4 ) ize
06014       integer ( kind = 4 ) j
06015       integer ( kind = 4 ) k
06016       integer ( kind = 4 ) m
06017       integer ( kind = 4 ) mplus1
06018       integer ( kind = 4 ) nb
06019       integer ( kind = 4 ) ncalc
06020       real ( kind = 8 ) one
06021       real ( kind = 8 ) p(8)
06022       real ( kind = 8 ) p0
06023       real ( kind = 8 ) q(7)
06024       real ( kind = 8 ) q0
06025       real ( kind = 8 ) r(5)
06026       real ( kind = 8 ) ratio
06027       real ( kind = 8 ) s(4)
06028       real ( kind = 8 ) sqxmin
06029       real ( kind = 8 ) t(6)
06030       real ( kind = 8 ) tinyx
06031       real ( kind = 8 ) two
06032       real ( kind = 8 ) twonu
06033       real ( kind = 8 ) twox
06034       real ( kind = 8 ) t1
06035       real ( kind = 8 ) t2
06036       real ( kind = 8 ) wminf
06037       real ( kind = 8 ) x
06038       real ( kind = 8 ) xinf
06039       real ( kind = 8 ) xmax
06040       real ( kind = 8 ) xmin
06041       real ( kind = 8 ) x2by4
06042       real ( kind = 8 ) zero
06043       !
06044       !  Mathematical constants
06045       !    A = LOG(2.D0) - Euler's constant
06046       !    D = SQRT(2.D0/PI)
06047       !
06048       data half / 0.5d0 /
06049       data one / 1.0d0 /
06050       data two / 2.0d0 /
06051       data zero /0.0d0/
06052       data four /4.0d0 /
06053       data tinyx / 1.0d-10/
06054       data a / 0.11593151565841244881d0/
06055       data d /0.797884560802865364d0/
06056       !
06057       !  Machine dependent parameters
06058       !
06059       data eps /2.22d-16/
06060       data sqxmin /1.49d-154/
06061       data xinf /1.79d+308/
06062       data xmin /2.23d-308/
06063       data xmax /705.342d0/
06064       !
06065       !  P, Q - Approximation for LOG(GAMMA(1+ALPHA))/ALPHA
06066       !                                         + Euler's constant
06067       !  Coefficients converted from hex to decimal and modified
06068       !  by W. J. Cody, 2/26/82
06069       !
06070       !  R, S - Approximation for (1-ALPHA*PI/SIN(ALPHA*PI))/(2.D0*ALPHA)
06071       !  T    - Approximation for SINH(Y)/Y
06072       !
06073       data p/ 0.805629875690432845d00,    0.204045500205365151d02, &
06074            0.157705605106676174d03,    0.536671116469207504d03, &
06075            0.900382759291288778d03,    0.730923886650660393d03, &
06076            0.229299301509425145d03,    0.822467033424113231d00/
06077       data q/ 0.294601986247850434d02,    0.277577868510221208d03, &
06078            0.120670325591027438d04,    0.276291444159791519d04, &
06079            0.344374050506564618d04,    0.221063190113378647d04, &
06080            0.572267338359892221d03/
06081       data r/-0.48672575865218401848d+0,  0.13079485869097804016d+2, &
06082            -0.10196490580880537526d+3,  0.34765409106507813131d+3, &
06083            0.34958981245219347820d-3/
06084       data s/-0.25579105509976461286d+2,  0.21257260432226544008d+3, &
06085            -0.61069018684944109624d+3,  0.42269668805777760407d+3/
06086       data t/ 0.16125990452916363814d-9, 0.25051878502858255354d-7, &
06087            0.27557319615147964774d-5, 0.19841269840928373686d-3, &
06088            0.83333333333334751799d-2, 0.16666666666666666446d+0/
06089       data estm/5.20583d1, 5.7607d0, 2.7782d0, 1.44303d1, 1.853004d2, &
06090            9.3715d0/
06091       data estf/4.18341d1, 7.1075d0, 6.4306d0, 4.25110d1, 1.35633d0, &
06092            8.45096d1, 2.0d1/
06093 
06094       ex = x
06095       enu = alpha
06096       ncalc = min ( nb, 0 ) - 2
06097 
06098       if ( 0 < nb .and. &
06099            ( zero <= enu .and. enu < one ) .and. &
06100            ( 1 <= ize .and. ize <= 2 ) .and. &
06101            ( ize /= 1 .or. ex <= xmax ) .and. &
06102            zero < ex )  then
06103 
06104          k = 0
06105          if ( enu < sqxmin ) then
06106             enu = zero
06107          end if
06108 
06109          if ( half < enu ) then
06110             k = 1
06111             enu = enu - one
06112          end if
06113 
06114          twonu = enu + enu
06115          iend = nb + k - 1
06116          c = enu * enu
06117          d3 = -c
06118 
06119          if ( ex <= one ) then
06120             !
06121             !  Calculation of P0 = GAMMA(1+ALPHA) * (2/X)**ALPHA,
06122             !                 Q0 = GAMMA(1-ALPHA) * (X/2)**ALPHA.
06123             !
06124             d1 = zero
06125             d2 = p(1)
06126             t1 = one
06127             t2 = q(1)
06128 
06129             do i = 2, 7, 2
06130                d1 = c * d1 + p(i)
06131                d2 = c * d2 + p(i+1)
06132                t1 = c * t1 + q(i)
06133                t2 = c * t2 + q(i+1)
06134             end do
06135 
06136             d1 = enu * d1
06137             t1 = enu * t1
06138             f1 = log ( ex )
06139             f0 = a + enu * ( p(8) &
06140                  - enu * ( d1 + d2 ) / ( t1 + t2 ) ) - f1
06141             q0 = exp ( -enu * ( a - enu * &
06142                  ( p(8) + enu * ( d1 - d2 ) / ( t1 - t2 ) ) - f1 ) )
06143             f1 = enu * f0
06144             p0 = exp ( f1 )
06145             !
06146             !  Calculation of F0.
06147             !
06148             d1 = r(5)
06149             t1 = one
06150             do i = 1, 4
06151                d1 = c * d1 + r(i)
06152                t1 = c * t1 + s(i)
06153             end do
06154 
06155             if ( abs ( f1 ) <= half ) then
06156                f1 = f1 * f1
06157                d2 = zero
06158                do i = 1, 6
06159                   d2 = f1 * d2 + t(i)
06160                end do
06161                d2 = f0 + f0 * f1 * d2
06162             else
06163                d2 = sinh ( f1 ) / enu
06164             end if
06165 
06166             f0 = d2 - enu * d1 / ( t1 * p0 )
06167             !
06168             !  X <= 1.0E-10.
06169             !
06170             !  Calculation of K(ALPHA,X) and X*K(ALPHA+1,X)/K(ALPHA,X).
06171             !
06172             if ( ex <= tinyx ) then
06173 
06174                bk(1) = f0 + ex * f0
06175 
06176                if ( ize == 1 ) then
06177                   bk(1) = bk(1) - ex * bk(1)
06178                end if
06179 
06180                ratio = p0 / f0
06181                c = ex * xinf
06182                !
06183                !  Calculation of K(ALPHA,X) and X*K(ALPHA+1,X)/K(ALPHA,X),
06184                !  1/2 <= ALPHA.
06185                !
06186                if ( k /= 0 ) then
06187 
06188                   ncalc = -1
06189 
06190                   if ( c / ratio <= bk(1) ) then
06191                      return
06192                   end if
06193 
06194                   bk(1) = ratio * bk(1) / ex
06195                   twonu = twonu + two
06196                   ratio = twonu
06197 
06198                end if
06199 
06200                ncalc = 1
06201 
06202                if ( nb == 1 ) then
06203                   return
06204                end if
06205                !
06206                !  Calculate  K(ALPHA+L,X)/K(ALPHA+L-1,X),  L  =  1, 2, ... , NB-1.
06207                !
06208                ncalc = -1
06209                do i = 2, nb
06210                   if ( c <= ratio ) then
06211                      return
06212                   end if
06213                   bk(i) = ratio / ex
06214                   twonu = twonu + two
06215                   ratio = twonu
06216                end do
06217 
06218                ncalc = 1
06219                j = ncalc + 1
06220 
06221                do i = j, nb
06222                   if ( xinf / bk(i) <= bk(ncalc) ) then
06223                      return
06224                   end if
06225                   bk(i) = bk(ncalc) * bk(i)
06226                   ncalc = i
06227                end do
06228 
06229                return
06230                !
06231                !  1.0E-10 < X <= 1.0.
06232                !
06233             else
06234 
06235                c = one
06236                x2by4 = ex * ex / four
06237                p0 = half * p0
06238                q0 = half * q0
06239                d1 = -one
06240                d2 = zero
06241                bk1 = zero
06242                bk2 = zero
06243                f1 = f0
06244                f2 = p0
06245 
06246 100            continue
06247 
06248                d1 = d1 + two
06249                d2 = d2 + one
06250                d3 = d1 + d3
06251                c = x2by4 * c / d2
06252                f0 = ( d2 * f0 + p0 + q0 ) / d3
06253                p0 = p0 / ( d2 - enu )
06254                q0 = q0 / ( d2 + enu )
06255                t1 = c * f0
06256                t2 = c * ( p0 - d2 * f0 )
06257                bk1 = bk1 + t1
06258                bk2 = bk2 + t2
06259 
06260                if ( eps < abs ( t1 / ( f1 + bk1 ) ) .or. &
06261                     eps < abs ( t2 / ( f2 + bk2 ) ) )  then
06262                   go to 100
06263                end if
06264 
06265                bk1 = f1 + bk1
06266                bk2 = two * ( f2 + bk2 ) / ex
06267 
06268                if ( ize == 2 ) then
06269                   d1 = exp ( ex )
06270                   bk1 = bk1 * d1
06271                   bk2 = bk2 * d1
06272                end if
06273 
06274                wminf = estf(1) * ex + estf(2)
06275 
06276             end if
06277             !
06278             !  1/EPS < X.
06279             !
06280          else if ( one < eps * ex ) then
06281 
06282             ncalc = nb
06283             bk1 = one / ( d * sqrt ( ex ) )
06284             do i = 1, nb
06285                bk(i) = bk1
06286             end do
06287 
06288             return
06289 
06290          else
06291             !
06292             !  1 < X.
06293             !
06294             twox = ex + ex
06295             blpha = zero
06296             ratio = zero
06297 
06298             if ( ex <= four ) then
06299                !
06300                !  Calculation of K(ALPHA+1,X)/K(ALPHA,X),  1.0 <= X <= 4.0.
06301                !
06302                d2 = aint ( estm(1) / ex + estm(2) )
06303                m = int ( d2 )
06304                d1 = d2 + d2
06305                d2 = d2 - half
06306                d2 = d2 * d2
06307                do i = 2, m
06308                   d1 = d1 - two
06309                   d2 = d2 - d1
06310                   ratio = ( d3 + d2 ) / ( twox + d1 - ratio )
06311                end do
06312                !
06313                !  Calculation of I(|ALPHA|,X) and I(|ALPHA|+1,X) by backward
06314                !  recurrence and K(ALPHA,X) from the Wronskian.
06315                !
06316                d2 = aint ( estm(3) * ex + estm(4) )
06317                m = int ( d2 )
06318                c = abs ( enu )
06319                d3 = c + c
06320                d1 = d3 - one
06321                f1 = xmin
06322                f0 = ( two * ( c + d2 ) / ex &
06323                     + half * ex / ( c + d2 + one ) ) * xmin
06324 
06325                do i = 3, m
06326                   d2 = d2 - one
06327                   f2 = ( d3 + d2 + d2 ) * f0
06328                   blpha = ( one + d1 / d2 ) * ( f2 + blpha )
06329                   f2 = f2 / ex + f1
06330                   f1 = f0
06331                   f0 = f2
06332                end do
06333 
06334                f1 = ( d3 + two ) * f0 / ex + f1
06335                d1 = zero
06336                t1 = one
06337                do i = 1, 7
06338                   d1 = c * d1 + p(i)
06339                   t1 = c * t1 + q(i)
06340                end do
06341 
06342                p0 = exp ( c * ( a + c * ( p(8) &
06343                     - c * d1 / t1 ) - log ( ex ) ) ) / ex
06344                f2 = ( c + half - ratio ) * f1 / ex
06345                bk1 = p0 + ( d3 * f0 - f2 + f0 + blpha ) &
06346                     / ( f2 + f1 + f0 ) * p0
06347 
06348                if ( ize == 1 ) then
06349                   bk1 = bk1 * exp ( - ex )
06350                end if
06351 
06352                wminf = estf(3) * ex + estf(4)
06353 
06354             else
06355                !
06356                !  Calculation of K(ALPHA,X) and K(ALPHA+1,X)/K(ALPHA,X), by backward
06357                !  recurrence, for 4 < X.
06358                !
06359                dm = aint ( estm(5) / ex + estm(6) )
06360                m = int ( dm )
06361                d2 = dm - half
06362                d2 = d2 * d2
06363                d1 = dm + dm
06364 
06365                do i = 2, m
06366                   dm = dm - one
06367                   d1 = d1 - two
06368                   d2 = d2 - d1
06369                   ratio = ( d3 + d2 ) / ( twox + d1 - ratio )
06370                   blpha = ( ratio + ratio * blpha ) / dm
06371                end do
06372 
06373                bk1 = one / ( ( d + d * blpha ) * sqrt ( ex ) )
06374 
06375                if ( ize == 1 ) then
06376                   bk1 = bk1 * exp ( - ex )
06377                end if
06378 
06379                wminf = estf(5) * ( ex - abs ( ex - estf(7) ) ) + estf(6)
06380 
06381             end if
06382             !
06383             !  Calculation of K(ALPHA+1,X) from K(ALPHA,X) and
06384             !  K(ALPHA+1,X)/K(ALPHA,X).
06385             !
06386             bk2 = bk1 + bk1 * ( enu + half - ratio ) / ex
06387 
06388          end if
06389          !
06390          !  Calculation of 'NCALC', K(ALPHA+I,X), I  =  0, 1, ... , NCALC-1,
06391          !  K(ALPHA+I,X)/K(ALPHA+I-1,X), I  =  NCALC, NCALC+1, ... , NB-1.
06392          !
06393          ncalc = nb
06394          bk(1) = bk1
06395 
06396          if ( iend == 0 ) then
06397             return
06398          end if
06399 
06400          j = 2 - k
06401 
06402          if ( 0 < j ) then
06403             bk(j) = bk2
06404          end if
06405 
06406          if ( iend == 1 ) then
06407             return
06408          end if
06409 
06410          m = min ( int ( wminf - enu ), iend )
06411 
06412          do i = 2, m
06413 
06414             t1 = bk1
06415             bk1 = bk2
06416             twonu = twonu + two
06417 
06418             if ( ex < one ) then
06419 
06420                if ( ( xinf / twonu ) * ex <= bk1 ) then
06421                   exit
06422                end if
06423 
06424             else
06425 
06426                if ( xinf / twonu <= bk1 / ex ) then
06427                   exit
06428                end if
06429 
06430             end if
06431 
06432             bk2 = twonu / ex * bk1 + t1
06433             itemp = i
06434             j = j + 1
06435 
06436             if ( 0 < j ) then
06437                bk(j) = bk2
06438             end if
06439 
06440          end do
06441 
06442          m = itemp
06443 
06444          if ( m == iend ) then
06445             return
06446          end if
06447 
06448          ratio = bk2 / bk1
06449          mplus1 = m + 1
06450          ncalc = -1
06451 
06452          do i = mplus1, iend
06453 
06454             twonu = twonu + two
06455             ratio = twonu / ex + one / ratio
06456             j = j + 1
06457 
06458             if ( 1 < j ) then
06459                bk(j) = ratio
06460             else
06461                if ( xinf / ratio <= bk2 ) then
06462                   return
06463                end if
06464                bk2 = ratio * bk2
06465             end if
06466 
06467          end do
06468 
06469          ncalc = max ( mplus1 - k, 1 )
06470 
06471          if ( ncalc == 1 ) then
06472             bk(1) = bk2
06473          end if
06474 
06475          if ( nb == 1 ) then
06476             return
06477          end if
06478 
06479          j = ncalc + 1
06480 
06481          do i = j, nb
06482             if ( xinf / bk(i) <= bk(ncalc) ) then
06483                return
06484             end if
06485             bk(i) = bk(ncalc) * bk(i)
06486             ncalc = i
06487          end do
06488 
06489       end if
06490 
06491       return
06492    end subroutine rkbesl
06493    subroutine rybesl ( x, alpha, nb, by, ncalc )
06494 
06495       !*****************************************************************************80
06496       !
06497       !! RYBESL calculates Y Bessel function with non-integer orders.
06498       !
06499       !  Discussion:
06500       !
06501       !    This routine calculates Bessel functions Y SUB(N+ALPHA) (X)
06502       !    for non-negative argument X, and non-negative order N+ALPHA.
06503       !
06504       !    This program draws heavily on Temme's Algol program for Y(a,x)
06505       !    and Y(a+1,x) and on Campbell's programs for Y_nu(x).  Temme's
06506       !    scheme is used for x < THRESH, and Campbell's scheme is used
06507       !    in the asymptotic region.  Segments of code from both sources
06508       !    have been translated into Fortran77, merged, and heavily modified.
06509       !    Modifications include parameterization of machine dependencies,
06510       !    use of a new approximation for ln(gamma(x)), and built-in
06511       !    protection against over/underflow.
06512       !
06513       !    In case of an error, NCALC .NE. NB, and not all Y's are
06514       !    calculated to the desired accuracy.
06515       !
06516       !    NCALC < -1:  An argument is out of range. For example,
06517       !    NB <= 0, IZE is not 1 or 2, or IZE=1 and ABS(X) .GE.
06518       !    XMAX.  In this case, BY(1) = 0.0, the remainder of the
06519       !    BY-vector is not calculated, and NCALC is set to
06520       !    MIN0(NB,0)-2  so that NCALC .NE. NB.
06521       !
06522       !    NCALC = -1:  Y(ALPHA,X) .GE. XINF.  The requested function
06523       !    values are set to 0.0.
06524       !
06525       !    1 < NCALC < NB: Not all requested function values could
06526       !    be calculated accurately.  BY(I) contains correct function
06527       !    values for I <= NCALC, and and the remaining NB-NCALC
06528       !    array elements contain 0.0.
06529       !
06530       !  Licensing:
06531       !
06532       !    This code is distributed under the GNU LGPL license.
06533       !
06534       !  Modified:
06535       !
06536       !    03 April 2007
06537       !
06538       !  Author:
06539       !
06540       !    Original FORTRAN77 version by William Cody.
06541       !    FORTRAN90 version by John Burkardt.
06542       !
06543       !  Reference:
06544       !
06545       !    JB Campbell,
06546       !    Bessel functions J_nu(x) and Y_nu(x) of real order and real argument,
06547       !    Computational Physics Communications,
06548       !    Volume 18, 1979, pages 133-142.
06549       !
06550       !    NM Temme,
06551       !    On the numerical evaluation of the ordinary Bessel function
06552       !    of the second kind,
06553       !    Journal of Computational Physics,
06554       !    Volume 21, 1976, pages 343-350.
06555       !
06556       !  Parameters:
06557       !
06558       !    Input, real ( kind = 8 ) X, the argument.  0 <= X.
06559       !
06560       !    Input, real ( kind = 8 ) ALPHA, the fractional part of the order
06561       !    for which the Y's are to be calculated.  0 <= ALPHA < 1.0.
06562       !
06563       !    Input, integer ( kind = 4 ) NB, the number of functions to be calculated, 
06564       !    NB .GT. 0.  The first function calculated is of order ALPHA, and the
06565       !    last is of order (NB - 1 + ALPHA).
06566       !
06567       !    Output, real ( kind = 8 ) BY(NB).  If the routine terminates normally
06568       !    (NCALC=NB), the vector BY contains the functions Y(ALPHA,X) through
06569       !    Y(NB-1+ALPHA,X),  If (0 < NCALC < NB), BY(I) contains correct
06570       !    function values for I <= NCALC, and contains the ratios
06571       !    Y(ALPHA+I-1,X)/Y(ALPHA+I-2,X) for the rest of the array.
06572       !
06573       !    Output, integer ( kind = 4 ) NCALC, error flag.  Before using the vector 
06574       !    BY, the user should check that NCALC=NB, i.e., all orders have been 
06575       !    calculated to the desired accuracy.
06576       !
06577       implicit none
06578 
06579       integer ( kind = 4 ) nb
06580 
06581       real ( kind = 8 ) alfa
06582       real ( kind = 8 ) alpha
06583       real ( kind = 8 ) aye
06584       real ( kind = 8 ) b
06585       real ( kind = 8 ) by(nb)
06586       real ( kind = 8 ) c
06587       real ( kind = 8 ) ch(21)
06588       real ( kind = 8 ) cosmu
06589       real ( kind = 8 ) d
06590       real ( kind = 8 ) del
06591       real ( kind = 8 ) den
06592       real ( kind = 8 ) ddiv
06593       real ( kind = 8 ) div
06594       real ( kind = 8 ) dmu
06595       real ( kind = 8 ) d1
06596       real ( kind = 8 ) d2
06597       real ( kind = 8 ) e
06598       real ( kind = 8 ) eight
06599       real ( kind = 8 ) en
06600       real ( kind = 8 ) enu
06601       real ( kind = 8 ) en1
06602       real ( kind = 8 ) eps
06603       real ( kind = 8 ) even
06604       real ( kind = 8 ) ex
06605       real ( kind = 8 ) f
06606       real ( kind = 8 ) fivpi
06607       real ( kind = 8 ) g
06608       real ( kind = 8 ) gamma
06609       real ( kind = 8 ) h
06610       real ( kind = 8 ) half
06611       integer ( kind = 4 ) i
06612       integer ( kind = 4 ) k
06613       integer ( kind = 4 ) na
06614       integer ( kind = 4 ) ncalc
06615       real ( kind = 8 ) odd
06616       real ( kind = 8 ) onbpi
06617       real ( kind = 8 ) one
06618       real ( kind = 8 ) one5
06619       real ( kind = 8 ) p
06620       real ( kind = 8 ) pa
06621       real ( kind = 8 ) pa1
06622       real ( kind = 8 ) pi
06623       real ( kind = 8 ) piby2
06624       real ( kind = 8 ) pim5
06625       real ( kind = 8 ) q
06626       real ( kind = 8 ) qa
06627       real ( kind = 8 ) qa1
06628       real ( kind = 8 ) q0
06629       real ( kind = 8 ) r
06630       real ( kind = 8 ) s
06631       real ( kind = 8 ) sinmu
06632       real ( kind = 8 ) sq2bpi
06633       real ( kind = 8 ) ten9
06634       real ( kind = 8 ) term
06635       real ( kind = 8 ) three
06636       real ( kind = 8 ) thresh
06637       real ( kind = 8 ) two
06638       real ( kind = 8 ) twobyx
06639       real ( kind = 8 ) x
06640       real ( kind = 8 ) xinf
06641       real ( kind = 8 ) xlarge
06642       real ( kind = 8 ) xmin
06643       real ( kind = 8 ) xna
06644       real ( kind = 8 ) x2
06645       real ( kind = 8 ) ya
06646       real ( kind = 8 ) ya1
06647       real ( kind = 8 ) zero
06648       !
06649       !  Mathematical constants
06650       !    FIVPI = 5*PI
06651       !    PIM5 = 5*PI - 15
06652       !    ONBPI = 1/PI
06653       !    PIBY2 = PI/2
06654       !    SQ2BPI = SQUARE ROOT OF 2/PI
06655       !
06656       data zero / 0.0d0 /
06657       data half / 0.5d0 /
06658       data one / 1.0d0 /
06659       data two / 2.0d0 /
06660       data three / 3.0d0 /
06661       data eight /8.0d0 /
06662       data one5 / 15.0d0 /
06663       data ten9 / 1.9d1/
06664       data fivpi /1.5707963267948966192d1 /
06665       data piby2 / 1.5707963267948966192d0/
06666       data pi /3.1415926535897932385d0 /
06667       data sq2bpi / 7.9788456080286535588d-1/
06668       data pim5 /7.0796326794896619231d-1/
06669       data onbpi / 3.1830988618379067154d-1/
06670       !
06671       !  Machine-dependent constants
06672       !
06673       data del / 1.0d-8 /
06674       data xmin / 4.46d-308 /
06675       data xinf / 1.79d308 /
06676       data eps / 1.11d-16 /
06677       data thresh / 16.0d0 /
06678       data xlarge / 1.0d8 /
06679       !
06680       !  Coefficients for Chebyshev polynomial expansion of
06681       !  1/gamma(1-x), abs(x) <= .5
06682       !
06683       data ch/-0.67735241822398840964d-23,-0.61455180116049879894d-22, &
06684            0.29017595056104745456d-20, 0.13639417919073099464d-18, &
06685            0.23826220476859635824d-17,-0.90642907957550702534d-17, &
06686            -0.14943667065169001769d-14,-0.33919078305362211264d-13, &
06687            -0.17023776642512729175d-12, 0.91609750938768647911d-11, &
06688            0.24230957900482704055d-09, 0.17451364971382984243d-08, &
06689            -0.33126119768180852711d-07,-0.86592079961391259661d-06, &
06690            -0.49717367041957398581d-05, 0.76309597585908126618d-04, &
06691            0.12719271366545622927d-02, 0.17063050710955562222d-02, &
06692            -0.76852840844786673690d-01,-0.28387654227602353814d+00, &
06693            0.92187029365045265648d+00/
06694 
06695       ex = x
06696       enu = alpha
06697 
06698       if ( 0 < nb .and. &
06699            xmin <= x .and. &
06700            ex < xlarge .and. &
06701            zero <= enu .and. &
06702            enu < one )  then
06703 
06704          xna = aint ( enu + half )
06705          na = int ( xna )
06706 
06707          if ( na == 1 ) then
06708             enu = enu - xna
06709          end if
06710 
06711          if ( enu == -half ) then
06712 
06713             p = sq2bpi / sqrt ( ex )
06714             ya = p * sin ( ex )
06715             ya1 = -p * cos ( ex )
06716             !
06717             !  Use Temme's scheme for small X.
06718             !
06719          else if ( ex < three ) then
06720 
06721             b = ex * half
06722             d = - log ( b )
06723             f = enu * d
06724             e = b**( -enu )
06725 
06726             if ( abs ( enu ) < del ) then
06727                c = onbpi
06728             else
06729                c = enu / sin ( enu * pi )
06730             end if
06731             !
06732             !  Computation of sinh(f)/f.
06733             !
06734             if ( abs ( f ) < one ) then
06735                x2 = f * f
06736                en = ten9
06737                s = one
06738                do i = 1, 9
06739                   s = s * x2 / en / ( en - one ) + one
06740                   en = en - two
06741                end do
06742             else
06743                s = ( e - one / e ) * half / f
06744             end if
06745             !
06746             !  Computation of 1/gamma(1-a) using Chebyshev polynomials.
06747             !
06748             x2 = enu * enu * eight
06749             aye = ch(1)
06750             even = zero
06751             alfa = ch(2)
06752             odd = zero
06753 
06754             do i = 3, 19, 2
06755                even = -( aye + aye + even )
06756                aye = -even * x2 - aye + ch(i)
06757                odd = -( alfa + alfa + odd )
06758                alfa = -odd * x2 - alfa + ch(i+1)
06759             end do
06760 
06761             even = ( even * half + aye ) * x2 - aye + ch(21)
06762             odd = ( odd + alfa ) * two
06763             gamma = odd * enu + even
06764             !
06765             !  End of computation of 1/gamma(1-a).
06766             !
06767             g = e * gamma
06768             e = ( e + one / e ) * half
06769             f = two * c * ( odd * e + even * s * d )
06770             e = enu * enu
06771             p = g * c
06772             q = onbpi / g
06773             c = enu * piby2
06774 
06775             if ( abs ( c ) < del ) then
06776                r = one
06777             else
06778                r = sin ( c ) / c
06779             end if
06780 
06781             r = pi * c * r * r
06782             c = one
06783             d = - b * b
06784             h = zero
06785             ya = f + r * q
06786             ya1 = p
06787             en = zero
06788 
06789 100         continue
06790 
06791             en = en + one
06792 
06793             if ( eps < abs ( g / ( one + abs ( ya ) ) ) &
06794                  + abs ( h / ( one + abs ( ya1 ) ) ) ) then
06795                f = ( f * en + p + q ) / ( en * en - e )
06796                c = c * d / en
06797                p = p / ( en - enu )
06798                q = q / ( en + enu )
06799                g = c * ( f + r * q )
06800                h = c * p - en * g
06801                ya = ya + g
06802                ya1 = ya1 + h
06803                go to 100
06804             end if
06805 
06806             ya = -ya
06807             ya1 = -ya1 / b
06808 
06809          else if ( ex < thresh ) then
06810             !
06811             !  Use Temme's scheme for moderate X.
06812             !
06813             c = ( half - enu ) * ( half + enu )
06814             b = ex + ex
06815             e = ( ex * onbpi * cos ( enu * pi ) / eps )
06816             e = e * e
06817             p = one
06818             q = -ex
06819             r = one + ex * ex
06820             s = r
06821             en = two
06822 
06823 200         continue
06824 
06825             if ( r * en * en < e ) then
06826                en1 = en + one
06827                d = ( en - one + c / en ) / s
06828                p = ( en + en - p * d ) / en1
06829                q = ( -b + q * d ) / en1
06830                s = p * p + q * q
06831                r = r * s
06832                en = en1
06833                go to 200
06834             end if
06835 
06836             f = p / s
06837             p = f
06838             g = -q / s
06839             q = g
06840 
06841 220         continue
06842 
06843             en = en - one
06844 
06845             if ( zero < en ) then
06846                r = en1 * ( two - p ) - two
06847                s = b + en1 * q
06848                d = ( en - one + c / en ) / ( r * r + s * s )
06849                p = d * r
06850                q = d * s
06851                e = f + one
06852                f = p * e - g * q
06853                g = q * e + p * g
06854                en1 = en
06855                go to 220
06856             end if
06857 
06858             f = one + f
06859             d = f * f + g * g
06860             pa = f / d
06861             qa = -g / d
06862             d = enu + half -p
06863             q = q + ex
06864             pa1 = ( pa * q - qa * d ) / ex
06865             qa1 = ( qa * q + pa * d ) / ex
06866             b = ex - piby2 * ( enu + half )
06867             c = cos ( b )
06868             s = sin ( b )
06869             d = sq2bpi / sqrt ( ex )
06870             ya = d * ( pa * s + qa * c )
06871             ya1 = d * ( qa1 * s - pa1 * c )
06872          else
06873             !
06874             !  Use Campbell's asymptotic scheme.
06875             !
06876             na = 0
06877             d1 = aint ( ex / fivpi )
06878             i = int ( d1 )
06879             dmu = (( ex - one5 * d1 ) - d1 * pim5 ) &
06880                  - ( alpha + half ) * piby2
06881 
06882             if ( i - 2 * ( i / 2 ) == 0 ) then
06883                cosmu = cos ( dmu )
06884                sinmu = sin ( dmu )
06885             else
06886                cosmu = -cos ( dmu )
06887                sinmu = -sin ( dmu )
06888             end if
06889 
06890             ddiv = eight * ex
06891             dmu = alpha
06892             den = sqrt ( ex )
06893 
06894             do k = 1, 2
06895 
06896                p = cosmu
06897                cosmu = sinmu
06898                sinmu = -p
06899                d1 = ( two * dmu - one ) * ( two * dmu + one )
06900                d2 = zero
06901                div = ddiv
06902                p = zero
06903                q = zero
06904                q0 = d1 / div
06905                term = q0
06906 
06907                do i = 2, 20
06908 
06909                   d2 = d2 + eight
06910                   d1 = d1 - d2
06911                   div = div + ddiv
06912                   term = -term * d1 / div
06913                   p = p + term
06914                   d2 = d2 + eight
06915                   d1 = d1 - d2
06916                   div = div + ddiv
06917                   term = term * d1 / div
06918                   q = q + term
06919 
06920                   if ( abs ( term ) <= eps ) then
06921                      exit
06922                   end if
06923 
06924                end do
06925 
06926                p = p + one
06927                q = q + q0
06928 
06929                if ( k == 1 ) then
06930                   ya = sq2bpi * ( p * cosmu - q * sinmu ) / den
06931                else
06932                   ya1 = sq2bpi * ( p * cosmu - q * sinmu ) / den
06933                end if
06934 
06935                dmu = dmu + one
06936 
06937             end do
06938 
06939          end if
06940 
06941          if ( na == 1 ) then
06942             h = two * ( enu + one ) / ex
06943             if ( one < h ) then
06944                if ( xinf / h < abs ( ya1 ) ) then
06945                   h = zero
06946                   ya = zero
06947                end if
06948             end if
06949             h = h * ya1 - ya
06950             ya = ya1
06951             ya1 = h
06952          end if
06953          !
06954          !  Now have first one or two Y's.
06955          !
06956          by(1) = ya
06957          by(2) = ya1
06958 
06959          if ( ya1 == zero ) then
06960 
06961             ncalc = 1
06962 
06963          else
06964 
06965             aye = one + alpha
06966             twobyx = two / ex
06967             ncalc = 2
06968 
06969             do i = 3, nb
06970 
06971                if ( twobyx < one ) then
06972                   if ( xinf / aye <= abs ( by(i-1) ) * twobyx ) then
06973                      go to 450
06974                   end if
06975                else
06976                   if ( xinf / aye / twobyx <= abs ( by(i-1) ) ) then
06977                      go to 450
06978                   end if
06979                end if
06980 
06981                by(i) = twobyx * aye * by(i-1) - by(i-2)
06982                aye = aye + one
06983                ncalc = ncalc + 1
06984 
06985             end do
06986 
06987          end if
06988 
06989 450      continue
06990 
06991          do i = ncalc + 1, nb
06992             by(i) = zero
06993          end do
06994 
06995       else
06996 
06997          by(1) = zero
06998          ncalc = min ( nb, 0 ) - 1
06999 
07000       end if
07001 
07002       return
07003    end subroutine rybesl
07004    subroutine timestamp ( )
07005 
07006       !*****************************************************************************80
07007       !
07008       !! TIMESTAMP prints the current YMDHMS date as a time stamp.
07009       !
07010       !  Example:
07011       !
07012       !    31 May 2001   9:45:54.872 AM
07013       !
07014       !  Licensing:
07015       !
07016       !    This code is distributed under the GNU LGPL license.
07017       !
07018       !  Modified:
07019       !
07020       !    06 August 2005
07021       !
07022       !  Author:
07023       !
07024       !    John Burkardt
07025       !
07026       !  Parameters:
07027       !
07028       !    None
07029       !
07030       implicit none
07031 
07032       character ( len = 8 ) ampm
07033       integer ( kind = 4 ) d
07034       integer ( kind = 4 ) h
07035       integer ( kind = 4 ) m
07036       integer ( kind = 4 ) mm
07037       character ( len = 9 ), parameter, dimension(12) :: month = (/ 
07038            'January  ', 'February ', 'March    ', 'April    ', 
07039            'May      ', 'June     ', 'July     ', 'August   ', 
07040            'September', 'October  ', 'November ', 'December ' /)
07041       integer ( kind = 4 ) n
07042       integer ( kind = 4 ) s
07043       integer ( kind = 4 ) values(8)
07044       integer ( kind = 4 ) y
07045 
07046       call date_and_time ( values = values )
07047 
07048       y = values(1)
07049       m = values(2)
07050       d = values(3)
07051       h = values(5)
07052       n = values(6)
07053       s = values(7)
07054       mm = values(8)
07055 
07056       if ( h < 12 ) then
07057          ampm = 'AM'
07058       else if ( h == 12 ) then
07059          if ( n == 0 .and. s == 0 ) then
07060             ampm = 'Noon'
07061          else
07062             ampm = 'PM'
07063          end if
07064       else
07065          h = h - 12
07066          if ( h < 12 ) then
07067             ampm = 'PM'
07068          else if ( h == 12 ) then
07069             if ( n == 0 .and. s == 0 ) then
07070                ampm = 'Midnight'
07071             else
07072                ampm = 'AM'
07073             end if
07074          end if
07075       end if
07076 
07077       write ( *, '(i2,1x,a,1x,i4,2x,i2,a1,i2.2,a1,i2.2,a1,i3.3,1x,a)' ) &
07078            d, trim ( month(m) ), y, h, ':', n, ':', s, '.', mm, trim ( ampm )
07079 
07080       return
07081    end subroutine timestamp
07082 
07083 END MODULE SpecFun
 All Classes Files Functions Variables