hstcyl.f 14 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461
  1. *DECK HSTCYL
  2. SUBROUTINE HSTCYL (A, B, M, MBDCND, BDA, BDB, C, D, N, NBDCND,
  3. + BDC, BDD, ELMBDA, F, IDIMF, PERTRB, IERROR, W)
  4. C***BEGIN PROLOGUE HSTCYL
  5. C***PURPOSE Solve the standard five-point finite difference
  6. C approximation on a staggered grid to the modified
  7. C Helmholtz equation in cylindrical coordinates.
  8. C***LIBRARY SLATEC (FISHPACK)
  9. C***CATEGORY I2B1A1A
  10. C***TYPE SINGLE PRECISION (HSTCYL-S)
  11. C***KEYWORDS CYLINDRICAL, ELLIPTIC, FISHPACK, HELMHOLTZ, PDE
  12. C***AUTHOR Adams, J., (NCAR)
  13. C Swarztrauber, P. N., (NCAR)
  14. C Sweet, R., (NCAR)
  15. C***DESCRIPTION
  16. C
  17. C HSTCYL solves the standard five-point finite difference
  18. C approximation on a staggered grid to the modified Helmholtz
  19. C equation in cylindrical coordinates
  20. C
  21. C (1/R)(d/dR)(R(dU/dR)) + (d/dZ)(dU/dZ)C
  22. C + LAMBDA*(1/R**2)*U = F(R,Z)
  23. C
  24. C This two-dimensional modified Helmholtz equation results
  25. C from the Fourier transform of a three-dimensional Poisson
  26. C equation.
  27. C
  28. C * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
  29. C
  30. C * * * * * * * * Parameter Description * * * * * * * * * *
  31. C
  32. C * * * * * * On Input * * * * * *
  33. C
  34. C A,B
  35. C The range of R, i.e. A .LE. R .LE. B. A must be less than B and
  36. C A must be non-negative.
  37. C
  38. C M
  39. C The number of grid points in the interval (A,B). The grid points
  40. C in the R-direction are given by R(I) = A + (I-0.5)DR for
  41. C I=1,2,...,M where DR =(B-A)/M. M must be greater than 2.
  42. C
  43. C MBDCND
  44. C Indicates the type of boundary conditions at R = A and R = B.
  45. C
  46. C = 1 If the solution is specified at R = A (see note below) and
  47. C R = B.
  48. C
  49. C = 2 If the solution is specified at R = A (see note below) and
  50. C the derivative of the solution with respect to R is
  51. C specified at R = B.
  52. C
  53. C = 3 If the derivative of the solution with respect to R is
  54. C specified at R = A (see note below) and R = B.
  55. C
  56. C = 4 If the derivative of the solution with respect to R is
  57. C specified at R = A (see note below) and the solution is
  58. C specified at R = B.
  59. C
  60. C = 5 If the solution is unspecified at R = A = 0 and the solution
  61. C is specified at R = B.
  62. C
  63. C = 6 If the solution is unspecified at R = A = 0 and the
  64. C derivative of the solution with respect to R is specified at
  65. C R = B.
  66. C
  67. C NOTE: If A = 0, do not use MBDCND = 1,2,3, or 4, but instead
  68. C use MBDCND = 5 or 6. The resulting approximation gives
  69. C the only meaningful boundary condition, i.e. dU/dR = 0.
  70. C (see D. Greenspan, 'Introductory Numerical Analysis Of
  71. C Elliptic Boundary Value Problems,' Harper and Row, 1965,
  72. C Chapter 5.)
  73. C
  74. C BDA
  75. C A one-dimensional array of length N that specifies the boundary
  76. C values (if any) of the solution at R = A. When MBDCND = 1 or 2,
  77. C
  78. C BDA(J) = U(A,Z(J)) , J=1,2,...,N.
  79. C
  80. C When MBDCND = 3 or 4,
  81. C
  82. C BDA(J) = (d/dR)U(A,Z(J)) , J=1,2,...,N.
  83. C
  84. C When MBDCND = 5 or 6, BDA is a dummy variable.
  85. C
  86. C BDB
  87. C A one-dimensional array of length N that specifies the boundary
  88. C values of the solution at R = B. When MBDCND = 1,4, or 5,
  89. C
  90. C BDB(J) = U(B,Z(J)) , J=1,2,...,N.
  91. C
  92. C When MBDCND = 2,3, or 6,
  93. C
  94. C BDB(J) = (d/dR)U(B,Z(J)) , J=1,2,...,N.
  95. C
  96. C C,D
  97. C The range of Z, i.e. C .LE. Z .LE. D. C must be less
  98. C than D.
  99. C
  100. C N
  101. C The number of unknowns in the interval (C,D). The unknowns in
  102. C the Z-direction are given by Z(J) = C + (J-0.5)DZ,
  103. C J=1,2,...,N, where DZ = (D-C)/N. N must be greater than 2.
  104. C
  105. C NBDCND
  106. C Indicates the type of boundary conditions at Z = C
  107. C and Z = D.
  108. C
  109. C = 0 If the solution is periodic in Z, i.e.
  110. C U(I,J) = U(I,N+J).
  111. C
  112. C = 1 If the solution is specified at Z = C and Z = D.
  113. C
  114. C = 2 If the solution is specified at Z = C and the derivative
  115. C of the solution with respect to Z is specified at
  116. C Z = D.
  117. C
  118. C = 3 If the derivative of the solution with respect to Z is
  119. C specified at Z = C and Z = D.
  120. C
  121. C = 4 If the derivative of the solution with respect to Z is
  122. C specified at Z = C and the solution is specified at
  123. C Z = D.
  124. C
  125. C BDC
  126. C A one dimensional array of length M that specifies the boundary
  127. C values of the solution at Z = C. When NBDCND = 1 or 2,
  128. C
  129. C BDC(I) = U(R(I),C) , I=1,2,...,M.
  130. C
  131. C When NBDCND = 3 or 4,
  132. C
  133. C BDC(I) = (d/dZ)U(R(I),C), I=1,2,...,M.
  134. C
  135. C When NBDCND = 0, BDC is a dummy variable.
  136. C
  137. C BDD
  138. C A one-dimensional array of length M that specifies the boundary
  139. C values of the solution at Z = D. when NBDCND = 1 or 4,
  140. C
  141. C BDD(I) = U(R(I),D) , I=1,2,...,M.
  142. C
  143. C When NBDCND = 2 or 3,
  144. C
  145. C BDD(I) = (d/dZ)U(R(I),D) , I=1,2,...,M.
  146. C
  147. C When NBDCND = 0, BDD is a dummy variable.
  148. C
  149. C ELMBDA
  150. C The constant LAMBDA in the modified Helmholtz equation. If
  151. C LAMBDA is greater than 0, a solution may not exist. However,
  152. C HSTCYL will attempt to find a solution. LAMBDA must be zero
  153. C when MBDCND = 5 or 6.
  154. C
  155. C F
  156. C A two-dimensional array that specifies the values of the right
  157. C side of the modified Helmholtz equation. For I=1,2,...,M
  158. C and J=1,2,...,N
  159. C
  160. C F(I,J) = F(R(I),Z(J)) .
  161. C
  162. C F must be dimensioned at least M X N.
  163. C
  164. C IDIMF
  165. C The row (or first) dimension of the array F as it appears in the
  166. C program calling HSTCYL. This parameter is used to specify the
  167. C variable dimension of F. IDIMF must be at least M.
  168. C
  169. C W
  170. C A one-dimensional array that must be provided by the user for
  171. C work space. W may require up to 13M + 4N + M*INT(log2(N))
  172. C locations. The actual number of locations used is computed by
  173. C HSTCYL and is returned in the location W(1).
  174. C
  175. C
  176. C * * * * * * On Output * * * * * *
  177. C
  178. C F
  179. C Contains the solution U(I,J) of the finite difference
  180. C approximation for the grid point (R(I),Z(J)) for
  181. C I=1,2,...,M, J=1,2,...,N.
  182. C
  183. C PERTRB
  184. C If a combination of periodic, derivative, or unspecified
  185. C boundary conditions is specified for a Poisson equation
  186. C (LAMBDA = 0), a solution may not exist. PERTRB is a con-
  187. C stant, calculated and subtracted from F, which ensures
  188. C that a solution exists. HSTCYL then computes this
  189. C solution, which is a least squares solution to the
  190. C original approximation. This solution plus any constant is also
  191. C a solution; hence, the solution is not unique. The value of
  192. C PERTRB should be small compared to the right side F.
  193. C Otherwise, a solution is obtained to an essentially different
  194. C problem. This comparison should always be made to insure that
  195. C a meaningful solution has been obtained.
  196. C
  197. C IERROR
  198. C An error flag that indicates invalid input parameters.
  199. C Except for numbers 0 and 11, a solution is not attempted.
  200. C
  201. C = 0 No error
  202. C
  203. C = 1 A .LT. 0
  204. C
  205. C = 2 A .GE. B
  206. C
  207. C = 3 MBDCND .LT. 1 or MBDCND .GT. 6
  208. C
  209. C = 4 C .GE. D
  210. C
  211. C = 5 N .LE. 2
  212. C
  213. C = 6 NBDCND .LT. 0 or NBDCND .GT. 4
  214. C
  215. C = 7 A = 0 and MBDCND = 1,2,3, or 4
  216. C
  217. C = 8 A .GT. 0 and MBDCND .GE. 5
  218. C
  219. C = 9 M .LE. 2
  220. C
  221. C = 10 IDIMF .LT. M
  222. C
  223. C = 11 LAMBDA .GT. 0
  224. C
  225. C = 12 A=0, MBDCND .GE. 5, ELMBDA .NE. 0
  226. C
  227. C Since this is the only means of indicating a possibly
  228. C incorrect call to HSTCYL, the user should test IERROR after
  229. C the call.
  230. C
  231. C W
  232. C W(1) contains the required length of W.
  233. C
  234. C *Long Description:
  235. C
  236. C * * * * * * * Program Specifications * * * * * * * * * * * *
  237. C
  238. C Dimension OF BDA(N),BDB(N),BDC(M),BDD(M),F(IDIMF,N),
  239. C Arguments W(see argument list)
  240. C
  241. C Latest June 1, 1977
  242. C Revision
  243. C
  244. C Subprograms HSTCYL,POISTG,POSTG2,GENBUN,POISD2,POISN2,POISP2,
  245. C Required COSGEN,MERGE,TRIX,TRI3,PIMACH
  246. C
  247. C Special NONE
  248. C Conditions
  249. C
  250. C Common NONE
  251. C Blocks
  252. C
  253. C I/O NONE
  254. C
  255. C Precision Single
  256. C
  257. C Specialist Roland Sweet
  258. C
  259. C Language FORTRAN
  260. C
  261. C History Written by Roland Sweet at NCAR in March, 1977
  262. C
  263. C Algorithm This subroutine defines the finite-difference
  264. C equations, incorporates boundary data, adjusts the
  265. C right side when the system is singular and calls
  266. C either POISTG or GENBUN which solves the linear
  267. C system of equations.
  268. C
  269. C Space 8228(decimal) = 20044(octal) locations on the
  270. C Required NCAR Control Data 7600
  271. C
  272. C Timing and The execution time T on the NCAR Control Data
  273. C Accuracy 7600 for subroutine HSTCYL is roughly proportional
  274. C to M*N*log2(N). Some typical values are listed in
  275. C the table below.
  276. C The solution process employed results in a loss
  277. C of no more than four significant digits for N and M
  278. C as large as 64. More detailed information about
  279. C accuracy can be found in the documentation for
  280. C subroutine POISTG which is the routine that
  281. C actually solves the finite difference equations.
  282. C
  283. C
  284. C M(=N) MBDCND NBDCND T(MSECS)
  285. C ----- ------ ------ --------
  286. C
  287. C 32 1-6 1-4 56
  288. C 64 1-6 1-4 230
  289. C
  290. C Portability American National Standards Institute Fortran.
  291. C The machine dependent constant PI is defined in
  292. C function PIMACH.
  293. C
  294. C Required COS
  295. C Resident
  296. C Routines
  297. C
  298. C Reference Schumann, U. and R. Sweet,'A Direct Method For
  299. C The Solution of Poisson's Equation With Neumann
  300. C Boundary Conditions On A Staggered Grid Of
  301. C Arbitrary Size,' J. Comp. Phys. 20(1976),
  302. C pp. 171-182.
  303. C
  304. C * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
  305. C
  306. C***REFERENCES U. Schumann and R. Sweet, A direct method for the
  307. C solution of Poisson's equation with Neumann boundary
  308. C conditions on a staggered grid of arbitrary size,
  309. C Journal of Computational Physics 20, (1976),
  310. C pp. 171-182.
  311. C***ROUTINES CALLED GENBUN, POISTG
  312. C***REVISION HISTORY (YYMMDD)
  313. C 801001 DATE WRITTEN
  314. C 890531 Changed all specific intrinsics to generic. (WRB)
  315. C 890531 REVISION DATE from Version 3.2
  316. C 891214 Prologue converted to Version 4.0 format. (BAB)
  317. C 920501 Reformatted the REFERENCES section. (WRB)
  318. C***END PROLOGUE HSTCYL
  319. C
  320. C
  321. DIMENSION F(IDIMF,*) ,BDA(*) ,BDB(*) ,BDC(*) ,
  322. 1 BDD(*) ,W(*)
  323. C***FIRST EXECUTABLE STATEMENT HSTCYL
  324. IERROR = 0
  325. IF (A .LT. 0.) IERROR = 1
  326. IF (A .GE. B) IERROR = 2
  327. IF (MBDCND.LE.0 .OR. MBDCND.GE.7) IERROR = 3
  328. IF (C .GE. D) IERROR = 4
  329. IF (N .LE. 2) IERROR = 5
  330. IF (NBDCND.LT.0 .OR. NBDCND.GE.5) IERROR = 6
  331. IF (A.EQ.0. .AND. MBDCND.NE.5 .AND. MBDCND.NE.6) IERROR = 7
  332. IF (A.GT.0. .AND. MBDCND.GE.5) IERROR = 8
  333. IF (IDIMF .LT. M) IERROR = 10
  334. IF (M .LE. 2) IERROR = 9
  335. IF (A.EQ.0. .AND. MBDCND.GE.5 .AND. ELMBDA.NE.0.) IERROR = 12
  336. IF (IERROR .NE. 0) RETURN
  337. DELTAR = (B-A)/M
  338. DLRSQ = DELTAR**2
  339. DELTHT = (D-C)/N
  340. DLTHSQ = DELTHT**2
  341. NP = NBDCND+1
  342. C
  343. C DEFINE A,B,C COEFFICIENTS IN W-ARRAY.
  344. C
  345. IWB = M
  346. IWC = IWB+M
  347. IWR = IWC+M
  348. DO 101 I=1,M
  349. J = IWR+I
  350. W(J) = A+(I-0.5)*DELTAR
  351. W(I) = (A+(I-1)*DELTAR)/(DLRSQ*W(J))
  352. K = IWC+I
  353. W(K) = (A+I*DELTAR)/(DLRSQ*W(J))
  354. K = IWB+I
  355. W(K) = ELMBDA/W(J)**2-2./DLRSQ
  356. 101 CONTINUE
  357. C
  358. C ENTER BOUNDARY DATA FOR R-BOUNDARIES.
  359. C
  360. GO TO (102,102,104,104,106,106),MBDCND
  361. 102 A1 = 2.*W(1)
  362. W(IWB+1) = W(IWB+1)-W(1)
  363. DO 103 J=1,N
  364. F(1,J) = F(1,J)-A1*BDA(J)
  365. 103 CONTINUE
  366. GO TO 106
  367. 104 A1 = DELTAR*W(1)
  368. W(IWB+1) = W(IWB+1)+W(1)
  369. DO 105 J=1,N
  370. F(1,J) = F(1,J)+A1*BDA(J)
  371. 105 CONTINUE
  372. 106 CONTINUE
  373. GO TO (107,109,109,107,107,109),MBDCND
  374. 107 W(IWC) = W(IWC)-W(IWR)
  375. A1 = 2.*W(IWR)
  376. DO 108 J=1,N
  377. F(M,J) = F(M,J)-A1*BDB(J)
  378. 108 CONTINUE
  379. GO TO 111
  380. 109 W(IWC) = W(IWC)+W(IWR)
  381. A1 = DELTAR*W(IWR)
  382. DO 110 J=1,N
  383. F(M,J) = F(M,J)-A1*BDB(J)
  384. 110 CONTINUE
  385. C
  386. C ENTER BOUNDARY DATA FOR THETA-BOUNDARIES.
  387. C
  388. 111 A1 = 2./DLTHSQ
  389. GO TO (121,112,112,114,114),NP
  390. 112 DO 113 I=1,M
  391. F(I,1) = F(I,1)-A1*BDC(I)
  392. 113 CONTINUE
  393. GO TO 116
  394. 114 A1 = 1./DELTHT
  395. DO 115 I=1,M
  396. F(I,1) = F(I,1)+A1*BDC(I)
  397. 115 CONTINUE
  398. 116 A1 = 2./DLTHSQ
  399. GO TO (121,117,119,119,117),NP
  400. 117 DO 118 I=1,M
  401. F(I,N) = F(I,N)-A1*BDD(I)
  402. 118 CONTINUE
  403. GO TO 121
  404. 119 A1 = 1./DELTHT
  405. DO 120 I=1,M
  406. F(I,N) = F(I,N)-A1*BDD(I)
  407. 120 CONTINUE
  408. 121 CONTINUE
  409. C
  410. C ADJUST RIGHT SIDE OF SINGULAR PROBLEMS TO INSURE EXISTENCE OF A
  411. C SOLUTION.
  412. C
  413. PERTRB = 0.
  414. IF (ELMBDA) 130,123,122
  415. 122 IERROR = 11
  416. GO TO 130
  417. 123 GO TO (130,130,124,130,130,124),MBDCND
  418. 124 GO TO (125,130,130,125,130),NP
  419. 125 CONTINUE
  420. DO 127 I=1,M
  421. A1 = 0.
  422. DO 126 J=1,N
  423. A1 = A1+F(I,J)
  424. 126 CONTINUE
  425. J = IWR+I
  426. PERTRB = PERTRB+A1*W(J)
  427. 127 CONTINUE
  428. PERTRB = PERTRB/(M*N*0.5*(A+B))
  429. DO 129 I=1,M
  430. DO 128 J=1,N
  431. F(I,J) = F(I,J)-PERTRB
  432. 128 CONTINUE
  433. 129 CONTINUE
  434. 130 CONTINUE
  435. C
  436. C MULTIPLY I-TH EQUATION THROUGH BY DELTHT**2
  437. C
  438. DO 132 I=1,M
  439. W(I) = W(I)*DLTHSQ
  440. J = IWC+I
  441. W(J) = W(J)*DLTHSQ
  442. J = IWB+I
  443. W(J) = W(J)*DLTHSQ
  444. DO 131 J=1,N
  445. F(I,J) = F(I,J)*DLTHSQ
  446. 131 CONTINUE
  447. 132 CONTINUE
  448. LP = NBDCND
  449. W(1) = 0.
  450. W(IWR) = 0.
  451. C
  452. C CALL GENBUN TO SOLVE THE SYSTEM OF EQUATIONS.
  453. C
  454. IF (NBDCND .EQ. 0) GO TO 133
  455. CALL POISTG (LP,N,1,M,W,W(IWB+1),W(IWC+1),IDIMF,F,IERR1,W(IWR+1))
  456. GO TO 134
  457. 133 CALL GENBUN (LP,N,1,M,W,W(IWB+1),W(IWC+1),IDIMF,F,IERR1,W(IWR+1))
  458. 134 CONTINUE
  459. W(1) = W(IWR+1)+3*M
  460. RETURN
  461. END