hstplr.f 15 KB

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