hwscsp.f 16 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405
  1. *DECK HWSCSP
  2. SUBROUTINE HWSCSP (INTL, TS, TF, M, MBDCND, BDTS, BDTF, RS, RF, N,
  3. + NBDCND, BDRS, BDRF, ELMBDA, F, IDIMF, PERTRB, IERROR, W)
  4. C***BEGIN PROLOGUE HWSCSP
  5. C***PURPOSE Solve a finite difference approximation to the modified
  6. C Helmholtz equation in spherical coordinates assuming
  7. C axisymmetry (no dependence on longitude).
  8. C***LIBRARY SLATEC (FISHPACK)
  9. C***CATEGORY I2B1A1A
  10. C***TYPE SINGLE PRECISION (HWSCSP-S)
  11. C***KEYWORDS ELLIPTIC, FISHPACK, HELMHOLTZ, PDE, SPHERICAL
  12. C***AUTHOR Adams, J., (NCAR)
  13. C Swarztrauber, P. N., (NCAR)
  14. C Sweet, R., (NCAR)
  15. C***DESCRIPTION
  16. C
  17. C Subroutine HWSCSP solves a finite difference approximation to the
  18. C modified Helmholtz equation in spherical coordinates assuming
  19. C axisymmetry (no dependence on longitude)
  20. C
  21. C (1/R**2)(d/dR)((R**2)(d/dR)U)
  22. C
  23. C + (1/(R**2)SIN(THETA))(d/dTHETA)(SIN(THETA)(d/dTHETA)U)
  24. C
  25. C + (LAMBDA/(RSIN(THETA))**2)U = F(THETA,R).
  26. C
  27. C This two dimensional modified Helmholtz equation results from
  28. C the Fourier transform of the three dimensional Poisson equation
  29. C
  30. C * * * * * * * * * * On Input * * * * * * * * * *
  31. C
  32. C INTL
  33. C = 0 On initial entry to HWSCSP or if any of the arguments
  34. C RS, RF, N, NBDCND are changed from a previous call.
  35. C = 1 If RS, RF, N, NBDCND are all unchanged from previous call
  36. C to HWSCSP.
  37. C
  38. C NOTE A call with INTL=0 takes approximately 1.5 times as
  39. C much time as a call with INTL = 1. Once a call with
  40. C INTL = 0 has been made then subsequent solutions
  41. C corresponding to different F, BDTS, BDTF, BDRS, BDRF can
  42. C be obtained faster with INTL = 1 since initialization is
  43. C not repeated.
  44. C
  45. C TS,TF
  46. C The range of THETA (colatitude), i.e., TS .LE. THETA .LE. TF.
  47. C TS must be less than TF. TS and TF are in radians. A TS of
  48. C zero corresponds to the north pole and a TF of PI corresponds
  49. C to the south pole.
  50. C
  51. C * * * * * * * * * * * * * * IMPORTANT * * * * * * * * * * * * * *
  52. C
  53. C If TF is equal to PI then it must be computed using the statement
  54. C TF = PIMACH(DUM). This insures that TF in the users program is
  55. C equal to PI in this program which permits several tests of the
  56. C input parameters that otherwise would not be possible.
  57. C
  58. C M
  59. C The number of panels into which the interval (TS,TF) is
  60. C subdivided. Hence, there will be M+1 grid points in the
  61. C THETA-direction given by THETA(K) = (I-1)DTHETA+TS for
  62. C I = 1,2,...,M+1, where DTHETA = (TF-TS)/M is the panel width.
  63. C
  64. C MBDCND
  65. C Indicates the type of boundary condition at THETA = TS and
  66. C THETA = TF.
  67. C
  68. C = 1 If the solution is specified at THETA = TS and THETA = TF.
  69. C = 2 If the solution is specified at THETA = TS and the
  70. C derivative of the solution with respect to THETA is
  71. C specified at THETA = TF (see note 2 below).
  72. C = 3 If the derivative of the solution with respect to THETA is
  73. C specified at THETA = TS and THETA = TF (see notes 1,2
  74. C below).
  75. C = 4 If the derivative of the solution with respect to THETA is
  76. C specified at THETA = TS (see note 1 below) and the
  77. C solution is specified at THETA = TF.
  78. C = 5 If the solution is unspecified at THETA = TS = 0 and the
  79. C solution is specified at THETA = TF.
  80. C = 6 If the solution is unspecified at THETA = TS = 0 and the
  81. C derivative of the solution with respect to THETA is
  82. C specified at THETA = TF (see note 2 below).
  83. C = 7 If the solution is specified at THETA = TS and the
  84. C solution is unspecified at THETA = TF = PI.
  85. C = 8 If the derivative of the solution with respect to THETA is
  86. C specified at THETA = TS (see note 1 below) and the solution
  87. C is unspecified at THETA = TF = PI.
  88. C = 9 If the solution is unspecified at THETA = TS = 0 and
  89. C THETA = TF = PI.
  90. C
  91. C NOTES: 1. If TS = 0, do not use MBDCND = 3,4, or 8, but
  92. C instead use MBDCND = 5,6, or 9 .
  93. C 2. If TF = PI, do not use MBDCND = 2,3, or 6, but
  94. C instead use MBDCND = 7,8, or 9 .
  95. C
  96. C BDTS
  97. C A one-dimensional array of length N+1 that specifies the values
  98. C of the derivative of the solution with respect to THETA at
  99. C THETA = TS. When MBDCND = 3,4, or 8,
  100. C
  101. C BDTS(J) = (d/dTHETA)U(TS,R(J)), J = 1,2,...,N+1 .
  102. C
  103. C When MBDCND has any other value, BDTS is a dummy variable.
  104. C
  105. C BDTF
  106. C A one-dimensional array of length N+1 that specifies the values
  107. C of the derivative of the solution with respect to THETA at
  108. C THETA = TF. When MBDCND = 2,3, or 6,
  109. C
  110. C BDTF(J) = (d/dTHETA)U(TF,R(J)), J = 1,2,...,N+1 .
  111. C
  112. C When MBDCND has any other value, BDTF is a dummy variable.
  113. C
  114. C RS,RF
  115. C The range of R, i.e., RS .LE. R .LT. RF. RS must be less than
  116. C RF. RS must be non-negative.
  117. C
  118. C N
  119. C The number of panels into which the interval (RS,RF) is
  120. C subdivided. Hence, there will be N+1 grid points in the
  121. C R-direction given by R(J) = (J-1)DR+RS for J = 1,2,...,N+1,
  122. C where DR = (RF-RS)/N is the panel width.
  123. C N must be greater than 2
  124. C
  125. C NBDCND
  126. C Indicates the type of boundary condition at R = RS and R = RF.
  127. C
  128. C = 1 If the solution is specified at R = RS and R = RF.
  129. C = 2 If the solution is specified at R = RS and the derivative
  130. C of the solution with respect to R is specified at R = RF.
  131. C = 3 If the derivative of the solution with respect to R is
  132. C specified at R = RS and R = RF.
  133. C = 4 If the derivative of the solution with respect to R is
  134. C specified at RS and the solution is specified at R = RF.
  135. C = 5 If the solution is unspecified at R = RS = 0 (see note
  136. C below) and the solution is specified at R = RF.
  137. C = 6 If the solution is unspecified at R = RS = 0 (see note
  138. C below) and the derivative of the solution with respect to
  139. C R is specified at R = RF.
  140. C
  141. C NOTE: NBDCND = 5 or 6 cannot be used with
  142. C MBDCND = 1,2,4,5, or 7 (the former indicates that the
  143. C solution is unspecified at R = 0, the latter
  144. C indicates that the solution is specified).
  145. C Use instead
  146. C NBDCND = 1 or 2 .
  147. C
  148. C BDRS
  149. C A one-dimensional array of length M+1 that specifies the values
  150. C of the derivative of the solution with respect to R at R = RS.
  151. C When NBDCND = 3 or 4,
  152. C
  153. C BDRS(I) = (d/dR)U(THETA(I),RS), I = 1,2,...,M+1 .
  154. C
  155. C When NBDCND has any other value, BDRS is a dummy variable.
  156. C
  157. C BDRF
  158. C A one-dimensional array of length M+1 that specifies the values
  159. C of the derivative of the solution with respect to R at R = RF.
  160. C When NBDCND = 2,3, or 6,
  161. C
  162. C BDRF(I) = (d/dR)U(THETA(I),RF), I = 1,2,...,M+1 .
  163. C
  164. C When NBDCND has any other value, BDRF is a dummy variable.
  165. C
  166. C ELMBDA
  167. C The constant LAMBDA in the Helmholtz equation. If
  168. C LAMBDA .GT. 0, a solution may not exist. However, HWSCSP will
  169. C attempt to find a solution. If NBDCND = 5 or 6 or
  170. C MBDCND = 5,6,7,8, or 9, ELMBDA must be zero.
  171. C
  172. C F
  173. C A two-dimensional array that specifies the value of the right
  174. C side of the Helmholtz equation and boundary values (if any).
  175. C for I = 2,3,...,M and J = 2,3,...,N
  176. C
  177. C F(I,J) = F(THETA(I),R(J)).
  178. C
  179. C On the boundaries F is defined by
  180. C
  181. C MBDCND F(1,J) F(M+1,J)
  182. C ------ ---------- ----------
  183. C
  184. C 1 U(TS,R(J)) U(TF,R(J))
  185. C 2 U(TS,R(J)) F(TF,R(J))
  186. C 3 F(TS,R(J)) F(TF,R(J))
  187. C 4 F(TS,R(J)) U(TF,R(J))
  188. C 5 F(0,R(J)) U(TF,R(J)) J = 1,2,...,N+1
  189. C 6 F(0,R(J)) F(TF,R(J))
  190. C 7 U(TS,R(J)) F(PI,R(J))
  191. C 8 F(TS,R(J)) F(PI,R(J))
  192. C 9 F(0,R(J)) F(PI,R(J))
  193. C
  194. C NBDCND F(I,1) F(I,N+1)
  195. C ------ -------------- --------------
  196. C
  197. C 1 U(THETA(I),RS) U(THETA(I),RF)
  198. C 2 U(THETA(I),RS) F(THETA(I),RF)
  199. C 3 F(THETA(I),RS) F(THETA(I),RF)
  200. C 4 F(THETA(I),RS) U(THETA(I),RF) I = 1,2,...,M+1
  201. C 5 F(TS,0) U(THETA(I),RF)
  202. C 6 F(TS,0) F(THETA(I),RF)
  203. C
  204. C F must be dimensioned at least (M+1)*(N+1).
  205. C
  206. C NOTE
  207. C
  208. C If the table calls for both the solution U and the right side F
  209. C at a corner then the solution must be specified.
  210. C
  211. C IDIMF
  212. C The row (or first) dimension of the array F as it appears in the
  213. C program calling HWSCSP. This parameter is used to specify the
  214. C variable dimension of F. IDIMF must be at least M+1 .
  215. C
  216. C W
  217. C A one-dimensional array that must be provided by the user for
  218. C work space. Its length can be computed from the formula below
  219. C which depends on the value of NBDCND.
  220. C
  221. C If NBDCND=2,4 or 6 define NUNK=N
  222. C If NBDCND=1 or 5 define NUNK=N-1
  223. C If NBDCND=3 define NUNK=N+1
  224. C
  225. C Now set K=INT(log2(NUNK))+1 and L=2**(K+1) then W must be
  226. C dimensioned at least (K-2)*L+K+5*(M+N)+MAX(2*N,6*M)+23
  227. C
  228. C **IMPORTANT** For purposes of checking, the required length
  229. C of W is computed by HWSCSP and stored in W(1)
  230. C in floating point format.
  231. C
  232. C
  233. C * * * * * * * * * * On Output * * * * * * * * * *
  234. C
  235. C F
  236. C Contains the solution U(I,J) of the finite difference
  237. C approximation for the grid point (THETA(I),R(J)),
  238. C I = 1,2,...,M+1, J = 1,2,...,N+1 .
  239. C
  240. C PERTRB
  241. C If a combination of periodic or derivative boundary conditions
  242. C is specified for a Poisson equation (LAMBDA = 0), a solution may
  243. C not exist. PERTRB is a constant, calculated and subtracted from
  244. C F, which ensures that a solution exists. HWSCSP then computes
  245. C this solution, which is a least squares solution to the original
  246. C approximation. This solution is not unique and is unnormalized.
  247. C The value of PERTRB should be small compared to the right side
  248. C F. Otherwise , a solution is obtained to an essentially
  249. C different problem. This comparison should always be made to
  250. C insure that a meaningful solution has been obtained.
  251. C
  252. C IERROR
  253. C An error flag that indicates invalid input parameters. Except
  254. C for numbers 0 and 10, a solution is not attempted.
  255. C
  256. C = 1 TS.LT.0. or TF.GT.PI
  257. C = 2 TS.GE.TF
  258. C = 3 M.LT.5
  259. C = 4 MBDCND.LT.1 or MBDCND.GT.9
  260. C = 5 RS.LT.0
  261. C = 6 RS.GE.RF
  262. C = 7 N.LT.5
  263. C = 8 NBDCND.LT.1 or NBDCND.GT.6
  264. C = 9 ELMBDA.GT.0
  265. C = 10 IDIMF.LT.M+1
  266. C = 11 ELMBDA.NE.0 and MBDCND.GE.5
  267. C = 12 ELMBDA.NE.0 and NBDCND equals 5 or 6
  268. C = 13 MBDCND equals 5,6 or 9 and TS.NE.0
  269. C = 14 MBDCND.GE.7 and TF.NE.PI
  270. C = 15 TS.EQ.0 and MBDCND equals 3,4 or 8
  271. C = 16 TF.EQ.PI and MBDCND equals 2,3 or 6
  272. C = 17 NBDCND.GE.5 and RS.NE.0
  273. C = 18 NBDCND.GE.5 and MBDCND equals 1,2,4,5 or 7
  274. C
  275. C Since this is the only means of indicating a possibly incorrect
  276. C call to HWSCSP, the user should test IERROR after a call.
  277. C
  278. C W
  279. C Contains intermediate values that must not be destroyed if
  280. C HWSCSP will be called again with INTL = 1. W(1) contains the
  281. C number of locations which W must have.
  282. C
  283. C *Long Description:
  284. C
  285. C * * * * * * * Program Specifications * * * * * * * * * * * *
  286. C
  287. C Dimension of BDTS(N+1),BDTF(N+1),BDRS(M+1),BDRF(M+1),
  288. C Arguments F(IDIMF,N+1),W(see argument list)
  289. C
  290. C Latest June 1979
  291. C Revision
  292. C
  293. C Subprograms HWSCSP,HWSCS1,BLKTRI,BLKTR1,PROD,PRODP,CPROD,CPRODP
  294. C Required ,COMBP,PPADD,PSGF,BSRH,PPSGF,PPSPF,TEVLS,INDXA,
  295. C ,INDXB,INDXC,R1MACH
  296. C
  297. C Special
  298. C Conditions
  299. C
  300. C Common CBLKT
  301. C Blocks
  302. C
  303. C I/O NONE
  304. C
  305. C Precision Single
  306. C
  307. C Specialist Paul N Swarztrauber
  308. C
  309. C Language FORTRAN
  310. C
  311. C History Version 1 September 1973
  312. C Version 2 April 1976
  313. C Version 3 June 1979
  314. C
  315. C Algorithm The routine defines the finite difference
  316. C equations, incorporates boundary data, and adjusts
  317. C the right side of singular systems and then calls
  318. C BLKTRI to solve the system.
  319. C
  320. C Space
  321. C Required
  322. C
  323. C Portability American National Standards Institute FORTRAN.
  324. C The machine accuracy is set using function R1MACH.
  325. C
  326. C Required NONE
  327. C Resident
  328. C Routines
  329. C
  330. C Reference Swarztrauber,P. and R. Sweet, 'Efficient FORTRAN
  331. C Subprograms for The Solution Of Elliptic Equations'
  332. C NCAR TN/IA-109, July, 1975, 138 pp.
  333. C
  334. C * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
  335. C
  336. C***REFERENCES P. N. Swarztrauber and R. Sweet, Efficient Fortran
  337. C subprograms for the solution of elliptic equations,
  338. C NCAR TN/IA-109, July 1975, 138 pp.
  339. C***ROUTINES CALLED HWSCS1, PIMACH
  340. C***REVISION HISTORY (YYMMDD)
  341. C 801001 DATE WRITTEN
  342. C 890531 Changed all specific intrinsics to generic. (WRB)
  343. C 890531 REVISION DATE from Version 3.2
  344. C 891214 Prologue converted to Version 4.0 format. (BAB)
  345. C 920501 Reformatted the REFERENCES section. (WRB)
  346. C***END PROLOGUE HWSCSP
  347. C
  348. DIMENSION F(IDIMF,*) ,BDTS(*) ,BDTF(*) ,BDRS(*) ,
  349. 1 BDRF(*) ,W(*)
  350. C***FIRST EXECUTABLE STATEMENT HWSCSP
  351. PI = PIMACH(DUM)
  352. IERROR = 0
  353. IF (TS.LT.0. .OR. TF.GT.PI) IERROR = 1
  354. IF (TS .GE. TF) IERROR = 2
  355. IF (M .LT. 5) IERROR = 3
  356. IF (MBDCND.LT.1 .OR. MBDCND.GT.9) IERROR = 4
  357. IF (RS .LT. 0.) IERROR = 5
  358. IF (RS .GE. RF) IERROR = 6
  359. IF (N .LT. 5) IERROR = 7
  360. IF (NBDCND.LT.1 .OR. NBDCND.GT.6) IERROR = 8
  361. IF (ELMBDA .GT. 0.) IERROR = 9
  362. IF (IDIMF .LT. M+1) IERROR = 10
  363. IF (ELMBDA.NE.0. .AND. MBDCND.GE.5) IERROR = 11
  364. IF (ELMBDA.NE.0. .AND. (NBDCND.EQ.5 .OR. NBDCND.EQ.6)) IERROR = 12
  365. IF ((MBDCND.EQ.5 .OR. MBDCND.EQ.6 .OR. MBDCND.EQ.9) .AND.
  366. 1 TS.NE.0.) IERROR = 13
  367. IF (MBDCND.GE.7 .AND. TF.NE.PI) IERROR = 14
  368. IF (TS.EQ.0. .AND.
  369. 1 (MBDCND.EQ.4 .OR. MBDCND.EQ.8 .OR. MBDCND.EQ.3)) IERROR = 15
  370. IF (TF.EQ.PI .AND.
  371. 1 (MBDCND.EQ.2 .OR. MBDCND.EQ.3 .OR. MBDCND.EQ.6)) IERROR = 16
  372. IF (NBDCND.GE.5 .AND. RS.NE.0.) IERROR = 17
  373. IF (NBDCND.GE.5 .AND. (MBDCND.EQ.1 .OR. MBDCND.EQ.2 .OR.
  374. 1 MBDCND.EQ.5 .OR. MBDCND.EQ.7))
  375. 2 IERROR = 18
  376. IF (IERROR.NE.0 .AND. IERROR.NE.9) RETURN
  377. NCK = N
  378. GO TO (101,103,102,103,101,103),NBDCND
  379. 101 NCK = NCK-1
  380. GO TO 103
  381. 102 NCK = NCK+1
  382. 103 L = 2
  383. K = 1
  384. 104 L = L+L
  385. K = K+1
  386. IF (NCK-L) 105,105,104
  387. 105 L = L+L
  388. NP1 = N+1
  389. MP1 = M+1
  390. I1 = (K-2)*L+K+MAX(2*N,6*M)+13
  391. I2 = I1+NP1
  392. I3 = I2+NP1
  393. I4 = I3+NP1
  394. I5 = I4+NP1
  395. I6 = I5+NP1
  396. I7 = I6+MP1
  397. I8 = I7+MP1
  398. I9 = I8+MP1
  399. I10 = I9+MP1
  400. W(1) = I10+M
  401. CALL HWSCS1 (INTL,TS,TF,M,MBDCND,BDTS,BDTF,RS,RF,N,NBDCND,BDRS,
  402. 1 BDRF,ELMBDA,F,IDIMF,PERTRB,W(2),W(I1),W(I2),W(I3),
  403. 2 W(I4),W(I5),W(I6),W(I7),W(I8),W(I9),W(I10))
  404. RETURN
  405. END