hstssp.f 18 KB

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