hstcsp.f 16 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446
  1. *DECK HSTCSP
  2. SUBROUTINE HSTCSP (INTL, A, B, M, MBDCND, BDA, BDB, C, D, N,
  3. + NBDCND, BDC, BDD, ELMBDA, F, IDIMF, PERTRB, IERROR, W)
  4. C***BEGIN PROLOGUE HSTCSP
  5. C***PURPOSE Solve the standard five-point finite difference
  6. C approximation on a staggered grid to the modified Helmholtz
  7. C equation in spherical coordinates assuming axisymmetry
  8. C (no dependence on longitude).
  9. C***LIBRARY SLATEC (FISHPACK)
  10. C***CATEGORY I2B1A1A
  11. C***TYPE SINGLE PRECISION (HSTCSP-S)
  12. C***KEYWORDS ELLIPTIC, FISHPACK, HELMHOLTZ, PDE, SPHERICAL
  13. C***AUTHOR Adams, J., (NCAR)
  14. C Swarztrauber, P. N., (NCAR)
  15. C Sweet, R., (NCAR)
  16. C***DESCRIPTION
  17. C
  18. C HSTCSP solves the standard five-point finite difference
  19. C approximation on a staggered grid to the modified Helmholtz
  20. C equation spherical coordinates assuming axisymmetry (no dependence
  21. C on longitude).
  22. C
  23. C (1/R**2)(d/dR)(R**2(dU/dR)) +
  24. C
  25. C 1/(R**2*SIN(THETA))(d/dTHETA)(SIN(THETA)(dU/dTHETA)) +
  26. C
  27. C (LAMBDA/(R*SIN(THETA))**2)U = F(THETA,R)
  28. C
  29. C where THETA is colatitude and R is the radial coordinate.
  30. C This two-dimensional modified Helmholtz equation results from
  31. C the Fourier transform of the three-dimensional Poisson equation.
  32. C
  33. C * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
  34. C
  35. C
  36. C * * * * * * * * Parameter Description * * * * * * * * * *
  37. C
  38. C
  39. C * * * * * * On Input * * * * * *
  40. C
  41. C INTL
  42. C = 0 On initial entry to HSTCSP or if any of the arguments
  43. C C, D, N, or NBDCND are changed from a previous call.
  44. C
  45. C = 1 If C, D, N, and NBDCND are all unchanged from previous
  46. C call to HSTCSP.
  47. C
  48. C NOTE: A call with INTL = 0 takes approximately 1.5 times as much
  49. C time as a call with INTL = 1. Once a call with INTL = 0
  50. C has been made then subsequent solutions corresponding to
  51. C different F, BDA, BDB, BDC, and BDD can be obtained
  52. C faster with INTL = 1 since initialization is not repeated.
  53. C
  54. C A,B
  55. C The range of THETA (colatitude), i.e. A .LE. THETA .LE. B. A
  56. C must be less than B and A must be non-negative. A and B are in
  57. C radians. A = 0 corresponds to the north pole and B = PI
  58. C corresponds to the south pole.
  59. C
  60. C * * * IMPORTANT * * *
  61. C
  62. C If B is equal to PI, then B must be computed using the statement
  63. C
  64. C B = PIMACH(DUM)
  65. C
  66. C This insures that B in the user's program is equal to PI in this
  67. C program which permits several tests of the input parameters that
  68. C otherwise would not be possible.
  69. C
  70. C * * * * * * * * * * * *
  71. C
  72. C M
  73. C The number of grid points in the interval (A,B). The grid points
  74. C in the THETA-direction are given by THETA(I) = A + (I-0.5)DTHETA
  75. C for I=1,2,...,M where DTHETA =(B-A)/M. M must be greater than 4.
  76. C
  77. C MBDCND
  78. C Indicates the type of boundary conditions at THETA = A and
  79. C THETA = B.
  80. C
  81. C = 1 If the solution is specified at THETA = A and THETA = B.
  82. C (See notes 1, 2 below)
  83. C
  84. C = 2 If the solution is specified at THETA = A and the derivative
  85. C of the solution with respect to THETA is specified at
  86. C THETA = B (See notes 1, 2 below).
  87. C
  88. C = 3 If the derivative of the solution with respect to THETA is
  89. C specified at THETA = A (See notes 1, 2 below) and THETA = B.
  90. C
  91. C = 4 If the derivative of the solution with respect to THETA is
  92. C specified at THETA = A (See notes 1, 2 below) and the
  93. C solution is specified at THETA = B.
  94. C
  95. C = 5 If the solution is unspecified at THETA = A = 0 and the
  96. C solution is specified at THETA = B. (See note 2 below)
  97. C
  98. C = 6 If the solution is unspecified at THETA = A = 0 and the
  99. C derivative of the solution with respect to THETA is
  100. C specified at THETA = B (See note 2 below).
  101. C
  102. C = 7 If the solution is specified at THETA = A and the
  103. C solution is unspecified at THETA = B = PI.
  104. C
  105. C = 8 If the derivative of the solution with respect to
  106. C THETA is specified at THETA = A (See note 1 below)
  107. C and the solution is unspecified at THETA = B = PI.
  108. C
  109. C = 9 If the solution is unspecified at THETA = A = 0 and
  110. C THETA = B = PI.
  111. C
  112. C NOTES: 1. If A = 0, do not use MBDCND = 1,2,3,4,7 or 8,
  113. C but instead use MBDCND = 5, 6, or 9.
  114. C
  115. C 2. if B = PI, do not use MBDCND = 1,2,3,4,5 or 6,
  116. C but instead use MBDCND = 7, 8, or 9.
  117. C
  118. C When A = 0 and/or B = PI the only meaningful boundary
  119. C condition is dU/dTHETA = 0. (See D. Greenspan, 'Numerical
  120. C Analysis of Elliptic Boundary Value Problems,' Harper and
  121. C Row, 1965, Chapter 5.)
  122. C
  123. C BDA
  124. C A one-dimensional array of length N that specifies the boundary
  125. C values (if any) of the solution at THETA = A. When
  126. C MBDCND = 1, 2, or 7,
  127. C
  128. C BDA(J) = U(A,R(J)) , J=1,2,...,N.
  129. C
  130. C When MBDCND = 3, 4, or 8,
  131. C
  132. C BDA(J) = (d/dTHETA)U(A,R(J)) , J=1,2,...,N.
  133. C
  134. C When MBDCND has any other value, BDA is a dummy variable.
  135. C
  136. C BDB
  137. C A one-dimensional array of length N that specifies the boundary
  138. C values of the solution at THETA = B. When MBDCND = 1, 4, or 5,
  139. C
  140. C BDB(J) = U(B,R(J)) , J=1,2,...,N.
  141. C
  142. C When MBDCND = 2,3, or 6,
  143. C
  144. C BDB(J) = (d/dTHETA)U(B,R(J)) , J=1,2,...,N.
  145. C
  146. C When MBDCND has any other value, BDB is a dummy variable.
  147. C
  148. C C,D
  149. C The range of R , i.e. C .LE. R .LE. D.
  150. C C must be less than D. C must be non-negative.
  151. C
  152. C N
  153. C The number of unknowns in the interval (C,D). The unknowns in
  154. C the R-direction are given by R(J) = C + (J-0.5)DR,
  155. C J=1,2,...,N, where DR = (D-C)/N. N must be greater than 4.
  156. C
  157. C NBDCND
  158. C Indicates the type of boundary conditions at R = C
  159. C and R = D.
  160. C
  161. C = 1 If the solution is specified at R = C and R = D.
  162. C
  163. C = 2 If the solution is specified at R = C and the derivative
  164. C of the solution with respect to R is specified at
  165. C R = D. (See note 1 below)
  166. C
  167. C = 3 If the derivative of the solution with respect to R is
  168. C specified at R = C and R = D.
  169. C
  170. C = 4 If the derivative of the solution with respect to R is
  171. C specified at R = C and the solution is specified at
  172. C R = D.
  173. C
  174. C = 5 If the solution is unspecified at R = C = 0 (See note 2
  175. C below) and the solution is specified at R = D.
  176. C
  177. C = 6 If the solution is unspecified at R = C = 0 (See note 2
  178. C below) and the derivative of the solution with respect to R
  179. C is specified at R = D.
  180. C
  181. C NOTE 1: If C = 0 and MBDCND = 3,6,8 or 9, the system of equations
  182. C to be solved is singular. The unique solution is
  183. C determined by extrapolation to the specification of
  184. C U(THETA(1),C). But in these cases the right side of the
  185. C system will be perturbed by the constant PERTRB.
  186. C
  187. C NOTE 2: NBDCND = 5 or 6 cannot be used with MBDCND = 1, 2, 4, 5,
  188. C or 7 (the former indicates that the solution is
  189. C unspecified at R = 0; the latter indicates that the
  190. C solution is specified). Use instead NBDCND = 1 or 2.
  191. C
  192. C BDC
  193. C A one dimensional array of length M that specifies the boundary
  194. C values of the solution at R = C. When NBDCND = 1 or 2,
  195. C
  196. C BDC(I) = U(THETA(I),C) , I=1,2,...,M.
  197. C
  198. C When NBDCND = 3 or 4,
  199. C
  200. C BDC(I) = (d/dR)U(THETA(I),C), I=1,2,...,M.
  201. C
  202. C When NBDCND has any other value, BDC is a dummy variable.
  203. C
  204. C BDD
  205. C A one-dimensional array of length M that specifies the boundary
  206. C values of the solution at R = D. When NBDCND = 1 or 4,
  207. C
  208. C BDD(I) = U(THETA(I),D) , I=1,2,...,M.
  209. C
  210. C When NBDCND = 2 or 3,
  211. C
  212. C BDD(I) = (d/dR)U(THETA(I),D) , I=1,2,...,M.
  213. C
  214. C When NBDCND has any other value, BDD is a dummy variable.
  215. C
  216. C ELMBDA
  217. C The constant LAMBDA in the modified Helmholtz equation. If
  218. C LAMBDA is greater than 0, a solution may not exist. However,
  219. C HSTCSP will attempt to find a solution.
  220. C
  221. C F
  222. C A two-dimensional array that specifies the values of the right
  223. C side of the modified Helmholtz equation. For I=1,2,...,M and
  224. C J=1,2,...,N
  225. C
  226. C F(I,J) = F(THETA(I),R(J)) .
  227. C
  228. C F must be dimensioned at least M X N.
  229. C
  230. C IDIMF
  231. C The row (or first) dimension of the array F as it appears in the
  232. C program calling HSTCSP. This parameter is used to specify the
  233. C variable dimension of F. IDIMF must be at least M.
  234. C
  235. C W
  236. C A one-dimensional array that must be provided by the user for
  237. C work space. With K = INT(log2(N))+1 and L = 2**(K+1), W may
  238. C require up to (K-2)*L+K+MAX(2N,6M)+4(N+M)+5 locations. The
  239. C actual number of locations used is computed by HSTCSP and is
  240. C returned in the location W(1).
  241. C
  242. C
  243. C * * * * * * On Output * * * * * *
  244. C
  245. C F
  246. C Contains the solution U(I,J) of the finite difference
  247. C approximation for the grid point (THETA(I),R(J)) for
  248. C I=1,2,...,M, J=1,2,...,N.
  249. C
  250. C PERTRB
  251. C If a combination of periodic, derivative, or unspecified
  252. C boundary conditions is specified for a Poisson equation
  253. C (LAMBDA = 0), a solution may not exist. PERTRB is a con-
  254. C stant, calculated and subtracted from F, which ensures
  255. C that a solution exists. HSTCSP then computes this
  256. C solution, which is a least squares solution to the
  257. C original approximation. This solution plus any constant is also
  258. C a solution; hence, the solution is not unique. The value of
  259. C PERTRB should be small compared to the right side F.
  260. C Otherwise, a solution is obtained to an essentially different
  261. C problem. This comparison should always be made to insure that
  262. C a meaningful solution has been obtained.
  263. C
  264. C IERROR
  265. C An error flag that indicates invalid input parameters.
  266. C Except for numbers 0 and 10, a solution is not attempted.
  267. C
  268. C = 0 No error
  269. C
  270. C = 1 A .LT. 0 or B .GT. PI
  271. C
  272. C = 2 A .GE. B
  273. C
  274. C = 3 MBDCND .LT. 1 or MBDCND .GT. 9
  275. C
  276. C = 4 C .LT. 0
  277. C
  278. C = 5 C .GE. D
  279. C
  280. C = 6 NBDCND .LT. 1 or NBDCND .GT. 6
  281. C
  282. C = 7 N .LT. 5
  283. C
  284. C = 8 NBDCND = 5 or 6 and MBDCND = 1, 2, 4, 5, or 7
  285. C
  286. C = 9 C .GT. 0 and NBDCND .GE. 5
  287. C
  288. C = 10 ELMBDA .GT. 0
  289. C
  290. C = 11 IDIMF .LT. M
  291. C
  292. C = 12 M .LT. 5
  293. C
  294. C = 13 A = 0 and MBDCND =1,2,3,4,7 or 8
  295. C
  296. C = 14 B = PI and MBDCND .LE. 6
  297. C
  298. C = 15 A .GT. 0 and MBDCND = 5, 6, or 9
  299. C
  300. C = 16 B .LT. PI and MBDCND .GE. 7
  301. C
  302. C = 17 LAMBDA .NE. 0 and NBDCND .GE. 5
  303. C
  304. C Since this is the only means of indicating a possibly
  305. C incorrect call to HSTCSP, the user should test IERROR after
  306. C the call.
  307. C
  308. C W
  309. C W(1) contains the required length of W. Also W contains
  310. C intermediate values that must not be destroyed if HSTCSP
  311. C will be called again with INTL = 1.
  312. C
  313. C *Long Description:
  314. C
  315. C * * * * * * * Program Specifications * * * * * * * * * * * *
  316. C
  317. C Dimension of BDA(N),BDB(N),BDC(M),BDD(M),F(IDIMF,N),
  318. C Arguments W(See argument list)
  319. C
  320. C Latest June 1979
  321. C Revision
  322. C
  323. C Subprograms HSTCSP,HSTCS1,BLKTRI,BLKTR1,INDXA,INDXB,INDXC,
  324. C Required PROD,PRODP,CPROD,CPRODP,PPADD,PSGF,BSRH,PPSGF,
  325. C PPSPF,COMPB,TEVLS,R1MACH
  326. C
  327. C Special NONE
  328. C Conditions
  329. C
  330. C Common CBLKT
  331. C Blocks
  332. C
  333. C I/O NONE
  334. C
  335. C Precision Single
  336. C
  337. C Specialist Roland Sweet
  338. C
  339. C Language FORTRAN
  340. C
  341. C History Written by Roland Sweet at NCAR in May, 1977
  342. C
  343. C Algorithm This subroutine defines the finite-difference
  344. C equations, incorporates boundary data, adjusts the
  345. C right side when the system is singular and calls
  346. C BLKTRI which solves the linear system of equations.
  347. C
  348. C Space 5269(decimal) = 12225(octal) locations on the
  349. C Required NCAR Control Data 7600
  350. C
  351. C Timing and The execution time T on the NCAR Control Data
  352. C Accuracy 7600 for subroutine HSTCSP is roughly proportional
  353. C to M*N*log2(N), but depends on the input parameter
  354. C INTL. Some values are listed in the table below.
  355. C The solution process employed results in a loss
  356. C of no more than FOUR significant digits for N and M
  357. C as large as 64. More detailed information about
  358. C accuracy can be found in the documentation for
  359. C subroutine BLKTRI which is the routine that
  360. C actually solves the finite difference equations.
  361. C
  362. C
  363. C M(=N) INTL MBDCND(=NBDCND) T(MSECS)
  364. C ----- ---- --------------- --------
  365. C
  366. C 32 0 1-6 132
  367. C 32 1 1-6 88
  368. C 64 0 1-6 546
  369. C 64 1 1-6 380
  370. C
  371. C Portability American National Standards Institute Fortran.
  372. C The machine accuracy is set using function R1MACH.
  373. C
  374. C Required COS,SIN,ABS,SQRT
  375. C Resident
  376. C Routines
  377. C
  378. C Reference Swarztrauber, P.N., 'A Direct Method For The
  379. C Discrete Solution Of Separable Elliptic Equations,'
  380. C SIAM J. Numer. Anal. 11(1974), pp. 1136-1150.
  381. C
  382. C * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
  383. C
  384. C***REFERENCES P. N. Swarztrauber and R. Sweet, Efficient Fortran
  385. C subprograms for the solution of elliptic equations,
  386. C NCAR TN/IA-109, July 1975, 138 pp.
  387. C P. N. Swarztrauber, A direct method for the discrete
  388. C solution of separable elliptic equations, SIAM Journal
  389. C on Numerical Analysis 11, (1974), pp. 1136-1150.
  390. C***ROUTINES CALLED HSTCS1, PIMACH
  391. C***REVISION HISTORY (YYMMDD)
  392. C 801001 DATE WRITTEN
  393. C 890531 Changed all specific intrinsics to generic. (WRB)
  394. C 890531 REVISION DATE from Version 3.2
  395. C 891214 Prologue converted to Version 4.0 format. (BAB)
  396. C 920501 Reformatted the REFERENCES section. (WRB)
  397. C***END PROLOGUE HSTCSP
  398. C
  399. C
  400. DIMENSION F(IDIMF,*) ,BDA(*) ,BDB(*) ,BDC(*) ,
  401. 1 BDD(*) ,W(*)
  402. C***FIRST EXECUTABLE STATEMENT HSTCSP
  403. PI = PIMACH(DUM)
  404. C
  405. C CHECK FOR INVALID INPUT PARAMETERS
  406. C
  407. IERROR = 0
  408. IF (A.LT.0. .OR. B.GT.PI) IERROR = 1
  409. IF (A .GE. B) IERROR = 2
  410. IF (MBDCND.LT.1 .OR. MBDCND.GT.9) IERROR = 3
  411. IF (C .LT. 0.) IERROR = 4
  412. IF (C .GE. D) IERROR = 5
  413. IF (NBDCND.LT.1 .OR. NBDCND.GT.6) IERROR = 6
  414. IF (N .LT. 5) IERROR = 7
  415. IF ((NBDCND.EQ.5 .OR. NBDCND.EQ.6) .AND. (MBDCND.EQ.1 .OR.
  416. 1 MBDCND.EQ.2 .OR. MBDCND.EQ.4 .OR. MBDCND.EQ.5 .OR.
  417. 2 MBDCND.EQ.7))
  418. 3 IERROR = 8
  419. IF (C.GT.0. .AND. NBDCND.GE.5) IERROR = 9
  420. IF (IDIMF .LT. M) IERROR = 11
  421. IF (M .LT. 5) IERROR = 12
  422. IF (A.EQ.0. .AND. MBDCND.NE.5 .AND. MBDCND.NE.6 .AND. MBDCND.NE.9)
  423. 1 IERROR = 13
  424. IF (B.EQ.PI .AND. MBDCND.LE.6) IERROR = 14
  425. IF (A.GT.0. .AND. (MBDCND.EQ.5 .OR. MBDCND.EQ.6 .OR. MBDCND.EQ.9))
  426. 1 IERROR = 15
  427. IF (B.LT.PI .AND. MBDCND.GE.7) IERROR = 16
  428. IF (ELMBDA.NE.0. .AND. NBDCND.GE.5) IERROR = 17
  429. IF (IERROR .NE. 0) GO TO 101
  430. IWBM = M+1
  431. IWCM = IWBM+M
  432. IWAN = IWCM+M
  433. IWBN = IWAN+N
  434. IWCN = IWBN+N
  435. IWSNTH = IWCN+N
  436. IWRSQ = IWSNTH+M
  437. IWWRK = IWRSQ+N
  438. IERR1 = 0
  439. CALL HSTCS1 (INTL,A,B,M,MBDCND,BDA,BDB,C,D,N,NBDCND,BDC,BDD,
  440. 1 ELMBDA,F,IDIMF,PERTRB,IERR1,W,W(IWBM),W(IWCM),
  441. 2 W(IWAN),W(IWBN),W(IWCN),W(IWSNTH),W(IWRSQ),W(IWWRK))
  442. W(1) = W(IWWRK)+IWWRK-1
  443. IERROR = IERR1
  444. 101 CONTINUE
  445. RETURN
  446. END