hwscrt.f 15 KB

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