-
Notifications
You must be signed in to change notification settings - Fork 337
/
greenspline.c
2558 lines (2345 loc) · 111 KB
/
greenspline.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-2020 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: Paul Wessel
* Date: 1-JAN-2010
* Version: 6 API
*
* Brief synopsis: greenspline grids data using Green's functions for a selected spline.
* The data may be Cartesian or geographical, and gridding can be done
* in a Cartesian 1-D, 2-D, 3-D space or on a sphere. The spline may be evaluated on
* a grid, on parts of a grid, or at specified arbitrary locations.
* Five classes of splines (for a total of 10 splines) are implemented:
* 1. Minimum curvature Cartesian spline [Sandwell, Geophys. Res. Lett, 1987]
* 2. Minimum curvature Cartesian spline in tension [Wessel & Bercovici, Math., Geol., 1998]
* 3. Regularized Cartesian spline in tension [Mitasova & Mitas, Math. Geol., 1993]
* 4. Minimum curvature spherical spline [Parker, "Geophysical Inverse Theory", 1994]
* 5. Minimum curvature spherical spline in tension [Wessel & Becker, Geophys. J. Int, 2008]
*
* Originally published as:
* "Wessel, P., 2009. A general-purpose Green's function-based interpolator,
* Computers & Geosciences, 35: 1247-1254".
*
* PW Update June 2013. The numerical implementation of the Green's function found by Wessel & Becker [2008]
* was unstable (it required the difference between k*P_v and log, and P_v is difficult to compute accurately).
* Bob Parker (Scripps) helped develop a series solution for canceling out the two singular behaviors,
* resulting in a more stable expression that converges reasonably rapidly. We now use this new series
* solution for -Sq, combined with a (new) cubic spline interpolation. This replaces the old -SQ machinery
* with linear interpolation which is now deprecated.
*
* PW Update July 2015. With help from Dong Ju Choi, San Diego Supercomputing Center, we have added Open MP
* support in greenspline and the matrix solvers in gmt_vector.c. Requires open MP support and use of -x.
*/
#include "gmt_dev.h"
#define THIS_MODULE_CLASSIC_NAME "greenspline"
#define THIS_MODULE_MODERN_NAME "greenspline"
#define THIS_MODULE_LIB "core"
#define THIS_MODULE_PURPOSE "Interpolate using Green's functions for splines in 1-3 dimensions"
#define THIS_MODULE_KEYS "<D{,AD(=,ED),ND(,TG(,CD)=f,G?},GDN"
#define THIS_MODULE_NEEDS "R"
#define THIS_MODULE_OPTIONS "-:>Vbdefghioqrs" GMT_OPT("FH") GMT_ADD_x_OPT
EXTERN_MSC int gmtlib_cspline (struct GMT_CTRL *GMT, double *x, double *y, uint64_t n, double *c);
/* Control structure for greenspline */
struct GREENSPLINE_CTRL {
struct GREENSPLINE_A { /* -A<gradientfile> */
bool active;
unsigned int mode; /* 0 = azimuths, 1 = directions, 2 = dx,dy components, 3 = dx, dy, dz components */
char *file;
} A ;
struct GREENSPLINE_C { /* -C[n]<cutoff>[+f<file>] */
bool active;
unsigned int movie; /* Undocumented and not-yet-working movie mode +m incremental grids, +M total grids vs eigenvalue */
unsigned int mode;
double value;
char *file;
} C;
struct GREENSPLINE_D { /* -D<distflag> */
bool active;
int mode; /* Can be negative */
} D;
struct GREENSPLINE_E { /* -E[<file>] */
bool active;
unsigned int mode;
char *file;
} E;
struct GREENSPLINE_G { /* -G<output_grdfile> */
bool active;
char *file;
} G;
struct GREENSPLINE_I { /* -Idx[/dy[/dz]] */
bool active;
double inc[3];
} I;
struct GREENSPLINE_L { /* -L */
bool active;
} L;
struct GREENSPLINE_M { /* -M<gfuncfile> */
bool active;
unsigned int mode; /* GMT_IN or GMT_OUT */
char *file;
} M;
struct GREENSPLINE_N { /* -N<outputnode_file> */
bool active;
char *file;
} N;
struct GREENSPLINE_Q { /* -Qdaz */
bool active;
double az;
double dir[3];
} Q;
struct GREENSPLINE_R3 { /* -Rxmin/xmax[/ymin/ymax[/zmin/zmaz]] | -Ggridfile */
bool active;
bool mode; /* true if settings came from a grid file */
unsigned int dimension; /* 1, 2, or 3 */
unsigned int offset; /* 0 or 1 */
double range[6]; /* Min/max for each dimension */
double inc[2]; /* xinc/yinc when -Rgridfile was given*/
} R3;
struct GREENSPLINE_S { /* -S<mode>[<tension][+<mod>[args]] */
bool active;
unsigned int mode;
double value[4];
double rval[2];
char *arg;
} S;
struct GREENSPLINE_T { /* -T<mask_grdfile> */
bool active;
char *file;
} T;
struct GREENSPLINE_W { /* -W[w] */
bool active;
unsigned int mode; /* 0 = got sigmas, 1 = got weights */
} W;
struct GREENSPLINE_Z { /* -Z undocumented debugging option */
bool active;
} Z;
};
enum Greenspline_modes { /* Various integer mode flags */
SANDWELL_1987_1D = 0,
SANDWELL_1987_2D = 1,
SANDWELL_1987_3D = 2,
WESSEL_BERCOVICI_1998_1D = 3,
WESSEL_BERCOVICI_1998_2D = 4,
WESSEL_BERCOVICI_1998_3D = 5,
MITASOVA_MITAS_1993_2D = 6,
MITASOVA_MITAS_1993_3D = 7,
PARKER_1994 = 8,
WESSEL_BECKER_2008 = 9,
LINEAR_1D = 10,
LINEAR_2D = 11,
N_METHODS = 12,
N_PARAMS = 11,
GREENSPLINE_TREND = 1, /* Remove/Restore linear trend */
GREENSPLINE_NORM = 2, /* Normalize residual data to 0-1 range */
SQ_N_NODES = 10001, /* Default number of nodes in the precalculated -Sq spline */
GSP_GOT_SIG = 0,
GSP_GOT_WEIGHTS = 1,
GREENSPLINE_NO_MOVIE = 0,
GREENSPLINE_INC_MOVIE = 1,
GREENSPLINE_CUM_MOVIE = 2,
GREENSPLINE_EIGEN_RATIO = 0,
GREENSPLINE_EIGEN_CUTOFF = 1
};
#ifndef M_LOG_2
#define M_LOG_2 0.69314718055994530942
#endif
#ifndef M_GAMMA
#define M_GAMMA 0.577215664901532860606512
#endif
#ifndef M_SQRT_PI
#define M_SQRT_PI 1.772453850905516027298167483341
#endif
#ifndef M_INV_SQRT_PI
#define M_INV_SQRT_PI (1.0 / M_SQRT_PI)
#endif
#define SQ_TRUNC_ERROR 1.0e-6 /* Max truncation error in Parker's simplified sum for WB'08 */
enum Greenspline_index { /* Indices for coeff array for normalization */
GSP_MEAN_X = 0,
GSP_MEAN_Y = 1,
GSP_MEAN_Z = 2,
GSP_SLP_X = 3,
GSP_SLP_Y = 4,
GSP_RANGE = 5,
GSP_LENGTH = 6};
struct GREENSPLINE_LOOKUP { /* Used to spline interpolation of precalculated function */
uint64_t n; /* Number of values in the spline setup */
double *y; /* Function values */
double *c; /* spline coefficients */
double *A, *B, *C; /* power/ratios of order l terms */
};
struct GREENSPLINE_ZGRID {
unsigned int nz;
double z_min, z_max, z_inc;
};
#ifdef DEBUG
static bool TEST = false; /* Global variable used for undocumented testing [under -DDEBUG only; see -+ hidden option] */
#endif
static void *New_Ctrl (struct GMT_CTRL *GMT) { /* Allocate and initialize a new control structure */
struct GREENSPLINE_CTRL *C;
C = gmt_M_memory (GMT, NULL, 1, struct GREENSPLINE_CTRL);
/* Initialize values whose defaults are not 0/false/NULL */
C->S.mode = SANDWELL_1987_2D;
C->S.rval[0] = -1.0; C->S.rval[1] = 1.0;
C->S.value[3] = (double)SQ_N_NODES; /* Default number of spline nodes */
C->S.value[2] = SQ_TRUNC_ERROR; /* Default truncation error for Legendre sum in -Sq */
return (C);
}
static void Free_Ctrl (struct GMT_CTRL *GMT, struct GREENSPLINE_CTRL *C) { /* Deallocate control structure */
if (!C) return;
gmt_M_str_free (C->A.file);
gmt_M_str_free (C->C.file);
gmt_M_str_free (C->G.file);
gmt_M_str_free (C->N.file);
gmt_M_str_free (C->T.file);
gmt_M_str_free (C->S.arg);
gmt_M_free (GMT, C);
}
static int usage (struct GMTAPI_CTRL *API, int level) {
const char *name = gmt_show_name_and_purpose (API, THIS_MODULE_LIB, THIS_MODULE_CLASSIC_NAME, THIS_MODULE_PURPOSE);
if (level == GMT_MODULE_PURPOSE) return (GMT_NOERROR);
GMT_Message (API, GMT_TIME_NONE, "usage: %s [<table>] -G<outfile> [-A<gradientfile>+f<format>] [-C[n]<val>[%%][+f<file>][+m|M]]\n", name);
GMT_Message (API, GMT_TIME_NONE, "\t[-D<mode>] [-E[<misfitfile>] [-I<dx>[/<dy>[/<dz>]] [-L] [-N<nodefile>] [-Q<az>]\n");
GMT_Message (API, GMT_TIME_NONE, "\t[-R<xmin>/<xmax[/<ymin>/<ymax>[/<zmin>/<zmax>]]] [-Sc|l|t|r|p|q[<pars>]] [-T<maskgrid>]\n");
GMT_Message (API, GMT_TIME_NONE, "\t[%s] [-W[w]] [%s] [%s] [%s]\n\t[%s] [%s]\n\t[%s] [%s]\n\t[%s] [%s] [%s]%s[%s] [%s]\n\n", GMT_V_OPT,
GMT_bi_OPT, GMT_d_OPT, GMT_e_OPT, GMT_g_OPT, GMT_h_OPT, GMT_i_OPT, GMT_o_OPT, GMT_q_OPT, GMT_r_OPT, GMT_s_OPT, GMT_x_OPT, GMT_colon_OPT, GMT_PAR_OPT);
if (level == GMT_SYNOPSIS) return (GMT_MODULE_SYNOPSIS);
GMT_Message (API, GMT_TIME_NONE, "\tChoose one of three ways to specify where to evaluate the spline:\n");
GMT_Message (API, GMT_TIME_NONE, "\t1. Specify a rectangular grid domain with options -R, -I [and optionally -r].\n");
GMT_Message (API, GMT_TIME_NONE, "\t2. Supply a mask file via -T whose values are NaN or 0. The spline will then\n");
GMT_Message (API, GMT_TIME_NONE, "\t only be evaluated at the nodes originally set to zero.\n");
GMT_Message (API, GMT_TIME_NONE, "\t3. Specify a set of output locations via the -N option.\n\n");
GMT_Message (API, GMT_TIME_NONE, "\t-G Output data. Give name of output file.\n");
GMT_Message (API, GMT_TIME_NONE, "\n\tOPTIONS:\n");
GMT_Option (API, "<");
GMT_Message (API, GMT_TIME_NONE, "\t-A ASCII file with surface gradients V to use in the modeling. Specify format:\n");
GMT_Message (API, GMT_TIME_NONE, "\t (0) For 1-D: x, slope, (1) X, Vmagnitude, Vazimuth(s), (2) X, Vazimuth(s), Vmagnitude,\n");
GMT_Message (API, GMT_TIME_NONE, "\t (3) X, Vmagnitude, Vangle(s), (4) X, Vcomponents, or (5) X, Vunit-vector, Vmagnitude.\n");
GMT_Message (API, GMT_TIME_NONE, "\t Here, X = (x, y[, z]) is the position vector, V is the gradient vector.\n");
GMT_Message (API, GMT_TIME_NONE, "\t-C Solve by SVD and eliminate eigenvalues whose ratio to largest eigenvalue is less than <val> [0].\n");
GMT_Message (API, GMT_TIME_NONE, "\t Optionally append +f<filename> to save the eigenvalues to this file.\n");
GMT_Message (API, GMT_TIME_NONE, "\t A negative cutoff will stop execution after saving the eigenvalues.\n");
GMT_Message (API, GMT_TIME_NONE, "\t Use -Cn to select only the largest <val> eigenvalues [all].\n");
GMT_Message (API, GMT_TIME_NONE, "\t Use <val>%% to select a percentage of the eigenvalues instead.\n");
GMT_Message (API, GMT_TIME_NONE, "\t [Default uses Gauss-Jordan elimination to solve the linear system]\n");
GMT_Message (API, GMT_TIME_NONE, "\t The +m|M modifiers are valid for 2-D gridding only and will create a series of intermediate\n");
GMT_Message (API, GMT_TIME_NONE, "\t grids for each eigenvalue holding the incremental (+m) or cumulative (+M) result.\n");
GMT_Message (API, GMT_TIME_NONE, "\t Requires -G with a filename template containing an integer C formatting specifier (e.g., %%d).\n");
GMT_Message (API, GMT_TIME_NONE, "\t-D Distance flag determines how we calculate distances between (x,y) points:\n");
GMT_Message (API, GMT_TIME_NONE, "\t Options 0 apples to Cartesian 1-D spline interpolation.\n");
GMT_Message (API, GMT_TIME_NONE, "\t -D0 x in user units, Cartesian distances.\n");
GMT_Message (API, GMT_TIME_NONE, "\t Options 1-3 apply to Cartesian 2-D surface spline interpolation.\n");
GMT_Message (API, GMT_TIME_NONE, "\t -D1 x,y in user units, Cartesian distances.\n");
GMT_Message (API, GMT_TIME_NONE, "\t -D2 x,y in degrees, flat Earth distances in meters.\n");
GMT_Message (API, GMT_TIME_NONE, "\t -D3 x,y in degrees, spherical distances in meters.\n");
GMT_Message (API, GMT_TIME_NONE, "\t Option 4 applies to 2-D spherical surface spline interpolation.\n");
GMT_Message (API, GMT_TIME_NONE, "\t -D4 x,y in degrees, use cosine of spherical distances in degrees.\n");
GMT_Message (API, GMT_TIME_NONE, "\t Option 5 applies to Cartesian 3-D volume interpolation.\n");
GMT_Message (API, GMT_TIME_NONE, "\t -D5 x,y,z in user units, Cartesian distances.\n");
GMT_Message (API, GMT_TIME_NONE, "\t For option 3-4, use PROJ_ELLIPSOID to select geodesic or great circle arcs.\n");
GMT_Message (API, GMT_TIME_NONE, "\t-E Evaluate solution at input locations and report misfit statistics.\n");
GMT_Message (API, GMT_TIME_NONE, "\t Append filename to save all data with two extra columns for model and misfit.\n");
GMT_Message (API, GMT_TIME_NONE, "\t-I Specify a regular set of output locations. Give equidistant increment for each dimension.\n");
GMT_Message (API, GMT_TIME_NONE, "\t Requires -R for specifying the output domain.\n");
GMT_Message (API, GMT_TIME_NONE, "\t-L Leave trend alone. Do not remove least squares plane from data before spline fit.\n");
GMT_Message (API, GMT_TIME_NONE, "\t Only applies to -D0-2. [Default removes linear trend, fits residuals, and restores trend].\n");
GMT_Message (API, GMT_TIME_NONE, "\t-N ASCII file with desired output locations.\n");
GMT_Message (API, GMT_TIME_NONE, "\t The resulting ASCII coordinates and interpolation are written to file given in -G\n");
GMT_Message (API, GMT_TIME_NONE, "\t or stdout if no file specified (see -bo for binary output).\n");
GMT_Message (API, GMT_TIME_NONE, "\t-Q Calculate the directional derivative in the <az> direction and return it instead of surface elevation.\n");
GMT_Option (API, "R");
if (gmt_M_showusage (API)) {
GMT_Message (API, GMT_TIME_NONE, "\t-R Specify a regular set of output locations. Give min and max coordinates for each dimension.\n");
GMT_Message (API, GMT_TIME_NONE, "\t Requires -I for specifying equidistant increments. For 2D-gridding a gridfile may be given;\n");
GMT_Message (API, GMT_TIME_NONE, "\t this then also sets -I (and perhaps -r); use those options to override the grid settings.\n");
}
GMT_Message (API, GMT_TIME_NONE, "\t-S Specify which spline to use; except for c|p, append normalized <tension> between 0 and 1:\n");
GMT_Message (API, GMT_TIME_NONE, "\t -Sc is minimum curvature spline (Sandwell, 1987) [Default].\n");
GMT_Message (API, GMT_TIME_NONE, "\t -Sl is a linear (1-D) or bilinear (2-D) spline.\n");
GMT_Message (API, GMT_TIME_NONE, "\t -St<tension>[/<scale>] is a Cartesian spline in tension (Wessel & Bercovici, 1998).\n");
GMT_Message (API, GMT_TIME_NONE, "\t Optionally, specify a length-scale [Default is the given output spacing].\n");
GMT_Message (API, GMT_TIME_NONE, "\t -Sr<tension> is a regularized spline in tension (Mitasova & Mitas, 1993).\n");
GMT_Message (API, GMT_TIME_NONE, "\t Optionally, specify a length-scale [Default is given output spacing].\n");
GMT_Message (API, GMT_TIME_NONE, "\t -Sp is a spherical surface spline (Parker, 1994); automatically sets -D4.\n");
GMT_Message (API, GMT_TIME_NONE, "\t -Sq is a spherical surface spline in tension (Wessel & Becker, 2008); automatically sets -D4.\n");
GMT_Message (API, GMT_TIME_NONE, "\t Append +e<error> to change maximum error in series truncation [%g].\n", SQ_TRUNC_ERROR);
GMT_Message (API, GMT_TIME_NONE, "\t Append +n<n> to change the (odd) number of precalculated nodes for spline interpolation [%d].\n", SQ_N_NODES);
GMT_Message (API, GMT_TIME_NONE, "\t-T Mask grid file whose values are NaN or 0; its header implicitly sets -R, -I (and -r).\n");
GMT_Message (API, GMT_TIME_NONE, "\t-W Expects one extra input column with data errors sigma_i.\n");
GMT_Message (API, GMT_TIME_NONE, "\t Append w to indicate this column carries weights instead.\n");
GMT_Message (API, GMT_TIME_NONE, "\t [Default makes weights via w_i = 1/sigma_i^2] for the least squares solution.\n");
GMT_Message (API, GMT_TIME_NONE, "\t Note this will only have an effect if -C is used.\n");
GMT_Option (API, "V,bi");
if (gmt_M_showusage (API)) GMT_Message (API, GMT_TIME_NONE, "\t Default is 2-5 input columns depending on dimensionality (see -D) and weights (see -W).\n");
GMT_Option (API, "d,e,g,h,i,o,q,r,s,x,:,.");
return (GMT_MODULE_USAGE);
}
static int parse (struct GMT_CTRL *GMT, struct GREENSPLINE_CTRL *Ctrl, struct GMT_OPTION *options) {
/* This parses the options provided to greenspline and sets parameters in CTRL.
* Any GMT common options will override values set previously by other commands.
* It also replaces any file names specified as input or output with the data ID
* returned when registering these sources/destinations with the API.
*/
int n_items;
unsigned int n_errors = 0, dimension, k, pos = 0;
char txt[6][GMT_LEN64], p[GMT_BUFSIZ] = {""}, *c = NULL;
struct GMT_OPTION *opt = NULL;
struct GMTAPI_CTRL *API = GMT->parent;
for (opt = options; opt; opt = opt->next) {
switch (opt->option) {
case '<': /* Skip input files */
if (GMT_Get_FilePath (GMT->parent, GMT_IS_DATASET, GMT_IN, GMT_FILE_REMOTE, &(opt->arg))) n_errors++;;
break;
case 'R': /* Normally processed internally but must be handled separately since it can take 1,2,3 dimensions */
GMT->common.R.active[RSET] = true;
Ctrl->R3.dimension = 1; /* At least */
if (opt->arg[0] == 'g' && opt->arg[1] == '\0') { /* Got -Rg */
Ctrl->R3.range[0] = 0.0; Ctrl->R3.range[1] = 360.0; Ctrl->R3.range[2] = -90.0; Ctrl->R3.range[3] = 90.0;
Ctrl->R3.dimension = 2;
break;
}
if (opt->arg[0] == 'g' && opt->arg[1] == '/') { /* Got -Rg/zmin/zmax */
Ctrl->R3.range[0] = 0.0; Ctrl->R3.range[1] = 360.0; Ctrl->R3.range[2] = -90.0; Ctrl->R3.range[3] = 90.0;
n_items = sscanf (&opt->arg[2], "%[^/]/%s", txt[4], txt[5]);
if (n_items != 2) {
GMT_Report (API, GMT_MSG_ERROR, "Option -Rg/z0/z1: Append the z-range\n");
n_errors++;
}
Ctrl->R3.dimension = 3;
break;
}
if (opt->arg[0] == 'd' && opt->arg[1] == '\0') { /* Got -Rd */
Ctrl->R3.range[0] = -180.0; Ctrl->R3.range[1] = 180.0; Ctrl->R3.range[2] = -90.0; Ctrl->R3.range[3] = 90.0;
Ctrl->R3.dimension = 2;
break;
}
if (opt->arg[0] == 'd' && opt->arg[1] == '/') { /* Got -Rd/zmin/zmax */
Ctrl->R3.range[0] = -180.0; Ctrl->R3.range[1] = 180.0; Ctrl->R3.range[2] = -90.0; Ctrl->R3.range[3] = 90.0;
n_items = sscanf (&opt->arg[2], "%[^/]/%s", txt[4], txt[5]);
if (n_items != 2) {
GMT_Report (API, GMT_MSG_ERROR, "Option -Rd/z0/z1: Append the z-range\n");
n_errors++;
}
Ctrl->R3.dimension = 3;
break;
}
if (!gmt_access (GMT, opt->arg, R_OK)) { /* Gave a readable file, presumably a grid */
struct GMT_GRID *G = NULL;
if ((G = GMT_Read_Data (API, GMT_IS_GRID, GMT_IS_FILE, GMT_IS_SURFACE, GMT_CONTAINER_ONLY, NULL, opt->arg, NULL)) == NULL) { /* Get header only */
return (API->error);
}
Ctrl->R3.range[0] = G->header->wesn[XLO]; Ctrl->R3.range[1] = G->header->wesn[XHI];
Ctrl->R3.range[2] = G->header->wesn[YLO]; Ctrl->R3.range[3] = G->header->wesn[YHI];
Ctrl->R3.inc[GMT_X] = G->header->inc[GMT_X]; Ctrl->R3.inc[GMT_Y] = G->header->inc[GMT_Y];
Ctrl->R3.offset = G->header->registration;
Ctrl->R3.dimension = 2;
Ctrl->R3.mode = true;
if (GMT_Destroy_Data (API, &G) != GMT_NOERROR) {
return (API->error);
}
break;
}
/* Only get here if the above cases did not trip */
n_items = sscanf (opt->arg, "%[^/]/%[^/]/%[^/]/%[^/]/%[^/]/%s", txt[0], txt[1], txt[2], txt[3], txt[4], txt[5]);
if (!(n_items == 2 || n_items == 4 || n_items == 6)) {
GMT_Report (API, GMT_MSG_ERROR, "Option -R: Give 2, 4, or 6 coordinates\n");
n_errors++;
}
n_errors += gmt_verify_expectations (GMT, gmt_M_type (GMT, GMT_IN, GMT_X), gmt_scanf_arg (GMT, txt[0], gmt_M_type (GMT, GMT_IN, GMT_X), true, &Ctrl->R3.range[0]), txt[0]);
n_errors += gmt_verify_expectations (GMT, gmt_M_type (GMT, GMT_IN, GMT_X), gmt_scanf_arg (GMT, txt[1], gmt_M_type (GMT, GMT_IN, GMT_X), true, &Ctrl->R3.range[1]), txt[1]);
if (n_items > 2) {
Ctrl->R3.dimension = 2;
n_errors += gmt_verify_expectations (GMT, gmt_M_type (GMT, GMT_IN, GMT_Y), gmt_scanf_arg (GMT, txt[2], gmt_M_type (GMT, GMT_IN, GMT_Y), true, &Ctrl->R3.range[2]), txt[2]);
n_errors += gmt_verify_expectations (GMT, gmt_M_type (GMT, GMT_IN, GMT_Y), gmt_scanf_arg (GMT, txt[3], gmt_M_type (GMT, GMT_IN, GMT_Y), true, &Ctrl->R3.range[3]), txt[3]);
}
if (n_items == 6) {
Ctrl->R3.dimension = 3;
n_errors += gmt_verify_expectations (GMT, gmt_M_type (GMT, GMT_IN, GMT_Z), gmt_scanf_arg (GMT, txt[4], gmt_M_type (GMT, GMT_IN, GMT_Z), false, &Ctrl->R3.range[4]), txt[4]);
n_errors += gmt_verify_expectations (GMT, gmt_M_type (GMT, GMT_IN, GMT_Z), gmt_scanf_arg (GMT, txt[5], gmt_M_type (GMT, GMT_IN, GMT_Z), false, &Ctrl->R3.range[5]), txt[5]);
}
if (Ctrl->R3.dimension > 1) gmt_M_memcpy (GMT->common.R.wesn, Ctrl->R3.range, 4, double);
break;
/* Processes program-specific parameters */
case 'A': /* Gradient data: -A<gradientfile>+f<format> */
Ctrl->A.active = true;
if (strchr (opt->arg, ',')) { /* Old syntax: Specified a particular format with -A<mode>,<file> */
if (gmt_M_compat_check (API->GMT, 5)) {
GMT_Report (API, GMT_MSG_COMPAT, "Option -A<format>,<gradientfile> is deprecated; use -A<gradientfile>+f<format> instead\n");
Ctrl->A.mode = (int)(opt->arg[0] - '0');
Ctrl->A.file = strdup (&opt->arg[2]);
}
else {
GMT_Report (API, GMT_MSG_ERROR, "Option -A: Expect -A>gradientfile>+f<format>\n");
n_errors++;
}
break;
}
/* New syntax */
if ((c = strstr (opt->arg, "+f")) == NULL) {
GMT_Report (API, GMT_MSG_ERROR, "Option -A: Expect -A>gradientfile>+f<format>\n");
n_errors++;
}
else {
Ctrl->A.mode = (int)(c[2] - '0');
c[0] = '\0'; /* Temporarily chop off the modifier */
if (opt->arg[0] == 0) {
GMT_Report (API, GMT_MSG_ERROR, "Option -A: No file given\n");
n_errors++;
}
else
Ctrl->A.file = strdup (opt->arg);
c[0] = '+'; /* Restore the modifier */
}
break;
case 'C': /* Solve by SVD */
Ctrl->C.active = true;
if (opt->arg[0] == 'n') Ctrl->C.mode = GREENSPLINE_EIGEN_CUTOFF;
k = (Ctrl->C.mode) ? 1 : 0;
if (gmt_get_modifier (opt->arg, 'f', p))
Ctrl->C.file = strdup (p);
if (gmt_get_modifier (opt->arg, 'm', p))
Ctrl->C.movie = GREENSPLINE_INC_MOVIE;
else if (gmt_get_modifier (opt->arg, 'M', p))
Ctrl->C.movie = GREENSPLINE_CUM_MOVIE;
if (strchr (opt->arg, '/')) { /* Old-style file specification */
if (gmt_M_compat_check (API->GMT, 5)) { /* OK */
sscanf (&opt->arg[k], "%lf/%s", &Ctrl->C.value, p);
Ctrl->C.file = strdup (p);
}
else {
GMT_Report (API, GMT_MSG_ERROR, "Option -C: Expected -C[n]<cut>[+<file>]\n");
n_errors++;
}
}
else {
if (opt->arg[k]) /* Check if we got percentages or a value */
Ctrl->C.value = (Ctrl->C.mode && strchr (opt->arg, '%')) ? 0.01 * atof (&opt->arg[k]) : atof (&opt->arg[k]);
else
Ctrl->C.value = 0.0;
}
break;
case 'D': /* Distance mode */
Ctrl->D.active = true;
Ctrl->D.mode = atoi (opt->arg); /* Since I added 0 to be 1-D later so now this is mode -1 */
break;
case 'E': /* Evaluate misfit -E[<file>]*/
Ctrl->E.active = true;
if (opt->arg) {
Ctrl->E.file = strdup (opt->arg);
Ctrl->E.mode = 1;
}
break;
case 'G': /* Output file */
Ctrl->G.active = true;
Ctrl->G.file = strdup (opt->arg);
break;
case 'I': /* Table or grid spacings */
Ctrl->I.active = true;
k = gmt_getincn (GMT, opt->arg, Ctrl->I.inc, 3);
if (k < 1) {
gmt_inc_syntax (GMT, 'I', 1);
n_errors++;
}
if (Ctrl->I.inc[GMT_Y] == 0.0) Ctrl->I.inc[GMT_Y] = Ctrl->I.inc[GMT_X];
if (Ctrl->I.inc[GMT_Z] == 0.0) Ctrl->I.inc[GMT_Z] = Ctrl->I.inc[GMT_X];
break;
case 'L': /* Leave trend alone */
Ctrl->L.active = true;
break;
case 'M': /* Read or write list of Green's function forces */
Ctrl->M.active = true;
Ctrl->M.file = strdup (opt->arg);
break;
case 'N': /* Output locations */
Ctrl->N.active = true;
if (opt->arg[0]) Ctrl->N.file = strdup (opt->arg);
if (GMT_Get_FilePath (GMT->parent, GMT_IS_DATASET, GMT_IN, GMT_FILE_REMOTE, &(Ctrl->N.file))) n_errors++;
break;
case 'Q': /* Directional derivative */
Ctrl->Q.active = true;
if (strchr (opt->arg, '/')) { /* Got 3-D vector components */
n_items = sscanf (opt->arg, "%lf/%lf/%lf", &Ctrl->Q.dir[0], &Ctrl->Q.dir[1], &Ctrl->Q.dir[2]);
if (n_items != 3) {
GMT_Report (API, GMT_MSG_ERROR, "Option -Q: Append azimuth (2-D) or x/y/z components (3-D)\n");
n_errors++;
}
gmt_normalize3v (GMT, Ctrl->Q.dir); /* Normalize to unit vector */
}
else if (opt->arg[0]) /* 2-D azimuth */
Ctrl->Q.az = atof(opt->arg);
else {
GMT_Report (API, GMT_MSG_ERROR, "Option -Q: Append azimuth (2-D) or x/y/z components (3-D)\n");
n_errors++;
}
break;
case 'S': /* Spline selection */
Ctrl->S.active = true;
Ctrl->S.arg = strdup (opt->arg);
switch (opt->arg[0]) {
case 'l': /* Cartesian linear spline in 1-D or 2-D (bilinear) */
Ctrl->S.mode = LINEAR_1D;
break;
case 'c': /* Cartesian minimum curvature spline */
Ctrl->S.mode = SANDWELL_1987_1D;
break;
case 't': /* Cartesian minimum curvature spline in tension */
Ctrl->S.mode = WESSEL_BERCOVICI_1998_1D;
if (strchr (opt->arg, '/'))
sscanf (&opt->arg[1], "%lf/%lf", &Ctrl->S.value[0], &Ctrl->S.value[1]);
else
Ctrl->S.value[0] = atof (&opt->arg[1]);
break;
case 'r': /* Regularized minimum curvature spline in tension in 2-D or 3-D */
Ctrl->S.mode = MITASOVA_MITAS_1993_2D;
if (strchr (opt->arg, '/'))
sscanf (&opt->arg[1], "%lf/%lf", &Ctrl->S.value[0], &Ctrl->S.value[1]);
else
Ctrl->S.value[0] = atof (&opt->arg[1]);
break;
case 'p': /* Spherical minimum curvature spline */
Ctrl->S.mode = PARKER_1994;
break;
case 'Q': /* Spherical minimum curvature spline in tension */
if (gmt_M_compat_check (GMT, 4)) {
GMT_Report (API, GMT_MSG_COMPAT, "Option -SQ is deprecated; see -Sq syntax instead.\n");
Ctrl->S.mode = WESSEL_BECKER_2008;
if (strchr (opt->arg, '/')) {
n_items = sscanf (&opt->arg[1], "%lf/%lf/%lf/%lf", &Ctrl->S.value[0], &Ctrl->S.value[1], &Ctrl->S.rval[0], &Ctrl->S.rval[1]);
if (n_items == 2) {
Ctrl->S.rval[0] = -1.0;
Ctrl->S.rval[1] = +1.0;
}
}
else
Ctrl->S.value[0] = atof (&opt->arg[1]);
}
else {
GMT_Report (API, GMT_MSG_ERROR, "Option -S: Append c|l|t|g|p|q\n");
n_errors++;
}
break;
case 'q': /* Spherical minimum curvature spline in tension */
Ctrl->S.mode = WESSEL_BECKER_2008;
Ctrl->S.value[0] = atof (&opt->arg[1]);
if (Ctrl->S.value[0] == 0.0) /* Switch to Parker_1994 since tension is zero */
Ctrl->S.mode = PARKER_1994;
if ((c = strchr (opt->arg, '+')) != NULL) {
while (gmt_strtok (c, "+", &pos, p)) {
switch (p[0]) {
case 'e': Ctrl->S.value[2] = atof (&p[1]); break; /* Change the truncation error limit */
case 'n': Ctrl->S.value[3] = atof (&p[1]); break; /* Change the number of nodes for the spline lookup */
case 'l': Ctrl->S.rval[0] = atof (&p[1]); break; /* Min value for spline, undocumented for testing only */
case 'u': Ctrl->S.rval[1] = atof (&p[1]); break; /* Max value for spline, undocumented for testing only */
default:
GMT_Report (API, GMT_MSG_ERROR, "Option -Sq: Unknown modifier %s\n", p);
n_errors++;
break;
}
}
}
break;
default:
GMT_Report (API, GMT_MSG_ERROR, "Option -S: Append c|l|t|g|p|q[<params>]\n");
n_errors++;
break;
}
break;
case 'T': /* Input mask grid */
Ctrl->T.active = true;
if (opt->arg[0]) Ctrl->T.file = strdup (opt->arg);
if (GMT_Get_FilePath (GMT->parent, GMT_IS_GRID, GMT_IN, GMT_FILE_REMOTE, &(Ctrl->T.file)))
n_errors++;
else { /* Obtain -R -I -r from file */
struct GMT_GRID *G = NULL;
if ((G = GMT_Read_Data (API, GMT_IS_GRID, GMT_IS_FILE, GMT_IS_SURFACE, GMT_CONTAINER_ONLY, NULL, Ctrl->T.file, NULL)) == NULL) { /* Get header only */
return (API->error);
}
Ctrl->R3.range[0] = G->header->wesn[XLO]; Ctrl->R3.range[1] = G->header->wesn[XHI];
Ctrl->R3.range[2] = G->header->wesn[YLO]; Ctrl->R3.range[3] = G->header->wesn[YHI];
Ctrl->R3.inc[GMT_X] = G->header->inc[GMT_X]; Ctrl->R3.inc[GMT_Y] = G->header->inc[GMT_Y];
Ctrl->R3.offset = G->header->registration;
Ctrl->R3.dimension = 2;
Ctrl->R3.mode = true;
if (GMT_Destroy_Data (API, &G) != GMT_NOERROR) {
return (API->error);
}
GMT->common.R.active[RSET] = true;
}
break;
case 'W': /* Expect data uncertainty or weights in last column */
Ctrl->W.active = true;
if (opt->arg[0] == 'w') Ctrl->W.mode = GSP_GOT_WEIGHTS; /* Got weights instead of sigmas */
break;
#ifdef DEBUG
case 'Z': /* Dump matrices */
Ctrl->Z.active = true;
break;
case '+': /* Turn on TEST mode */
TEST = true;
break;
#endif
default: /* Report bad options */
n_errors += gmt_default_error (GMT, opt->option);
break;
}
}
if (Ctrl->M.active) { /* Determine if this is for reading or writing Green's function forces */
if (Ctrl->S.active) /* Writing, since -S was given */
Ctrl->M.mode = GMT_OUT;
else if (gmt_access (GMT, Ctrl->M.file, F_OK)) {
GMT_Report (API, GMT_MSG_ERROR, "-M option given but file %s not found\n", Ctrl->M.file);
n_errors++;
}
else /* Read in previous Green's function forces */
Ctrl->M.mode = GMT_IN;
}
if (Ctrl->S.mode == WESSEL_BECKER_2008) { /* Check that nodes is an odd integer */
double fn = rint (Ctrl->S.value[3]);
int64_t n = lrint (fn);
if (!doubleAlmostEqual (Ctrl->S.value[3], fn) || ((n%2) == 0)) {
GMT_Report (API, GMT_MSG_ERROR, "Option -Sq option +n<N> modifier: <N> must be an odd integer\n");
n_errors++;
}
if (Ctrl->S.value[2] < 0.0 || Ctrl->S.value[2] > 1.0e-4) {
GMT_Report (API, GMT_MSG_ERROR, "Option -Sq option +e<err> modifier: <err> must be positive and < 1.0e-4\n");
n_errors++;
}
}
if (Ctrl->S.mode == PARKER_1994 || Ctrl->S.mode == WESSEL_BECKER_2008) Ctrl->D.mode = 4; /* Automatically set */
dimension = (Ctrl->D.mode == 0) ? 1 : ((Ctrl->D.mode == 5) ? 3 : 2);
if (dimension == 2 && Ctrl->R3.mode) { /* Set -R via a gridfile */
/* Here, -R<grdfile> was used and we will use the settings supplied by the grid file (unless overridden) */
if (!Ctrl->I.active) { /* -I was not set separately; set indirectly */
Ctrl->I.inc[GMT_X] = Ctrl->R3.inc[GMT_X];
Ctrl->I.inc[GMT_Y] = Ctrl->R3.inc[GMT_Y];
Ctrl->I.active = true;
}
/* Here, -r means toggle the grids registration */
if (GMT->common.R.active[GSET]) {
GMT->common.R.active[GSET] = !Ctrl->R3.offset;
GMT->common.R.registration = !Ctrl->R3.offset;
}
else {
GMT->common.R.active[GSET] = Ctrl->R3.offset;
GMT->common.R.registration = Ctrl->R3.offset;
}
}
n_errors += gmt_M_check_condition (GMT, Ctrl->A.active && gmt_access (GMT, Ctrl->A.file, R_OK), "Option -A: Cannot read file %s!\n", Ctrl->A.file);
n_errors += gmt_M_check_condition (GMT, Ctrl->A.active && Ctrl->A.mode > 5, "Option -A: format must be in 0-5 range\n");
n_errors += gmt_M_check_condition (GMT, !(GMT->common.R.active[RSET] || Ctrl->N.active || Ctrl->T.active), "No output locations specified (use either [-R -I], -N, or -T)\n");
n_errors += gmt_M_check_condition (GMT, Ctrl->R3.mode && dimension != 2, "The -R<gridfile> or -T<gridfile> option only applies to 2-D gridding\n");
n_errors += gmt_M_check_condition (GMT, Ctrl->C.movie && dimension != 2, "The -C +m|M modifier only applies to 2-D gridding\n");
#ifdef DEBUG
n_errors += gmt_M_check_condition (GMT, !TEST && !Ctrl->N.active && Ctrl->R3.dimension != dimension, "The -R and -D options disagree on the dimension\n");
#else
n_errors += gmt_M_check_condition (GMT, !Ctrl->N.active && Ctrl->R3.dimension != dimension, "The -R and -D options disagree on the dimension\n");
#endif
n_errors += gmt_check_binary_io (GMT, dimension + 1);
n_errors += gmt_M_check_condition (GMT, Ctrl->S.value[0] < 0.0 || Ctrl->S.value[0] >= 1.0, "Option -S: Tension must be in range 0 <= t < 1\n");
n_errors += gmt_M_check_condition (GMT, !(Ctrl->S.mode == PARKER_1994 || Ctrl->S.mode == WESSEL_BECKER_2008) && Ctrl->D.mode == 3, "Option -Sc|t|r: Cannot select -D3\n");
n_errors += gmt_M_check_condition (GMT, Ctrl->S.mode == LINEAR_1D && Ctrl->D.mode > 3, "Option -Sl: Cannot select -D4 or higher\n");
n_errors += gmt_M_check_condition (GMT, Ctrl->I.active && (Ctrl->I.inc[GMT_X] <= 0.0 || (dimension > 1 && Ctrl->I.inc[GMT_Y] <= 0.0) || (dimension == 3 && Ctrl->I.inc[GMT_Z] <= 0.0)), "Option -I: Must specify positive increment(s)\n");
n_errors += gmt_M_check_condition (GMT, dimension == 2 && !Ctrl->N.active && !(Ctrl->G.active || Ctrl->G.file), "Option -G: Must specify output grid file name\n");
n_errors += gmt_M_check_condition (GMT, Ctrl->C.active && Ctrl->C.value < 0.0 && !Ctrl->C.file, "Option -C: Must specify file name for eigenvalues if cut < 0\n");
n_errors += gmt_M_check_condition (GMT, Ctrl->T.active && !Ctrl->T.file, "Option -T: Must specify mask grid file name\n");
n_errors += gmt_M_check_condition (GMT, Ctrl->T.active && dimension != 2, "Option -T: Only applies to 2-D gridding\n");
n_errors += gmt_M_check_condition (GMT, Ctrl->N.active && !Ctrl->N.file, "Option -N: Must specify node file name\n");
n_errors += gmt_M_check_condition (GMT, Ctrl->N.active && Ctrl->N.file && gmt_access (GMT, Ctrl->N.file, R_OK), "Option -N: Cannot read file %s!\n", Ctrl->N.file);
n_errors += gmt_M_check_condition (GMT, (Ctrl->I.active + GMT->common.R.active[RSET]) == 1 && dimension == 2, "Must specify -R, -I, [-r], -G for gridding\n");
n_errors += gmt_M_check_condition (GMT, Ctrl->C.movie && strchr (Ctrl->G.file, '%') == NULL, "Must give a filename template via -G when -C..+m|M is selected\n");
n_errors += gmt_M_check_condition (GMT, Ctrl->C.movie && dimension != 2, "The -C..+m|M modifier is only available for 2-D gridding\n");
return (n_errors ? GMT_PARSE_ERROR : GMT_NOERROR);
}
#ifdef DEBUG
/* Dump a table of x, G, dGdx for test purposes [requires option -+ and compilation with -DDEBUG] */
GMT_LOCAL void greenspline_dump_green (struct GMT_CTRL *GMT, double (*G) (struct GMT_CTRL *, double, double *, struct GREENSPLINE_LOOKUP *), double (*D) (struct GMT_CTRL *, double, double *, struct GREENSPLINE_LOOKUP *), double par[], double x0, double x1, int N, struct GREENSPLINE_LOOKUP *Lz, struct GREENSPLINE_LOOKUP *Lg) {
int i;
double x, dx, dy, y, t, ry, rdy;
double min_y, max_y, min_dy, max_dy;
min_y = min_dy = DBL_MAX;
max_y = max_dy = -DBL_MAX;
dx = (x1 - x0) / (N - 1);
for (i = 0; i < N; i++) {
x = x0 + i * dx;
y = G (GMT, x, par, Lz);
dy = D (GMT, x, par, Lg);
if (y < min_y) min_y = y;
if (y > max_y) max_y = y;
if (dy < min_dy) min_dy = dy;
if (dy > max_dy) max_dy = dy;
}
ry = max_y - min_y;
rdy = max_dy - min_dy;
for (i = 0; i < N; i++) {
x = x0 + i * dx;
t = (x0 < 0.0) ? acosd (x) : x;
y = G (GMT, x, par, Lz);
dy = D (GMT, x, par, Lg);
dy = (rdy > 0.0) ? (dy - min_dy)/rdy : 1.0;
printf ("%g\t%g\t%g\t%g\n", x, (y - min_y) / ry, dy, t);
}
}
#endif
/* Below are all the individual Green's functions. Note that most of them take an argument
* that is unused except in the spline lookup version of WB_08. */
/*---------------------- ONE DIMENSION ---------------------- */
/* Basic linear spline (bilinear in 2-D) */
GMT_LOCAL double greenspline_spline1d_linear (struct GMT_CTRL *GMT, double r, double par[], struct GREENSPLINE_LOOKUP *unused) {
/* Dumb linear spline */
gmt_M_unused(GMT); gmt_M_unused(par); gmt_M_unused(unused);
return (r); /* Just regular spline; par not used */
}
GMT_LOCAL double greenspline_grad_spline1d_linear (struct GMT_CTRL *GMT, double r, double par[], struct GREENSPLINE_LOOKUP *unused) {
/* d/dr of r is 1 */
gmt_M_unused(GMT); gmt_M_unused(r); gmt_M_unused(par); gmt_M_unused(unused);
return (1.0);
}
/* greenspline_spline1d_sandwell computes the Green function for a 1-d spline
* as per Sandwell [1987], G(r) = r^3. All r must be >= 0.
*/
GMT_LOCAL double greenspline_spline1d_sandwell (struct GMT_CTRL *GMT, double r, double par[], struct GREENSPLINE_LOOKUP *unused) {
gmt_M_unused(GMT); gmt_M_unused(par); gmt_M_unused(unused);
if (r == 0.0) return (0.0);
return (pow (r, 3.0)); /* Just regular spline; par not used */
}
GMT_LOCAL double greenspline_grad_spline1d_sandwell (struct GMT_CTRL *GMT, double r, double par[], struct GREENSPLINE_LOOKUP *unused) {
gmt_M_unused(GMT); gmt_M_unused(par); gmt_M_unused(unused);
return (r); /* Just regular spline; par not used */
}
GMT_LOCAL double greenspline_spline1d_Wessel_Bercovici (struct GMT_CTRL *GMT, double r, double par[], struct GREENSPLINE_LOOKUP *unused) {
/* greenspline_spline1d_Wessel_Bercovici computes the Green function for a 1-d spline
* in tension as per Wessel and Bercovici [1988], G(u) = exp(-u) + u - 1,
* where u = par[0] * r and par[0] = sqrt (t/(1-t)).
* All r must be >= 0. par[0] = c
*/
double cx;
gmt_M_unused(GMT); gmt_M_unused(unused);
if (r == 0.0) return (0.0);
cx = par[0] * r;
return (exp (-cx) + cx - 1.0);
}
GMT_LOCAL double greenspline_grad_spline1d_Wessel_Bercovici (struct GMT_CTRL *GMT, double r, double par[], struct GREENSPLINE_LOOKUP *unused) {
double cx;
gmt_M_unused(GMT); gmt_M_unused(unused);
if (r == 0.0) return (0.0);
cx = par[0] * r;
return (1.0 - exp (-cx));
}
/*---------------------- TWO DIMENSIONS ---------------------- */
GMT_LOCAL double greenspline_spline2d_sandwell (struct GMT_CTRL *GMT, double r, double par[], struct GREENSPLINE_LOOKUP *unused) {
gmt_M_unused(GMT); gmt_M_unused(par); gmt_M_unused(unused);
if (r == 0.0) return (0.0);
return (r * r * (log (r) - 1.0)); /* Just regular spline; par not used */
}
GMT_LOCAL double greenspline_grad_spline2d_sandwell (struct GMT_CTRL *GMT, double r, double par[], struct GREENSPLINE_LOOKUP *unused) {
gmt_M_unused(GMT); gmt_M_unused(par); gmt_M_unused(unused);
if (r == 0.0) return (0.0);
return (r * (2.0 * log (r) - 1.0)); /* Just regular spline; par not used */
}
/* greenspline_spline2d_Wessel_Bercovici computes the Green function for a 2-d spline
* in tension as per Wessel and Bercovici [1988], G(u) = K(u) - log(u),
* where u = par[0] * r and par[0] = sqrt (t/(1-t)).
* K is the modified Bessel function K of order zero.
* All r must be >= 0.
* par[0] = c
* par[1] = 2/c
*/
GMT_LOCAL double greenspline_spline2d_Wessel_Bercovici (struct GMT_CTRL *GMT, double r, double par[], struct GREENSPLINE_LOOKUP *unused) {
double y, z, cx, t, g;
gmt_M_unused(GMT); gmt_M_unused(unused);
if (r == 0.0) return (0.0);
cx = par[0] * r;
if (r <= par[1]) {
y = 0.25 * (t = cx * cx);
z = t / 14.0625;
g = (-log(0.5*cx) * (z * (3.5156229 + z * (3.0899424 + z * (1.2067492 + z * (0.2659732 +
z * (0.360768e-1 + z * 0.45813e-2))))))) + (y * (0.42278420 + y * (0.23069756 +
y * (0.3488590e-1 + y * (0.262698e-2 + y * (0.10750e-3 + y * 0.74e-5))))));
}
else {
y = par[1] / r;
g = (exp (-cx) / sqrt (cx)) * (1.25331414 + y * (-0.7832358e-1 + y * (0.2189568e-1 +
y * (-0.1062446e-1 + y * (0.587872e-2 + y * (-0.251540e-2 + y * 0.53208e-3))))))
+ log (cx) - M_LOG_2 + M_GAMMA;
}
return (g);
}
GMT_LOCAL double greenspline_grad_spline2d_Wessel_Bercovici (struct GMT_CTRL *GMT, double r, double par[], struct GREENSPLINE_LOOKUP *unused) {
double y, z, cx, t, dgdr;
gmt_M_unused(GMT); gmt_M_unused(unused);
if (r == 0.0) return (0.0);
cx = par[0] * r;
if (r <= par[1]) {
y = 0.25 * (t = cx * cx);
z = t / 14.0625;
dgdr = -((log(0.5*cx) * (cx * (0.5 + z * (0.87890594 + z * (0.51498869 + z * (0.15084934 +
z * (0.2658733e-1 + z * (0.301532e-2 + z * 0.32411e-3)))))))) + (1.0/cx) * (y * (0.15443144 +
y * (-0.67278579 + y * (-0.18156897 + y * (-0.1919402e-1 + y * (-0.110404e-2 + y * (-0.4686e-4))))))));
}
else {
y = par[1] / r;
dgdr = 0.5*y - ((exp (-cx) / sqrt (cx)) * (1.25331414 + y * (0.23498619 + y * (-0.3655620e-1 +
y * (0.1504268e-1 + y * (-0.780353e-2 + y * (0.325614e-2 + y * (-0.68245e-3))))))));
}
return (dgdr*par[0]);
}
/* greenspline_spline2d_Mitasova_Mitas computes the regularized Green function for a 2-d spline
* in tension as per Mitasova and Mitas [1993], G(u) = log (u) + Ei(u),
* where u = par[1] * r^2. All r must be >= 0.
* par[0] = phi (par[0] unused by function)
* par[1] = phi^2/4
*/
GMT_LOCAL double greenspline_spline2d_Mitasova_Mitas (struct GMT_CTRL *GMT, double r, double par[], struct GREENSPLINE_LOOKUP *unused) {
double x, z, g, En, Ed;
gmt_M_unused(GMT); gmt_M_unused(unused);
if (r == 0.0) return (0.0);
x = z = par[1] * r * r;
if (z <= 1.0) {
g = 0.99999193 * x;
x *= x;
g -= 0.24991055 * x;
x *= x;
g += 0.05519968 * x;
x *= x;
g -= 0.00976004 * x;
x *= x;
g += 0.00107857 * x;
}
else {
g = log (x) + M_GAMMA;
En = 0.2677737343 + 8.6347608925 * x;
Ed = 3.9584869228 + 21.0996530827 * x;
x *= x;
En += 18.0590169730 * x;
Ed += 25.6329561486 * x;
x *= x;
En += 8.5733287401 * x;
Ed += 9.5733223454 * x;
x *= x;
En += x;
Ed += x;
g += (En / Ed) / (z * exp(z));
}
return (g);
}
GMT_LOCAL double greenspline_grad_spline2d_Mitasova_Mitas (struct GMT_CTRL *GMT, double r, double par[], struct GREENSPLINE_LOOKUP *unused) {
double u, dgdr;
gmt_M_unused(GMT); gmt_M_unused(unused);
if (r == 0.0) return (0.0);
u = par[1] * r * r;
dgdr = 2.0 * (1.0 - exp (-u))/r;
return (dgdr);
}
/*---------------------- TWO DIMENSIONS (SPHERE) ---------------------- */
/* greenspline_spline2d_Parker computes the Green function for a 2-d surface spline
* as per Parker [1994], G(x) = dilog(),
* where x is cosine of distances. All x must be -1 <= x <= +1.
* Parameters passed are:
* par[0] = 6/M_PI^2 (to normalize results)
*/
GMT_LOCAL double greenspline_spline2d_Parker (struct GMT_CTRL *GMT, double x, double par[], struct GREENSPLINE_LOOKUP *unused) {
/* Normalized to 0-1 */
gmt_M_unused(unused);
if (x == +1.0) return (1.0);
if (x == -1.0) return (0.0);
return (gmt_dilog (GMT, 0.5 - 0.5 * x) * par[0]);
}
GMT_LOCAL double greenspline_grad_spline2d_Parker (struct GMT_CTRL *GMT, double x, double par[], struct GREENSPLINE_LOOKUP *unused) {
/* Normalized to 0-1 */
gmt_M_unused(GMT); gmt_M_unused(par); gmt_M_unused(unused);
if (x == +1.0 || x == -1.0) return (0.0);
return (log (0.5 - 0.5 * x) * sqrt ((1.0 - x) / (1.0 + x)));
}
GMT_LOCAL void greenspline_series_prepare (struct GMT_CTRL *GMT, double p, unsigned int L, struct GREENSPLINE_LOOKUP *Lz, struct GREENSPLINE_LOOKUP *Lg) {
/* Precalculate Legendre series terms involving various powers/ratios of l and p */
unsigned int l;
double pp, t1, t2;
GMT_Report (GMT->parent, GMT_MSG_INFORMATION, "Precalculate max %u terms for Legendre summation\n", L+1);
Lz->A = gmt_M_memory (GMT, NULL, L+1, double);
Lz->B = gmt_M_memory (GMT, NULL, L+1, double);
Lz->C = gmt_M_memory (GMT, NULL, L+1, double);
if (Lg) Lg->A = gmt_M_memory (GMT, NULL, L+1, double);
pp = p * p;
for (l = 1; l <= L; l++) {
t1 = l + 1.0;
t2 = 2.0 * l + 1.0;
Lz->A[l] = t2 * pp / (l*t1*(l*t1+pp));
Lz->B[l] = t2 / t1;
Lz->C[l] = l / t1;
if (Lg) Lg->A[l] = t2 * pp / (t1*(l*t1+pp));
}
if (Lg) Lg->B = Lz->B; /* Lg needs the same B,C coefficients as Lz, so just share them via the link */
if (Lg) Lg->C = Lz->C;
}
GMT_LOCAL unsigned int greenspline_get_max_L (struct GMT_CTRL *GMT, double p, double err) {
/* Return max L needed given p and truncation err for Parker's simplified loop expression */
gmt_M_unused(GMT);
return ((unsigned int)lrint (p / sqrt (err) + 10.0));
}
GMT_LOCAL void greenspline_free_lookup (struct GMT_CTRL *GMT, struct GREENSPLINE_LOOKUP **Lptr, unsigned int mode) {
/* Free all items allocated under the lookup structures.
* mode = 0 means Lz and mode = 1 means Lg; the latter has no B & C arrays */
struct GREENSPLINE_LOOKUP *L = *Lptr;
if (L == NULL) return; /* Nothing to free */
gmt_M_free (GMT, L->y);
gmt_M_free (GMT, L->c);
gmt_M_free (GMT, L->A);
if (mode == 0) { /* Only Lz has A,B,C while Lg borrows B,C */
gmt_M_free (GMT, L->B);
gmt_M_free (GMT, L->C);
}
gmt_M_free (GMT, L);
*Lptr = NULL;
}
GMT_LOCAL unsigned int greenspline_get_L (double x, double p, double err) {
/* Determines the truncation order L_max given x, p, and err.
* See ParkerNotesJan2013.pdf in gurudocs for explanations and derivation */
unsigned int L_max;
double s, c, pp, Lf, gam, lam;
/* Get all the trig functions needed; x is cos(theta), but we need */
/* s = sin (theta/2); c=cos(theta/2); */
/* cos(theta) = x = 1.0 - 2.0 *s^2, so */
s = sqrt (0.5 * (1 - x));
c = sqrt (0.5 * (1 + x));
pp = p * p;
/* Approximate the highest order L_max we need to sum given specified err limit */
if (s == 0.0)
Lf = p / sqrt (err);
else if (c == 0.0)