-
Notifications
You must be signed in to change notification settings - Fork 375
/
Copy pathgmt_vector.c
1929 lines (1712 loc) · 72.6 KB
/
gmt_vector.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
/*--------------------------------------------------------------------
*
* Copyright (c) 1991-2025 by the GMT Team (https://www.generic-mapping-tools.org/team.html)
* See LICENSE.TXT file for copying and redistribution conditions.
*
* This program is free software; you can redistribute it and/or modify
* it under the terms of the GNU Lesser General Public License as published by
* the Free Software Foundation; version 3 or any later version.
*
* This program is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU Lesser General Public License for more details.
*
* Contact info: www.generic-mapping-tools.org
*--------------------------------------------------------------------*/
/*
* Author: Walter H.F. Smith
* Date: 1-JAN-2010
* Version: 5.x
*/
/*
* A) List of exported gmt_* functions available to modules and libraries via gmt_dev.h:
*
* gmt_add3v
* gmt_cart_to_geo
* gmt_cart_to_polar
* gmt_chol_dcmp
* gmt_chol_recover
* gmt_chol_solv
* gmt_cross3v
* gmt_dot2v
* gmt_dot3v
* gmt_fix_up_path
* gmt_gauss
* gmt_gaussjordan
* gmt_geo_to_cart
* gmt_jacobi
* gmt_mag3v
* gmt_make_rot_matrix
* gmt_make_rot_matrix2
* gmt_matrix_matrix_add
* gmt_matrix_matrix_mult
* gmt_matrix_vect_mult
* gmt_matrix_vector_mult
* gmt_n_cart_to_geo
* gmt_normalize2v
* gmt_normalize3v
* gmt_polar_to_cart
* gmt_resample_path
* gmt_set_line_resampling
* gmt_solve_svd
* gmt_sort_svd_values
* gmt_sub3v
* gmt_svdcmp
*/
#include "gmt_dev.h"
#include "gmt_internals.h"
#define MAX_SWEEPS 50
/* Local functions */
GMT_LOCAL void gmtvector_switchrows (double *a, double *b, unsigned int n1, unsigned int n2, unsigned int n) {
double *oa = (double *)malloc (sizeof(double)*n);
memcpy(oa, a+n*n1, sizeof(double)*n);
memcpy(a+n*n1, a+n*n2, sizeof(double)*n);
memcpy(a+n*n2, oa, sizeof(double)*n);
gmt_M_double_swap (b[n1], b[n2]);
gmt_M_str_free (oa);
}
/* Given a matrix a[0..m-1][0...n-1], this routine computes its singular
value decomposition, A=UWVt. The matrix U replaces a on output.
The diagonal matrix of singular values W is output as a vector
w[0...n-1]. The matrix V (Not V transpose) is output as
v[0...n-1][0....n-1]. m must be greater than or equal to n; if it is
smaller, then a should be filled up to square with zero rows.
Modified from Numerical Recipes -> page 68.
*/
#define SIGN(a,b) ((b) >= 0.0 ? fabs(a) : -fabs(a))
#ifndef HAVE_LAPACK
/* Use version provided by DJ 2015-07-06 [SDSC] */
GMT_LOCAL int gmtvector_svdcmp_nr (struct GMT_CTRL *GMT, double *a, unsigned int m_in, unsigned int n_in, double *w, double *v) {
/* void svdcmp(double *a,int m,int n,double *w,double *v) */
int flag,i,its,j,jj,k,l=0,nm = 0, n = n_in, m = m_in;
double c,f,h,s,x,y,z;
double anorm=0.0,tnorm, g=0.0,scale=0.0;
double *rv1 = NULL;
if (m < n) {
GMT_Report (GMT->parent, GMT_MSG_ERROR, "Failure in gmt_svdcmp: m < n augment A with additional rows\n");
return (GMT_DIM_TOO_SMALL);
}
/* allocate work space */
if ((rv1 = gmt_M_memory (GMT, NULL, n, double)) == NULL) return (GMT_MEMORY_ERROR);
/* do householder reduction to bidiagonal form */
for (i=0;i<n;i++) {
l=i+1;
rv1[i]=scale*g;
g=s=scale=0.0;
if (i < m) {
for (k=i;k<m;k++) scale += fabs (a[k*n+i]); /* a[k][i] */
if (scale) {
for (k=i;k<m;k++) {
a[k*n+i] /= scale; /* a[k][i] */
s += a[k*n+i]*a[k*n+i]; /* a[k][i] */
}
f=a[i*n+i]; /* a[i][i] */
g= -1.0*SIGN(sqrt(s),f);
h=f*g-s;
a[i*n+i]=f-g; /* a[i][i] */
if (i != n-1) {
for (j=l;j<n;j++) {
for (s=0.0,k=i;k<m;k++) s += a[k*n+i]*a[k*n+j]; /* a[k][i] a[k][j] */
f=s/h;
for (k=i;k<m;k++) a[k*n+j] += f*a[k*n+i]; /* a[k][j] a[k][i] */
}
}
for (k=i;k<m;k++) a[k*n+i] *= scale; /* a[k][i] */
}
}
w[i]=scale*g;
g=s=scale=0.0;
if (i <= m-1 && i != n-1) {
for (k=l;k<n;k++) scale += fabs (a[i*n+k]); /* a[i][k] */
if (scale) {
for (k=l;k<n;k++) {
a[i*n+k] /= scale; /* a[i][k] */
s += a[i*n+k]*a[i*n+k]; /* a[i][k] */
}
f=a[i*n+l]; /* a[i][l] */
g = -1.0*SIGN(sqrt(s),f);
h=f*g-s;
a[i*n+l]=f-g; /* a[i][l] */
for (k=l;k<n;k++) rv1[k]=a[i*n+k]/h; /* a[i][k] */
if (i != m-1) {
for (j=l;j<m;j++) {
for (s=0.0,k=l;k<n;k++) s += a[j*n+k]*a[i*n+k]; /*a[j][k] a[i][k] */
for (k=l;k<n;k++) a[j*n+k] += s*rv1[k]; /* a[j][k] */
}
}
for (k=l;k<n;k++) a[i*n+k] *= scale; /* a[i][k] */
}
}
tnorm=fabs (w[i])+fabs (rv1[i]);
anorm=MAX(anorm,tnorm);
}
/* accumulation of right-hand transforms */
for (i=n-1;i>=0;i--) {
if (i < n-1) {
if (g) {
for (j=l;j<n;j++) v[j*n+i]=(a[i*n+j]/a[i*n+l])/g; /* v[j][i] a[i][j] a[i][l] */
for (j=l;j<n;j++) {
for (s=0.0,k=l;k<n;k++) s += a[i*n+k]*v[k*n+j]; /* a[i][k] v[k][j] */
for (k=l;k<n;k++) v[k*n+j] += s*v[k*n+i]; /* v[k][j] v[k][i] */
}
}
for (j=l;j<n;j++) v[i*n+j]=v[j*n+i]=0.0; /* v[i][j] v[j][i] */
}
v[i*n+i]=1.0; /* v[i][i] */
g=rv1[i];
l=i;
}
/* accumulation of left-hand transforms */
for (i=n-1;i>=0;i--) {
l=i+1;
g=w[i];
if (i < n-1) for (j=l;j<n;j++) a[i*n+j]=0.0; /* a[i][j] */
if (g) {
g=1.0/g;
if (i != n-1) {
#ifdef _OPENMP
#pragma omp parallel for private(j,k,f,s) shared(a,n,m,g,l)
#endif
for (j=l;j<n;j++) {
for (s=0.0,k=l;k<m;k++) s += a[k*n+i]*a[k*n+j]; /* a[k][i] a[k][j] */
f=(s/a[i*n+i])*g; /* a[i][i] */
for (k=i;k<m;k++) a[k*n+j] += f*a[k*n+i]; /* a[k][j] a[k][i] */
}
}
for (j=i;j<m;j++) a[j*n+i] *= g; /* a[j][i] */
}
else {
for (j=i;j<m;j++) a[j*n+i]=0.0; /* a[j][i] */
}
++a[i*n+i]; /* a[i][i] */
}
/* diagonalization of the bidiagonal form */
for (k=n-1;k>=0;k--) { /* loop over singular values */
for (its=1;its<=30;its++) { /* loop over allowed iterations */
flag=1;
for (l=k;l>=0;l--) { /* test for splitting */
nm=l-1;
if (fabs(rv1[l])+anorm == anorm) {
flag=0;
break;
}
if (fabs (w[nm])+anorm == anorm) break;
}
if (flag) {
c=0.0; /* cancellation of rv1[l] if l > 1 */
s=1.0;
for (i=l;i<=k;i++) {
f=s*rv1[i];
if (fabs (f)+anorm != anorm) {
g=w[i];
h=hypot (f,g);
w[i]=h;
h=1.0/h;
c=g*h;
s=(-1.0*f*h);
for (j=0;j<m;j++) {
y=a[j*n+nm]; /* a[j][nm] */
z=a[j*n+i]; /* a[j][i] */
a[j*n+nm]=(y*c)+(z*s); /* a[j][nm] */
a[j*n+i]=(z*c)-(y*s); /* a[j][i] */
}
}
}
}
z=w[k];
if (l == k) { /* convergence */
if (z < 0.0) { /* singular value is made positive */
w[k]= -1.0*z;
for (j=0;j<n;j++) v[j*n+k] *= (-1.0); /* v[j][k] */
}
break;
}
if (its == 30) {
GMT_Report (GMT->parent, GMT_MSG_WARNING, "Failure in gmt_svdcmp: No convergence in 30 iterations\n");
#ifndef _OPENMP
return (GMT_RUNTIME_ERROR);
#endif
}
x=w[l]; /* shift from bottom 2-by-2 minor */
nm=k-1;
y=w[nm];
g=rv1[nm];
h=rv1[k];
f=((y-z)*(y+z)+(g-h)*(g+h))/(2.0*h*y);
g=hypot (f,1.0);
f=((x-z)*(x+z)+h*((y/(f+SIGN(g,f)))-h))/x;
/* next QR transformation */
c=s=1.0;
for (j=l;j<=nm;j++) {
i=j+1;
g=rv1[i];
y=w[i];
h=s*g;
g=c*g;
z=hypot(f,h);
rv1[j]=z;
c=f/z;
s=h/z;
f=(x*c)+(g*s);
g=(g*c)-(x*s);
h=y*s;
y=y*c;
for (jj=0;jj<n;jj++) {
x=v[jj*n+j]; /* v[jj][j] */
z=v[jj*n+i]; /* v[jj][i] */
v[jj*n+j]=(x*c)+(z*s); /* v[jj][j] */
v[jj*n+i]=(z*c)-(x*s); /* v[jj][i] */
}
z=hypot(f,h);
w[j]=z; /* rotation can be arbitrary if z=0 */
if (z) {
z=1.0/z;
c=f*z;
s=h*z;
}
f=(c*g)+(s*y);
x=(c*y)-(s*g);
for (jj=0;jj<m;jj++) {
y=a[jj*n+j]; /* a[jj][j] */
z=a[jj*n+i]; /* a[jj][i] */
a[jj*n+j]=(y*c)+(z*s); /* a[jj][j] */
a[jj*n+i]=(z*c)-(y*s); /* a[jj][i] */
}
}
rv1[l]=0.0;
rv1[k]=f;
w[k]=x;
}
}
gmt_M_free (GMT, rv1);
return (GMT_NOERROR);
}
#endif
GMT_LOCAL int gmtvector_compare_singular_values (const void *point_1v, const void *point_2v) {
/* Routine for qsort to sort struct GMT_SINGULAR_VALUE on decreasing eigenvalues
* keeping track of the original order before sorting.
*/
const struct GMT_SINGULAR_VALUE *E_1 = point_1v, *E_2 = point_2v;
bool bad_1, bad_2;
bad_1 = gmt_M_is_dnan (E_1->value);
bad_2 = gmt_M_is_dnan (E_2->value);
if (bad_1 && bad_2) return (0);
if (bad_1) return (+1);
if (bad_2) return (-1);
if (fabs(E_1->value) < fabs(E_2->value)) return (+1);
if (fabs(E_1->value) > fabs(E_2->value)) return (-1);
return (0);
}
GMT_LOCAL uint64_t gmtvector_fix_up_path_cartonly (struct GMT_CTRL *GMT, double **a_x, double **a_y, uint64_t n, unsigned int mode) {
/* Takes pointers to a list of <n> Cartesian x/y pairs and adds
* auxiliary points to build a staircase curve.
* If mode=1: staircase; first follows y, then x
* If mode=2: staircase; first follows x, then y
* Returns the new number of points (original plus auxiliary).
*/
uint64_t i, k, n_new = 2*n - 1;
double *x = NULL, *y = NULL;
x = *a_x; y = *a_y; /* Default is to return the input unchanged */
if (n < 2 || mode == 0) return n; /* Nothing to do */
gmt_prep_tmp_arrays (GMT, GMT_NOTSET, 1, 2); /* Init or reallocate two tmp vectors */
/* Start at first point */
GMT->hidden.mem_coord[GMT_X][0] = x[0]; GMT->hidden.mem_coord[GMT_Y][0] = y[0];
for (i = k = 1; i < n; i++) { /* For remaining points we must insert an intermediate node */
gmt_prep_tmp_arrays (GMT, GMT_NOTSET, k+1, 2); /* Init or reallocate two tmp vectors */
if (mode == GMT_STAIRS_X) { /* First follow x, then y */
GMT->hidden.mem_coord[GMT_X][k] = x[i];
GMT->hidden.mem_coord[GMT_Y][k] = y[i-1];
}
else if (mode == GMT_STAIRS_Y) { /* First follow y, then x */
GMT->hidden.mem_coord[GMT_X][k] = x[i-1];
GMT->hidden.mem_coord[GMT_Y][k] = y[i];
}
k++;
/* Then add original point */
gmt_prep_tmp_arrays (GMT, GMT_NOTSET, k+1, 2); /* Init or reallocate two tmp vectors */
GMT->hidden.mem_coord[GMT_X][k] = x[i]; GMT->hidden.mem_coord[GMT_Y][k] = y[i];
k++;
}
/* Destroy old allocated memory and put the new one in place */
gmt_M_free (GMT, x); gmt_M_free (GMT, y);
*a_x = gmtlib_assign_vector (GMT, n_new, GMT_X);
*a_y = gmtlib_assign_vector (GMT, n_new, GMT_Y);
return (n_new);
}
GMT_LOCAL uint64_t gmtvector_fix_up_path_polar (struct GMT_CTRL *GMT, double **a_x, double **a_y, uint64_t n, double step, unsigned int mode) {
/* Takes pointers to a list of <n> theta/r pairs (in user units) and adds
* auxiliary points if the distance between two given points exceeds
* <step> units.
* If mode=0: returns points along a straight line
* If mode=1: staircase; first follows r, then theta
* If mode=2: staircase; first follows theta, then r
* Returns the new number of points (original plus auxiliary).
*/
uint64_t i, j, n_new, n_step = 0;
double *x = NULL, *y = NULL, c;
if (mode == GMT_STAIRS_OFF) return n; /* Nothing to do */
x = *a_x; y = *a_y;
gmt_prep_tmp_arrays (GMT, GMT_NOTSET, 1, 2); /* Init or reallocate tmp vectors */
GMT->hidden.mem_coord[GMT_X][0] = x[0]; GMT->hidden.mem_coord[GMT_Y][0] = y[0]; n_new = 1;
if (step <= 0.0) step = 1.0; /* Sanity valve; if step not given we set it to 1 */
for (i = 1; i < n; i++) {
if (mode == GMT_STAIRS_X) { /* First follow theta, then r */
n_step = lrint (fabs (x[i] - x[i-1]) / step);
for (j = 1; j <= n_step; j++) {
c = j / (double)n_step;
gmt_prep_tmp_arrays (GMT, GMT_NOTSET, n_new, 2); /* Init or reallocate tmp read vectors */
GMT->hidden.mem_coord[GMT_X][n_new] = x[i-1] * (1 - c) + x[i] * c;
GMT->hidden.mem_coord[GMT_Y][n_new] = y[i-1];
n_new++;
}
GMT->hidden.mem_coord[GMT_X][n_new] = x[i]; GMT->hidden.mem_coord[GMT_Y][n_new] = y[i]; n_new++;
}
else if (mode == GMT_STAIRS_Y) { /* First follow r, then theta */
GMT->hidden.mem_coord[GMT_X][n_new] = x[i-1]; GMT->hidden.mem_coord[GMT_Y][n_new] = y[i]; n_new++;
n_step = lrint (fabs (x[i] - x[i-1]) / step);
for (j = 1; j <= n_step; j++) {
c = j / (double)n_step;
gmt_prep_tmp_arrays (GMT, GMT_NOTSET, n_new, 2); /* Init or reallocate tmp read vectors */
GMT->hidden.mem_coord[GMT_X][n_new] = x[i-1] * (1 - c) + x[i] * c;
GMT->hidden.mem_coord[GMT_Y][n_new] = y[i];
n_new++;
}
}
gmt_prep_tmp_arrays (GMT, GMT_NOTSET, n_new, 2); /* Init or reallocate tmp read vectors */
}
/* Destroy old allocated memory and put the new one in place */
gmt_M_free (GMT, x); gmt_M_free (GMT, y);
*a_x = gmtlib_assign_vector (GMT, n_new, GMT_X);
*a_y = gmtlib_assign_vector (GMT, n_new, GMT_Y);
return (n_new);
}
GMT_LOCAL uint64_t gmtvector_fix_up_path_cartesian (struct GMT_CTRL *GMT, double **a_x, double **a_y, uint64_t n, double step, unsigned int mode) {
/* Takes pointers to a list of <n> x/y pairs (in user units) and adds
* auxiliary points if the distance between two given points exceeds
* <step> units.
* If mode=0: returns points along a straight line
* If mode=1: staircase; first follows y, then x
* If mode=2: staircase; first follows x, then y
* Returns the new number of points (original plus auxiliary).
*/
unsigned int k = 1;
uint64_t i, j, n_new, n_step = 0;
double *x = NULL, *y = NULL, c;
x = *a_x; y = *a_y;
gmt_prep_tmp_arrays (GMT, GMT_NOTSET, 1, 2); /* Init or reallocate tmp vectors */
GMT->hidden.mem_coord[GMT_X][0] = x[0]; GMT->hidden.mem_coord[GMT_Y][0] = y[0]; n_new = 1;
if (step <= 0.0) step = 1.0; /* Sanity valve; if step not given we set it to 1 */
for (i = 1; i < n; i++) {
if (mode == GMT_STAIRS_X) { /* First follow x, then y */
n_step = lrint (fabs (x[i] - x[i-1]) / step);
for (j = 1; j <= n_step; j++) {
c = j / (double)n_step;
gmt_prep_tmp_arrays (GMT, GMT_NOTSET, n_new, 2); /* Init or reallocate tmp read vectors */
GMT->hidden.mem_coord[GMT_X][n_new] = x[i-1] * (1 - c) + x[i] * c;
GMT->hidden.mem_coord[GMT_Y][n_new] = y[i-1];
n_new++;
}
n_step = lrint (fabs (y[i]-y[i-1]) / step);
for (j = k; j < n_step; j++) { /* Start at 0 to make sure corner point is saved */
c = j / (double)n_step;
gmt_prep_tmp_arrays (GMT, GMT_NOTSET, n_new, 2); /* Init or reallocate tmp read vectors */
GMT->hidden.mem_coord[GMT_X][n_new] = x[i];
GMT->hidden.mem_coord[GMT_Y][n_new] = y[i-1] * (1 - c) + y[i] * c;
n_new++;
}
k = 0;
}
else if (mode == GMT_STAIRS_Y) { /* First follow y, then x */
n_step = lrint (fabs (y[i]-y[i-1]) / step);
for (j = 1; j <= n_step; j++) {
c = j / (double)n_step;
gmt_prep_tmp_arrays (GMT, GMT_NOTSET, n_new, 2); /* Init or reallocate tmp read vectors */
GMT->hidden.mem_coord[GMT_X][n_new] = x[i-1];
GMT->hidden.mem_coord[GMT_Y][n_new] = y[i-1] * (1 - c) + y[i] * c;
n_new++;
}
n_step = lrint (fabs (x[i]-x[i-1]) / step);
for (j = k; j < n_step; j++) { /* Start at 0 to make sure corner point is saved */
c = j / (double)n_step;
gmt_prep_tmp_arrays (GMT, GMT_NOTSET, n_new, 2); /* Init or reallocate tmp read vectors */
GMT->hidden.mem_coord[GMT_X][n_new] = x[i-1] * (1 - c) + x[i] * c;
GMT->hidden.mem_coord[GMT_Y][n_new] = y[i];
n_new++;
}
k = 0;
}
/* Follow straight line */
else if ((n_step = lrint (hypot (x[i]-x[i-1], y[i]-y[i-1]) / step)) > 1) { /* Must insert (n_step - 1) points, i.e. create n_step intervals */
for (j = 1; j < n_step; j++) {
c = j / (double)n_step;
gmt_prep_tmp_arrays (GMT, GMT_NOTSET, n_new, 2); /* Init or reallocate tmp read vectors */
GMT->hidden.mem_coord[GMT_X][n_new] = x[i-1] * (1 - c) + x[i] * c;
GMT->hidden.mem_coord[GMT_Y][n_new] = y[i-1] * (1 - c) + y[i] * c;
n_new++;
}
}
gmt_prep_tmp_arrays (GMT, GMT_NOTSET, n_new, 2); /* Init or reallocate tmp read vectors */
GMT->hidden.mem_coord[GMT_X][n_new] = x[i]; GMT->hidden.mem_coord[GMT_Y][n_new] = y[i]; n_new++;
}
/* Destroy old allocated memory and put the new one in place */
gmt_M_free (GMT, x); gmt_M_free (GMT, y);
*a_x = gmtlib_assign_vector (GMT, n_new, GMT_X);
*a_y = gmtlib_assign_vector (GMT, n_new, GMT_Y);
return (n_new);
}
GMT_LOCAL uint64_t gmtvector_resample_path_spherical (struct GMT_CTRL *GMT, double **lon, double **lat, uint64_t n_in, double step_out, enum GMT_enum_track mode) {
/* See gmt_resample_path below for details. */
bool meridian, new_pair;
uint64_t last_row_in = 0, row_in, row_out, n_out;
unsigned int k;
double dist_out, gap, d_lon = 0.0, L = 0.0, frac_to_a, frac_to_b, minlon, maxlon, a[3], b[3], c[3];
double P[3], Rot0[3][3], Rot[3][3], total_angle_rad = 0.0, angle_rad, ya = 0.0, yb = 0.0;
double *dist_in = NULL, *lon_out = NULL, *lat_out = NULL, *lon_in = *lon, *lat_in = *lat;
if (step_out < 0.0) { /* Safety valve */
GMT_Report (GMT->parent, GMT_MSG_ERROR, "Internal error: gmtvector_resample_path_spherical given negative step-size\n");
return (GMT_RUNTIME_ERROR);
}
if (mode > GMT_TRACK_SAMPLE_ADJ) { /* Bad mode*/
GMT_Report (GMT->parent, GMT_MSG_ERROR, "Internal error: gmtvector_resample_path_spherical given bad mode %d\n", mode);
return (GMT_RUNTIME_ERROR);
}
if (mode < GMT_TRACK_SAMPLE_FIX) {
if (GMT->current.map.dist[GMT_MAP_DIST].arc) /* Gave an increment in arc length (degree, min, sec) */
step_out /= GMT->current.map.dist[GMT_MAP_DIST].scale; /* Get degrees */
else /* Gave increment in spatial distance (km, meter, etc.) */
step_out = (step_out / GMT->current.map.dist[GMT_MAP_DIST].scale) / GMT->current.proj.DIST_M_PR_DEG; /* Get degrees */
return (gmt_fix_up_path (GMT, lon, lat, n_in, step_out, mode)); /* Insert extra points only */
}
dist_in = gmt_dist_array (GMT, lon_in, lat_in, n_in, true); /* Compute cumulative distances along line */
if (step_out == 0.0) step_out = (dist_in[n_in-1] - dist_in[0])/100.0; /* If nothing is selected we get 101 points */
/* Determine n_out, the number of output points */
if (mode == GMT_TRACK_SAMPLE_ADJ) { /* Round to nearest multiple of step_out, then adjust step to match exactly */
n_out = lrint (dist_in[n_in-1] / step_out);
step_out = dist_in[n_in-1] / n_out; /* Ensure exact fit */
}
else { /* Stop when last multiple is reached */
n_out = lrint (floor (dist_in[n_in-1] / step_out));
}
n_out++; /* Since number of points = number of segments + 1 */
if ((lon_out = gmt_M_memory (GMT, NULL, n_out, double)) == NULL) return (GMT_MEMORY_ERROR);
if ((lat_out = gmt_M_memory (GMT, NULL, n_out, double)) == NULL) return (GMT_MEMORY_ERROR);
lon_out[0] = lon_in[0]; lat_out[0] = lat_in[0]; /* Start at same origin */
for (row_in = row_out = 1; row_out < n_out; row_out++) { /* For remaining output points */
dist_out = row_out * step_out; /* Rhe desired output distance */
while (row_in < n_in && dist_in[row_in] < dist_out) row_in++; /* Wind so row points to the next input point with a distance >= dist_out */
gap = dist_in[row_in] - dist_out; /* Distance to next input datapoint */
new_pair = (row_in > last_row_in);
if (gmt_M_is_zero (gap)) { /* Use the input point as is */
lon_out[row_out] = lon_in[row_in]; lat_out[row_out] = lat_in[row_in];
}
else { /* Must interpolate along great-circle arc from a (point row_in-1) to b (point row_in) */
if (new_pair) { /* Must perform calculations on the two points */
L = dist_in[row_in] - dist_in[row_in-1]; /* Distance between points a and b */
ya = gmt_lat_swap (GMT, lat_in[row_in-1], GMT_LATSWAP_G2O); /* Convert to geocentric */
yb = gmt_lat_swap (GMT, lat_in[row_in], GMT_LATSWAP_G2O); /* Convert to geocentric */
}
frac_to_a = gap / L;
frac_to_b = 1.0 - frac_to_a;
if (GMT->current.map.loxodrome) { /* Linear resampling along Mercator straight line */
if (new_pair) gmt_M_set_delta_lon (lon_in[row_in-1], lon_in[row_in], d_lon);
a[0] = D2R * lon_in[row_in-1]; a[1] = d_log (GMT, tand (45.0 + 0.5 * ya));
b[0] = D2R * (lon_in[row_in-1] + d_lon); b[1] = d_log (GMT, tand (45.0 + 0.5 * yb));
for (k = 0; k < 2; k++) c[k] = a[k] * frac_to_a + b[k] * frac_to_b; /* Linear interpolation to find output point c */
lon_out[row_out] = c[0] * R2D;
lat_out[row_out] = atand (sinh (c[1]));
lat_out[row_out] = gmt_lat_swap (GMT, lat_out[row_out], GMT_LATSWAP_O2G); /* Convert back to geodetic */
}
else { /* Spherical resampling along segment */
if (new_pair) { /* Must perform calculations on the two points */
gmt_geo_to_cart (GMT, ya, lon_in[row_in-1], a, true);
gmt_geo_to_cart (GMT, yb, lon_in[row_in], b, true);
total_angle_rad = d_acos (gmt_dot3v (GMT, a, b)); /* Get spherical angle from a to b in radians; this is same distance as L */
gmt_cross3v (GMT, a, b, P); /* Get pole P of plane trough a and b (and center of Earth) */
gmt_normalize3v (GMT, P); /* Make sure P has unit length */
gmtlib_init_rot_matrix (Rot0, P); /* Get partial rotation matrix since no actual angle is applied yet */
}
gmt_M_memcpy (Rot, Rot0, 9, double); /* Get a copy of the "0-angle" rotation matrix */
angle_rad = total_angle_rad * frac_to_b; /* Angle we need to rotate from a to c */
gmtlib_load_rot_matrix (angle_rad, Rot, P); /* Build the actual rotation matrix for this angle */
gmt_matrix_vect_mult (GMT, 3U, Rot, a, c); /* Rotate from a to get c */
gmt_cart_to_geo (GMT, &lat_out[row_out], &lon_out[row_out], c, true);
lat_out[row_out] = gmt_lat_swap (GMT, lat_out[row_out], GMT_LATSWAP_O2G); /* Convert back to geodetic */
}
minlon = MIN (lon_in[row_in-1], lon_in[row_in]);
maxlon = MAX (lon_in[row_in-1], lon_in[row_in]);
meridian = doubleAlmostEqualZero (maxlon, minlon); /* A meridian; make sure we get right lon value */
if (meridian)
lon_out[row_out] = minlon;
else if (lon_out[row_out] < minlon)
lon_out[row_out] += 360.0;
else if (lon_out[row_out] > maxlon)
lon_out[row_out] -= 360.0;
}
last_row_in = row_in;
}
/* Destroy old allocated memory and put the new one in place */
gmt_M_free (GMT, lon_in);
gmt_M_free (GMT, lat_in);
gmt_M_free (GMT, dist_in);
*lon = lon_out;
*lat = lat_out;
return (n_out);
}
GMT_LOCAL uint64_t gmtvector_resample_path_cartesian (struct GMT_CTRL *GMT, double **x, double **y, uint64_t n_in, double step_out, enum GMT_enum_track mode) {
/* See gmt_resample_path below for details. */
uint64_t last_row_in = 0, row_in, row_out, n_out;
bool new_pair;
double dist_out, gap, L = 0.0, frac_to_a, frac_to_b;
double *dist_in = NULL, *x_out = NULL, *y_out = NULL, *x_in = *x, *y_in = *y;
if (step_out < 0.0) { /* Safety valve */
GMT_Report (GMT->parent, GMT_MSG_ERROR, "Internal error: gmtvector_resample_path_cartesian given negative step-size\n");
return (GMT_RUNTIME_ERROR);
}
if (mode > GMT_TRACK_SAMPLE_ADJ) { /* Bad mode*/
GMT_Report (GMT->parent, GMT_MSG_ERROR, "Internal error: gmtvector_resample_path_cartesian given bad mode %d\n", mode);
return (GMT_RUNTIME_ERROR);
}
if (mode < GMT_TRACK_SAMPLE_FIX) return (gmtvector_fix_up_path_cartesian (GMT, x, y, n_in, step_out, mode)); /* Insert extra points only */
dist_in = gmt_dist_array (GMT, x_in, y_in, n_in, true); /* Compute cumulative distances along line */
if (step_out == 0.0) step_out = (dist_in[n_in-1] - dist_in[0])/100.0; /* If nothing is selected we get 101 points */
/* Determine n_out, the number of output points */
if (mode == GMT_TRACK_SAMPLE_ADJ) { /* Round to nearest multiples, then adjust step to match exactly */
n_out = lrint (dist_in[n_in-1] / step_out);
step_out = dist_in[n_in-1] / n_out;
}
else { /* Stop when last multiple is reached */
n_out = lrint (floor (dist_in[n_in-1] / step_out));
}
n_out++; /* Since number of points = number of segments + 1 */
if ((x_out = gmt_M_memory (GMT, NULL, n_out, double)) == NULL) return (GMT_MEMORY_ERROR);
if ((y_out = gmt_M_memory (GMT, NULL, n_out, double)) == NULL) return (GMT_MEMORY_ERROR);
x_out[0] = x_in[0]; y_out[0] = y_in[0]; /* Start at same origin */
for (row_in = row_out = 1; row_out < n_out; row_out++) { /* For remaining output points */
dist_out = row_out * step_out; /* The desired output distance */
while (row_in < n_in && dist_in[row_in] < dist_out) row_in++; /* Wind so row points to the next input point with a distance >= dist_out */
gap = dist_in[row_in] - dist_out; /* Distance to next input datapoint */
new_pair = (row_in > last_row_in);
if (gmt_M_is_zero (gap)) { /* Use the input point as is */
x_out[row_out] = x_in[row_in]; y_out[row_out] = y_in[row_in];
}
else { /* Must interpolate along great-circle arc from a (point row_in-1) to b (point row_in) */
if (new_pair) L = dist_in[row_in] - dist_in[row_in-1]; /* Distance between points a and b */
frac_to_a = gap / L;
frac_to_b = 1.0 - frac_to_a;
x_out[row_out] = x_in[row_in-1] * frac_to_a + x_in[row_in] * frac_to_b; /* Linear interpolation to find output point */
y_out[row_out] = y_in[row_in-1] * frac_to_a + y_in[row_in] * frac_to_b; /* Linear interpolation to find output point */
}
last_row_in = row_in;
}
/* Destroy old allocated memory and put the new one in place */
gmt_M_free (GMT, x_in);
gmt_M_free (GMT, y_in);
gmt_M_free (GMT, dist_in);
*x = x_out;
*y = y_out;
return (n_out);
}
int gmt_jacobi (struct GMT_CTRL *GMT, double *a, unsigned int n, unsigned int m, double *d, double *v, double *b, double *z, unsigned int *nrots) {
/*
*
* Find eigenvalues & eigenvectors of a square symmetric matrix by Jacobi's
* method. Given A, find V and D such that A = V * D * V-transpose, with
* V an orthogonal matrix and D a diagonal matrix. The eigenvalues of A
* are on diag(D), and the j-th column of V is the eigenvector corresponding
* to the j-th diagonal element of D. Returns 0 if OK, -1 if it fails to
* converge in MAX_SWEEPS.
*
* a is sent as a square symmetric matrix, of size n, and row dimension m.
* Only the diagonal and super-diagonal elements of a will be used, so the
* sub-diagonal elements could be used to preserve a, or could have been
* destroyed by an earlier attempt to form the Cholesky decomposition of a.
* On return, the super-diagonal elements are destroyed. The diagonal and
* sub-diagonal elements are unchanged.
* d is returned as an n-vector containing the eigenvalues of a, sorted
* so that d[i] >= d[j] when i < j. d = diag(D).
* v is returned as an n by n matrix, V, with row dimension m, and the
* columns of v are the eigenvectors corresponding to the values in d.
* b is an n-vector of workspace, used to keep a copy of the diagonal
* elements which is updated only after a full sweep.
* z is an n-vector of workspace, used to accumulate the updates to
* the diagonal values of a during each sweep. This reduces round-
* off problems.
* nrots is the number of rotations performed. Bounds on round-off error
* can be estimated from this if desired.
*
* Numerical Details:
* The basic algorithms is in many textbooks. The idea is to make an
* infinite series (which turns out to be at quadratically convergent)
* of steps, in each of which A_new = P-transpose * A_old * P, where P is
* a plane-rotation matrix in the p,q plane, through an angle chosen to
* zero A_new(p,q) and A_new(q,p). The sum of the diagonal elements
* of A is unchanged by these operations, but the sum of squares of
* diagonal elements of a is increased by 2 * |A_old(p,q)| at each step.
* Although later steps make non-zero again the previously zeroed entries,
* the sum of squares of diagonal elements increases with each rotation,
* while the sum of squares of off-diagonals keeps decreasing, so that
* eventually A_new is diagonal to machine precision. This should
* happen in a few (3 to 7) sweeps.
*
* If only the eigenvalues are wanted then there are faster methods, but
* if all eigenvalues and eigenvectors are needed, then this method is
* only somewhat slower than the fastest method (Householder tri-
* diagonalization followed by symmetric QR iterations), and this method
* is numerically extremely stable.
*
* C G J Jacobi ("Ueber ein leichtes Vefahren, die in der Theorie der
* Saekularstoerungen vorkommenden Gelichungen numerisch aufzuloesen",
* Crelle's Journal, v. 30, pp. 51--94, 1846) originally searched the
* entire (half) matrix for the largest |A(p,q)| to select each step.
* When the method was developed for machine computation (R T Gregory,
* "Computing eigenvalues and eigenvectors of a symmetric matrix on
* the ILLIAC", Math. Tab. and other Aids to Comp., v. 7, pp. 215--220,
* 1953) it was done with a series of "sweeps" through the upper triangle,
* visiting all p,q in turn. Later, D A Pope and C Tompkins ("Maximizing
* functions of rotations - experiments concerning speed of diagonalization
* of symmetric matrices using Jacobi's method", J Assoc. Comput. Mach.
* v. 4, pp. 459--466, 1957) introduced a variant that skips small
* elements on the first few sweeps. The algorithm here was given by
* Heinz Rutishauser (1918--1970) and published in Numer. Math. v. 9,
* pp 1--10, 1966, and in Linear Algebra (the Handbook for Automatic
* Computation, v. II), by James Hardy Wilkinson and C. Reinsch (Springer-
* Verlag, 1971). It also appears in Numerical Recipes.
*
* This algorithm takes care to avoid round-off error in several ways.
* First, although there are four values of theta in (-pi, pi] that
* would zero A(p,q), there is only one with magnitude <= pi/4.
* This one is used. This is most stable, and also has the effect
* that, if A_old(p,p) >= A_old(q,q) then A_new(p,p) > A_new(q,q).
* Two copies of the diagonal elements are maintained in d[] and b[].
* d[] is updated immediately in each rotation, and each new rotation
* is computed based on d[], so that each rotation gets the benefit
* of the previous ones. However, z[] is also used to accumulate
* the sum of all the changes in the diagonal elements during one sweep,
* and z[] is used to update b[] after each sweep. Then b is copied
* to d. In this way, at the end of each sweep, d is reset to avoid
* accumulating round-off.
*
* This routine determines whether y is small compared to x by testing
* if (fabs(y) + fabs(x) == fabs(x) ). It is assumed that the
* underflow which may occur here is nevertheless going to allow this
* expression to be evaluated as true or false and execution to
* continue. If the run environment doesn't allow this, the routine
* won't work properly.
*
* programmer: W. H. F. Smith, 7 June, 1991.
* Revised: PW: 12-MAR-1998 for GMT 3.1
* Revision by WHF Smith, March 03, 2000, to speed up loop indexes.
*/
unsigned int p, q, pp, pq, mp1, pm, qm, nsweeps, j, jm, i, k;
double sum, threshold, g, h, t, theta, c, s, tau;
/* Begin by initializing v, b, d, and z. v = identity matrix,
b = d = diag(a), and z = 0: */
gmt_M_memset (v, m*n, double);
gmt_M_memset (z, n, double);
mp1 = m + 1;
for (p = 0, pp = 0; p < n; p++, pp+=mp1) {
v[pp] = 1.0;
b[p] = a[pp];
d[p] = b[p];
}
/* End of initializations. Set counters and begin: */
(*nrots) = 0;
nsweeps = 0;
while (nsweeps < MAX_SWEEPS) {
/* Sum off-diagonal elements of upper triangle. */
sum = 0.0;
for (q = 1, qm = m; q < n; q++, qm += m ) {
for (p = 0, pq = qm; p < q; p++, pq++) sum += fabs(a[pq]);
}
/* Exit this loop (converged) when sum == 0.0 */
if (sum == 0.0) break;
/* If (nsweeps < 3) do only bigger elements; else all */
threshold = (nsweeps < 3) ? 0.2 * sum / ( n * n ) : 0.0;
/* Now sweep whole upper triangle doing Givens rotations: */
for (q = 1, qm = m; q < n; q++, qm += m ) {
for (p = 0, pm = 0, pq = qm; p < q; p++, pm += m, pq++) {
/* In 3/2000 I swapped order of these loops,
to allow simple incrementing of pq */
if (a[pq] == 0.0) continue; /* New 3/2000 */
g = 100.0 * fabs (a[pq]);
/* After four sweeps, if g is small relative
to a(p,p) and a(q,q), skip the
rotation and set a(p,q) to zero. */
if ((nsweeps > 3) && ((fabs (d[p])+g) == fabs (d[p])) && ((fabs (d[q])+g) == fabs (d[q]))) {
a[pq] = 0.0;
}
else if (fabs (a[pq]) > threshold) {
h = d[q] - d[p];
if (h == 0.0)
t = 1.0; /* This if block is new 3/2000 */
else if (fabs (h) + g == fabs (h))
t = a[pq] / h;
else {
theta = 0.5 * h / a[pq];
t = 1.0 / (fabs (theta) + sqrt (1.0 + theta*theta) );
if (theta < 0.0) t = -t;
}
c = 1.0 / sqrt (1.0 + t*t);
s = t * c;
tau = s / (1.0 + c);
h = t * a[pq];
z[p] -= h;
z[q] += h;
d[p] -= h;
d[q] += h;
a[pq] = 0.0;
for (j = 0; j < p; j++) {
g = a[j + pm];
h = a[j + qm];
a[j + pm] = g - s * (h + g * tau);
a[j + qm] = h + s * (g - h * tau);
}
for (j = p+1, jm = m*(p+1); j < q; j++, jm += m ) {
g = a[p + jm];
h = a[j + qm];
a[p + jm] = g - s * (h + g * tau);
a[j + qm] = h + s * (g - h * tau);
}
for (j = q+1, jm = m*(q+1); j < n; j++, jm += m ) {
g = a[p + jm];
h = a[q + jm];
a[p + jm] = g - s * (h + g * tau);
a[q + jm] = h + s * (g - h * tau);
}
for (j = 0; j < n; j++) {
g = v[j + pm];
h = v[j + qm];
v[j + pm] = g - s * (h + g * tau);
v[j + qm] = h + s * (g - h * tau);
}
(*nrots)++;
}
}
}
/* End of one sweep of the upper triangle. */
nsweeps++;
for (p = 0; p < n; p++) {
b[p] += z[p]; /* Update the b copy of diagonal */
d[p] = b[p]; /* Replace d with b to reduce round-off error */
z[p] = 0.0; /* Clear z. */
}
}
/* Get here via break when converged, or when nsweeps == MAX_SWEEPS.
Sort eigenvalues by insertion: */
for (i = 0; i < n-1; i++) {
k = i;
g = d[i];
for (j = i+1; j < n; j++) { /* Find max location */
if (d[j] >= g) {
k = j;
g = d[j];
}
}
if (k != i) { /* Need to swap value and vector */
d[k] = d[i];
d[i] = g;
p = i * m;
q = k * m;
for (j = 0; j < n; j++) {
g = v[j + p];
v[j + p] = v[j + q];
v[j + q] = g;
}
}
}
/* Return 0 if converged; else print warning and return -1: */
if (nsweeps == MAX_SWEEPS) {
GMT_Report (GMT->parent, GMT_MSG_WARNING, "gmt_jacobi failed to converge in %d sweeps\n", nsweeps);
return(-1);
}
return(0);
}
int gmt_gauss (struct GMT_CTRL *GMT, double *a, double *vec, unsigned int n, unsigned int nstore, bool itriag) {
/* subroutine gauss, by william menke */
/* july 1978 (modified feb 1983, nov 85) */
/* a subroutine to solve a system of n linear equations in n unknowns
* gaussian reduction with partial pivoting is used
* a (sent, destroyed) n by n matrix
* vec (sent, overwritten) n vector, replaced w/ solution
* nstore (sent) dimension of a
* ierror (returned) zero on no error
* itriag (sent) matrix triangularized only
* on true useful when solving
* multiple systems with same a
*/
static unsigned int l1;
unsigned int *line = NULL, i = 0, j, k, l, j1, j2, *isub = NULL;
int iet, ieb;
size_t n_alloc = 0;
double big, testa, b, sum;
iet = 0; /* initial error flags, one for triagularization*/
ieb = 0; /* one for backsolving */
gmt_M_malloc2 (GMT, line, isub, n, &n_alloc, unsigned int);
/* triangularize the matrix a */
/* replacing the zero elements of the triangularized matrix */
/* with the coefficients needed to transform the vector vec */
if (itriag) { /* triangularize matrix */
for (j = 0; j < n; j++) { /*line is an array of flags*/
line[j] = 0;
/* elements of a are not moved during pivoting*/
/* line=0 flags unused lines */
}
for (j=0; j<n-1; j++) {
/* triangularize matrix by partial pivoting */
big = 0.0; /* find biggest element in j-th column*/
/* of unused portion of matrix*/
for (l1=0; l1<n; l1++) {
if (line[l1]==0) {
testa = fabs ((*(a+l1*nstore+j)));
if (testa>big) {
i=l1;
big=testa;
}
}
}
if (big<=DBL_EPSILON) iet=1; /* test for div by 0 */
line[i]=1; /* selected unused line becomes used line */
isub[j]=i; /* isub points to j-th row of tri. matrix */
sum=1.0/(*(a+i*nstore+j));
/*reduce matrix towards triangle */
for (k=0; k<n; k++) {
if (line[k]==0) {
b=(*(a+k*nstore+j))*sum;
for (l=j+1; l<n; l++) *(a+k*nstore+l) -= (b*(*(a+i*nstore+l)));
*(a+k*nstore+j)=b;
}
}
}
for( j=0; j<n; j++ ) {
/* find last unused row and set its pointer */
/* this row contains the apex of the triangle */
if (line[j]==0) {
l1=j; /* apex of triangle */
isub[n-1]=j;
break;
}
}