-
Notifications
You must be signed in to change notification settings - Fork 2
/
shana.c
2153 lines (2001 loc) · 67 KB
/
shana.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
/* part of the shansyn spherical harmonics package,
see COPYRIGHT for license
$Id: shana.c,v 1.30 2009/04/16 00:55:29 becker Exp becker $
*/
#include "shana.h"
int finite(double ); /* why doesn't this get included with math.h? */
#define myfinite(x) (finite((double)(x)))
void spline(COMP_PRECISION *,COMP_PRECISION *,
int ,COMP_PRECISION ,COMP_PRECISION ,
COMP_PRECISION *);
/*
shana
calculates a real spherical harmonic expansion using theoretical
physics-type fully normalized spherical harmonics as specified in
Dahlen and Tromp (Theoretical Global Seismology, 1998, Appendix B)
uses
- binary or ASCII block gridded data on the surface of the Earth
(can be in GMT/netcdf format) and solves by straightforward
integration
- irregularly distributed data and performs least squares fit
- contour lines along which Delta functions are integrated, or
- individual points whose Delta function contributions are added up
type shana -h for information on the parameter settings
TODO:
clean up coding
implement inverse FFT in longitude
(c) Thorsten Becker 1999 - 2009 twb@usc.edu
*/
// check for floating exceptions in velocity routines
//
//#define CHECK_FLOATING_EXCEPTION
int main(int argc, char *argv[] )
{
int lmax,lmc,lmax1,lmsize,nlon,nlat,i,j,opmode,nrp=0,intmode,verbose=BE_VERBOSE,flip=0,pixelreg,
calculate_derivatives,vectors,nlontimesnlat,os1,os2,rc;
DATA_PRECISION *func=NULL,avg,min,max,*cloc,minphi,maxphi,
mintheta,maxtheta,dtheta=0.,dphi=0.,*tmp;
COMP_PRECISION damping,theta,phi;
double tmpd;
char *filename;
#ifdef USE_GMT4
int ret_code;
struct GRD_HEADER *header;
GMT_LONG dummy[4]={0,0,0,0};
header=(struct GRD_HEADER *)calloc(1,sizeof(struct GRD_HEADER));
double wesn[6];
#else /* new GMT API */
void *API; /* The API control structure */
struct GMT_GRID *G = NULL; /* Structure to hold output grid */
struct GMT_GRID_HEADER *header = NULL;
#endif
filename=(char *)malloc(sizeof(char)*STRING_LENGTH*2);
#ifdef USE_GMT4
/* */
GMT_begin (argc, argv);
GMT_grd_init (header, argc, argv, FALSE);
GMT_program = argv[0];
#else
API = GMT_Create_Session (argv[0], 2U, 0, NULL);
#endif
if(argc > 1 && strcmp(argv[1],"-h") == 0 )argc=-9;
calculate_derivatives=vectors=FALSE;
lmax=DEF_LMAX;
opmode=DEF_OPMODE;
intmode=DEF_INT_MODE;
damping = DEF_DAMPING;
/*
argument assignment , check some of the flags
*/
switch(argc){
case 1:{
break;
}
case 2:{
sscanf(argv[1],"%i",&lmax);
break;}
case 3:{
sscanf(argv[1],"%i",&lmax);
check_opmode(argv[2],&opmode,filename,
&vectors,&calculate_derivatives);
break;
}
case 4:{
sscanf(argv[1],"%i",&lmax);
check_opmode(argv[2],&opmode,filename,
&vectors,&calculate_derivatives);
sscanf(argv[3],"%i",&intmode);
break;
}
case 5:{
sscanf(argv[1],"%i",&lmax);
check_opmode(argv[2],&opmode,filename,
&vectors,&calculate_derivatives);
sscanf(argv[3],"%i",&intmode);
sscanf(argv[4],"%lf",&tmpd);damping = (COMP_PRECISION) tmpd;
break;
}
default:{
phelp(argv[0]);
exit(-1);}}
lmsize= (int)((((COMP_PRECISION)lmax)+1.0)*(((COMP_PRECISION)lmax)+2)/2.0);
lmax1 = lmax+1;
lmc = lmsize*2 - lmax1;
if(lmax >= 0){
if(verbose)fprintf(stderr,"%s: expansion to l_max=%i, DOF %i-1\n",argv[0],lmax,lmc);
}else {
fprintf(stderr,"%s: expansion to l_max=%i doesn't make sense\n",argv[0],lmax);
phelp(argv[0]);exit(-1);
}
if(opmode < 0){
flip=1;
opmode*= -1;
if(verbose)fprintf(stderr,"%s: latitudes are flipped\n",argv[0]);
}
if((lmax < 30) && (intmode == GAUSSIAN))
fprintf(stderr,"%s: WARNING: Gaussian integration for low order expansion\n",argv[0]);
/*
DATA INPUT
*/
// flags and default settings for input
pixelreg=FALSE;
switch(opmode){
/*
read in binary block data of a global scalar function
*/
case BLOCK:{
if(verbose)fprintf(stderr,"%s: expecting binary block on offset (grid) coordinates\n",argv[0]);
rc =fread(&nlat,sizeof(int),1,stdin);
rc+=fread(&nlon,sizeof(int),1,stdin);
if(verbose)fprintf(stderr,"%s: read dimensions expecting nlon:%i times nlat:%i values\n",
argv[0],nlon,nlat);
myrealloc_dp(&func,(nlontimesnlat=nlon*nlat));zero_dp(func,nlontimesnlat);
rc+=fread(func,sizeof(DATA_PRECISION)*nlontimesnlat,1,stdin);
minmax(&min,&max,&avg,func,nlontimesnlat);
//
// this are the expected input data coordinates
// dx/2 <= x <= 360-dx/2 and -90+dy/2 <= y <= 90-dy/2
//
dphi= TWOPI/((COMP_PRECISION)nlon);
minphi= 0.5*dphi;
maxphi= minphi+(COMP_PRECISION)(nlon-1)*dphi;
//
dtheta= PI/((COMP_PRECISION)nlat);
mintheta= PI - 0.5*dtheta;
maxtheta= mintheta-(COMP_PRECISION)(nlat-1)*dtheta;
if(verbose)
field_message(minphi,dphi,maxphi,mintheta,dtheta,maxtheta,
nlon,nlat,min,max,avg,argv[0]);
break;
}
/*
lon lat value data
CONTOUR:
read in ASCII contours for integration of Delta functions along
contours
POINTS_FOR_AB:
POINTS_FOR_A_MATRIX:
read in ASCII lon lat coordinates for creation of an A matrix for
least squares solution of fitting problem OR for spherical
harmonic analysis of spotted delta function type points that have
a scalar value
*/
case POINTS_FOR_AB:
case POINTS_FOR_A_MATRIX:
case POINTS_FOR_A_MATRIX_ASCII:
case POINTS_FOR_AVEC_MATRIX:
case POINTS_FOR_AVEC_MATRIX_ASCII:
case CONTOUR:{
if(verbose){
if(intmode == LEAST_SQUARES)
fprintf(stderr,"%s: reading lon lat z for least squares fit of spherical harmonics\n",
argv[0]);
else{
if(opmode == CONTOUR)
fprintf(stderr,"%s: reading contour lon lat z tripels along which to integrate\n",argv[0]);
else if(opmode == POINTS_FOR_A_MATRIX)
fprintf(stderr,"%s: expecting lon lat z for A matrix (binary, single precision) output (z unused)\n",
argv[0]);
else if(opmode == POINTS_FOR_AVEC_MATRIX)
fprintf(stderr,"%s: expecting lon lat z for A matrix vector field (binary, single precision) output (z unused)\n",
argv[0]);
else if(opmode == POINTS_FOR_A_MATRIX_ASCII)
fprintf(stderr,"%s: expecting lon lat z for A matrix (ascii) output (z unused)\n",
argv[0]);
else if(opmode == POINTS_FOR_AVEC_MATRIX_ASCII)
fprintf(stderr,"%s: expecting lon lat z for A matrix vector field (ascii) output (z unused)\n",
argv[0]);
else if(opmode == POINTS_FOR_AB)
fprintf(stderr,"%s: expecting lon lat z for sum over Delta points\n",argv[0]);
}
}
nrp=1;
mymalloc_dp(&cloc,2*nrp);
myrealloc_dp(&func,nrp);
os1 = 0;
while(fscanf(stdin,THREE_GMTDATA_FSCAN_FORMAT,
(cloc+os1),(cloc+os1+1),(func+(nrp-1))) == 3){
/* reformat coordinates */
*(cloc+os1) =LONGITUDE2PHI( *(cloc+os1));
*(cloc+os1+1)=LATITUDE2THETA(*(cloc+os1+1));
nrp++;
os1 += 2;
myrealloc_dp(&cloc,2*nrp);myrealloc_dp(&func,nrp);
}
nrp--;
if(verbose)fprintf(stderr,"%s: read %i points\n",argv[0],nrp);
if(!nrp){
fprintf(stderr,"%s: Exiting, no data.\n",argv[0]);exit(-1);}
minmax(&min,&max,&avg,func,nrp);
if((opmode == CONTOUR)&&(intmode != LEAST_SQUARES)) /* contour input */
if(nrp!=1)
for(i=0;i<nrp;i++)
if(dist_rad(cloc,i, (i<nrp-1)?(i+1):(i-1)) > 0.5*PI/(COMP_PRECISION)lmax){
fprintf(stderr,"%s: pts. %i and %i have low spacing of %g degrees\n",
argv[0],i,i+1,dist_rad(cloc,i, (i<nrp-1)?(i+1):(i-1))/PIOVERONEEIGHTY);
fprintf(stderr,"%s: maximum should be around %g for lmax=%i\n",
argv[0],90.0/(COMP_PRECISION)lmax,lmax);
}
if(verbose)fprintf(stderr,"%s: min=%10g avg=%10g max=%10g\n",
argv[0],min,avg/((COMP_PRECISION)(nrp)),max);
break;
}
/*
read in ASCII block data of a global scalar function
*/
case ASCII_BLOCK:{
if(verbose)
fprintf(stderr,"%s: reading ascii block on grid registration coordinates\n",argv[0]);
nlon=1;
myrealloc_dp(&func,nlon);
while((fscanf(stdin,"%f ",(func+nlon-1)))==1){
nlon++;
myrealloc_dp(&func,nlon);
}
nlon--;
nlat=(int)sqrt((COMP_PRECISION)nlon/2.);
nlon=2*nlat;
// swap ordering
mymalloc_dp(&tmp,nlon*nlat);
for(os1=i=0;i<nlat;i++,os1+=nlon)
for(os2=j=0;j<nlon;j++,os2+=nlat)
*(tmp+os2+i)= *(func+os1+j);
j=nlon*nlat;
for(i=0;i<j;i++)
func[i]= tmp[i];
free(tmp);
// dx/2 <= x <= 360-dx/2 and -90+dy/2 <= y <= 90-dy/2
if(verbose)
fprintf(stderr,"%s: nlon=%i, nlat=%i\n",argv[0],nlon,nlat);
minmax(&min,&max,&avg,func,nlat*nlon);
dphi= TWOPI/((COMP_PRECISION)nlon);
minphi=0.5*dphi;
maxphi=minphi+(COMP_PRECISION)(nlon-1)*dphi;
dtheta= PI/((COMP_PRECISION)nlat);
mintheta=PI - 0.5*dtheta;
maxtheta=mintheta - (COMP_PRECISION)(nlat-1)*dtheta;
if(verbose)
field_message(minphi,dphi,maxphi,mintheta,dtheta,maxtheta,nlon,nlat,min,max,avg,argv[0]);
break;
}
/*
read in ASCII block data of a global scalar function with header
*/
case ASCII_BLOCK_HEADER:{
if(verbose)
fprintf(stderr,"%s: reading ascii block with header\n",argv[0]);
if(fscanf(stdin,"%i %i %f %f %f %f\n",&nlon,&nlat,&minphi,&maxphi,&mintheta,&maxtheta)!=6){
fprintf(stderr,"%s: need header line with nlon nlat minphi maxphi mintheta maxtheta\n\n",argv[0]);
exit(-1);
}
minphi= LONGITUDE2PHI(minphi);
maxphi= LONGITUDE2PHI(maxphi);
mintheta=LATITUDE2THETA(mintheta);
maxtheta=LATITUDE2THETA(maxtheta);
dphi= (maxphi-minphi) /((COMP_PRECISION)(nlon-1));
dtheta=fabs((maxtheta-mintheta)/((COMP_PRECISION)(nlat-1)));
myrealloc_dp(&func,(nlontimesnlat=nlon*nlat));
/* read in data */
for(i=0;i<nlat;i++)
for(os1=j=0;j<nlon;j++,os1+=nlat)
if(fscanf(stdin,"%f ",(func+os1+i))!=1){
fprintf(stderr,"%s: read error at %i %i\n",argv[0],i,j);
exit(-1);
}
if(flip)/* flip the latitudes */
gmt2myconvention_rotate4(func,nlon,nlat,1);
minmax(&min,&max,&avg,func,nlontimesnlat);
if(verbose)
field_message(minphi,dphi,maxphi,mintheta,dtheta,maxtheta,nlon,nlat,min,max,avg,argv[0]);
break;
}
/*
read in GMT grd file of scalar data
*/
case GMT_BLOCK:{
#ifdef USE_GMT4
if((ret_code = GMT_read_grd_info (filename,header))< 0){
/* return code has changed for GMT4, i presume? */
fprintf(stderr,"%s: error opening %s for grdinfo\n",argv[0],filename);
exit(-1);
}else{
pixelreg=(header->node_offset?TRUE:FALSE);
minphi=LONGITUDE2PHI(header->x_min+(pixelreg?header->x_inc/2.0:0.0));
maxphi=LONGITUDE2PHI(header->x_max-(pixelreg?header->x_inc/2.0:0.0));
mintheta=LATITUDE2THETA(header->y_min+(pixelreg?header->y_inc/2.0:0.0));
maxtheta=LATITUDE2THETA(header->y_max-(pixelreg?header->y_inc/2.0:0.0));
dphi= DEG2RAD( header->x_inc);
dtheta=DEG2RAD( header->y_inc);
// for grid line registered nx_i = (x_i^max-x_i^min)/dx_i + 1
// for pixel reg nx_i = (x_i^max-x_i^min)/dx_i
nlon=header->nx;
nlat=header->ny;
}
#else
/* read header info */
if((G = GMT_Read_Data (API, GMT_IS_GRID, GMT_IS_FILE,
GMT_IS_SURFACE, GMT_CONTAINER_ONLY, NULL, filename, NULL))==NULL)
return (-1);
header = G->header;
pixelreg=(header->registration == GMT_GRID_PIXEL_REG)?TRUE:FALSE;
minphi=LONGITUDE2PHI(header->wesn[XLO]+(pixelreg?header->inc[GMT_X]/2.0:0.0));
maxphi=LONGITUDE2PHI(header->wesn[XHI]-(pixelreg?header->inc[GMT_X]/2.0:0.0));
mintheta=LATITUDE2THETA(header->wesn[YLO]+(pixelreg?header->inc[GMT_Y]/2.0:0.0));
maxtheta=LATITUDE2THETA(header->wesn[YHI]-(pixelreg?header->inc[GMT_Y]/2.0:0.0));
dphi= DEG2RAD( header->inc[GMT_X]);
dtheta=DEG2RAD( header->inc[GMT_Y]);
// for grid line registered nx_i = (x_i^max-x_i^min)/dx_i + 1
// for pixel reg nx_i = (x_i^max-x_i^min)/dx_i
nlon=header->n_columns;
nlat=header->n_rows;
#endif
if(!vectors){
/*
read in scalar data field from GMT grd file
*/
fprintf(stderr,"%s: reading from grd-file %s\n",argv[0],filename);
myrealloc_dp(&func,nlon*nlat);
#ifdef USE_GMT4
GMT_read_grd(filename,header, func, 0.0, 0.0, 0.0, 0.0, dummy, 0);
gmt2myconvention_rotate4(func,nlon,nlat,1.0);
#else
/* read data */
if((G = GMT_Read_Data (API, GMT_IS_GRID, GMT_IS_FILE,
GMT_IS_SURFACE, GMT_DATA_ONLY, NULL, filename, G))==NULL)
return (-1);
gmt2myconvention_rotate(func,nlon,nlat,1,G);
#endif
minmax(&min,&max,&avg,func,(nlontimesnlat=nlat*nlon));
if(verbose)
field_message(minphi,dphi,maxphi,mintheta,dtheta,maxtheta,nlon,nlat,min,max,avg,argv[0]);
if((minphi == maxphi) || (mintheta == maxtheta)){
fprintf(stderr,"%s: use arrays with non-overlapping borders: %g <= x <= %g, %g <= y <= %g\n",
argv[0],PHI2LONGITUDE(minphi),PHI2LONGITUDE(maxphi),
THETA2LATITUDE(mintheta),THETA2LATITUDE(maxtheta));
exit(-1);
}
}else{
/*
read in vectors from two grd files
use dimensions specs from first one and check if second has identical
bounds
*/
if(verbose)
fprintf(stderr,"%s: reading phi component from grd-file %s\n",argv[0],filename);
myrealloc_dp(&func,(nlontimesnlat=nlon*nlat)*2);
/* read from grid */
#ifdef USE_GMT4
GMT_read_grd (filename,header, func, 0.0, 0.0, 0.0, 0.0, dummy, 0);
gmt2myconvention_rotate4(func,nlon,nlat,1.0);
#else /* read datA */
if((G = GMT_Read_Data (API, GMT_IS_GRID, GMT_IS_FILE,
GMT_IS_SURFACE, GMT_DATA_ONLY, NULL, filename, G))==NULL)
return (-1);
gmt2myconvention_rotate(func,nlon,nlat,1.0,G);
#endif
minmax(&min,&max,&avg,func,nlontimesnlat);
if(verbose)
field_message(minphi,dphi,maxphi,mintheta,dtheta,maxtheta,nlon,nlat,min,max,avg,argv[0]);
/* second velocity file */
if(verbose)
fprintf(stderr,"%s: and theta component from grd-file %s\n",argv[0],(filename+STRING_LENGTH));
#ifdef USE_GMT4
if(GMT_read_grd_info ((filename+STRING_LENGTH),header)== -1){
fprintf(stderr,"%s: error opening %s for grdinfo\n",argv[0],filename);
exit(-1);
}
if(SIGNIFICANTLY_DIFFERENT(header->x_min, PHI2LONGITUDE(minphi)) ||
SIGNIFICANTLY_DIFFERENT(header->y_min, THETA2LATITUDE(mintheta)) ||
SIGNIFICANTLY_DIFFERENT(header->x_max, PHI2LONGITUDE(maxphi)) ||
SIGNIFICANTLY_DIFFERENT(header->y_max, THETA2LATITUDE(maxtheta)) ||
SIGNIFICANTLY_DIFFERENT(header->nx, nlon) ||
SIGNIFICANTLY_DIFFERENT(header->ny, nlat)){
fprintf(stderr,"%s: %s and %s have to have identical size\n",
argv[0],(filename),(filename+STRING_LENGTH));
fprintf(stderr,"%s: delta xmin: %g\n",
argv[0],header->x_min-PHI2LONGITUDE(minphi));
fprintf(stderr,"%s: delta ymin: %g\n",
argv[0],header->y_min-THETA2LATITUDE(mintheta));
fprintf(stderr,"%s: delta xmax: %g\n",
argv[0],header->x_max-PHI2LONGITUDE(maxphi));
fprintf(stderr,"%s: delta ymax: %g\n",
argv[0],header->y_max-THETA2LATITUDE(maxtheta));
fprintf(stderr,"%s: delta nx: %i\n",
argv[0],header->nx -nlon );
fprintf(stderr,"%s: delta ny: %i\n",
argv[0],header->ny -nlat );
exit(-1);
}
GMT_read_grd ((filename+STRING_LENGTH),header,
(func+nlontimesnlat), 0.0, 0.0, 0.0, 0.0, dummy, 0);
gmt2myconvention_rotate4((func+nlontimesnlat),nlon,nlat,1.0);
#else
if((G = GMT_Read_Data(API, GMT_IS_GRID, GMT_IS_FILE, GMT_IS_SURFACE, GMT_CONTAINER_ONLY, NULL, (filename + STRING_LENGTH), NULL))==NULL)
return (-1);
header = G->header;
if(SIGNIFICANTLY_DIFFERENT(header->wesn[XLO], PHI2LONGITUDE(minphi)) ||
SIGNIFICANTLY_DIFFERENT(header->wesn[YLO], THETA2LATITUDE(mintheta)) ||
SIGNIFICANTLY_DIFFERENT(header->wesn[XHI], PHI2LONGITUDE(maxphi)) ||
SIGNIFICANTLY_DIFFERENT(header->wesn[YHI], THETA2LATITUDE(maxtheta)) ||
SIGNIFICANTLY_DIFFERENT(header->n_columns, nlon) ||
SIGNIFICANTLY_DIFFERENT(header->n_rows, nlat)){
fprintf(stderr,"%s: %s and %s have to have identical size\n",
argv[0],(filename),(filename+STRING_LENGTH));
fprintf(stderr,"%s: delta xmin: %g\n",
argv[0],header->wesn[XLO]-PHI2LONGITUDE(minphi));
fprintf(stderr,"%s: delta ymin: %g\n",
argv[0],header->wesn[YLO]-THETA2LATITUDE(mintheta));
fprintf(stderr,"%s: delta xmax: %g\n",
argv[0],header->wesn[XHI]-PHI2LONGITUDE(maxphi));
fprintf(stderr,"%s: delta ymax: %g\n",
argv[0],header->wesn[YHI]-THETA2LATITUDE(maxtheta));
fprintf(stderr,"%s: delta nx: %i\n",
argv[0],header->n_columns -nlon );
fprintf(stderr,"%s: delta ny: %i\n",
argv[0],header->n_rows -nlat );
exit(-1);
}
if((G = GMT_Read_Data (API, GMT_IS_GRID, GMT_IS_FILE,
GMT_IS_SURFACE, GMT_DATA_ONLY, NULL, (filename+STRING_LENGTH), G))==NULL)
return (-1);
gmt2myconvention_rotate((func+nlontimesnlat),nlon,nlat,1.0,G);
#endif
minmax(&min,&max,&avg,(func+nlontimesnlat),nlontimesnlat);
if(verbose)
field_message(minphi,dphi,maxphi,mintheta,
dtheta,maxtheta,nlon,nlat,min,max,avg,argv[0]);
} /* end vectors branch */
if(intmode == LEAST_SQUARES){
/*
need to initialize cloc array, my convention
*/
mymalloc_dp(&cloc,nlontimesnlat*2);
for(j=0,phi=minphi;j < nlon;j++,phi += dphi){ /* */
for(i=0,theta=mintheta;i<nlat;i++,theta -= dtheta){ /* theta is the fast index */
*(cloc+(j*nlat+i)*2) =phi;
*(cloc+(j*nlat+i)*2+1)=theta;
}
}
nrp = nlontimesnlat;
}
break;
}
default:{
fprintf(stderr,"%s: opmode %i is not supported\n",argv[0],opmode);exit(-1);
}}
if((opmode!=CONTOUR) && (opmode!=POINTS_FOR_A_MATRIX) && (opmode!=POINTS_FOR_AVEC_MATRIX)
&& (opmode!=POINTS_FOR_A_MATRIX_ASCII)&& (opmode!=POINTS_FOR_AVEC_MATRIX_ASCII) \
&& (opmode!=POINTS_FOR_AB)
&& (nlat < lmax+1)){
fprintf(stderr,"%s: nlat (%i) is lower than l_{max}+1 (%i).\n",argv[0],nlat,lmax+1);
fprintf(stderr,"%s: you might want to use lmax=%i instead\n",argv[0],nlat-1);
if(intmode == GAUSSIAN){
fprintf(stderr,"%s: for Gaussian integration mode, we need at least lmax+1 points\n",
argv[0]);
exit(-1);
}
}
if(verbose)fprintf(stderr,"%s: input done, computing coefficients...\n",argv[0]);
calc_coeff(&func,lmax,nlon,nlat,argv[0],opmode,cloc,nrp,intmode,verbose,
minphi,maxphi,dphi,mintheta,maxtheta,dtheta,
vectors, calculate_derivatives,pixelreg,damping);
if(verbose)
fprintf(stderr,"%s: Done.\n",argv[0]);
#ifdef USE_GMT4
// this would be nice, but it saometimes gave a double free error, need to check
//GMT_end (argc, argv);
#else
GMT_Destroy_Session (API);
#endif
return 0;
}
/*************************************************************/
/*
routine to do the SPHERICAL HARMONIC EXPANSION
*/
/*************************************************************/
void calc_coeff(DATA_PRECISION **func,int lmax,
int nlon,int nlat, char *program,
int opmode,DATA_PRECISION *cloc, int nrp,
int intmode,int verbose,DATA_PRECISION minphi,
DATA_PRECISION maxphi,
DATA_PRECISION dphi,DATA_PRECISION mintheta,
DATA_PRECISION maxtheta,
DATA_PRECISION dtheta, int vectors,
int calculate_derivatives,
int pixelreg, COMP_PRECISION damping)
{
COMP_PRECISION *y=NULL,*x=NULL,
*p=NULL,theta,phi,f[3][2],*sinfac=NULL,*dptheta=NULL,tmp,tmp2,
*cosarr=NULL,*sinarr=NULL,atmp[2],btmp[2],ctmp[2],dist,
*oneoverl=NULL,sinphi,cosphi,
*tmp_func,*absc=NULL,sintheta,costheta,dummy,*sndder,phirange,
thetarange,norm,chi2,sum,dnorm;
static COMP_PRECISION normfac,sqrt2normfac;
int i=0,j,k,l,m,n,lmsize,lmax1,memneeded,lookup_p=1,coeff,nlontimesnlat,klim,
lmc,nr_gauss_pts,lmax_max,os1,os2,ic1,ic2,os3,nsize[3],nrp2,
jc1,jc2,phisumlim,thetasumlim,ysumn,warned,damping_mode;
float *aflt = NULL,*sol=NULL, *sigma = NULL,*aflt_copy=NULL;
DATA_PRECISION *tmp_dfunc,*tmp_cloc;
/*
y is cos(theta)
x is phi
cosarr and sinarr hold sin(m phi) and cos(m phi)
p is the legendre function arrar
DPTHETA and DP the derivatives of the legendre functions (factors)
*/
y=x=sinarr=cosarr=NULL;
#ifdef REDUCE_LMAX
/*
this is some sort of check if we have enough data points
for trapeziodal integration: degree l would normally need
2l data points to be exact;
*/
lmax_max=(int)(nlat-1)/2;
#else
lmax_max=lmax;
#endif
// some size and length variables
lmsize= (int)((((COMP_PRECISION)lmax)+1.0)*(((COMP_PRECISION)lmax)+2)/2.0);
lmax1=lmax+1;
// number of Gauss integration points, should be lmax+1
nr_gauss_pts=lmax1;
nlontimesnlat=nlon*nlat;
if(intmode == LEAST_SQUARES){
/*
LEAST SQUARES ESTIMATE
*/
/* switch damping mode */
if(fabs(damping) > 0)
if(damping < 0){
damping_mode = 2; /* gradient */
damping = -damping;
}
else
damping_mode = 1; /* norm */
else
damping_mode = 0;
if(verbose){
fprintf(stderr,"%s: least squares inversion, ",program);
if(damping_mode)
fprintf(stderr,"%s damping with %g\n",(damping_mode == 1)?("norm"):("roughness"),damping);
else
fprintf(stderr,"no damping\n");
}
/*
start inversion part
*/
switch(opmode){
case POINTS_FOR_A_MATRIX:
case POINTS_FOR_A_MATRIX_ASCII:
case POINTS_FOR_AVEC_MATRIX:
case POINTS_FOR_AVEC_MATRIX_ASCII:
case POINTS_FOR_AB:
case CONTOUR:
case GMT_BLOCK:
/*
weights for every data point
*/
nrp2 = vectors ? nrp * 2 : nrp;
myrealloc_dp(&sigma,nrp2);
/* assign sigma = 1 / weight */
if(opmode != GMT_BLOCK){
for(i=0;i < nrp2;i++)
sigma[i] = 1.0;
}else{ /* GMT block */
/*
block data, weigh by latitude
*/
k = n = 0;
tmp_dfunc = tmp_cloc = NULL;
for(i=0;i < nrp;i++){
sum = sin(cloc[i*2+1]); /* sin(theta) */
if(fabs(sum) > EPS_DATA_PREC){
/* use this point */
sigma[n] = 1/sum;
/* location */
myrealloc_dp(&tmp_cloc,2*(k+1));
tmp_cloc[k*2] = cloc[i*2];
tmp_cloc[k*2+1]=cloc[i*2+1];
/* first function value */
myrealloc_dp(&tmp_dfunc,(n+1));
tmp_dfunc[n] = *(*func+i);
n++;
if(vectors){ /* */
myrealloc_dp(&tmp_dfunc,(n+1));
tmp_dfunc[n] = *(*func+nrp+i);
sigma[n] = 1/sum;
n++;
}
k++;
}
}
if(n != nrp2){fprintf(stderr,"logic error 001, %i %i\n",n,nrp2);exit(-1);}
fprintf(stderr,"%s: used %i out of %i total points (%i total values), rejected %i polar latitude points in GMT block\n",
program,k,nrp,n,nrp-k);
nrp = k; /* shrink */
nrp2 = vectors ? nrp * 2 : nrp;
myrealloc_dp(&sigma,nrp2);
myrealloc_dp(&cloc,nrp*2);
for(i=0;i < nrp;i++){ /* assign */
cloc[i*2] = tmp_cloc[i*2];
cloc[i*2+1] = tmp_cloc[i*2+1];
}
if(vectors){
if(n != nrp2){
fprintf(stderr,"logic error 12\n");exit(-1);}
}
myrealloc_dp(func,nrp2);
for(i=0;i < nrp2;i++) /* assign */
*(*func+i) = tmp_dfunc[i];
free(tmp_dfunc);free(tmp_cloc);
}
/*
get the A matrix, fortran style, and use sigma
*/
make_a(lmax,lmsize,nrp,&aflt,cloc,verbose,
program,damping_mode,nsize,damping,TRUE,TRUE,sigma,vectors);
/*
save, for variance reduction (choose speed over mem usage)
*/
myrealloc_dp(&aflt_copy,nsize[2]);
memcpy(aflt_copy,aflt,(size_t)(nsize[2])*sizeof(float));
/*
assign data to solution vector which will get overwritten
*/
myrealloc_dp(&sol,nsize[0]);
for(i=0;i < nrp2;i++) /* function values */
sol[i] = *(*func+i)/sigma[i];
for(;i < nsize[0];i++) /* for damping */
sol[i] = 0.0;
/* solve and overwrite function values with solution
coefficients and a matrix with something else */
#ifdef USE_LAPACK
solver_ab_lls(aflt,nsize[0],nsize[1],sol,program);
#else
fprintf(stderr,"%s: recompile with LAPACK support, need solver_ab_lls for least squares \n",program);
exit(-1);
#endif
/* don't need that output */
free(aflt);
/*
assess solution
*/
/* solution norm */
norm = chi2 = 0.0;
for(i=0;i < nsize[1];i++)
norm += sol[i]*sol[i];
norm = sqrt(norm);
/*
misfit
*/
chi2 = dnorm = 0.0;
klim = vectors ? 2 : 1;
n = 0;
for(i=0;i < nrp;i++){
for(k = 0; k < klim;k++){
sum = 0.0;
for(j=0;j < nsize[1];j++) /* aflt is fortran style and
weighted, need to undo */
sum += aflt_copy[j*nsize[0]+n] * sol[j];
sum *= sigma[n]; /* undo 1/sigma scaling for A(i,j) */
/* misfit/sigma */
sum = (sum - (*(*func+k*nrp+i)))/sigma[n];
sum = sum * sum;
chi2 += sum;
dnorm += (*(*func+k*nrp+i)) * (*(*func+k*nrp+i));
n++;
}
}
dnorm = sqrt(dnorm);
fprintf(stderr,"%s: solution DOF: %i solution norm: %.5e VR(1-chi/|x|): %.1f%%\n",
program,nsize[1],norm,(1-sqrt(chi2)/dnorm)*100);
free(aflt_copy);free(sigma); /* no use anymore */
/*
output of solution
*/
fprintf(stdout,"%i\n",lmax);
if(vectors){
fprintf(stdout,COEFF_DATA_FORMAT,0.0,0.0);fprintf(stdout,"\t");fprintf(stdout,COEFF_DATA_FORMAT,0.0,0.0); /* no degree zero term */
fprintf(stdout,"\n");
for(j=0,l=1;l <= lmax;l++){
for(m=0;m <= l;m++){
if(m == 0){
fprintf(stdout,COEFF_DATA_FORMAT,sol[j],0.0); /* A */
j++;
fprintf(stdout,"\t"); /* C */
fprintf(stdout,COEFF_DATA_FORMAT,sol[j],0.0);
j++;
}else{
fprintf(stdout,COEFF_DATA_FORMAT,sol[j],sol[j+1]); /* AB */
j += 2;
fprintf(stdout,"\t");
fprintf(stdout,COEFF_DATA_FORMAT,sol[j],sol[j+1]); /* CD */
j += 2;
}
fprintf(stdout,"\n");
}
}
}else{ /* scalar */
for(j=l=0;l <= lmax;l++){
for(m=0;m <= l;m++){
/* scalar */
if(m == 0){
fprintf(stdout,COEFF_DATA_FORMAT,sol[j],0.0);
j++;
}else{
fprintf(stdout,COEFF_DATA_FORMAT,sol[j],sol[j+1]);
j += 2;
}
fprintf(stdout,"\n");
}
}
}
free(sol);
break;
default:
fprintf(stderr,"opmode %i not implemented for least squares yet\n",opmode);
exit(-1);
}
/*
END LEAST SQUARES
*/
}else{
/*
ALL OTHER STRAIGHTFORWARD INTEGRATION MODES
*/
//
// decide on integration mode depending on which opmode is chosen
// and determine summing bounds and alike
//
switch(opmode){
case GMT_BLOCK:
case ASCII_BLOCK_HEADER:
case ASCII_BLOCK:
case BLOCK:{
/* use block data */
phirange=fabs(maxphi-minphi);
thetarange=fabs(maxtheta-mintheta);
// check the boundaries in longitude to see if full geographical coverage
ic1=(fabs(phirange-TWOPI)<1e-6);
ic2=(fabs(phirange-TWOPI+dphi)<1e-6);
if(!ic1 && !ic2){
fprintf(stderr,"%s: expecting full longitudinal coverage, ics: %i %i, range: %g dx: %g\n",
program,ic1,ic2,PHI2LONGITUDE(phirange),PHI2LONGITUDE(dphi));
exit(-1);
}else{
if(ic1){
/*
full 360 range with redundant repetition
of phi=0 value. limit summation to 0...360-dx
*/
phisumlim=nlon-1;
// check if values at Greenwich are same
for(os1=(nlon-1)*nlat,i=warned=0;i<nlat;i++)
if((*func)[i] != ((*func)[os1+i])){
if(!warned){
fprintf(stderr,
"%s: at least one mismatch at the Greenwich, averaging\n",
program);
warned=1;
}
// average in this case
(*func)[i]=((*func)[i] + ((*func)[os1+i]))/2.0;
}
}else{ // Greenwich only in field once, use all points in x direction
phisumlim=nlon;
}
}
// now check latitudes
jc1=(fabs(thetarange-PI)<1e-2);
jc2=(fabs(thetarange-PI+dtheta)<1e-2);
if(!jc1 && !jc2){
fprintf(stderr,"%s: expecting full latitudinal coverage, range: %g (%g), %g(%g) %g(%g)\n",
program,RAD2DEG(thetarange),RAD2DEG(PI-dtheta),mintheta,THETA2LATITUDE(mintheta),
maxtheta,THETA2LATITUDE(maxtheta));
exit(-1);
}else{
thetasumlim=nlat;
if(jc1){
/*
full 180 range with both North and South pole values
use only half values at poles if we use sum rule, for Gauss, it's OK
*/
if(intmode == TRAPEZOIDAL)
for(os1=nlat-1,j=i=0;i<nlon;i++,j+=nlat){
((*func)[j]) /= 2.0;
((*func)[j+os1]) /= 2.0;
}
// we should count both polar entries only as a half box, thus
ysumn=nlat-1;
// NOTE: this is different from the integration loop end which is
// thetasumlim and always equal to nlat
}else{
// check if whole lat range covered, ie. if offset like -90+dy/2 and 90-dy/2
if((fabs(mintheta-PI+dtheta/2.0)>1e-6)||
(fabs(maxtheta-dtheta/2.0)>1e-6)){
fprintf(stderr,"%s: error: use latitude range from %g to %g for %g spacing\n",
program,THETA2LATITUDE(PI-dtheta/2.0),
THETA2LATITUDE(dtheta/2.0),RAD2DEG(dtheta));
exit(-1);
}
// all boxes count in this case to divide PI
ysumn=nlat;
}
}
/*
x coordinates (0 to 2 pi) if pixel registration, shift from original values
to mid of cell
*/
myrealloc_cp(&x,nlon);zero_cp(x,nlon);
for(j=0,phi=minphi; j < nlon;j++,phi += dphi)
x[j]= phi;
/*
y coordinates, form 0 (north) to pi (south), y is cos(theta)
*/
myrealloc_cp(&y,nlat);zero_cp(y,nlat);
// shouldn't use dtheta since this is positive by definition
tmp=(COMP_PRECISION)(maxtheta-mintheta)/(COMP_PRECISION)(nlat-1);
// y coords from now on are in cos(theta)
for(i=0,theta=mintheta; i<nlat;i++,theta += tmp){
y[i] = cos(theta);
/* fprintf(stderr,"i: %i nlat: %i y[i]: %12g\n",i,nlat,y[i]); */
}
#ifdef DEBUG
fprintf(stderr,"%s: integration coord range: x: %g to %g y: %g(%g) to %g(%g)\n",
program,
PHI2LONGITUDE(x[0]), PHI2LONGITUDE(x[phisumlim-1]),
THETA2LATITUDE(acos(y[0])),acos(y[0]),
THETA2LATITUDE(acos(y[thetasumlim-1])),
acos(y[thetasumlim-1]));
#endif
if((intmode == TRAPEZOIDAL) || (intmode == GAUSSIAN)){
/*
cos and sin factors for R l m and S l m,
cos(m phi) and sin(m phi)
would not need that for FFT method
*/
myrealloc_cp(&cosarr,nlon*lmax1);zero_cp(cosarr,nlon*lmax1);
myrealloc_cp(&sinarr,nlon*lmax1);zero_cp(sinarr,nlon*lmax1);
// m=0 terms
for(j=0; j < nlon ;j++){
cosarr[j]=1.0; sinarr[j]=0.0;
}
// 1 <= m <= lmax
for(os1=nlon,m=1;m<=lmax;m++,os1+=nlon){
for(j=0; j < nlon ;j++){
tmp = ((COMP_PRECISION)m)*x[j];
cosarr[os1+j]=cos(tmp);
sinarr[os1+j]=sin(tmp);
}
}
}
/*
Legendre functions, assign P_lm(theta) to array P
using y=cos(theta)[j]
*/
if(verbose)
fprintf(stderr,
"%s: constructing Assoc. legendre functions...\n",program);
myrealloc_cp(&p,lmsize*nlat);zero_cp(p,lmsize*nlat);
if(verbose)
fprintf(stderr,"%s: using %6g MBs for P array\n",
program,
lmsize*sizeof(COMP_PRECISION)*nlat/ONE_MEGABYTE);
if(intmode != GAUSSIAN){
//
// obtain Legendre functions at the latitudinal points of the input data
//
plgndr(y,nlat,p,lmax);
// set these such that derivatives and sinfac etc.
// work for Gauss points and original input which is given on nlat points
//
nr_gauss_pts=nlat;
absc=y;// absc values are the y coordinates
}else{
/*
obtain Legendre functions at Gauss integration points
(should be lmax+1), the P_lm will then also be already
multiplied by the weight factor for quadrature
we furthermore obtain the abscissae values from this
routine, to replace the y[j] later (still need y for interpolation later)
*/
myrealloc_cp(&absc,nr_gauss_pts);
plgndr_g(p,lmax,nr_gauss_pts,absc);
}
/*
integretation area sin (theta) d theta d phi, sine factor,
sin(arccos(y))=sqrt(1-y^2)
also used for derivatives, therefore calculate here
for real abscissae values, not necessary y
*/
myrealloc_cp(&sinfac,nr_gauss_pts);zero_cp(sinfac,nr_gauss_pts);
for(i=0;i<nr_gauss_pts;i++)
sinfac[i] = (((tmp=(1.0-absc[i]*absc[i]))>=0.0)?(sqrt(tmp)):(0.0));
if(calculate_derivatives){
/*
to expand vector fields calculate derivative of P with respect to theta
--> DPTHETA,