1 | ! This module encapsulates misc. math functions |
---|
2 | |
---|
3 | MODULE lmdz_libmath_pch |
---|
4 | USE slatec_xer, ONLY: xermsg |
---|
5 | IMPLICIT NONE; PRIVATE |
---|
6 | PUBLIC pchfe_95, pchsp_95 |
---|
7 | CONTAINS |
---|
8 | REAL FUNCTION cbrt(x) |
---|
9 | ! Return the cubic root for positive or negative x |
---|
10 | REAL :: x |
---|
11 | |
---|
12 | cbrt = sign(1., x) * (abs(x)**(1. / 3.)) |
---|
13 | END FUNCTION cbrt |
---|
14 | |
---|
15 | SUBROUTINE CHFEV(X1, X2, F1, F2, D1, D2, NE, XE, FE, NEXT, IERR) |
---|
16 | !***BEGIN PROLOGUE CHFEV |
---|
17 | !***PURPOSE Evaluate a cubic polynomial given in Hermite form at an |
---|
18 | ! array of points. While designed for use by PCHFE, it may |
---|
19 | ! be useful directly as an evaluator for a piecewise cubic |
---|
20 | ! Hermite function in applications, such as graphing, where |
---|
21 | ! the interval is known in advance. |
---|
22 | !***LIBRARY SLATEC (PCHIP) |
---|
23 | !***CATEGORY E3 |
---|
24 | !***TYPE SINGLE PRECISION (CHFEV-S, DCHFEV-D) |
---|
25 | !***KEYWORDS CUBIC HERMITE EVALUATION, CUBIC POLYNOMIAL EVALUATION, |
---|
26 | ! PCHIP |
---|
27 | !***AUTHOR Fritsch, F. N., (LLNL) |
---|
28 | ! Lawrence Livermore National Laboratory |
---|
29 | ! P.O. Box 808 (L-316) |
---|
30 | ! Livermore, CA 94550 |
---|
31 | ! FTS 532-4275, (510) 422-4275 |
---|
32 | !***DESCRIPTION |
---|
33 | |
---|
34 | ! CHFEV: Cubic Hermite Function EValuator |
---|
35 | |
---|
36 | ! Evaluates the cubic polynomial determined by function values |
---|
37 | ! F1,F2 and derivatives D1,D2 on interval (X1,X2) at the points |
---|
38 | ! XE(J), J=1(1)NE. |
---|
39 | |
---|
40 | ! ---------------------------------------------------------------------- |
---|
41 | |
---|
42 | ! Calling sequence: |
---|
43 | |
---|
44 | ! INTEGER NE, NEXT(2), IERR |
---|
45 | ! REAL X1, X2, F1, F2, D1, D2, XE(NE), FE(NE) |
---|
46 | |
---|
47 | ! CALL CHFEV (X1,X2, F1,F2, D1,D2, NE, XE, FE, NEXT, IERR) |
---|
48 | |
---|
49 | ! Parameters: |
---|
50 | |
---|
51 | ! X1,X2 -- (input) endpoints of interval of definition of cubic. |
---|
52 | ! (Error return if X1.EQ.X2 .) |
---|
53 | |
---|
54 | ! F1,F2 -- (input) values of function at X1 and X2, respectively. |
---|
55 | |
---|
56 | ! D1,D2 -- (input) values of derivative at X1 and X2, respectively. |
---|
57 | |
---|
58 | ! NE -- (input) number of evaluation points. (Error return if |
---|
59 | ! NE.LT.1 .) |
---|
60 | |
---|
61 | ! XE -- (input) real array of points at which the function is to be |
---|
62 | ! evaluated. If any of the XE are outside the interval |
---|
63 | ! [X1,X2], a warning error is returned in NEXT. |
---|
64 | |
---|
65 | ! FE -- (output) real array of values of the cubic function defined |
---|
66 | ! by X1,X2, F1,F2, D1,D2 at the points XE. |
---|
67 | |
---|
68 | ! NEXT -- (output) integer array indicating number of extrapolation |
---|
69 | ! points: |
---|
70 | ! NEXT(1) = number of evaluation points to left of interval. |
---|
71 | ! NEXT(2) = number of evaluation points to right of interval. |
---|
72 | |
---|
73 | ! IERR -- (output) error flag. |
---|
74 | ! Normal return: |
---|
75 | ! IERR = 0 (no errors). |
---|
76 | ! "Recoverable" errors: |
---|
77 | ! IERR = -1 if NE.LT.1 . |
---|
78 | ! IERR = -2 if X1.EQ.X2 . |
---|
79 | ! (The FE-array has not been changed in either case.) |
---|
80 | |
---|
81 | !***REFERENCES (NONE) |
---|
82 | !***ROUTINES CALLED XERMSG |
---|
83 | !***REVISION HISTORY (YYMMDD) |
---|
84 | ! 811019 DATE WRITTEN |
---|
85 | ! 820803 Minor cosmetic changes for release 1. |
---|
86 | ! 890411 Added SAVE statements (Vers. 3.2). |
---|
87 | ! 890531 Changed all specific intrinsics to generic. (WRB) |
---|
88 | ! 890703 Corrected category record. (WRB) |
---|
89 | ! 890703 REVISION DATE from Version 3.2 |
---|
90 | ! 891214 Prologue converted to Version 4.0 format. (BAB) |
---|
91 | ! 900315 CALLs to XERROR changed to CALLs to XERMSG. (THJ) |
---|
92 | !***END PROLOGUE CHFEV |
---|
93 | ! Programming notes: |
---|
94 | |
---|
95 | ! To produce a double precision version, simply: |
---|
96 | ! a. Change CHFEV to DCHFEV wherever it occurs, |
---|
97 | ! b. Change the real declaration to double precision, and |
---|
98 | ! c. Change the constant ZERO to double precision. |
---|
99 | |
---|
100 | ! DECLARE ARGUMENTS. |
---|
101 | |
---|
102 | INTEGER :: NE, NEXT(2), IERR |
---|
103 | REAL :: X1, X2, F1, F2, D1, D2, XE(*), FE(*) |
---|
104 | |
---|
105 | ! DECLARE LOCAL VARIABLES. |
---|
106 | |
---|
107 | INTEGER :: I |
---|
108 | REAL :: C2, C3, DEL1, DEL2, DELTA, H, X, XMI, XMA, ZERO |
---|
109 | SAVE ZERO |
---|
110 | DATA ZERO /0./ |
---|
111 | |
---|
112 | ! VALIDITY-CHECK ARGUMENTS. |
---|
113 | |
---|
114 | !***FIRST EXECUTABLE STATEMENT CHFEV |
---|
115 | IF (NE < 1) GO TO 5001 |
---|
116 | H = X2 - X1 |
---|
117 | IF (H == ZERO) GO TO 5002 |
---|
118 | |
---|
119 | ! INITIALIZE. |
---|
120 | |
---|
121 | IERR = 0 |
---|
122 | NEXT(1) = 0 |
---|
123 | NEXT(2) = 0 |
---|
124 | XMI = MIN(ZERO, H) |
---|
125 | XMA = MAX(ZERO, H) |
---|
126 | |
---|
127 | ! COMPUTE CUBIC COEFFICIENTS (EXPANDED ABOUT X1). |
---|
128 | |
---|
129 | DELTA = (F2 - F1) / H |
---|
130 | DEL1 = (D1 - DELTA) / H |
---|
131 | DEL2 = (D2 - DELTA) / H |
---|
132 | ! (DELTA IS NO LONGER NEEDED.) |
---|
133 | C2 = -(DEL1 + DEL1 + DEL2) |
---|
134 | C3 = (DEL1 + DEL2) / H |
---|
135 | ! (H, DEL1 AND DEL2 ARE NO LONGER NEEDED.) |
---|
136 | |
---|
137 | ! EVALUATION LOOP. |
---|
138 | |
---|
139 | DO I = 1, NE |
---|
140 | X = XE(I) - X1 |
---|
141 | FE(I) = F1 + X * (D1 + X * (C2 + X * C3)) |
---|
142 | ! COUNT EXTRAPOLATION POINTS. |
---|
143 | IF (X<XMI) NEXT(1) = NEXT(1) + 1 |
---|
144 | IF (X>XMA) NEXT(2) = NEXT(2) + 1 |
---|
145 | ! (NOTE REDUNDANCY--IF EITHER CONDITION IS TRUE, OTHER IS FALSE.) |
---|
146 | END DO |
---|
147 | |
---|
148 | ! NORMAL RETURN. |
---|
149 | |
---|
150 | RETURN |
---|
151 | |
---|
152 | ! ERROR RETURNS. |
---|
153 | |
---|
154 | 5001 CONTINUE |
---|
155 | ! NE.LT.1 RETURN. |
---|
156 | IERR = -1 |
---|
157 | CALL XERMSG ('SLATEC', 'CHFEV', & |
---|
158 | 'NUMBER OF EVALUATION POINTS LESS THAN ONE', IERR, 1) |
---|
159 | RETURN |
---|
160 | |
---|
161 | 5002 CONTINUE |
---|
162 | ! X1.EQ.X2 RETURN. |
---|
163 | IERR = -2 |
---|
164 | CALL XERMSG ('SLATEC', 'CHFEV', 'INTERVAL ENDPOINTS EQUAL', IERR, & |
---|
165 | 1) |
---|
166 | RETURN |
---|
167 | !------------- LAST LINE OF CHFEV FOLLOWS ------------------------------ |
---|
168 | END SUBROUTINE CHFEV |
---|
169 | |
---|
170 | SUBROUTINE PCHFE(N, X, F, D, INCFD, SKIP, NE, XE, FE, IERR) |
---|
171 | !***BEGIN PROLOGUE PCHFE |
---|
172 | !***PURPOSE Evaluate a piecewise cubic Hermite function at an array of |
---|
173 | ! points. May be used by itself for Hermite interpolation, |
---|
174 | ! or as an evaluator for PCHIM or PCHIC. |
---|
175 | !***LIBRARY SLATEC (PCHIP) |
---|
176 | !***CATEGORY E3 |
---|
177 | !***TYPE SINGLE PRECISION (PCHFE-S, DPCHFE-D) |
---|
178 | !***KEYWORDS CUBIC HERMITE EVALUATION, HERMITE INTERPOLATION, PCHIP, |
---|
179 | ! PIECEWISE CUBIC EVALUATION |
---|
180 | !***AUTHOR Fritsch, F. N., (LLNL) |
---|
181 | ! Lawrence Livermore National Laboratory |
---|
182 | ! P.O. Box 808 (L-316) |
---|
183 | ! Livermore, CA 94550 |
---|
184 | ! FTS 532-4275, (510) 422-4275 |
---|
185 | !***DESCRIPTION |
---|
186 | |
---|
187 | ! PCHFE: Piecewise Cubic Hermite Function Evaluator |
---|
188 | |
---|
189 | ! Evaluates the cubic Hermite function defined by N, X, F, D at |
---|
190 | ! the points XE(J), J=1(1)NE. |
---|
191 | |
---|
192 | ! To provide compatibility with PCHIM and PCHIC, includes an |
---|
193 | ! increment between successive values of the F- and D-arrays. |
---|
194 | |
---|
195 | ! ---------------------------------------------------------------------- |
---|
196 | |
---|
197 | ! Calling sequence: |
---|
198 | |
---|
199 | ! PARAMETER (INCFD = ...) |
---|
200 | ! INTEGER N, NE, IERR |
---|
201 | ! REAL X(N), F(INCFD,N), D(INCFD,N), XE(NE), FE(NE) |
---|
202 | ! LOGICAL SKIP |
---|
203 | |
---|
204 | ! CALL PCHFE (N, X, F, D, INCFD, SKIP, NE, XE, FE, IERR) |
---|
205 | |
---|
206 | ! Parameters: |
---|
207 | |
---|
208 | ! N -- (input) number of data points. (Error return if N.LT.2 .) |
---|
209 | |
---|
210 | ! X -- (input) real array of independent variable values. The |
---|
211 | ! elements of X must be strictly increasing: |
---|
212 | ! X(I-1) .LT. X(I), I = 2(1)N. |
---|
213 | ! (Error return if not.) |
---|
214 | |
---|
215 | ! F -- (input) real array of function values. F(1+(I-1)*INCFD) is |
---|
216 | ! the value corresponding to X(I). |
---|
217 | |
---|
218 | ! D -- (input) real array of derivative values. D(1+(I-1)*INCFD) is |
---|
219 | ! the value corresponding to X(I). |
---|
220 | |
---|
221 | ! INCFD -- (input) increment between successive values in F and D. |
---|
222 | ! (Error return if INCFD.LT.1 .) |
---|
223 | |
---|
224 | ! SKIP -- (input/output) LOGICAL variable which should be set to |
---|
225 | ! .TRUE. if the user wishes to skip checks for validity of |
---|
226 | ! preceding parameters, or to .FALSE. otherwise. |
---|
227 | ! This will save time in case these checks have already |
---|
228 | ! been performed (say, in PCHIM or PCHIC). |
---|
229 | ! SKIP will be set to .TRUE. on normal return. |
---|
230 | |
---|
231 | ! NE -- (input) number of evaluation points. (Error return if |
---|
232 | ! NE.LT.1 .) |
---|
233 | |
---|
234 | ! XE -- (input) real array of points at which the function is to be |
---|
235 | ! evaluated. |
---|
236 | |
---|
237 | ! NOTES: |
---|
238 | ! 1. The evaluation will be most efficient if the elements |
---|
239 | ! of XE are increasing relative to X; |
---|
240 | ! that is, XE(J) .GE. X(I) |
---|
241 | ! implies XE(K) .GE. X(I), all K.GE.J . |
---|
242 | ! 2. If any of the XE are outside the interval [X(1),X(N)], |
---|
243 | ! values are extrapolated from the nearest extreme cubic, |
---|
244 | ! and a warning error is returned. |
---|
245 | |
---|
246 | ! FE -- (output) real array of values of the cubic Hermite function |
---|
247 | ! defined by N, X, F, D at the points XE. |
---|
248 | |
---|
249 | ! IERR -- (output) error flag. |
---|
250 | ! Normal return: |
---|
251 | ! IERR = 0 (no errors). |
---|
252 | ! Warning error: |
---|
253 | ! IERR.GT.0 means that extrapolation was performed at |
---|
254 | ! IERR points. |
---|
255 | ! "Recoverable" errors: |
---|
256 | ! IERR = -1 if N.LT.2 . |
---|
257 | ! IERR = -2 if INCFD.LT.1 . |
---|
258 | ! IERR = -3 if the X-array is not strictly increasing. |
---|
259 | ! IERR = -4 if NE.LT.1 . |
---|
260 | ! (The FE-array has not been changed in any of these cases.) |
---|
261 | ! NOTE: The above errors are checked in the order listed, |
---|
262 | ! and following arguments have **NOT** been validated. |
---|
263 | |
---|
264 | !***REFERENCES (NONE) |
---|
265 | !***ROUTINES CALLED CHFEV, XERMSG |
---|
266 | !***REVISION HISTORY (YYMMDD) |
---|
267 | ! 811020 DATE WRITTEN |
---|
268 | ! 820803 Minor cosmetic changes for release 1. |
---|
269 | ! 870707 Minor cosmetic changes to prologue. |
---|
270 | ! 890531 Changed all specific intrinsics to generic. (WRB) |
---|
271 | ! 890831 Modified array declarations. (WRB) |
---|
272 | ! 890831 REVISION DATE from Version 3.2 |
---|
273 | ! 891214 Prologue converted to Version 4.0 format. (BAB) |
---|
274 | ! 900315 CALLs to XERROR changed to CALLs to XERMSG. (THJ) |
---|
275 | !***END PROLOGUE PCHFE |
---|
276 | ! Programming notes: |
---|
277 | |
---|
278 | ! 1. To produce a double precision version, simply: |
---|
279 | ! a. Change PCHFE to DPCHFE, and CHFEV to DCHFEV, wherever they |
---|
280 | ! occur, |
---|
281 | ! b. Change the real declaration to double precision, |
---|
282 | |
---|
283 | ! 2. Most of the coding between the CALL to CHFEV and the end of |
---|
284 | ! the IR-loop could be eliminated if it were permissible to |
---|
285 | ! assume that XE is ordered relative to X. |
---|
286 | |
---|
287 | ! 3. CHFEV does not assume that X1 is less than X2. thus, it would |
---|
288 | ! be possible to write a version of PCHFE that assumes a strict- |
---|
289 | ! ly decreasing X-array by simply running the IR-loop backwards |
---|
290 | ! (and reversing the order of appropriate tests). |
---|
291 | |
---|
292 | ! 4. The present code has a minor bug, which I have decided is not |
---|
293 | ! worth the effort that would be required to fix it. |
---|
294 | ! If XE contains points in [X(N-1),X(N)], followed by points .LT. |
---|
295 | ! X(N-1), followed by points .GT.X(N), the extrapolation points |
---|
296 | ! will be counted (at least) twice in the total returned in IERR. |
---|
297 | |
---|
298 | ! DECLARE ARGUMENTS. |
---|
299 | |
---|
300 | INTEGER :: N, INCFD, NE, IERR |
---|
301 | REAL :: X(*), F(INCFD, *), D(INCFD, *), XE(*), FE(*) |
---|
302 | LOGICAL :: SKIP |
---|
303 | |
---|
304 | ! DECLARE LOCAL VARIABLES. |
---|
305 | |
---|
306 | INTEGER :: I, IERC, IR, J, JFIRST, NEXT(2), NJ |
---|
307 | |
---|
308 | ! VALIDITY-CHECK ARGUMENTS. |
---|
309 | |
---|
310 | !***FIRST EXECUTABLE STATEMENT PCHFE |
---|
311 | IF (SKIP) GO TO 5 |
---|
312 | |
---|
313 | IF (N<2) GO TO 5001 |
---|
314 | IF (INCFD<1) GO TO 5002 |
---|
315 | DO I = 2, N |
---|
316 | IF (X(I)<=X(I - 1)) GO TO 5003 |
---|
317 | END DO |
---|
318 | |
---|
319 | ! FUNCTION DEFINITION IS OK, GO ON. |
---|
320 | |
---|
321 | 5 CONTINUE |
---|
322 | IF (NE<1) GO TO 5004 |
---|
323 | IERR = 0 |
---|
324 | SKIP = .TRUE. |
---|
325 | |
---|
326 | ! LOOP OVER INTERVALS. ( INTERVAL INDEX IS IL = IR-1 . ) |
---|
327 | ! ( INTERVAL IS X(IL).LE.X.LT.X(IR) . ) |
---|
328 | JFIRST = 1 |
---|
329 | IR = 2 |
---|
330 | 10 CONTINUE |
---|
331 | |
---|
332 | ! SKIP OUT OF LOOP IF HAVE PROCESSED ALL EVALUATION POINTS. |
---|
333 | |
---|
334 | IF (JFIRST > NE) GO TO 5000 |
---|
335 | |
---|
336 | ! LOCATE ALL POINTS IN INTERVAL. |
---|
337 | |
---|
338 | DO J = JFIRST, NE |
---|
339 | IF (XE(J) >= X(IR)) GO TO 30 |
---|
340 | END DO |
---|
341 | J = NE + 1 |
---|
342 | GO TO 40 |
---|
343 | |
---|
344 | ! HAVE LOCATED FIRST POINT BEYOND INTERVAL. |
---|
345 | |
---|
346 | 30 CONTINUE |
---|
347 | IF (IR == N) J = NE + 1 |
---|
348 | |
---|
349 | 40 CONTINUE |
---|
350 | NJ = J - JFIRST |
---|
351 | |
---|
352 | ! SKIP EVALUATION IF NO POINTS IN INTERVAL. |
---|
353 | |
---|
354 | IF (NJ == 0) GO TO 50 |
---|
355 | |
---|
356 | ! EVALUATE CUBIC AT XE(I), I = JFIRST (1) J-1 . |
---|
357 | |
---|
358 | ! ---------------------------------------------------------------- |
---|
359 | CALL CHFEV (X(IR - 1), X(IR), F(1, IR - 1), F(1, IR), D(1, IR - 1), D(1, IR), & |
---|
360 | NJ, XE(JFIRST), FE(JFIRST), NEXT, IERC) |
---|
361 | ! ---------------------------------------------------------------- |
---|
362 | IF (IERC < 0) GO TO 5005 |
---|
363 | |
---|
364 | IF (NEXT(2) == 0) GO TO 42 |
---|
365 | ! IF (NEXT(2) .GT. 0) THEN |
---|
366 | ! IN THE CURRENT SET OF XE-POINTS, THERE ARE NEXT(2) TO THE |
---|
367 | ! RIGHT OF X(IR). |
---|
368 | |
---|
369 | IF (IR < N) GO TO 41 |
---|
370 | ! IF (IR .EQ. N) THEN |
---|
371 | ! THESE ARE ACTUALLY EXTRAPOLATION POINTS. |
---|
372 | IERR = IERR + NEXT(2) |
---|
373 | GO TO 42 |
---|
374 | 41 CONTINUE |
---|
375 | ! ELSE |
---|
376 | ! WE SHOULD NEVER HAVE GOTTEN HERE. |
---|
377 | GO TO 5005 |
---|
378 | ! ENDIF |
---|
379 | ! ENDIF |
---|
380 | 42 CONTINUE |
---|
381 | |
---|
382 | IF (NEXT(1) == 0) GO TO 49 |
---|
383 | ! IF (NEXT(1) .GT. 0) THEN |
---|
384 | ! IN THE CURRENT SET OF XE-POINTS, THERE ARE NEXT(1) TO THE |
---|
385 | ! LEFT OF X(IR-1). |
---|
386 | |
---|
387 | IF (IR > 2) GO TO 43 |
---|
388 | ! IF (IR .EQ. 2) THEN |
---|
389 | ! THESE ARE ACTUALLY EXTRAPOLATION POINTS. |
---|
390 | IERR = IERR + NEXT(1) |
---|
391 | GO TO 49 |
---|
392 | 43 CONTINUE |
---|
393 | ! ELSE |
---|
394 | ! XE IS NOT ORDERED RELATIVE TO X, SO MUST ADJUST |
---|
395 | ! EVALUATION INTERVAL. |
---|
396 | |
---|
397 | ! FIRST, LOCATE FIRST POINT TO LEFT OF X(IR-1). |
---|
398 | DO I = JFIRST, J - 1 |
---|
399 | IF (XE(I) < X(IR - 1)) GO TO 45 |
---|
400 | END DO |
---|
401 | ! NOTE-- CANNOT DROP THROUGH HERE UNLESS THERE IS AN ERROR |
---|
402 | ! IN CHFEV. |
---|
403 | GO TO 5005 |
---|
404 | |
---|
405 | 45 CONTINUE |
---|
406 | ! RESET J. (THIS WILL BE THE NEW JFIRST.) |
---|
407 | J = I |
---|
408 | |
---|
409 | ! NOW FIND OUT HOW FAR TO BACK UP IN THE X-ARRAY. |
---|
410 | DO I = 1, IR - 1 |
---|
411 | IF (XE(J) < X(I)) GO TO 47 |
---|
412 | END DO |
---|
413 | ! NB-- CAN NEVER DROP THROUGH HERE, SINCE XE(J).LT.X(IR-1). |
---|
414 | |
---|
415 | 47 CONTINUE |
---|
416 | ! AT THIS POINT, EITHER XE(J) .LT. X(1) |
---|
417 | ! OR X(I-1) .LE. XE(J) .LT. X(I) . |
---|
418 | ! RESET IR, RECOGNIZING THAT IT WILL BE INCREMENTED BEFORE |
---|
419 | ! CYCLING. |
---|
420 | IR = MAX(1, I - 1) |
---|
421 | ! ENDIF |
---|
422 | ! ENDIF |
---|
423 | 49 CONTINUE |
---|
424 | |
---|
425 | JFIRST = J |
---|
426 | |
---|
427 | ! END OF IR-LOOP. |
---|
428 | |
---|
429 | 50 CONTINUE |
---|
430 | IR = IR + 1 |
---|
431 | IF (IR <= N) GO TO 10 |
---|
432 | |
---|
433 | ! NORMAL RETURN. |
---|
434 | |
---|
435 | 5000 CONTINUE |
---|
436 | RETURN |
---|
437 | |
---|
438 | ! ERROR RETURNS. |
---|
439 | |
---|
440 | 5001 CONTINUE |
---|
441 | ! N.LT.2 RETURN. |
---|
442 | IERR = -1 |
---|
443 | CALL XERMSG ('SLATEC', 'PCHFE', & |
---|
444 | 'NUMBER OF DATA POINTS LESS THAN TWO', IERR, 1) |
---|
445 | RETURN |
---|
446 | |
---|
447 | 5002 CONTINUE |
---|
448 | ! INCFD.LT.1 RETURN. |
---|
449 | IERR = -2 |
---|
450 | CALL XERMSG ('SLATEC', 'PCHFE', 'INCREMENT LESS THAN ONE', IERR, & |
---|
451 | 1) |
---|
452 | RETURN |
---|
453 | |
---|
454 | 5003 CONTINUE |
---|
455 | ! X-ARRAY NOT STRICTLY INCREASING. |
---|
456 | IERR = -3 |
---|
457 | CALL XERMSG ('SLATEC', 'PCHFE', 'X-ARRAY NOT STRICTLY INCREASING' & |
---|
458 | , IERR, 1) |
---|
459 | RETURN |
---|
460 | |
---|
461 | 5004 CONTINUE |
---|
462 | ! NE.LT.1 RETURN. |
---|
463 | IERR = -4 |
---|
464 | CALL XERMSG ('SLATEC', 'PCHFE', & |
---|
465 | 'NUMBER OF EVALUATION POINTS LESS THAN ONE', IERR, 1) |
---|
466 | RETURN |
---|
467 | |
---|
468 | 5005 CONTINUE |
---|
469 | ! ERROR RETURN FROM CHFEV. |
---|
470 | ! *** THIS CASE SHOULD NEVER OCCUR *** |
---|
471 | IERR = -5 |
---|
472 | CALL XERMSG ('SLATEC', 'PCHFE', & |
---|
473 | 'ERROR RETURN FROM CHFEV -- FATAL', IERR, 2) |
---|
474 | RETURN |
---|
475 | !------------- LAST LINE OF PCHFE FOLLOWS ------------------------------ |
---|
476 | END SUBROUTINE PCHFE |
---|
477 | |
---|
478 | SUBROUTINE PCHFE_95(X, F, D, SKIP, XE, FE, IERR) |
---|
479 | |
---|
480 | ! PURPOSE Evaluate a piecewise cubic Hermite function at an array of |
---|
481 | ! points. May be used by itself for Hermite interpolation, |
---|
482 | ! or as an evaluator for PCHIM or PCHIC. |
---|
483 | ! CATEGORY E3 |
---|
484 | ! KEYWORDS CUBIC HERMITE EVALUATION, HERMITE INTERPOLATION, PCHIP, |
---|
485 | ! PIECEWISE CUBIC EVALUATION |
---|
486 | |
---|
487 | ! PCHFE: Piecewise Cubic Hermite Function Evaluator |
---|
488 | ! Evaluates the cubic Hermite function defined by X, F, D at |
---|
489 | ! the points XE. |
---|
490 | |
---|
491 | USE lmdz_assert_eq, ONLY: assert_eq |
---|
492 | |
---|
493 | REAL, INTENT(IN) :: X(:) ! real array of independent variable values |
---|
494 | ! The elements of X must be strictly increasing. |
---|
495 | |
---|
496 | REAL, INTENT(IN) :: F(:) ! real array of function values |
---|
497 | ! F(I) is the value corresponding to X(I). |
---|
498 | |
---|
499 | REAL, INTENT(IN) :: D(:) ! real array of derivative values |
---|
500 | ! D(I) is the value corresponding to X(I). |
---|
501 | |
---|
502 | LOGICAL, INTENT(INOUT) :: SKIP |
---|
503 | ! request to skip checks for validity of "x" |
---|
504 | ! If "skip" is false then "pchfe" will check that size(x) >= 2 and |
---|
505 | ! "x" is in strictly ascending order. |
---|
506 | ! Setting "skip" to true will save time in case these checks have |
---|
507 | ! already been performed (say, in "PCHIM" or "PCHIC"). |
---|
508 | ! "SKIP" will be set to TRUE on normal return. |
---|
509 | |
---|
510 | REAL, INTENT(IN) :: XE(:) ! points at which the function is to be evaluated |
---|
511 | ! NOTES: |
---|
512 | ! 1. The evaluation will be most efficient if the elements of XE |
---|
513 | ! are increasing relative to X. |
---|
514 | ! That is, XE(J) .GE. X(I) |
---|
515 | ! implies XE(K) .GE. X(I), all K.GE.J |
---|
516 | ! 2. If any of the XE are outside the interval [X(1),X(N)], values |
---|
517 | ! are extrapolated from the nearest extreme cubic, and a warning |
---|
518 | ! error is returned. |
---|
519 | |
---|
520 | REAL, INTENT(OUT) :: FE(:) ! values of the cubic Hermite function |
---|
521 | ! defined by X, F, D at the points XE |
---|
522 | |
---|
523 | INTEGER, INTENT(OUT) :: IERR ! error flag |
---|
524 | ! Normal return: |
---|
525 | ! IERR = 0 no error |
---|
526 | ! Warning error: |
---|
527 | ! IERR > 0 means that extrapolation was performed at IERR points |
---|
528 | ! "Recoverable" errors: |
---|
529 | ! IERR = -1 if N < 2 |
---|
530 | ! IERR = -3 if the X-array is not strictly increasing |
---|
531 | ! IERR = -4 if NE < 1 |
---|
532 | ! NOTE: The above errors are checked in the order listed, and |
---|
533 | ! following arguments have **NOT** been validated. |
---|
534 | |
---|
535 | ! Variables local to the procedure: |
---|
536 | |
---|
537 | INTEGER N, NE |
---|
538 | |
---|
539 | !--------------------------------------- |
---|
540 | |
---|
541 | n = assert_eq(size(x), size(f), size(d), "PCHFE_95 n") |
---|
542 | ne = assert_eq(size(xe), size(fe), "PCHFE_95 ne") |
---|
543 | CALL PCHFE(N, X, F, D, 1, SKIP, NE, XE, FE, IERR) |
---|
544 | |
---|
545 | END SUBROUTINE PCHFE_95 |
---|
546 | |
---|
547 | FUNCTION pchsp_95(x, f, ibeg, iend, vc_beg, vc_end) |
---|
548 | |
---|
549 | ! PURPOSE: Set derivatives needed to determine the Hermite |
---|
550 | ! representation of the cubic spline interpolant to given data, |
---|
551 | ! with specified boundary conditions. |
---|
552 | |
---|
553 | ! Part of the "pchip" package. |
---|
554 | |
---|
555 | ! CATEGORY: E1A |
---|
556 | |
---|
557 | ! KEYWORDS: cubic hermite interpolation, piecewise cubic |
---|
558 | ! interpolation, spline interpolation |
---|
559 | |
---|
560 | ! DESCRIPTION: "pchsp" stands for "Piecewise Cubic Hermite Spline" |
---|
561 | ! Computes the Hermite representation of the cubic spline |
---|
562 | ! interpolant to the data given in X and F satisfying the boundary |
---|
563 | ! conditions specified by Ibeg, iend, vc_beg and VC_end. |
---|
564 | |
---|
565 | ! The resulting piecewise cubic Hermite function may be evaluated |
---|
566 | ! by "pchfe" or "pchfd". |
---|
567 | |
---|
568 | ! NOTE: This is a modified version of C. de Boor's cubic spline |
---|
569 | ! routine "cubspl". |
---|
570 | |
---|
571 | ! REFERENCE: Carl de Boor, A Practical Guide to Splines, Springer, |
---|
572 | ! 2001, pages 43-47 |
---|
573 | |
---|
574 | USE lmdz_assert_eq, ONLY: assert_eq |
---|
575 | |
---|
576 | REAL, INTENT(IN) :: x(:) |
---|
577 | ! independent variable values |
---|
578 | ! The elements of X must be strictly increasing: |
---|
579 | ! X(I-1) < X(I), I = 2...N. |
---|
580 | ! (Error return if not.) |
---|
581 | ! (error if size(x) < 2) |
---|
582 | |
---|
583 | REAL, INTENT(IN) :: f(:) |
---|
584 | ! dependent variable values to be interpolated |
---|
585 | ! F(I) is value corresponding to X(I). |
---|
586 | |
---|
587 | INTEGER, INTENT(IN) :: ibeg |
---|
588 | ! desired boundary condition at beginning of data |
---|
589 | |
---|
590 | ! IBEG = 0 to set pchsp_95(1) so that the third derivative is con- |
---|
591 | ! tinuous at X(2). This is the "not a knot" condition |
---|
592 | ! provided by de Boor's cubic spline routine CUBSPL. |
---|
593 | ! This is the default boundary condition. |
---|
594 | ! IBEG = 1 if first derivative at X(1) is given in VC_BEG. |
---|
595 | ! IBEG = 2 if second derivative at X(1) is given in VC_BEG. |
---|
596 | ! IBEG = 3 to use the 3-point difference formula for pchsp_95(1). |
---|
597 | ! (Reverts to the default boundary condition if size(x) < 3 .) |
---|
598 | ! IBEG = 4 to use the 4-point difference formula for pchsp_95(1). |
---|
599 | ! (Reverts to the default boundary condition if size(x) < 4 .) |
---|
600 | |
---|
601 | ! NOTES: |
---|
602 | ! 1. An error return is taken if IBEG is out of range. |
---|
603 | ! 2. For the "natural" boundary condition, use IBEG=2 and |
---|
604 | ! VC_BEG=0. |
---|
605 | |
---|
606 | INTEGER, INTENT(IN) :: iend |
---|
607 | ! IC(2) = IEND, desired condition at end of data. |
---|
608 | ! IEND may take on the same values as IBEG, but applied to |
---|
609 | ! derivative at X(N). In case IEND = 1 or 2, The value is given in VC_END. |
---|
610 | |
---|
611 | ! NOTES: |
---|
612 | ! 1. An error return is taken if IEND is out of range. |
---|
613 | ! 2. For the "natural" boundary condition, use IEND=2 and |
---|
614 | ! VC_END=0. |
---|
615 | |
---|
616 | REAL, INTENT(IN), optional :: vc_beg |
---|
617 | ! desired boundary value, as indicated above. |
---|
618 | ! VC_BEG need be set only if IBEG = 1 or 2 . |
---|
619 | |
---|
620 | REAL, INTENT(IN), optional :: vc_end |
---|
621 | ! desired boundary value, as indicated above. |
---|
622 | ! VC_END need be set only if Iend = 1 or 2 . |
---|
623 | |
---|
624 | REAL pchsp_95(size(x)) |
---|
625 | ! derivative values at the data points |
---|
626 | ! These values will determine the cubic spline interpolant |
---|
627 | ! with the requested boundary conditions. |
---|
628 | ! The value corresponding to X(I) is stored in |
---|
629 | ! PCHSP_95(I), I=1...N. |
---|
630 | |
---|
631 | ! LOCAL VARIABLES: |
---|
632 | REAL wk(2, size(x)) ! real array of working storage |
---|
633 | INTEGER n ! number of data points |
---|
634 | INTEGER ierr, ic(2) |
---|
635 | REAL vc(2) |
---|
636 | |
---|
637 | !------------------------------------------------------------------- |
---|
638 | |
---|
639 | n = assert_eq(size(x), size(f), "pchsp_95 n") |
---|
640 | IF ((ibeg == 1 .OR. ibeg == 2) .AND. .NOT. present(vc_beg)) THEN |
---|
641 | PRINT *, "vc_beg required for IBEG = 1 or 2" |
---|
642 | stop 1 |
---|
643 | end if |
---|
644 | IF ((iend == 1 .OR. iend == 2) .AND. .NOT. present(vc_end)) THEN |
---|
645 | PRINT *, "vc_end required for IEND = 1 or 2" |
---|
646 | stop 1 |
---|
647 | end if |
---|
648 | ic = (/ibeg, iend/) |
---|
649 | IF (present(vc_beg)) vc(1) = vc_beg |
---|
650 | IF (present(vc_end)) vc(2) = vc_end |
---|
651 | CALL PCHSP(IC, VC, N, X, F, pchsp_95, 1, WK, size(WK), IERR) |
---|
652 | IF (ierr /= 0) stop 1 |
---|
653 | |
---|
654 | END FUNCTION pchsp_95 |
---|
655 | |
---|
656 | REAL FUNCTION PCHDF(K, X, S, IERR) |
---|
657 | !***BEGIN PROLOGUE PCHDF |
---|
658 | !***SUBSIDIARY |
---|
659 | !***PURPOSE Computes divided differences for PCHCE and PCHSP |
---|
660 | !***LIBRARY SLATEC (PCHIP) |
---|
661 | !***TYPE SINGLE PRECISION (PCHDF-S, DPCHDF-D) |
---|
662 | !***AUTHOR Fritsch, F. N., (LLNL) |
---|
663 | !***DESCRIPTION |
---|
664 | |
---|
665 | ! PCHDF: PCHIP Finite Difference Formula |
---|
666 | |
---|
667 | ! Uses a divided difference formulation to compute a K-point approx- |
---|
668 | ! imation to the derivative at X(K) based on the data in X and S. |
---|
669 | |
---|
670 | ! Called by PCHCE and PCHSP to compute 3- and 4-point boundary |
---|
671 | ! derivative approximations. |
---|
672 | |
---|
673 | ! ---------------------------------------------------------------------- |
---|
674 | |
---|
675 | ! On input: |
---|
676 | ! K is the order of the desired derivative approximation. |
---|
677 | ! K must be at least 3 (error return if not). |
---|
678 | ! X contains the K values of the independent variable. |
---|
679 | ! X need not be ordered, but the values **MUST** be |
---|
680 | ! distinct. (Not checked here.) |
---|
681 | ! S contains the associated slope values: |
---|
682 | ! S(I) = (F(I+1)-F(I))/(X(I+1)-X(I)), I=1(1)K-1. |
---|
683 | ! (Note that S need only be of length K-1.) |
---|
684 | |
---|
685 | ! On return: |
---|
686 | ! S will be destroyed. |
---|
687 | ! IERR will be set to -1 if K.LT.2 . |
---|
688 | ! PCHDF will be set to the desired derivative approximation if |
---|
689 | ! IERR=0 or to zero if IERR=-1. |
---|
690 | |
---|
691 | ! ---------------------------------------------------------------------- |
---|
692 | |
---|
693 | !***SEE ALSO PCHCE, PCHSP |
---|
694 | !***REFERENCES Carl de Boor, A Practical Guide to Splines, Springer- |
---|
695 | ! Verlag, New York, 1978, pp. 10-16. |
---|
696 | !***ROUTINES CALLED XERMSG |
---|
697 | !***REVISION HISTORY (YYMMDD) |
---|
698 | ! 820503 DATE WRITTEN |
---|
699 | ! 820805 Converted to SLATEC library version. |
---|
700 | ! 870813 Minor cosmetic changes. |
---|
701 | ! 890411 Added SAVE statements (Vers. 3.2). |
---|
702 | ! 890411 REVISION DATE from Version 3.2 |
---|
703 | ! 891214 Prologue converted to Version 4.0 format. (BAB) |
---|
704 | ! 900315 CALLs to XERROR changed to CALLs to XERMSG. (THJ) |
---|
705 | ! 900328 Added TYPE section. (WRB) |
---|
706 | ! 910408 Updated AUTHOR and DATE WRITTEN sections in prologue. (WRB) |
---|
707 | ! 920429 Revised format and order of references. (WRB,FNF) |
---|
708 | ! 930503 Improved purpose. (FNF) |
---|
709 | !***END PROLOGUE PCHDF |
---|
710 | |
---|
711 | !**End |
---|
712 | |
---|
713 | ! DECLARE ARGUMENTS. |
---|
714 | |
---|
715 | INTEGER :: K, IERR |
---|
716 | REAL :: X(K), S(K) |
---|
717 | |
---|
718 | ! DECLARE LOCAL VARIABLES. |
---|
719 | |
---|
720 | INTEGER :: I, J |
---|
721 | REAL :: VALUE, ZERO |
---|
722 | SAVE ZERO |
---|
723 | DATA ZERO /0./ |
---|
724 | |
---|
725 | ! CHECK FOR LEGAL VALUE OF K. |
---|
726 | |
---|
727 | !***FIRST EXECUTABLE STATEMENT PCHDF |
---|
728 | IF (K < 3) GO TO 5001 |
---|
729 | |
---|
730 | ! COMPUTE COEFFICIENTS OF INTERPOLATING POLYNOMIAL. |
---|
731 | |
---|
732 | DO J = 2, K - 1 |
---|
733 | DO I = 1, K - J |
---|
734 | S(I) = (S(I + 1) - S(I)) / (X(I + J) - X(I)) |
---|
735 | END DO |
---|
736 | END DO |
---|
737 | |
---|
738 | ! EVALUATE DERIVATIVE AT X(K). |
---|
739 | |
---|
740 | VALUE = S(1) |
---|
741 | DO I = 2, K - 1 |
---|
742 | VALUE = S(I) + VALUE * (X(K) - X(I)) |
---|
743 | END DO |
---|
744 | |
---|
745 | ! NORMAL RETURN. |
---|
746 | |
---|
747 | IERR = 0 |
---|
748 | PCHDF = VALUE |
---|
749 | RETURN |
---|
750 | |
---|
751 | ! ERROR RETURN. |
---|
752 | |
---|
753 | 5001 CONTINUE |
---|
754 | ! K.LT.3 RETURN. |
---|
755 | IERR = -1 |
---|
756 | CALL XERMSG ('SLATEC', 'PCHDF', 'K LESS THAN THREE', IERR, 1) |
---|
757 | PCHDF = ZERO |
---|
758 | RETURN |
---|
759 | !------------- LAST LINE OF PCHDF FOLLOWS ------------------------------ |
---|
760 | END FUNCTION PCHDF |
---|
761 | |
---|
762 | SUBROUTINE PCHSP(IC, VC, N, X, F, D, INCFD, WK, NWK, IERR) |
---|
763 | !***BEGIN PROLOGUE PCHSP |
---|
764 | !***PURPOSE Set derivatives needed to determine the Hermite represen- |
---|
765 | ! tation of the cubic spline interpolant to given data, with |
---|
766 | ! specified boundary conditions. |
---|
767 | !***LIBRARY SLATEC (PCHIP) |
---|
768 | !***CATEGORY E1A |
---|
769 | !***TYPE SINGLE PRECISION (PCHSP-S, DPCHSP-D) |
---|
770 | !***KEYWORDS CUBIC HERMITE INTERPOLATION, PCHIP, |
---|
771 | ! PIECEWISE CUBIC INTERPOLATION, SPLINE INTERPOLATION |
---|
772 | !***AUTHOR Fritsch, F. N., (LLNL) |
---|
773 | ! Lawrence Livermore National Laboratory |
---|
774 | ! P.O. Box 808 (L-316) |
---|
775 | ! Livermore, CA 94550 |
---|
776 | ! FTS 532-4275, (510) 422-4275 |
---|
777 | !***DESCRIPTION |
---|
778 | |
---|
779 | ! PCHSP: Piecewise Cubic Hermite Spline |
---|
780 | |
---|
781 | ! Computes the Hermite representation of the cubic spline inter- |
---|
782 | ! polant to the data given in X and F satisfying the boundary |
---|
783 | ! conditions specified by IC and VC. |
---|
784 | |
---|
785 | ! To facilitate two-dimensional applications, includes an increment |
---|
786 | ! between successive values of the F- and D-arrays. |
---|
787 | |
---|
788 | ! The resulting piecewise cubic Hermite function may be evaluated |
---|
789 | ! by PCHFE or PCHFD. |
---|
790 | |
---|
791 | ! NOTE: This is a modified version of C. de Boor's cubic spline |
---|
792 | ! routine CUBSPL. |
---|
793 | |
---|
794 | ! ---------------------------------------------------------------------- |
---|
795 | |
---|
796 | ! Calling sequence: |
---|
797 | |
---|
798 | ! PARAMETER (INCFD = ...) |
---|
799 | ! INTEGER IC(2), N, NWK, IERR |
---|
800 | ! REAL VC(2), X(N), F(INCFD,N), D(INCFD,N), WK(NWK) |
---|
801 | |
---|
802 | ! CALL PCHSP (IC, VC, N, X, F, D, INCFD, WK, NWK, IERR) |
---|
803 | |
---|
804 | ! Parameters: |
---|
805 | |
---|
806 | ! IC -- (input) integer array of length 2 specifying desired |
---|
807 | ! boundary conditions: |
---|
808 | ! IC(1) = IBEG, desired condition at beginning of data. |
---|
809 | ! IC(2) = IEND, desired condition at end of data. |
---|
810 | |
---|
811 | ! IBEG = 0 to set D(1) so that the third derivative is con- |
---|
812 | ! tinuous at X(2). This is the "not a knot" condition |
---|
813 | ! provided by de Boor's cubic spline routine CUBSPL. |
---|
814 | ! < This is the default boundary condition. > |
---|
815 | ! IBEG = 1 if first derivative at X(1) is given in VC(1). |
---|
816 | ! IBEG = 2 if second derivative at X(1) is given in VC(1). |
---|
817 | ! IBEG = 3 to use the 3-point difference formula for D(1). |
---|
818 | ! (Reverts to the default b.c. if N.LT.3 .) |
---|
819 | ! IBEG = 4 to use the 4-point difference formula for D(1). |
---|
820 | ! (Reverts to the default b.c. if N.LT.4 .) |
---|
821 | ! NOTES: |
---|
822 | ! 1. An error return is taken if IBEG is out of range. |
---|
823 | ! 2. For the "natural" boundary condition, use IBEG=2 and |
---|
824 | ! VC(1)=0. |
---|
825 | |
---|
826 | ! IEND may take on the same values as IBEG, but applied to |
---|
827 | ! derivative at X(N). In case IEND = 1 or 2, the value is |
---|
828 | ! given in VC(2). |
---|
829 | |
---|
830 | ! NOTES: |
---|
831 | ! 1. An error return is taken if IEND is out of range. |
---|
832 | ! 2. For the "natural" boundary condition, use IEND=2 and |
---|
833 | ! VC(2)=0. |
---|
834 | |
---|
835 | ! VC -- (input) real array of length 2 specifying desired boundary |
---|
836 | ! values, as indicated above. |
---|
837 | ! VC(1) need be set only if IC(1) = 1 or 2 . |
---|
838 | ! VC(2) need be set only if IC(2) = 1 or 2 . |
---|
839 | |
---|
840 | ! N -- (input) number of data points. (Error return if N.LT.2 .) |
---|
841 | |
---|
842 | ! X -- (input) real array of independent variable values. The |
---|
843 | ! elements of X must be strictly increasing: |
---|
844 | ! X(I-1) .LT. X(I), I = 2(1)N. |
---|
845 | ! (Error return if not.) |
---|
846 | |
---|
847 | ! F -- (input) real array of dependent variable values to be inter- |
---|
848 | ! polated. F(1+(I-1)*INCFD) is value corresponding to X(I). |
---|
849 | |
---|
850 | ! D -- (output) real array of derivative values at the data points. |
---|
851 | ! These values will determine the cubic spline interpolant |
---|
852 | ! with the requested boundary conditions. |
---|
853 | ! The value corresponding to X(I) is stored in |
---|
854 | ! D(1+(I-1)*INCFD), I=1(1)N. |
---|
855 | ! No other entries in D are changed. |
---|
856 | |
---|
857 | ! INCFD -- (input) increment between successive values in F and D. |
---|
858 | ! This argument is provided primarily for 2-D applications. |
---|
859 | ! (Error return if INCFD.LT.1 .) |
---|
860 | |
---|
861 | ! WK -- (scratch) real array of working storage. |
---|
862 | |
---|
863 | ! NWK -- (input) length of work array. |
---|
864 | ! (Error return if NWK.LT.2*N .) |
---|
865 | |
---|
866 | ! IERR -- (output) error flag. |
---|
867 | ! Normal return: |
---|
868 | ! IERR = 0 (no errors). |
---|
869 | ! "Recoverable" errors: |
---|
870 | ! IERR = -1 if N.LT.2 . |
---|
871 | ! IERR = -2 if INCFD.LT.1 . |
---|
872 | ! IERR = -3 if the X-array is not strictly increasing. |
---|
873 | ! IERR = -4 if IBEG.LT.0 or IBEG.GT.4 . |
---|
874 | ! IERR = -5 if IEND.LT.0 of IEND.GT.4 . |
---|
875 | ! IERR = -6 if both of the above are true. |
---|
876 | ! IERR = -7 if NWK is too small. |
---|
877 | ! NOTE: The above errors are checked in the order listed, |
---|
878 | ! and following arguments have **NOT** been validated. |
---|
879 | ! (The D-array has not been changed in any of these cases.) |
---|
880 | ! IERR = -8 in case of trouble solving the linear system |
---|
881 | ! for the interior derivative values. |
---|
882 | ! (The D-array may have been changed in this case.) |
---|
883 | ! ( Do **NOT** use it! ) |
---|
884 | |
---|
885 | !***REFERENCES Carl de Boor, A Practical Guide to Splines, Springer- |
---|
886 | ! Verlag, New York, 1978, pp. 53-59. |
---|
887 | !***ROUTINES CALLED PCHDF, XERMSG |
---|
888 | !***REVISION HISTORY (YYMMDD) |
---|
889 | ! 820503 DATE WRITTEN |
---|
890 | ! 820804 Converted to SLATEC library version. |
---|
891 | ! 870707 Minor cosmetic changes to prologue. |
---|
892 | ! 890411 Added SAVE statements (Vers. 3.2). |
---|
893 | ! 890703 Corrected category record. (WRB) |
---|
894 | ! 890831 Modified array declarations. (WRB) |
---|
895 | ! 890831 REVISION DATE from Version 3.2 |
---|
896 | ! 891214 Prologue converted to Version 4.0 format. (BAB) |
---|
897 | ! 900315 CALLs to XERROR changed to CALLs to XERMSG. (THJ) |
---|
898 | ! 920429 Revised format and order of references. (WRB,FNF) |
---|
899 | !***END PROLOGUE PCHSP |
---|
900 | ! Programming notes: |
---|
901 | |
---|
902 | ! To produce a double precision version, simply: |
---|
903 | ! a. Change PCHSP to DPCHSP wherever it occurs, |
---|
904 | ! b. Change the real declarations to double precision, and |
---|
905 | ! c. Change the constants ZERO, HALF, ... to double precision. |
---|
906 | |
---|
907 | ! DECLARE ARGUMENTS. |
---|
908 | |
---|
909 | INTEGER :: IC(2), N, INCFD, NWK, IERR |
---|
910 | REAL :: VC(2), X(*), F(INCFD, *), D(INCFD, *), WK(2, *) |
---|
911 | |
---|
912 | ! DECLARE LOCAL VARIABLES. |
---|
913 | |
---|
914 | INTEGER :: IBEG, IEND, INDEX, J, NM1 |
---|
915 | REAL :: G, HALF, ONE, STEMP(3), THREE, TWO, XTEMP(4), ZERO |
---|
916 | SAVE ZERO, HALF, ONE, TWO, THREE |
---|
917 | |
---|
918 | DATA ZERO /0./, HALF /0.5/, ONE /1./, TWO /2./, THREE /3./ |
---|
919 | |
---|
920 | ! VALIDITY-CHECK ARGUMENTS. |
---|
921 | |
---|
922 | !***FIRST EXECUTABLE STATEMENT PCHSP |
---|
923 | IF (N<2) GO TO 5001 |
---|
924 | IF (INCFD<1) GO TO 5002 |
---|
925 | DO J = 2, N |
---|
926 | IF (X(J)<=X(J - 1)) GO TO 5003 |
---|
927 | END DO |
---|
928 | |
---|
929 | IBEG = IC(1) |
---|
930 | IEND = IC(2) |
---|
931 | IERR = 0 |
---|
932 | IF ((IBEG<0).OR.(IBEG>4)) IERR = IERR - 1 |
---|
933 | IF ((IEND<0).OR.(IEND>4)) IERR = IERR - 2 |
---|
934 | IF (IERR<0) GO TO 5004 |
---|
935 | |
---|
936 | ! FUNCTION DEFINITION IS OK -- GO ON. |
---|
937 | |
---|
938 | IF (NWK < 2 * N) GO TO 5007 |
---|
939 | |
---|
940 | ! COMPUTE FIRST DIFFERENCES OF X SEQUENCE AND STORE IN WK(1,.). ALSO, |
---|
941 | ! COMPUTE FIRST DIVIDED DIFFERENCE OF DATA AND STORE IN WK(2,.). |
---|
942 | DO J = 2, N |
---|
943 | WK(1, J) = X(J) - X(J - 1) |
---|
944 | WK(2, J) = (F(1, J) - F(1, J - 1)) / WK(1, J) |
---|
945 | END DO |
---|
946 | |
---|
947 | ! SET TO DEFAULT BOUNDARY CONDITIONS IF N IS TOO SMALL. |
---|
948 | |
---|
949 | IF (IBEG>N) IBEG = 0 |
---|
950 | IF (IEND>N) IEND = 0 |
---|
951 | |
---|
952 | ! SET UP FOR BOUNDARY CONDITIONS. |
---|
953 | |
---|
954 | IF ((IBEG==1).OR.(IBEG==2)) THEN |
---|
955 | D(1, 1) = VC(1) |
---|
956 | ELSE IF (IBEG > 2) THEN |
---|
957 | ! PICK UP FIRST IBEG POINTS, IN REVERSE ORDER. |
---|
958 | DO J = 1, IBEG |
---|
959 | INDEX = IBEG - J + 1 |
---|
960 | ! INDEX RUNS FROM IBEG DOWN TO 1. |
---|
961 | XTEMP(J) = X(INDEX) |
---|
962 | IF (J < IBEG) STEMP(J) = WK(2, INDEX) |
---|
963 | END DO |
---|
964 | ! -------------------------------- |
---|
965 | D(1, 1) = PCHDF (IBEG, XTEMP, STEMP, IERR) |
---|
966 | ! -------------------------------- |
---|
967 | IF (IERR /= 0) GO TO 5009 |
---|
968 | IBEG = 1 |
---|
969 | ENDIF |
---|
970 | |
---|
971 | IF ((IEND==1).OR.(IEND==2)) THEN |
---|
972 | D(1, N) = VC(2) |
---|
973 | ELSE IF (IEND > 2) THEN |
---|
974 | ! PICK UP LAST IEND POINTS. |
---|
975 | DO J = 1, IEND |
---|
976 | INDEX = N - IEND + J |
---|
977 | ! INDEX RUNS FROM N+1-IEND UP TO N. |
---|
978 | XTEMP(J) = X(INDEX) |
---|
979 | IF (J < IEND) STEMP(J) = WK(2, INDEX + 1) |
---|
980 | END DO |
---|
981 | ! -------------------------------- |
---|
982 | D(1, N) = PCHDF (IEND, XTEMP, STEMP, IERR) |
---|
983 | ! -------------------------------- |
---|
984 | IF (IERR /= 0) GO TO 5009 |
---|
985 | IEND = 1 |
---|
986 | ENDIF |
---|
987 | |
---|
988 | ! --------------------( BEGIN CODING FROM CUBSPL )-------------------- |
---|
989 | |
---|
990 | ! **** A TRIDIAGONAL LINEAR SYSTEM FOR THE UNKNOWN SLOPES S(J) OF |
---|
991 | ! F AT X(J), J=1,...,N, IS GENERATED AND THEN SOLVED BY GAUSS ELIM- |
---|
992 | ! INATION, WITH S(J) ENDING UP IN D(1,J), ALL J. |
---|
993 | ! WK(1,.) AND WK(2,.) ARE USED FOR TEMPORARY STORAGE. |
---|
994 | |
---|
995 | ! CONSTRUCT FIRST EQUATION FROM FIRST BOUNDARY CONDITION, OF THE FORM |
---|
996 | ! WK(2,1)*S(1) + WK(1,1)*S(2) = D(1,1) |
---|
997 | |
---|
998 | IF (IBEG == 0) THEN |
---|
999 | IF (N == 2) THEN |
---|
1000 | ! NO CONDITION AT LEFT END AND N = 2. |
---|
1001 | WK(2, 1) = ONE |
---|
1002 | WK(1, 1) = ONE |
---|
1003 | D(1, 1) = TWO * WK(2, 2) |
---|
1004 | ELSE |
---|
1005 | ! NOT-A-KNOT CONDITION AT LEFT END AND N .GT. 2. |
---|
1006 | WK(2, 1) = WK(1, 3) |
---|
1007 | WK(1, 1) = WK(1, 2) + WK(1, 3) |
---|
1008 | D(1, 1) = ((WK(1, 2) + TWO * WK(1, 1)) * WK(2, 2) * WK(1, 3) & |
---|
1009 | + WK(1, 2)**2 * WK(2, 3)) / WK(1, 1) |
---|
1010 | ENDIF |
---|
1011 | ELSE IF (IBEG == 1) THEN |
---|
1012 | ! SLOPE PRESCRIBED AT LEFT END. |
---|
1013 | WK(2, 1) = ONE |
---|
1014 | WK(1, 1) = ZERO |
---|
1015 | ELSE |
---|
1016 | ! SECOND DERIVATIVE PRESCRIBED AT LEFT END. |
---|
1017 | WK(2, 1) = TWO |
---|
1018 | WK(1, 1) = ONE |
---|
1019 | D(1, 1) = THREE * WK(2, 2) - HALF * WK(1, 2) * D(1, 1) |
---|
1020 | ENDIF |
---|
1021 | |
---|
1022 | ! IF THERE ARE INTERIOR KNOTS, GENERATE THE CORRESPONDING EQUATIONS AND |
---|
1023 | ! CARRY OUT THE FORWARD PASS OF GAUSS ELIMINATION, AFTER WHICH THE J-TH |
---|
1024 | ! EQUATION READS WK(2,J)*S(J) + WK(1,J)*S(J+1) = D(1,J). |
---|
1025 | |
---|
1026 | NM1 = N - 1 |
---|
1027 | IF (NM1 > 1) THEN |
---|
1028 | DO J = 2, NM1 |
---|
1029 | IF (WK(2, J - 1) == ZERO) GO TO 5008 |
---|
1030 | G = -WK(1, J + 1) / WK(2, J - 1) |
---|
1031 | D(1, J) = G * D(1, J - 1) & |
---|
1032 | + THREE * (WK(1, J) * WK(2, J + 1) + WK(1, J + 1) * WK(2, J)) |
---|
1033 | WK(2, J) = G * WK(1, J - 1) + TWO * (WK(1, J) + WK(1, J + 1)) |
---|
1034 | END DO |
---|
1035 | ENDIF |
---|
1036 | |
---|
1037 | ! CONSTRUCT LAST EQUATION FROM SECOND BOUNDARY CONDITION, OF THE FORM |
---|
1038 | ! (-G*WK(2,N-1))*S(N-1) + WK(2,N)*S(N) = D(1,N) |
---|
1039 | |
---|
1040 | ! IF SLOPE IS PRESCRIBED AT RIGHT END, ONE CAN GO DIRECTLY TO BACK- |
---|
1041 | ! SUBSTITUTION, SINCE ARRAYS HAPPEN TO BE SET UP JUST RIGHT FOR IT |
---|
1042 | ! AT THIS POINT. |
---|
1043 | IF (IEND == 1) GO TO 30 |
---|
1044 | |
---|
1045 | IF (IEND == 0) THEN |
---|
1046 | IF (N==2 .AND. IBEG==0) THEN |
---|
1047 | ! NOT-A-KNOT AT RIGHT ENDPOINT AND AT LEFT ENDPOINT AND N = 2. |
---|
1048 | D(1, 2) = WK(2, 2) |
---|
1049 | GO TO 30 |
---|
1050 | ELSE IF ((N==2) .OR. (N==3 .AND. IBEG==0)) THEN |
---|
1051 | ! EITHER (N=3 AND NOT-A-KNOT ALSO AT LEFT) OR (N=2 AND *NOT* |
---|
1052 | ! NOT-A-KNOT AT LEFT END POINT). |
---|
1053 | D(1, N) = TWO * WK(2, N) |
---|
1054 | WK(2, N) = ONE |
---|
1055 | IF (WK(2, N - 1) == ZERO) GO TO 5008 |
---|
1056 | G = -ONE / WK(2, N - 1) |
---|
1057 | ELSE |
---|
1058 | ! NOT-A-KNOT AND N .GE. 3, AND EITHER N.GT.3 OR ALSO NOT-A- |
---|
1059 | ! KNOT AT LEFT END POINT. |
---|
1060 | G = WK(1, N - 1) + WK(1, N) |
---|
1061 | ! DO NOT NEED TO CHECK FOLLOWING DENOMINATORS (X-DIFFERENCES). |
---|
1062 | D(1, N) = ((WK(1, N) + TWO * G) * WK(2, N) * WK(1, N - 1) & |
---|
1063 | + WK(1, N)**2 * (F(1, N - 1) - F(1, N - 2)) / WK(1, N - 1)) / G |
---|
1064 | IF (WK(2, N - 1) == ZERO) GO TO 5008 |
---|
1065 | G = -G / WK(2, N - 1) |
---|
1066 | WK(2, N) = WK(1, N - 1) |
---|
1067 | ENDIF |
---|
1068 | ELSE |
---|
1069 | ! SECOND DERIVATIVE PRESCRIBED AT RIGHT ENDPOINT. |
---|
1070 | D(1, N) = THREE * WK(2, N) + HALF * WK(1, N) * D(1, N) |
---|
1071 | WK(2, N) = TWO |
---|
1072 | IF (WK(2, N - 1) == ZERO) GO TO 5008 |
---|
1073 | G = -ONE / WK(2, N - 1) |
---|
1074 | ENDIF |
---|
1075 | |
---|
1076 | ! COMPLETE FORWARD PASS OF GAUSS ELIMINATION. |
---|
1077 | |
---|
1078 | WK(2, N) = G * WK(1, N - 1) + WK(2, N) |
---|
1079 | IF (WK(2, N) == ZERO) GO TO 5008 |
---|
1080 | D(1, N) = (G * D(1, N - 1) + D(1, N)) / WK(2, N) |
---|
1081 | |
---|
1082 | ! CARRY OUT BACK SUBSTITUTION |
---|
1083 | |
---|
1084 | 30 CONTINUE |
---|
1085 | DO J = NM1, 1, -1 |
---|
1086 | IF (WK(2, J) == ZERO) GO TO 5008 |
---|
1087 | D(1, J) = (D(1, J) - WK(1, J) * D(1, J + 1)) / WK(2, J) |
---|
1088 | END DO |
---|
1089 | ! --------------------( END CODING FROM CUBSPL )-------------------- |
---|
1090 | |
---|
1091 | ! NORMAL RETURN. |
---|
1092 | |
---|
1093 | RETURN |
---|
1094 | |
---|
1095 | ! ERROR RETURNS. |
---|
1096 | |
---|
1097 | 5001 CONTINUE |
---|
1098 | ! N.LT.2 RETURN. |
---|
1099 | IERR = -1 |
---|
1100 | CALL XERMSG ('SLATEC', 'PCHSP', & |
---|
1101 | 'NUMBER OF DATA POINTS LESS THAN TWO', IERR, 1) |
---|
1102 | RETURN |
---|
1103 | |
---|
1104 | 5002 CONTINUE |
---|
1105 | ! INCFD.LT.1 RETURN. |
---|
1106 | IERR = -2 |
---|
1107 | CALL XERMSG ('SLATEC', 'PCHSP', 'INCREMENT LESS THAN ONE', IERR, & |
---|
1108 | 1) |
---|
1109 | RETURN |
---|
1110 | |
---|
1111 | 5003 CONTINUE |
---|
1112 | ! X-ARRAY NOT STRICTLY INCREASING. |
---|
1113 | IERR = -3 |
---|
1114 | CALL XERMSG ('SLATEC', 'PCHSP', 'X-ARRAY NOT STRICTLY INCREASING' & |
---|
1115 | , IERR, 1) |
---|
1116 | RETURN |
---|
1117 | |
---|
1118 | 5004 CONTINUE |
---|
1119 | ! IC OUT OF RANGE RETURN. |
---|
1120 | IERR = IERR - 3 |
---|
1121 | CALL XERMSG ('SLATEC', 'PCHSP', 'IC OUT OF RANGE', IERR, 1) |
---|
1122 | RETURN |
---|
1123 | |
---|
1124 | 5007 CONTINUE |
---|
1125 | ! NWK TOO SMALL RETURN. |
---|
1126 | IERR = -7 |
---|
1127 | CALL XERMSG ('SLATEC', 'PCHSP', 'WORK ARRAY TOO SMALL', IERR, 1) |
---|
1128 | RETURN |
---|
1129 | |
---|
1130 | 5008 CONTINUE |
---|
1131 | ! SINGULAR SYSTEM. |
---|
1132 | ! *** THEORETICALLY, THIS CAN ONLY OCCUR IF SUCCESSIVE X-VALUES *** |
---|
1133 | ! *** ARE EQUAL, WHICH SHOULD ALREADY HAVE BEEN CAUGHT (IERR=-3). *** |
---|
1134 | IERR = -8 |
---|
1135 | CALL XERMSG ('SLATEC', 'PCHSP', 'SINGULAR LINEAR SYSTEM', IERR, & |
---|
1136 | 1) |
---|
1137 | RETURN |
---|
1138 | |
---|
1139 | 5009 CONTINUE |
---|
1140 | ! ERROR RETURN FROM PCHDF. |
---|
1141 | ! *** THIS CASE SHOULD NEVER OCCUR *** |
---|
1142 | IERR = -9 |
---|
1143 | CALL XERMSG ('SLATEC', 'PCHSP', 'ERROR RETURN FROM PCHDF', IERR, & |
---|
1144 | 1) |
---|
1145 | RETURN |
---|
1146 | !------------- LAST LINE OF PCHSP FOLLOWS ------------------------------ |
---|
1147 | END SUBROUTINE PCHSP |
---|
1148 | |
---|
1149 | |
---|
1150 | END MODULE lmdz_libmath_pch |
---|