-
Notifications
You must be signed in to change notification settings - Fork 11
/
Copy pathtoms760.f90
1350 lines (1203 loc) · 41.4 KB
/
toms760.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
MODULE Grid_Interpolation
! ALGORITHM 760, COLLECTED ALGORITHMS FROM ACM.
! THIS WORK PUBLISHED IN TRANSACTIONS ON MATHEMATICAL SOFTWARE,
! VOL. 22, NO. 3, September, 1996, P. 357--361.
!AL Jun 2014: fixed f90 reshape
IMPLICIT NONE
Type Grid_Interpolation_Data
REAL, pointer :: wk(:,:,:) => NULL()
REAL, pointer :: x(:), y(:)
integer nx, ny
end Type Grid_Interpolation_Data
integer, parameter :: GI = KIND(1.0)
CONTAINS
!Utility wrappers by AL, 2012
subroutine Grid_Interpolation_Init(W, x,y, z)
Type(Grid_Interpolation_Data):: W
REAL, INTENT(IN) :: x(:)
REAL, INTENT(IN) :: y(:)
REAL, INTENT(IN) :: z(:,:)
W%nx = size(x)
W%ny = size(y)
allocate(W%Wk(3,W%NX,W%NY))
allocate(W%x(W%nx))
allocate(W%y(W%ny))
W%x = x
W%y = y
CALL rgpd3p(W%nx, W%ny, x, y, z, W%wk)
end subroutine Grid_Interpolation_Init
subroutine Grid_Interpolation_Free(W)
Type(Grid_Interpolation_Data):: W
deallocate(W%x)
deallocate(W%y)
deallocate(W%Wk)
end subroutine Grid_Interpolation_Free
function Grid_Interpolate(W,zin, x,y) result (res)
Type(Grid_Interpolation_Data) :: W
real res, z(1), xx(1),yy(1)
real, intent(in) :: x,y, zin(:,:)
integer md,ier
md=2
xx(1)=x
yy(1)=y
call rgbi3p(W%Wk,md, W%nx, W%ny, W%x, W%y, zin, 1, xx, yy, z, ier)
if (ier/=0) stop 'Grid_Interpolate error'
res = z(1)
end function Grid_Interpolate
subroutine Grid_InterpolatePoints(W, zin, nip, x,y,z)
Type(Grid_Interpolation_Data) :: W
integer, intent(in) :: nip
real, intent(out):: z(*)
real, intent(in) :: x(*),y(*), zin(:,:)
integer md,ier
md=2
call rgbi3p(W%Wk,md, W%nx, W%ny, W%x, W%y, zin, nip, x, y, z, ier)
if (ier/=0) stop 'Grid_Interpolate error'
end subroutine Grid_InterpolatePoints
SUBROUTINE rgbi3p(Wk,md, nxd, nyd, xd, yd, zd, nip, xi, yi, zi, ier)
! Code converted using TO_F90 by Alan Miller
! Date: 2003-06-11 Time: 10:11:03
! Rectangular-grid bivariate interpolation
! (a master subroutine of the RGBI3P/RGSF3P subroutine package)
! Hiroshi Akima
! U.S. Department of Commerce, NTIA/ITS
! Version of 1995/08
! This subroutine performs interpolation of a bivariate function, z(x,y), on a
! rectangular grid in the x-y plane. It is based on the revised Akima method.
! In this subroutine, the interpolating function is a piecewise function
! composed of a set of bicubic (bivariate third-degree) polynomials, each
! applicable to a rectangle of the input grid in the x-y plane.
! Each polynomial is determined locally.
! This subroutine has the accuracy of a bicubic polynomial, i.e., it
! interpolates accurately when all data points lie on a surface of a
! bicubic polynomial.
! The grid lines can be unevenly spaced.
! The input arguments are
! MD = mode of computation
! = 1 for new XD, YD, or ZD data (default)
! = 2 for old XD, YD, and ZD data,
! NXD = number of the input-grid data points in the x coordinate
! (must be 2 or greater),
! NYD = number of the input-grid data points in the y coordinate
! (must be 2 or greater),
! XD = array of dimension NXD containing the x coordinates of the
! input-grid data points (must be in a monotonic increasing order),
! YD = array of dimension NYD containing the y coordinates of the
! input-grid data points (must be in a monotonic increasing order),
! ZD = two-dimensional array of dimension NXD*NYD
! containing the z(x,y) values at the input-grid data points,
! NIP = number of the output points at which interpolation
! of the z value is desired (must be 1 or greater),
! XI = array of dimension NIP containing the x coordinates
! of the output points,
! YI = array of dimension NIP containing the y coordinates
! of the output points.
! The output arguments are
! ZI = array of dimension NIP where the interpolated z
! values at the output points are to be stored,
! IER = error flag
! = 0 for no errors
! = 1 for NXD = 1 or less
! = 2 for NYD = 1 or less
! = 3 for identical XD values or XD values out of sequence
! = 4 for identical YD values or YD values out of sequence
! = 5 for NIP = 0 or less.
! N.B. The workspace has been removed from the argument list.
! WK = three dimensional array of dimension 3*NXD*NYD used internally
! as a work area.
! The very fisrt call to this subroutine and the call with a new XD, YD, and
! ZD array must be made with MD=1. The call with MD=2 must be preceded by
! another call with the same XD, YD, and ZD arrays. Between the call with
! MD=2 and its preceding call, the WK array must not be disturbed.
! The constant in the PARAMETER statement below is
! NIPIMX = maximum number of output points to be processed at a time.
! The constant value has been selected empirically.
! This subroutine calls the RGPD3P, RGLCTN, and RGPLNL subroutines.
! Specification statements
! .. Parameters ..
INTEGER, INTENT(IN) :: md
INTEGER, INTENT(IN) :: nxd
INTEGER, INTENT(IN) :: nyd
REAL, INTENT(IN) :: xd(nxd)
REAL, INTENT(IN) :: yd(nyd)
REAL, INTENT(IN) :: zd(nxd,nyd)
INTEGER, INTENT(IN) :: nip
REAL, INTENT(IN) :: xi(nip)
REAL, INTENT(IN) :: yi(nip)
REAL, INTENT(OUT) :: zi(nip)
INTEGER, INTENT(OUT) :: ier
REAL, INTENT(INOUT) :: wk(3,nxd,nyd)
! ..
! .. Local Scalars ..
INTEGER, PARAMETER :: nipimx=51
INTEGER :: iip, ix, iy, nipi
! ..
! .. Local Arrays ..
INTEGER :: inxi(nipimx), inyi(nipimx)
! ..
! .. External Subroutines ..
! EXTERNAL rglctn, rgpd3p, rgplnl
! ..
! .. Intrinsic Functions ..
! INTRINSIC MIN
! ..
! Preliminary processing
! Error check
IF (nxd <= 1) GO TO 40
IF (nyd <= 1) GO TO 50
DO ix = 2,nxd
IF (xd(ix) <= xd(ix-1)) GO TO 60
END DO
DO iy = 2,nyd
IF (yd(iy) <= yd(iy-1)) GO TO 70
END DO
IF (nip <= 0) GO TO 80
ier = 0
! Calculation
! Estimates partial derivatives at all input-grid data points (for MD=1).
IF (md /= 2) THEN
CALL rgpd3p(nxd, nyd, xd, yd, zd, wk)
END IF
! DO-loop with respect to the output point
! Processes NIPIMX output points, at most, at a time.
DO iip = 1,nip,nipimx
nipi = MIN(nip- (iip-1),nipimx)
! Locates the output points.
CALL rglctn(nxd, nyd, xd, yd, nipi, xi(iip), yi(iip), inxi, inyi)
! Calculates the z values at the output points.
CALL rgplnl(nxd, nyd, xd, yd, zd, wk, nipi, xi(iip), yi(iip), inxi, inyi, &
zi(iip))
END DO
RETURN
! Error exit
40 WRITE (*,FMT=9000)
ier = 1
GO TO 90
50 WRITE (*,FMT=9010)
ier = 2
GO TO 90
60 WRITE (*,FMT=9020) ix,xd(ix)
ier = 3
GO TO 90
70 WRITE (*,FMT=9030) iy,yd(iy)
ier = 4
GO TO 90
80 WRITE (*,FMT=9040)
ier = 5
90 WRITE (*,FMT=9050) nxd,nyd,nip
RETURN
! Format statements for error messages
9000 FORMAT (/' *** RGBI3P Error 1: NXD = 1 or less')
9010 FORMAT (/' *** RGBI3P Error 2: NYD = 1 or less')
9020 FORMAT (/' *** RGBI3P Error 3: Identical XD values or', &
' XD values out of sequence'/ ' IX =', i6, ', XD(IX) =', e11.3)
9030 FORMAT (/' *** RGBI3P Error 4: Identical YD values or', &
' YD values out of sequence',/,' IY =',i6,', YD(IY) =', e11.3)
9040 FORMAT (/' *** RGBI3P Error 5: NIP = 0 or less')
9050 FORMAT (' NXD =', i5,', NYD =', i5,', NIP =', i5/)
END SUBROUTINE rgbi3p
SUBROUTINE rgsf3p(wk, md, nxd, nyd, xd, yd, zd, nxi, xi, nyi, yi, zi, ier)
! Rectangular-grid surface fitting
! (a master subroutine of the RGBI3P/RGSF3P subroutine package)
! Hiroshi Akima
! U.S. Department of Commerce, NTIA/ITS
! Version of 1995/08
! This subroutine performs surface fitting by interpolating values of a
! bivariate function, z(x,y), on a rectangular grid in the x-y plane.
! It is based on the revised Akima method.
! In this subroutine, the interpolating function is a piecewise function
! composed of a set of bicubic (bivariate third-degree) polynomials, each
! applicable to a rectangle of the input grid in the x-y plane.
! Each polynomial is determined locally.
! This subroutine has the accuracy of a bicubic polynomial, i.e., it fits the
! surface accurately when all data points lie on a surface of a bicubic
! polynomial.
! The grid lines of the input and output data can be unevenly spaced.
! The input arguments are
! MD = mode of computation
! = 1 for new XD, YD, or ZD data (default)
! = 2 for old XD, YD, and ZD data,
! NXD = number of the input-grid data points in the x
! coordinate (must be 2 or greater),
! NYD = number of the input-grid data points in the y
! coordinate (must be 2 or greater),
! XD = array of dimension NXD containing the x coordinates
! of the input-grid data points (must be in a
! monotonic increasing order),
! YD = array of dimension NYD containing the y coordinates
! of the input-grid data points (must be in a
! monotonic increasing order),
! ZD = two-dimensional array of dimension NXD*NYD
! containing the z(x,y) values at the input-grid data points,
! NXI = number of output grid points in the x coordinate
! (must be 1 or greater),
! XI = array of dimension NXI containing the x coordinates
! of the output grid points,
! NYI = number of output grid points in the y coordinate
! (must be 1 or greater),
! YI = array of dimension NYI containing the y coordinates
! of the output grid points.
! The output arguments are
! ZI = two-dimensional array of dimension NXI*NYI, where the interpolated
! z values at the output grid points are to be stored,
! IER = error flag
! = 0 for no error
! = 1 for NXD = 1 or less
! = 2 for NYD = 1 or less
! = 3 for identical XD values or XD values out of sequence
! = 4 for identical YD values or YD values out of sequence
! = 5 for NXI = 0 or less
! = 6 for NYI = 0 or less.
! N.B. The workspace has been removed from the argument list.
! WK = three-dimensional array of dimension 3*NXD*NYD used internally
! as a work area.
! The very first call to this subroutine and the call with a new XD, YD, or
! ZD array must be made with MD=1. The call with MD=2 must be preceded by
! another call with the same XD, YD, and ZD arrays. Between the call with
! MD=2 and its preceding call, the WK array must not be disturbed.
! The constant in the PARAMETER statement below is
! NIPIMX = maximum number of output points to be processed at a time.
! The constant value has been selected empirically.
! This subroutine calls the RGPD3P, RGLCTN, and RGPLNL subroutines.
! Specification statements
! .. Parameters ..
INTEGER, INTENT(IN) :: md
INTEGER, INTENT(IN) :: nxd
INTEGER, INTENT(IN) :: nyd
REAL, INTENT(IN) :: xd(nxd)
REAL, INTENT(IN) :: yd(nyd)
REAL, INTENT(IN OUT) :: zd(nxd,nyd)
INTEGER, INTENT(IN) :: nxi
REAL, INTENT(IN OUT) :: xi(nxi)
INTEGER, INTENT(IN) :: nyi
REAL, INTENT(IN) :: yi(nyi)
REAL, INTENT(IN OUT) :: zi(nxi,nyi)
INTEGER, INTENT(OUT) :: ier
REAL, INTENT(INOUT) :: wk(3,nxd,nyd)
! ..
! .. Local Scalars ..
INTEGER, PARAMETER :: nipimx=51
INTEGER :: ix, ixi, iy, iyi, nipi
! ..
! .. Local Arrays ..
REAL :: yii(nipimx)
INTEGER :: inxi(nipimx), inyi(nipimx)
! ..
! .. External Subroutines ..
! EXTERNAL rglctn,rgpd3p,rgplnl
! ..
! .. Intrinsic Functions ..
! INTRINSIC MIN
! ..
! Preliminary processing
! Error check
IF (nxd <= 1) GO TO 60
IF (nyd <= 1) GO TO 70
DO ix = 2,nxd
IF (xd(ix) <= xd(ix-1)) GO TO 80
END DO
DO iy = 2,nyd
IF (yd(iy) <= yd(iy-1)) GO TO 90
END DO
IF (nxi <= 0) GO TO 100
IF (nyi <= 0) GO TO 110
ier = 0
! Calculation
! Estimates partial derivatives at all input-grid data points
! (for MD=1).
IF (md /= 2) THEN
CALL rgpd3p(nxd, nyd, xd, yd, zd, wk)
END IF
! Outermost DO-loop with respect to the y coordinate of the output grid points
DO iyi = 1,nyi
DO ixi = 1,nipimx
yii(ixi) = yi(iyi)
END DO
! Second DO-loop with respect to the x coordinate of the output grid points
! Processes NIPIMX output-grid points, at most, at a time.
DO ixi = 1,nxi,nipimx
nipi = MIN(nxi- (ixi-1), nipimx)
! Locates the output-grid points.
CALL rglctn(nxd, nyd, xd, yd, nipi, xi(ixi), yii, inxi, inyi)
! Calculates the z values at the output-grid points.
CALL rgplnl(nxd, nyd, xd, yd, zd, wk, nipi, xi(ixi), yii, inxi, inyi, &
zi(ixi,iyi))
END DO
END DO
RETURN
! Error exit
60 WRITE (*,FMT=9000)
ier = 1
GO TO 120
70 WRITE (*,FMT=9010)
ier = 2
GO TO 120
80 WRITE (*,FMT=9020) ix,xd(ix)
ier = 3
GO TO 120
90 WRITE (*,FMT=9030) iy,yd(iy)
ier = 4
GO TO 120
100 WRITE (*,FMT=9040)
ier = 5
GO TO 120
110 WRITE (*,FMT=9050)
ier = 6
120 WRITE (*,FMT=9060) nxd,nyd,nxi,nyi
RETURN
! Format statements for error messages
9000 FORMAT (/' *** RGSF3P Error 1: NXD = 1 or less')
9010 FORMAT (/' *** RGSF3P Error 2: NYD = 1 or less')
9020 FORMAT (/' *** RGSF3P Error 3: Identical XD values or', &
' XD values out of sequence',/,' IX =',i6,', XD(IX) =', e11.3)
9030 FORMAT (/' *** RGSF3P Error 4: Identical YD values or', &
' YD values out of sequence',/,' IY =',i6,', YD(IY) =', e11.3)
9040 FORMAT (/' *** RGSF3P Error 5: NXI = 0 or less')
9050 FORMAT (/' *** RGSF3P Error 6: NYI = 0 or less')
9060 FORMAT (' NXD =', i5, ', NYD =', i5, ', NXI =', i5,', NYI =', i5 /)
END SUBROUTINE rgsf3p
! ..
! Statement Function definitions
! z2f(xx1,xx2,zz0,zz1) = (zz1-zz0)*xx2/xx1 + zz0
! z3f(xx1,xx2,xx3,zz0,zz1,zz2) = ((zz2-zz0)* (xx3-xx1)/xx2 - &
! (zz1-zz0)* (xx3-xx2)/xx1)* (xx3/ (xx2-xx1)) + zz0
FUNCTION z2f(xx1, xx2, zz0, zz1) RESULT(fn_val)
REAL(GI), INTENT(IN) :: xx1, xx2, zz0, zz1
REAL(GI) :: fn_val
fn_val = (zz1-zz0)*xx2/xx1 + zz0
RETURN
END FUNCTION z2f
FUNCTION z3f(xx1, xx2, xx3, zz0, zz1, zz2) RESULT(fn_val)
REAL(GI), INTENT(IN) :: xx1, xx2, xx3, zz0, zz1, zz2
REAL(GI) :: fn_val
fn_val = ((zz2-zz0)*(xx3-xx1)/xx2 - (zz1-zz0)*(xx3-xx2)/xx1) * &
(xx3/(xx2-xx1)) + zz0
RETURN
END FUNCTION z3f
SUBROUTINE rgpd3p(nxd, nyd, xd, yd, zd, pdd)
! Partial derivatives of a bivariate function on a rectangular grid
! (a supporting subroutine of the RGBI3P/RGSF3P subroutine package)
! Hiroshi Akima
! U.S. Department of Commerce, NTIA/ITS
! Version of 1995/08
! This subroutine estimates three partial derivatives, zx, zy, and
! zxy, of a bivariate function, z(x,y), on a rectangular grid in
! the x-y plane. It is based on the revised Akima method that has
! the accuracy of a bicubic polynomial.
! The input arguments are
! NXD = number of the input-grid data points in the x
! coordinate (must be 2 or greater),
! NYD = number of the input-grid data points in the y
! coordinate (must be 2 or greater),
! XD = array of dimension NXD containing the x coordinates of the
! input-grid data points (must be in a monotonic increasing order),
! YD = array of dimension NYD containing the y coordinates of the
! input-grid data points (must be in a monotonic increasing order),
! ZD = two-dimensional array of dimension NXD*NYD
! containing the z(x,y) values at the input-grid data points.
! The output argument is
! PDD = three-dimensional array of dimension 3*NXD*NYD,
! where the estimated zx, zy, and zxy values at the
! input-grid data points are to be stored.
! Specification statements
! .. Scalar Arguments ..
INTEGER, INTENT(IN) :: nxd
INTEGER, INTENT(IN) :: nyd
REAL, INTENT(IN) :: xd(nxd)
REAL, INTENT(IN) :: yd(nyd)
REAL, INTENT(IN) :: zd(nxd,nyd)
REAL, INTENT(OUT) :: pdd(3,nxd,nyd)
! ..
! .. Local Scalars ..
REAL(GI) :: b00, b00x, b00y, b01, b10, b11, cx1, cx2, cx3, cy1, cy2, &
cy3, disf, dnm, dz00, dz01, dz02, dz03, dz10, dz11, dz12, &
dz13, dz20, dz21, dz22, dz23, dz30, dz31, dz32, dz33, &
dzx10, dzx20, dzx30, dzxy11, dzxy12, dzxy13, dzxy21, &
dzxy22, dzxy23, dzxy31, dzxy32, dzxy33, dzy01, dzy02, &
dzy03, epsln, pezx, pezxy, pezy, smpef, smpei, smwtf, &
smwti, sx, sxx, sxxy, sxxyy, sxy, sxyy, sxyz, sxz, sy, syy, &
syz, sz, volf, wt, x0, x1, x2, x3, y0, y1, y2, &
y3, z00, z01, z02, z03, z10, z11, z12, z13, z20, z21, z22, &
z23, z30, z31, z32, z33, zxdi, zxydi, zydi
INTEGER :: ipex, ipey, ix0, ix1, ix2, ix3, iy0, iy1, iy2, iy3, nx0, ny0
! ..
! .. Local Arrays ..
REAL(GI) :: b00xa(4), b00ya(4), b01a(4), b10a(4), cxa(3,4), cya(3,4), &
sxa(4), sxxa(4), sya(4), syya(4), xa(3,4), ya(3,4), &
z0ia(3,4), zi0a(3,4)
! INTEGER :: idlt(3,4)
! ..
! .. Intrinsic Functions ..
! INTRINSIC MAX
! ..
! .. Statement Functions ..
! REAL :: z2f,z3f
! ..
INTEGER, save :: idlt(3,4) = RESHAPE([-3, -2, -1, -2, -1, 1,-1, 1,2, 1 ,2, 3 ], [ 3, 4 ])
!AL Jun 14: Fixed F90 translation
! Calculation
! Initial setting of some local variables
nx0 = MAX(4,nxd)
ny0 = MAX(4,nyd)
! Double DO-loop with respect to the input grid points
DO iy0 = 1,nyd
DO ix0 = 1,nxd
x0 = xd(ix0)
y0 = yd(iy0)
z00 = zd(ix0,iy0)
! Part 1. Estimation of ZXDI
! Initial setting
smpef = 0.0
smwtf = 0.0
smpei = 0.0
smwti = 0.0
! DO-loop with respect to the primary estimate
DO ipex = 1,4
! Selects necessary grid points in the x direction.
ix1 = ix0 + idlt(1,ipex)
ix2 = ix0 + idlt(2,ipex)
ix3 = ix0 + idlt(3,ipex)
IF ((ix1 < 1) .OR. (ix2 < 1) .OR. (ix3 < 1) .OR. &
(ix1 > nx0) .OR. (ix2 > nx0) .OR. (ix3 > nx0)) CYCLE
! Selects and/or supplements the x and z values.
x1 = xd(ix1) - x0
z10 = zd(ix1,iy0)
IF (nxd >= 4) THEN
x2 = xd(ix2) - x0
x3 = xd(ix3) - x0
z20 = zd(ix2,iy0)
z30 = zd(ix3,iy0)
ELSE IF (nxd == 3) THEN
x2 = xd(ix2) - x0
z20 = zd(ix2,iy0)
x3 = 2*xd(3) - xd(2) - x0
z30 = z3f(x1,x2,x3,z00,z10,z20)
ELSE IF (nxd == 2) THEN
x2 = 2*xd(2) - xd(1) - x0
z20 = z2f(x1,x2,z00,z10)
x3 = 2*xd(1) - xd(2) - x0
z30 = z2f(x1,x3,z00,z10)
END IF
dzx10 = (z10-z00)/x1
dzx20 = (z20-z00)/x2
dzx30 = (z30-z00)/x3
! Calculates the primary estimate of partial derivative zx as
! the coefficient of the bicubic polynomial.
cx1 = x2*x3/ ((x1-x2)* (x1-x3))
cx2 = x3*x1/ ((x2-x3)* (x2-x1))
cx3 = x1*x2/ ((x3-x1)* (x3-x2))
pezx = cx1*dzx10 + cx2*dzx20 + cx3*dzx30
! Calculates the volatility factor and distance factor in the x
! direction for the primary estimate of zx.
sx = x1 + x2 + x3
sz = z00 + z10 + z20 + z30
sxx = x1*x1 + x2*x2 + x3*x3
sxz = x1*z10 + x2*z20 + x3*z30
dnm = 4.0*sxx - sx*sx
b00 = (sxx*sz-sx*sxz)/dnm
b10 = (4.0*sxz-sx*sz)/dnm
dz00 = z00 - b00
dz10 = z10 - (b00+b10*x1)
dz20 = z20 - (b00+b10*x2)
dz30 = z30 - (b00+b10*x3)
volf = dz00**2 + dz10**2 + dz20**2 + dz30**2
disf = sxx
! Calculates the EPSLN value, which is used to decide whether or
! not the volatility factor is essentially zero.
epsln = (z00**2+z10**2+z20**2+z30**2)*1.0E-12
! Accumulates the weighted primary estimates of zx and their weights.
IF (volf > epsln) THEN
! - For a finite weight.
wt = 1.0/ (volf*disf)
smpef = smpef + wt*pezx
smwtf = smwtf + wt
ELSE
! - For an infinite weight.
smpei = smpei + pezx
smwti = smwti + 1.0
END IF
! Saves the necessary values for estimating zxy
xa(1,ipex) = x1
xa(2,ipex) = x2
xa(3,ipex) = x3
zi0a(1,ipex) = z10
zi0a(2,ipex) = z20
zi0a(3,ipex) = z30
cxa(1,ipex) = cx1
cxa(2,ipex) = cx2
cxa(3,ipex) = cx3
sxa(ipex) = sx
sxxa(ipex) = sxx
b00xa(ipex) = b00
b10a(ipex) = b10
END DO
! Calculates the final estimate of zx.
IF (smwti < 0.5) THEN
! - When no infinite weights exist.
zxdi = smpef/smwtf
ELSE
! - When infinite weights exist.
zxdi = smpei/smwti
END IF
! End of Part 1.
! Part 2. Estimation of ZYDI
! Initial setting
smpef = 0.0
smwtf = 0.0
smpei = 0.0
smwti = 0.0
! DO-loop with respect to the primary estimate
DO ipey = 1,4
! Selects necessary grid points in the y direction.
iy1 = iy0 + idlt(1,ipey)
iy2 = iy0 + idlt(2,ipey)
iy3 = iy0 + idlt(3,ipey)
IF ((iy1 < 1) .OR. (iy2 < 1) .OR. (iy3 < 1) .OR. &
(iy1 > ny0) .OR. (iy2 > ny0) .OR. (iy3 > ny0)) CYCLE
! Selects and/or supplements the y and z values.
y1 = yd(iy1) - y0
z01 = zd(ix0,iy1)
IF (nyd >= 4) THEN
y2 = yd(iy2) - y0
y3 = yd(iy3) - y0
z02 = zd(ix0,iy2)
z03 = zd(ix0,iy3)
ELSE IF (nyd == 3) THEN
y2 = yd(iy2) - y0
z02 = zd(ix0,iy2)
y3 = 2*yd(3) - yd(2) - y0
z03 = z3f(y1,y2,y3,z00,z01,z02)
ELSE IF (nyd == 2) THEN
y2 = 2*yd(2) - yd(1) - y0
z02 = z2f(y1,y2,z00,z01)
y3 = 2*yd(1) - yd(2) - y0
z03 = z2f(y1,y3,z00,z01)
END IF
dzy01 = (z01-z00)/y1
dzy02 = (z02-z00)/y2
dzy03 = (z03-z00)/y3
! Calculates the primary estimate of partial derivative zy as
! the coefficient of the bicubic polynomial.
cy1 = y2*y3/ ((y1-y2)* (y1-y3))
cy2 = y3*y1/ ((y2-y3)* (y2-y1))
cy3 = y1*y2/ ((y3-y1)* (y3-y2))
pezy = cy1*dzy01 + cy2*dzy02 + cy3*dzy03
! Calculates the volatility factor and distance factor in the y
! direction for the primary estimate of zy.
sy = y1 + y2 + y3
sz = z00 + z01 + z02 + z03
syy = y1*y1 + y2*y2 + y3*y3
syz = y1*z01 + y2*z02 + y3*z03
dnm = 4.0*syy - sy*sy
b00 = (syy*sz-sy*syz)/dnm
b01 = (4.0*syz-sy*sz)/dnm
dz00 = z00 - b00
dz01 = z01 - (b00+b01*y1)
dz02 = z02 - (b00+b01*y2)
dz03 = z03 - (b00+b01*y3)
volf = dz00**2 + dz01**2 + dz02**2 + dz03**2
disf = syy
! Calculates the EPSLN value, which is used to decide whether or
! not the volatility factor is essentially zero.
epsln = (z00**2+z01**2+z02**2+z03**2)*1.0E-12
! Accumulates the weighted primary estimates of zy and their weights.
IF (volf > epsln) THEN
! - For a finite weight.
wt = 1.0/ (volf*disf)
smpef = smpef + wt*pezy
smwtf = smwtf + wt
ELSE
! - For an infinite weight.
smpei = smpei + pezy
smwti = smwti + 1.0
END IF
! Saves the necessary values for estimating zxy
ya(1,ipey) = y1
ya(2,ipey) = y2
ya(3,ipey) = y3
z0ia(1,ipey) = z01
z0ia(2,ipey) = z02
z0ia(3,ipey) = z03
cya(1,ipey) = cy1
cya(2,ipey) = cy2
cya(3,ipey) = cy3
sya(ipey) = sy
syya(ipey) = syy
b00ya(ipey) = b00
b01a(ipey) = b01
END DO
! Calculates the final estimate of zy.
IF (smwti < 0.5) THEN
! - When no infinite weights exist.
zydi = smpef/smwtf
ELSE
! - When infinite weights exist.
zydi = smpei/smwti
END IF
! End of Part 2.
! Part 3. Estimation of ZXYDI
! Initial setting
smpef = 0.0
smwtf = 0.0
smpei = 0.0
smwti = 0.0
! Outer DO-loops with respect to the primary estimates in the x direction
DO ipex = 1,4
ix1 = ix0 + idlt(1,ipex)
ix2 = ix0 + idlt(2,ipex)
ix3 = ix0 + idlt(3,ipex)
IF ((ix1 < 1) .OR. (ix2 < 1) .OR. (ix3 < 1) .OR. &
(ix1 > nx0) .OR. (ix2 > nx0) .OR. (ix3 > nx0)) CYCLE
! Retrieves the necessary values for estimating zxy in the x direction.
x1 = xa(1,ipex)
x2 = xa(2,ipex)
x3 = xa(3,ipex)
z10 = zi0a(1,ipex)
z20 = zi0a(2,ipex)
z30 = zi0a(3,ipex)
cx1 = cxa(1,ipex)
cx2 = cxa(2,ipex)
cx3 = cxa(3,ipex)
sx = sxa(ipex)
sxx = sxxa(ipex)
b00x = b00xa(ipex)
b10 = b10a(ipex)
! Inner DO-loops with respect to the primary estimates in the y direction
DO ipey = 1,4
iy1 = iy0 + idlt(1,ipey)
iy2 = iy0 + idlt(2,ipey)
iy3 = iy0 + idlt(3,ipey)
IF ((iy1 < 1) .OR. (iy2 < 1) .OR. (iy3 < 1) .OR. (iy1 > ny0) .OR. &
(iy2 > ny0) .OR. (iy3 > ny0)) CYCLE
! Retrieves the necessary values for estimating zxy in the y direction.
y1 = ya(1,ipey)
y2 = ya(2,ipey)
y3 = ya(3,ipey)
z01 = z0ia(1,ipey)
z02 = z0ia(2,ipey)
z03 = z0ia(3,ipey)
cy1 = cya(1,ipey)
cy2 = cya(2,ipey)
cy3 = cya(3,ipey)
sy = sya(ipey)
syy = syya(ipey)
b00y = b00ya(ipey)
b01 = b01a(ipey)
! Selects and/or supplements the z values.
IF (nyd >= 4) THEN
z11 = zd(ix1,iy1)
z12 = zd(ix1,iy2)
z13 = zd(ix1,iy3)
IF (nxd >= 4) THEN
z21 = zd(ix2,iy1)
z22 = zd(ix2,iy2)
z23 = zd(ix2,iy3)
z31 = zd(ix3,iy1)
z32 = zd(ix3,iy2)
z33 = zd(ix3,iy3)
ELSE IF (nxd == 3) THEN
z21 = zd(ix2,iy1)
z22 = zd(ix2,iy2)
z23 = zd(ix2,iy3)
z31 = z3f(x1,x2,x3,z01,z11,z21)
z32 = z3f(x1,x2,x3,z02,z12,z22)
z33 = z3f(x1,x2,x3,z03,z13,z23)
ELSE IF (nxd == 2) THEN
z21 = z2f(x1,x2,z01,z11)
z22 = z2f(x1,x2,z02,z12)
z23 = z2f(x1,x2,z03,z13)
z31 = z2f(x1,x3,z01,z11)
z32 = z2f(x1,x3,z02,z12)
z33 = z2f(x1,x3,z03,z13)
END IF
ELSE IF (nyd == 3) THEN
z11 = zd(ix1,iy1)
z12 = zd(ix1,iy2)
z13 = z3f(y1,y2,y3,z10,z11,z12)
IF (nxd >= 4) THEN
z21 = zd(ix2,iy1)
z22 = zd(ix2,iy2)
z31 = zd(ix3,iy1)
z32 = zd(ix3,iy2)
ELSE IF (nxd == 3) THEN
z21 = zd(ix2,iy1)
z22 = zd(ix2,iy2)
z31 = z3f(x1,x2,x3,z01,z11,z21)
z32 = z3f(x1,x2,x3,z02,z12,z22)
ELSE IF (nxd == 2) THEN
z21 = z2f(x1,x2,z01,z11)
z22 = z2f(x1,x2,z02,z12)
z31 = z2f(x1,x3,z01,z11)
z32 = z2f(x1,x3,z02,z12)
END IF
z23 = z3f(y1,y2,y3,z20,z21,z22)
z33 = z3f(y1,y2,y3,z30,z31,z32)
ELSE IF (nyd == 2) THEN
z11 = zd(ix1,iy1)
z12 = z2f(y1,y2,z10,z11)
z13 = z2f(y1,y3,z10,z11)
IF (nxd >= 4) THEN
z21 = zd(ix2,iy1)
z31 = zd(ix3,iy1)
ELSE IF (nxd == 3) THEN
z21 = zd(ix2,iy1)
z31 = z3f(x1,x2,x3,z01,z11,z21)
ELSE IF (nxd == 2) THEN
z21 = z2f(x1,x2,z01,z11)
z31 = z2f(x1,x3,z01,z11)
END IF
z22 = z2f(y1,y2,z20,z21)
z23 = z2f(y1,y3,z20,z21)
z32 = z2f(y1,y2,z30,z31)
z33 = z2f(y1,y3,z30,z31)
END IF
! Calculates the primary estimate of partial derivative zxy as
! the coefficient of the bicubic polynomial.
dzxy11 = (z11-z10-z01+z00)/ (x1*y1)
dzxy12 = (z12-z10-z02+z00)/ (x1*y2)
dzxy13 = (z13-z10-z03+z00)/ (x1*y3)
dzxy21 = (z21-z20-z01+z00)/ (x2*y1)
dzxy22 = (z22-z20-z02+z00)/ (x2*y2)
dzxy23 = (z23-z20-z03+z00)/ (x2*y3)
dzxy31 = (z31-z30-z01+z00)/ (x3*y1)
dzxy32 = (z32-z30-z02+z00)/ (x3*y2)
dzxy33 = (z33-z30-z03+z00)/ (x3*y3)
pezxy = cx1* (cy1*dzxy11+cy2*dzxy12+cy3*dzxy13) + &
cx2* (cy1*dzxy21+cy2*dzxy22+cy3*dzxy23) + &
cx3* (cy1*dzxy31+cy2*dzxy32+cy3*dzxy33)
! Calculates the volatility factor and distance factor in the x
! and y directions for the primary estimate of zxy.
b00 = (b00x+b00y)/2.0
sxy = sx*sy
sxxy = sxx*sy
sxyy = sx*syy
sxxyy = sxx*syy
sxyz = x1* (y1*z11+y2*z12+y3*z13) + x2* (y1*z21+y2*z22+y3*z23) + &
x3* (y1*z31+y2*z32+y3*z33)
b11 = (sxyz-b00*sxy-b10*sxxy-b01*sxyy)/sxxyy
dz00 = z00 - b00
dz01 = z01 - (b00+b01*y1)
dz02 = z02 - (b00+b01*y2)
dz03 = z03 - (b00+b01*y3)
dz10 = z10 - (b00+b10*x1)
dz11 = z11 - (b00+b01*y1+x1* (b10+b11*y1))
dz12 = z12 - (b00+b01*y2+x1* (b10+b11*y2))
dz13 = z13 - (b00+b01*y3+x1* (b10+b11*y3))
dz20 = z20 - (b00+b10*x2)
dz21 = z21 - (b00+b01*y1+x2* (b10+b11*y1))
dz22 = z22 - (b00+b01*y2+x2* (b10+b11*y2))
dz23 = z23 - (b00+b01*y3+x2* (b10+b11*y3))
dz30 = z30 - (b00+b10*x3)
dz31 = z31 - (b00+b01*y1+x3* (b10+b11*y1))
dz32 = z32 - (b00+b01*y2+x3* (b10+b11*y2))
dz33 = z33 - (b00+b01*y3+x3* (b10+b11*y3))
volf = dz00**2 + dz01**2 + dz02**2 + dz03**2 + &
dz10**2 + dz11**2 + dz12**2 + dz13**2 + &
dz20**2 + dz21**2 + dz22**2 + dz23**2 + &
dz30**2 + dz31**2 + dz32**2 + dz33**2
disf = sxx*syy
! Calculates EPSLN.
epsln = (z00**2 + z01**2 + z02**2 + z03**2 + z10**2 + &
z11**2 + z12**2 + z13**2 + z20**2 + z21**2 + z22**2 + &
z23**2 + z30**2 + z31**2 + z32**2 + z33**2)* 1.0E-12
! Accumulates the weighted primary estimates of zxy and their weights.
IF (volf > epsln) THEN
! - For a finite weight.
wt = 1.0/ (volf*disf)
smpef = smpef + wt*pezxy
smwtf = smwtf + wt
ELSE
! - For an infinite weight.
smpei = smpei + pezxy
smwti = smwti + 1.0
END IF
END DO
END DO
! Calculates the final estimate of zxy.
IF (smwti < 0.5) THEN
! - When no infinite weights exist.
zxydi = smpef/smwtf
ELSE
! - When infinite weights exist.
zxydi = smpei/smwti
END IF
! End of Part 3
pdd(1,ix0,iy0) = zxdi
pdd(2,ix0,iy0) = zydi
pdd(3,ix0,iy0) = zxydi
END DO
END DO
RETURN
END SUBROUTINE rgpd3p
SUBROUTINE rglctn(nxd, nyd, xd, yd, nip, xi, yi, inxi, inyi)
! Location of the desired points in a rectangular grid
! (a supporting subroutine of the RGBI3P/RGSF3P subroutine package)
! Hiroshi Akima
! U.S. Department of Commerce, NTIA/ITS
! Version of 1995/08
! This subroutine locates the desired points in a rectangular grid
! in the x-y plane.
! The grid lines can be unevenly spaced.
! The input arguments are
! NXD = number of the input-grid data points in the x
! coordinate (must be 2 or greater),
! NYD = number of the input-grid data points in the y
! coordinate (must be 2 or greater),
! XD = array of dimension NXD containing the x coordinates of the
! input-grid data points (must be in a monotonic increasing order),
! YD = array of dimension NYD containing the y coordinates of the
! input-grid data points (must be in a monotonic increasing order),
! NIP = number of the output points to be located (must be 1 or greater),
! XI = array of dimension NIP containing the x coordinates
! of the output points to be located,
! YI = array of dimension NIP containing the y coordinates
! of the output points to be located.
! The output arguments are
! INXI = integer array of dimension NIP where the interval
! numbers of the XI array elements are to be stored,
! INYI = integer array of dimension NIP where the interval
! numbers of the YI array elements are to be stored.
! The interval numbers are between 0 and NXD and between 0 and NYD,
! respectively.
! Specification statements
! .. Scalar Arguments ..
INTEGER, INTENT(IN) :: nxd
INTEGER, INTENT(IN) :: nyd
REAL, INTENT(IN) :: xd(nxd)
REAL, INTENT(IN) :: yd(nyd)
INTEGER, INTENT(IN) :: nip
REAL, INTENT(IN) :: xi(nip)
REAL, INTENT(IN) :: yi(nip)
INTEGER, INTENT(OUT) :: inxi(nip)
INTEGER, INTENT(OUT) :: inyi(nip)
! ..
! .. Local Scalars ..
REAL :: xii, yii
INTEGER :: iip, imd, imn, imx, ixd, iyd, nintx, ninty
! ..
! DO-loop with respect to IIP, which is the point number of the output point