-
Notifications
You must be signed in to change notification settings - Fork 376
/
Copy pathssrfpack.c
4403 lines (3455 loc) · 134 KB
/
ssrfpack.c
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
/*
* ssrfpack.c: Translated via f2c then massaged so that f2c include and lib
* are not required to compile and link the sph supplement.
*
* All these functions are local (static) and included into gmt_sph.c where they are used.
*/
GMT_LOCAL double d_sign (doublereal *a, doublereal *b) {
double x;
x = (*a >= 0 ? *a : - *a);
return (*b >= 0 ? x : -x);
}
/* Table of constant values */
static doublereal c_b23 = 1.0;
GMT_LOCAL integer aplyr_(doublereal *x, doublereal *y, doublereal *z__, doublereal *cx, doublereal *sx, doublereal *cy, doublereal *sy, doublereal *xp, doublereal *yp, doublereal *zp) {
/* Builtin functions */
double sqrt(doublereal);
/* Local variables */
static doublereal t;
/* *********************************************************** */
/* From SSRFPACK */
/* Robert J. Renka */
/* Dept. of Computer Science */
/* Univ. of North Texas */
/* [email protected] */
/* 05/09/92 */
/* This subroutine applies the rotation R defined by Sub- */
/* routine CONSTR to the unit vector (X Y Z)**T, i,e. (X,Y,Z) */
/* is rotated to (XP,YP,ZP). If (XP,YP,ZP) lies in the */
/* southern hemisphere (ZP < 0), (XP,YP) are set to the */
/* coordinates of the nearest point of the equator, ZP re- */
/* maining unchanged. */
/* On input: */
/* X,Y,Z = Coordinates of a point on the unit sphere. */
/* CX,SX,CY,SY = Elements of the rotation defined by */
/* Subroutine CONSTR. */
/* Input parameters are not altered except as noted below. */
/* On output: */
/* XP,YP,ZP = Coordinates of the rotated point on the */
/* sphere unless ZP < 0, in which case */
/* (XP,YP,0) is the closest point of the */
/* equator to the rotated point. Storage */
/* for XP, YP, and ZP may coincide with */
/* storage for X, Y, and Z, respectively, */
/* if the latter need not be saved. */
/* Modules required by APLYR: None */
/* Intrinsic function called by APLYR: SQRT */
/* *********************************************************** */
/* Local parameter: */
/* T = Temporary variable */
t = *sx * *y + *cx * *z__;
*yp = *cx * *y - *sx * *z__;
*zp = *sy * *x + *cy * t;
*xp = *cy * *x - *sy * t;
if (*zp >= 0.) return 0;
/* Move (XP,YP,ZP) to the equator. */
t = sqrt(*xp * *xp + *yp * *yp);
if (t == 0.) {
goto L1;
}
*xp /= t;
*yp /= t;
return 0;
/* Move the south pole to an arbitrary point of the equator. */
L1:
*xp = 1.;
*yp = 0.;
return 0;
} /* aplyr_ */
GMT_LOCAL integer aplyrt_(doublereal *g1p, doublereal *g2p, doublereal *cx, doublereal *sx, doublereal *cy, doublereal *sy, doublereal *g) {
static doublereal t;
/* *********************************************************** */
/* From SSRFPACK */
/* Robert J. Renka */
/* Dept. of Computer Science */
/* Univ. of North Texas */
/* [email protected] */
/* 05/09/92 */
/* This subroutine applies the inverse (transpose) of the */
/* rotation defined by Subroutine CONSTR to the vector */
/* (G1P G2P 0)**T, i.e., the gradient (G1P,G2P,0) in the rot- */
/* ated coordinate system is mapped to (G1,G2,G3) in the */
/* original coordinate system. */
/* On input: */
/* G1P,G2P = X and Y components, respectively, of the */
/* gradient in the rotated coordinate system. */
/* CX,SX,CY,SY = Elements of the rotation R constructed */
/* by Subroutine CONSTR. */
/* Input parameters are not altered by this routine. */
/* On output: */
/* G = X, Y, and Z components (in that order) of the */
/* inverse rotation applied to (G1P,G2P,0) -- */
/* gradient in the original coordinate system. */
/* Modules required by APLYRT: None */
/* *********************************************************** */
/* Local parameters: */
/* T = Temporary variable */
/* Parameter adjustments */
--g;
/* Function Body */
t = *sy * *g1p;
g[1] = *cy * *g1p;
g[2] = *cx * *g2p - *sx * t;
g[3] = -(*sx) * *g2p - *cx * t;
return 0;
} /* aplyrt_ */
GMT_LOCAL doublereal arclen_(doublereal *p, doublereal *q) {
/* System generated locals */
doublereal ret_val, d__1;
/* Builtin functions */
double atan(doublereal), sqrt(doublereal);
/* Local variables */
static doublereal d__;
static integer i__;
/* *********************************************************** */
/* From SSRFPACK */
/* Robert J. Renka */
/* Dept. of Computer Science */
/* Univ. of North Texas */
/* [email protected] */
/* 05/09/92 */
/* This function computes the arc-length (angle in radians) */
/* between a pair of points on the unit sphere. */
/* On input: */
/* P,Q = Arrays of length 3 containing the X, Y, and Z */
/* coordinates (in that order) of points on the */
/* unit sphere. */
/* Input parameters are not altered by this function. */
/* On output: */
/* ARCLEN = Angle in radians between the unit vectors */
/* P and Q. 0 .LE. ARCLEN .LE. PI. */
/* Modules required by ARCLEN: None */
/* Intrinsic functions called by ARCLEN: ATAN, SQRT */
/* *********************************************************** */
/* Local parameters: */
/* D = Euclidean norm squared of P+Q */
/* I = DO-loop index */
/* Parameter adjustments */
--q;
--p;
/* Function Body */
d__ = 0.;
for (i__ = 1; i__ <= 3; ++i__) {
/* Computing 2nd power */
d__1 = p[i__] + q[i__];
d__ += d__1 * d__1;
/* L1: */
}
if (d__ == 0.) {
/* P and Q are separated by 180 degrees. */
ret_val = atan(1.) * 4.;
} else if (d__ >= 4.) {
/* P and Q coincide. */
ret_val = 0.;
} else {
ret_val = atan(sqrt((4. - d__) / d__)) * 2.;
}
return ret_val;
} /* arclen_ */
GMT_LOCAL integer snhcsh_(doublereal *x, doublereal *sinhm, doublereal *coshm, doublereal *coshmm) {
/* Initialized data */
static doublereal c1 = .1666666666659;
static doublereal c2 = .008333333431546;
static doublereal c3 = 1.984107350948e-4;
static doublereal c4 = 2.768286868175e-6;
/* Builtin functions */
double exp(doublereal);
/* Local variables */
static doublereal f, ax, xc, xs, xsd2, xsd4, expx;
/* *********************************************************** */
/* From SSRFPACK */
/* Robert J. Renka */
/* Dept. of Computer Science */
/* Univ. of North Texas */
/* [email protected] */
/* 03/18/90 */
/* This subroutine computes approximations to the modified */
/* hyperbolic functions defined below with relative error */
/* bounded by 4.7E-12 for a floating point number system with */
/* sufficient precision. For IEEE standard single precision, */
/* the relative error is less than 1.E-5 for all x. */
/* Note that the 13-digit constants in the data statements */
/* below may not be acceptable to all compilers. */
/* On input: */
/* X = Point at which the functions are to be */
/* evaluated. */
/* X is not altered by this routine. */
/* On output: */
/* SINHM = sinh(X) - X. */
/* COSHM = cosh(X) - 1. */
/* COSHMM = cosh(X) - 1 - X*X/2. */
/* Modules required by SNHCSH: None */
/* Intrinsic functions called by SNHCSH: ABS, EXP */
/* *********************************************************** */
ax = fabs(*x);
xs = ax * ax;
if (ax <= .5) {
/* Approximations for small X: */
xc = *x * xs;
*sinhm = xc * (((c4 * xs + c3) * xs + c2) * xs + c1);
xsd4 = xs * .25;
xsd2 = xsd4 + xsd4;
f = (((c4 * xsd4 + c3) * xsd4 + c2) * xsd4 + c1) * xsd4;
*coshmm = xsd2 * f * (f + 2.);
*coshm = *coshmm + xsd2;
} else {
/* Approximations for large X: */
expx = exp(ax);
*sinhm = -(1. / expx + ax + ax - expx) / 2.;
if (*x < 0.) {
*sinhm = -(*sinhm);
}
*coshm = (1. / expx - 2. + expx) / 2.;
*coshmm = *coshm - xs / 2.;
}
return 0;
} /* snhcsh_ */
GMT_LOCAL integer arcint_(doublereal *p, doublereal *p1, doublereal *p2, doublereal *f1, doublereal *f2, doublereal *g1, doublereal *g2, doublereal *sigma, doublereal *f, doublereal *g, doublereal *gn) {
/* Builtin functions */
double sqrt(doublereal), exp(doublereal);
/* Local variables */
static doublereal a, e;
static integer i__;
static doublereal s, b1, b2, d1, d2, e1, e2, al, cm, gt, sm, tm, un[3],
ts, cm2, sb1, sb2, sm2, tm1, tm2, tp1, tp2, cmm, sig, ems, tau1,
tau2, sinh__, sinh2, dummy, unorm;
/* *********************************************************** */
/* From SSRFPACK */
/* Robert J. Renka */
/* Dept. of Computer Science */
/* Univ. of North Texas */
/* [email protected] */
/* 11/21/96 */
/* Given 3 points P, P1, and P2 lying on a common geodesic */
/* of the unit sphere with P between P1 and P2, along with */
/* data values and gradients at P1 and P2, this subroutine */
/* computes an interpolated value F and a gradient vector G */
/* AT P. F and the tangential component of G are taken to be */
/* the value and derivative (with respect to arc-length) of */
/* a Hermite interpolatory tension spline defined by the end- */
/* point values and tangential gradient components. The nor- */
/* mal component of G is obtained by linear interpolation of */
/* the normal components of the gradients at P1 and P2. */
/* On input: */
/* P = Cartesian coordinates of a point lying on the */
/* arc defined by P1 and P2. P(1)**2 + P(2)**2 + */
/* P(3)**2 = 1. */
/* P1,P2 = Coordinates of distinct points on the unit */
/* sphere defining an arc with length less than */
/* 180 degrees. */
/* F1,F2 = Data values associated with P1 and P2, */
/* respectively. */
/* G1,G2 = Gradient vectors associated with P1 and P2. */
/* G1 and G2 are orthogonal to P1 and P2, */
/* respectively. */
/* SIGMA = Tension factor associated with P1-P2. */
/* The above parameters are not altered by this routine. */
/* G = Array of length 3. */
/* On output: */
/* F = Interpolated value at P. */
/* G = Interpolated gradient at P. */
/* GN = Normal component of G with the direction */
/* P1 X P2 taken to be positive. The extrapola- */
/* tion procedure requires this component. */
/* For each vector V, V(1), V(2), and V(3) contain X, Y, */
/* and Z components, respectively. */
/* SSRFPACK modules required by ARCINT: ARCLEN, SNHCSH */
/* Intrinsic functions called by ARCINT: ABS, EXP, SQRT */
/* *********************************************************** */
/* Parameter adjustments */
--g;
--g2;
--g1;
--p2;
--p1;
--p;
/* Function Body */
/* Local parameters: */
/* A = Angle in radians (arc-length) between P1 and */
/* P2 */
/* AL = Arc-length between P1 and P */
/* B1,B2 = Local coordinates of P with respect to P1-P2 */
/* CM,CMM = Coshm(SIG) and Coshmm(SIG) -- refer to SNHCSH */
/* CM2 = Coshm(SB2) */
/* DUMMY = Dummy parameter for SNHCSH */
/* D1,D2 = Scaled second differences */
/* E = CM**2 - SM*Sinh = SIG*SM - 2*CMM (scaled by */
/* 2*EMS if SIG > .5) */
/* EMS = Exp(-SIG) */
/* E1,E2 = Exp(-SB1), Exp(-SB2) */
/* GT = Tangential component of G -- component in the */
/* direction UN X P */
/* I = DO-loop index */
/* LUN = Logical unit for error messages */
/* S = Slope: (F2-F1)/A */
/* SB1,SB2 = SIG*B1, SIG*B2 */
/* SIG = Abs(SIGMA) */
/* SINH = Sinh(SIGMA) */
/* SINH2 = Sinh(SB2) */
/* SM,SM2 = Sinhm(SIG), Sinhm(SB2) */
/* TAU1,TAU2 = Tangential derivatives (components of G1,G2) */
/* at P1 and P2 */
/* TM = 1-EMS */
/* TM1,TM2 = 1-E1, 1-E2 */
/* TP1,TP2 = 1+E1, 1+E2 */
/* TS = TM**2 */
/* UN = Unit normal to the plane of P, P1, and P2 */
/* UNORM = Euclidean norm of P1 X P2 -- used to normalize */
/* UN */
/* Compute unit normal UN. */
un[0] = p1[2] * p2[3] - p1[3] * p2[2];
un[1] = p1[3] * p2[1] - p1[1] * p2[3];
un[2] = p1[1] * p2[2] - p1[2] * p2[1];
unorm = sqrt(un[0] * un[0] + un[1] * un[1] + un[2] * un[2]);
if (unorm == 0.) {
goto L2;
}
/* Normalize UN. */
for (i__ = 1; i__ <= 3; ++i__) {
un[i__ - 1] /= unorm;
/* L1: */
}
/* Compute tangential derivatives at the endpoints: */
/* TAU1 = (G1,UN X P1) = (G1,P2)/UNORM and */
/* TAU2 = (G2,UN X P2) = -(G2,P1)/UNORM. */
tau1 = (g1[1] * p2[1] + g1[2] * p2[2] + g1[3] * p2[3]) / unorm;
tau2 = -(g2[1] * p1[1] + g2[2] * p1[2] + g2[3] * p1[3]) / unorm;
/* Compute arc-lengths A, AL. */
a = arclen_(&p1[1], &p2[1]);
if (a == 0.) {
goto L2;
}
al = arclen_(&p1[1], &p[1]);
/* Compute local coordinates, slope, and second differences. */
b2 = al / a;
b1 = 1. - b2;
s = (*f2 - *f1) / a;
d1 = s - tau1;
d2 = tau2 - s;
/* Test the range of SIGMA. */
sig = fabs(*sigma);
if (sig < 1e-9) {
/* Hermite cubic interpolation. */
*f = *f1 + al * (tau1 + b2 * (d1 + b1 * (d1 - d2)));
gt = tau1 + b2 * (d1 + d2 + b1 * 3. * (d1 - d2));
} else if (sig <= .5) {
/* 0 < SIG .LE. .5. Use approximations designed to avoid */
/* cancellation error in the hyperbolic functions. */
sb2 = sig * b2;
snhcsh_(&sig, &sm, &cm, &cmm);
snhcsh_(&sb2, &sm2, &cm2, &dummy);
sinh__ = sm + sig;
sinh2 = sm2 + sb2;
e = sig * sm - cmm - cmm;
*f = *f1 + al * tau1 + a * ((cm * sm2 - sm * cm2) * (d1 + d2) + sig *
(cm * cm2 - sinh__ * sm2) * d1) / (sig * e);
gt = tau1 + ((cm * cm2 - sm * sinh2) * (d1 + d2) + sig * (cm * sinh2
- sinh__ * cm2) * d1) / e;
} else {
/* SIG > .5. Use negative exponentials in order to avoid */
/* overflow. Note that EMS = EXP(-SIG). */
sb1 = sig * b1;
sb2 = sig - sb1;
e1 = exp(-sb1);
e2 = exp(-sb2);
ems = e1 * e2;
tm = 1. - ems;
ts = tm * tm;
tm1 = 1. - e1;
tm2 = 1. - e2;
e = tm * (sig * (ems + 1.) - tm - tm);
*f = *f1 + al * s + a * (tm * tm1 * tm2 * (d1 + d2) + sig * ((e2 *
tm1 * tm1 - b1 * ts) * d1 + (e1 * tm2 * tm2 - b2 * ts) * d2))
/ (sig * e);
tp1 = e1 + 1.;
tp2 = e2 + 1.;
gt = s + (tm1 * (tm * tp2 - sig * e2 * tp1) * d1 - tm2 * (tm * tp1 -
sig * e1 * tp2) * d2) / e;
}
/* Compute GN. */
*gn = b1 * (un[0] * g1[1] + un[1] * g1[2] + un[2] * g1[3]) + b2 * (un[0] *
g2[1] + un[1] * g2[2] + un[2] * g2[3]);
/* Compute G = GT*(UN X P) + GN*UN. */
g[1] = gt * (un[1] * p[3] - un[2] * p[2]) + *gn * un[0];
g[2] = gt * (un[2] * p[1] - un[0] * p[3]) + *gn * un[1];
g[3] = gt * (un[0] * p[2] - un[1] * p[1]) + *gn * un[2];
return 0;
/* P1 X P2 = 0. Print an error message and terminate */
/* processing. */
/* 2 WRITE (LUN,100) (P1(I),I=1,3), (P2(I),I=1,3) */
/* 100 FORMAT ('1','ERROR IN ARCINT -- P1 = ',2(F9.6,', '), */
/* . F9.6/1X,19X,'P2 = ',2(F9.6,', '),F9.6) */
L2:
#ifdef SPH_DEBUG
fprintf (stderr, "ERROR IN ARCINT -- P1 = %9.6f %9.6f %9.6f P2 = %9.6f %9.6f %9.6f\n", p1[1], p1[2], p1[3], p2[1], p2[2], p2[3]);
#endif
i__ = 1;
return 0;
} /* arcint_ */
GMT_LOCAL integer constr_(doublereal *xk, doublereal *yk, doublereal *zk, doublereal *cx, doublereal *sx, doublereal *cy, doublereal *sy) {
/* Builtin functions */
double sqrt(doublereal);
/* *********************************************************** */
/* From SSRFPACK */
/* Robert J. Renka */
/* Dept. of Computer Science */
/* Univ. of North Texas */
/* [email protected] */
/* 05/09/92 */
/* This subroutine constructs the elements of a 3 by 3 */
/* orthogonal matrix R which rotates a point (XK,YK,ZK) on */
/* the unit sphere to the north pole, i.e., */
/* (XK) (CY 0 -SY) (1 0 0) (XK) (0) */
/* R * (YK) = ( 0 1 0) * (0 CX -SX) * (YK) = (0) */
/* (ZK) (SY 0 CY) (0 SX CX) (ZK) (1) */
/* On input: */
/* XK,YK,ZK = Components of a unit vector to be */
/* rotated to (0,0,1). */
/* Input parameters are not altered by this routine. */
/* On output: */
/* CX,SX,CY,SY = Elements of R: CX,SX define a rota- */
/* tion about the X-axis and CY,SY define */
/* a rotation about the Y-axis. */
/* Modules required by CONSTR: None */
/* Intrinsic function called by CONSTR: SQRT */
/* *********************************************************** */
*cy = sqrt(*yk * *yk + *zk * *zk);
*sy = *xk;
if (*cy != 0.) {
*cx = *zk / *cy;
*sx = *yk / *cy;
} else {
/* (XK,YK,ZK) lies on the X-axis. */
*cx = 1.;
*sx = 0.;
}
return 0;
} /* constr_ */
GMT_LOCAL doublereal hval_(doublereal *b, doublereal *h1, doublereal *h2, doublereal *hp1, doublereal *hp2, doublereal *sigma) {
/* System generated locals */
doublereal ret_val;
/* Builtin functions */
double exp(doublereal);
/* Local variables */
static doublereal e, s, b1, b2, d1, d2, e1, e2, cm, sm, tm, ts, cm2, sb1,
sb2, sm2, tm1, tm2, cmm, sig, ems, dummy;
/* *********************************************************** */
/* From SSRFPACK */
/* Robert J. Renka */
/* Dept. of Computer Science */
/* Univ. of North Texas */
/* [email protected] */
/* 11/21/96 */
/* Given a line segment P1-P2 containing a point P, along */
/* with values and derivatives at the endpoints, this func- */
/* tion returns the value H(P), where H is the Hermite inter- */
/* polatory tension spline defined by the endpoint data. */
/* On input: */
/* B = Local coordinate of P with respect to P1-P2: */
/* P = B*P1 + (1-B)*P2, and thus B = d(P,P2)/ */
/* d(P1,P2), where d(P1,P2) is the distance between */
/* P1 and P2. B < 0 or B > 1 results in extrapola- */
/* tion. */
/* H1,H2 = Values interpolated at P1 and P2, respec- */
/* tively. */
/* HP1,HP2 = Products of d(P1,P2) with first order der- */
/* ivatives at P1 and P2, respectively. HP1 */
/* may, for example, be the scalar product of */
/* P2-P1 with a gradient at P1. */
/* SIGMA = Nonnegative tension factor associated with */
/* the spline. SIGMA = 0 corresponds to a */
/* cubic spline, and H approaches the linear */
/* interpolant of H1 and H2 as SIGMA increases. */
/* Input parameters are not altered by this function. */
/* On output: */
/* HVAL = Interpolated value H(P). */
/* SSRFPACK module required by HVAL: SNHCSH */
/* Intrinsic functions called by HVAL: ABS, EXP */
/* *********************************************************** */
b1 = *b;
b2 = 1. - b1;
/* Compute slope S and second differences D1 and D2 scaled */
/* by the separation between P1 and P2. */
s = *h2 - *h1;
d1 = s - *hp1;
d2 = *hp2 - s;
/* Test the range of SIGMA. */
sig = fabs(*sigma);
if (sig < 1e-9) {
/* Hermite cubic interpolation: */
ret_val = *h1 + b2 * (*hp1 + b2 * (d1 + b1 * (d1 - d2)));
} else if (sig <= .5) {
/* 0 < SIG .LE. .5. Use approximations designed to avoid */
/* cancellation error in the hyperbolic functions. */
sb2 = sig * b2;
snhcsh_(&sig, &sm, &cm, &cmm);
snhcsh_(&sb2, &sm2, &cm2, &dummy);
e = sig * sm - cmm - cmm;
ret_val = *h1 + b2 * *hp1 + ((cm * sm2 - sm * cm2) * (d1 + d2) + sig *
(cm * cm2 - (sm + sig) * sm2) * d1) / (sig * e);
} else {
/* SIG > .5. Use negative exponentials in order to avoid */
/* overflow. Note that EMS = EXP(-SIG). */
sb1 = sig * b1;
sb2 = sig - sb1;
e1 = exp(-sb1);
e2 = exp(-sb2);
ems = e1 * e2;
tm = 1. - ems;
ts = tm * tm;
tm1 = 1. - e1;
tm2 = 1. - e2;
e = tm * (sig * (ems + 1.) - tm - tm);
ret_val = *h1 + b2 * s + (tm * tm1 * tm2 * (d1 + d2) + sig * ((e2 *
tm1 * tm1 - b1 * ts) * d1 + (e1 * tm2 * tm2 - b2 * ts) * d2))
/ (sig * e);
}
return ret_val;
} /* hval_ */
GMT_LOCAL doublereal fval_(doublereal *b1, doublereal *b2, doublereal *b3, doublereal *v1, doublereal *v2, doublereal *v3, doublereal *f1, doublereal *f2, doublereal *f3, doublereal *g1, doublereal *g2, doublereal *g3, doublereal *sig1, doublereal *sig2, doublereal *sig3) {
/* System generated locals */
doublereal ret_val;
/* Builtin functions */
double sqrt(doublereal);
/* Local variables */
static doublereal f, g[3];
static integer i__;
static doublereal c1, c2, c3, q1[3], q2[3], q3[3], s1, s2, s3, u1[3], u2[
3], u3[3], ds, dv, u1n, u2n, u3n, sig, val, dum, sum;
/* *********************************************************** */
/* From SSRFPACK */
/* Robert J. Renka */
/* Dept. of Computer Science */
/* Univ. of North Texas */
/* [email protected] */
/* 05/09/92 */
/* Given data values and gradients at the three vertices of */
/* a spherical triangle containing a point P, this routine */
/* computes the value of F at P where F interpolates the ver- */
/* tex data. Along the triangle sides, the interpolatory */
/* function F is the Hermite interpolatory tension spline */
/* defined by the values and tangential gradient components */
/* at the endpoints, and the gradient component normal to the */
/* triangle side varies linearly with respect to arc-length */
/* between the normal gradient components at the endpoints. */
/* A first-order C-1 blending method is used on the underly- */
/* ing planar triangle. Since values and gradients on an arc */
/* depend only on the vertex data, the method results in C-1 */
/* continuity when used to interpolate over a triangulation. */
/* The blending method consists of taking F(P) to be a */
/* weighted sum of the values at PP of three univariate Her- */
/* mite interpolatory tension splines defined on the line */
/* segments which join the vertices to the opposite sides and */
/* pass through PP: the central projection of P onto the */
/* underlying planar triangle. The tension factors for these */
/* splines are obtained by linear interpolation between the */
/* pair of tension factors associated with the triangle sides */
/* which join at the appropriate vertex. */
/* A tension factor SIGMA associated with a Hermite interp- */
/* olatory tension spline is a nonnegative parameter which */
/* determines the curviness of the spline. SIGMA = 0 results */
/* in a cubic spline, and the spline approaches the linear */
/* interpolant as SIGMA increases. */
/* On input: */
/* B1,B2,B3 = Barycentric coordinates of PP with re- */
/* spect to the (planar) underlying triangle */
/* (V1,V2,V3), where PP is the central */
/* projection of P onto this triangle. */
/* V1,V2,V3 = Cartesian coordinates of the vertices of */
/* a spherical triangle containing P. V3 */
/* Left V1->V2. */
/* F1,F2,F3 = Data values associated with the vertices. */
/* G1,G2,G3 = Gradients associated with the vertices. */
/* Gi is orthogonal to Vi for i = 1,2,3. */
/* SIG1,SIG2,SIG3 = Tension factors associated with the */
/* triangle sides opposite V1, V2, and */
/* V3, respectively. */
/* Input parameters are not altered by this function. */
/* On output: */
/* FVAL = Interpolated value at P. */
/* Each vector V above contains X, Y, and Z components in */
/* V(1), V(2), and V(3), respectively. */
/* SSRFPACK modules required by FVAL: ARCINT, ARCLEN, HVAL */
/* Intrinsic function called by FVAL: SQRT */
/* *********************************************************** */
/* Local parameters: */
/* C1,C2,C3 = Coefficients (weight functions) of partial */
/* interpolants. C1 = 1 on the edge opposite */
/* V1 and C1 = 0 on the other edges. Simi- */
/* larly for C2 and C3. C1+C2+C3 = 1. */
/* DS = Directional derivative (scaled by distance) */
/* at U1, U2, or U3: DS = (G,U1-V1)/U1N = */
/* -(G,V1)/U1N on side opposite V1, where G/ */
/* U1N (plus an orthogonal component) is the */
/* projection of G onto the planar triangle */
/* DUM = Dummy variable for calls to ARCINT */
/* DV = Directional derivatives (scaled by distance) */
/* at a vertex: D1 = (G1,U1-V1) = (G1,U1) */
/* F,G = Value and gradient at Q1 Q2, or Q3 obtained */
/* by interpolation along one of the arcs of */
/* the spherical triangle */
/* I = DO-loop index */
/* Q1,Q2,Q3 = Central projections of U1, U2, and U3 onto */
/* the sphere and thus lying on an arc of the */
/* spherical triangle */
/* SIG = Tension factor for a side-vertex (partial) */
/* interpolant: obtained by linear interpo- */
/* lation applied to triangle side tensions */
/* SUM = Quantity used to normalize C1, C2, and C3 */
/* S1,S2,S3 = Sums of pairs of barycentric coordinates: */
/* used to compute U1, U2, U3, and SIG */
/* U1,U2,U3 = Points on the boundary of the planar trian- */
/* gle and lying on the lines containing PP */
/* and the vertices. U1 is opposite V1, etc. */
/* U1N,U2N,U3N = Quantities used to compute Q1, Q2, and Q3 */
/* (magnitudes of U1, U2, and U3) */
/* VAL = Local variable used to accumulate the con- */
/* tributions to FVAL */
/* Compute weight functions C1, C2, and C3. */
/* Parameter adjustments */
--g3;
--g2;
--g1;
--v3;
--v2;
--v1;
/* Function Body */
c1 = *b2 * *b3;
c2 = *b3 * *b1;
c3 = *b1 * *b2;
sum = c1 + c2 + c3;
if (sum <= 0.) {
/* P coincides with a vertex. */
ret_val = *b1 * *f1 + *b2 * *f2 + *b3 * *f3;
return ret_val;
}
/* Normalize C1, C2, and C3. */
c1 /= sum;
c2 /= sum;
c3 /= sum;
/* Compute (S1,S2,S3), (U1,U2,U3) and (U1N,U2N,U3N). */
s1 = *b2 + *b3;
s2 = *b3 + *b1;
s3 = *b1 + *b2;
u1n = 0.;
u2n = 0.;
u3n = 0.;
for (i__ = 1; i__ <= 3; ++i__) {
u1[i__ - 1] = (*b2 * v2[i__] + *b3 * v3[i__]) / s1;
u2[i__ - 1] = (*b3 * v3[i__] + *b1 * v1[i__]) / s2;
u3[i__ - 1] = (*b1 * v1[i__] + *b2 * v2[i__]) / s3;
u1n += u1[i__ - 1] * u1[i__ - 1];
u2n += u2[i__ - 1] * u2[i__ - 1];
u3n += u3[i__ - 1] * u3[i__ - 1];
/* L1: */
}
/* Compute Q1, Q2, and Q3. */
u1n = sqrt(u1n);
u2n = sqrt(u2n);
u3n = sqrt(u3n);
for (i__ = 1; i__ <= 3; ++i__) {
q1[i__ - 1] = u1[i__ - 1] / u1n;
q2[i__ - 1] = u2[i__ - 1] / u2n;
q3[i__ - 1] = u3[i__ - 1] / u3n;
/* L2: */
}
/* Compute interpolated value (VAL) at P by looping on */
/* triangle sides. */
val = 0.;
/* Contribution from side opposite V1: */
/* Compute value and gradient at Q1 by interpolating */
/* between V2 and V3. */
arcint_(q1, &v2[1], &v3[1], f2, f3, &g2[1], &g3[1], sig1, &f, g, &dum);
/* Add in the contribution. */
dv = g1[1] * u1[0] + g1[2] * u1[1] + g1[3] * u1[2];
ds = -(g[0] * v1[1] + g[1] * v1[2] + g[2] * v1[3]) / u1n;
sig = (*b2 * *sig3 + *b3 * *sig2) / s1;
val += c1 * hval_(b1, f1, &f, &dv, &ds, &sig);
/* Contribution from side opposite V2: */
/* Compute value and gradient at Q2 by interpolating */
/* between V3 and V1. */
arcint_(q2, &v3[1], &v1[1], f3, f1, &g3[1], &g1[1], sig2, &f, g, &dum);
/* Add in the contribution. */
dv = g2[1] * u2[0] + g2[2] * u2[1] + g2[3] * u2[2];
ds = -(g[0] * v2[1] + g[1] * v2[2] + g[2] * v2[3]) / u2n;
sig = (*b3 * *sig1 + *b1 * *sig3) / s2;
val += c2 * hval_(b2, f2, &f, &dv, &ds, &sig);
/* Contribution from side opposite V3: */
/* Compute interpolated value and gradient at Q3 */
/* by interpolating between V1 and V2. */
arcint_(q3, &v1[1], &v2[1], f1, f2, &g1[1], &g2[1], sig3, &f, g, &dum);
/* Add in the final contribution. */
dv = g3[1] * u3[0] + g3[2] * u3[1] + g3[3] * u3[2];
ds = -(g[0] * v3[1] + g[1] * v3[2] + g[2] * v3[3]) / u3n;
sig = (*b1 * *sig2 + *b2 * *sig1) / s3;
ret_val = val + c3 * hval_(b3, f3, &f, &dv, &ds, &sig);
return ret_val;
} /* fval_ */
GMT_LOCAL integer getsig_(integer *n, doublereal *x, doublereal *y, doublereal *z__, doublereal *h__, integer *list, integer *lptr, integer *lend, doublereal *grad, doublereal *tol, doublereal *sigma, doublereal *dsmax, integer *ier) {
static doublereal sbig = 85.;
/* System generated locals */
integer i__1, i__2;
doublereal d__1, d__2;
/* Builtin functions */
double sqrt(doublereal), exp(doublereal), d_sign(doublereal *, doublereal
*);
/* Local variables */
static doublereal a, e, f, s, t, c1, c2, d0, d1, d2, f0;
static integer n1, n2;
static doublereal p1[3], p2[3], s1, s2, t0, t1, t2, al, fp, tm, un[3];
static integer nm1, lp1, lp2;
static doublereal tp1, d1d2, scm, dsm, ems, sig;
static integer lpl;
static doublereal sgn;
static integer nit;
static doublereal ssm, ems2, d1pd2, fneg, dsig, dmax__, fmax;
static integer icnt;
static doublereal ftol, rtol, stol, coshm, sigin, sinhm, ssinh;
static doublereal unorm;
static doublereal coshmm;
/* *********************************************************** */
/* From SSRFPACK */
/* Robert J. Renka */
/* Dept. of Computer Science */
/* Univ. of North Texas */
/* [email protected] */
/* 11/21/96 */
/* Given a triangulation of a set of nodes on the unit */
/* sphere, along with data values H and gradients GRAD at the */
/* nodes, this subroutine determines, for each triangulation */
/* arc, the smallest (nonnegative) tension factor SIGMA such */
/* that the Hermite interpolatory tension spline H(A), de- */
/* fined by SIGMA and the endpoint values and directional */
/* derivatives, preserves local shape properties of the data. */
/* In order to define the shape properties on an arc, it is */
/* convenient to map the arc to an interval (A1,A2). Then, */
/* denoting the endpoint data values by H1,H2 and the deriva- */
/* tives (tangential gradient components) by HP1,HP2, and */
/* letting S = (H2-H1)/(A2-A1), the data properties are */
/* Monotonicity: S, HP1, and HP2 are nonnegative or */
/* nonpositive, */
/* and */
/* Convexity: HP1 .LE. S .LE. HP2 or HP1 .GE. S */
/* .GE. HP2. */