hwsplr.f 18 KB

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