hwscyl.f 16 KB

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