forked from tomojitakasu/RTKLIB
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathrtkpos.c
1858 lines (1640 loc) · 68.1 KB
/
rtkpos.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
/*------------------------------------------------------------------------------
* rtkpos.c : precise positioning
*
* Copyright (C) 2007-2018 by T.TAKASU, All rights reserved.
*
* version : $Revision: 1.1 $ $Date: 2008/07/17 21:48:06 $
* history : 2007/01/12 1.0 new
* 2007/03/13 1.1 add slip detection by LLI flag
* 2007/04/18 1.2 add antenna pcv correction
* change rtkpos argin
* 2008/07/18 1.3 refactored
* 2009/01/02 1.4 modify rtk positioning api
* 2009/03/09 1.5 support glonass, gallileo and qzs
* 2009/08/27 1.6 fix bug on numerical exception
* 2009/09/03 1.7 add check of valid satellite number
* add check time sync for moving-base
* 2009/11/23 1.8 add api rtkopenstat(),rtkclosestat()
* add receiver h/w bias estimation
* add solution status output
* 2010/04/04 1.9 support ppp-kinematic and ppp-static modes
* support earth tide correction
* changed api:
* rtkpos()
* 2010/09/07 1.10 add elevation mask to hold ambiguity
* 2012/02/01 1.11 add extended receiver error model
* add glonass interchannel bias correction
* add slip detectior by L1-L5 gf jump
* output snr of rover receiver in residuals
* 2013/03/10 1.12 add otl and pole tides corrections
* 2014/05/26 1.13 support beidou and galileo
* add output of gal-gps and bds-gps time offset
* 2014/05/28 1.14 fix bug on memory exception with many sys and freq
* 2014/08/26 1.15 add functino to swap sol-stat file with keywords
* 2014/10/21 1.16 fix bug on beidou amb-res with pos2-bdsarmode=0
* 2014/11/08 1.17 fix bug on ar-degradation by unhealthy satellites
* 2015/03/23 1.18 residuals referenced to reference satellite
* 2018/01/29 1.19 unfix ambiguity between gps and qzss
*-----------------------------------------------------------------------------*/
#include <stdarg.h>
#include "rtklib.h"
static const char rcsid[]="$Id:$";
/* constants/macros ----------------------------------------------------------*/
#define SQR(x) ((x)*(x))
#define SQRT(x) ((x)<=0.0?0.0:sqrt(x))
#define MIN(x,y) ((x)<=(y)?(x):(y))
#define ROUND(x) (int)floor((x)+0.5)
#define VAR_POS SQR(30.0) /* initial variance of receiver pos (m^2) */
#define VAR_VEL SQR(10.0) /* initial variance of receiver vel ((m/s)^2) */
#define VAR_ACC SQR(10.0) /* initial variance of receiver acc ((m/ss)^2) */
#define VAR_HWBIAS SQR(1.0) /* initial variance of h/w bias ((m/MHz)^2) */
#define VAR_GRA SQR(0.001) /* initial variance of gradient (m^2) */
#define INIT_ZWD 0.15 /* initial zwd (m) */
#define PRN_HWBIAS 1E-6 /* process noise of h/w bias (m/MHz/sqrt(s)) */
#define GAP_RESION 120 /* gap to reset ionosphere parameters (epochs) */
#define MAXACC 30.0 /* max accel for doppler slip detection (m/s^2) */
#define VAR_HOLDAMB 0.001 /* constraint to hold ambiguity (cycle^2) */
#define TTOL_MOVEB (1.0+2*DTTOL)
/* time sync tolerance for moving-baseline (s) */
/* number of parameters (pos,ionos,tropos,hw-bias,phase-bias,real,estimated) */
#define NF(opt) ((opt)->ionoopt==IONOOPT_IFLC?1:(opt)->nf)
#define NP(opt) ((opt)->dynamics==0?3:9)
#define NI(opt) ((opt)->ionoopt!=IONOOPT_EST?0:MAXSAT)
#define NT(opt) ((opt)->tropopt<TROPOPT_EST?0:((opt)->tropopt<TROPOPT_ESTG?2:6))
#define NL(opt) ((opt)->glomodear!=2?0:NFREQGLO)
#define NB(opt) ((opt)->mode<=PMODE_DGPS?0:MAXSAT*NF(opt))
#define NR(opt) (NP(opt)+NI(opt)+NT(opt)+NL(opt))
#define NX(opt) (NR(opt)+NB(opt))
/* state variable index */
#define II(s,opt) (NP(opt)+(s)-1) /* ionos (s:satellite no) */
#define IT(r,opt) (NP(opt)+NI(opt)+NT(opt)/2*(r)) /* tropos (r:0=rov,1:ref) */
#define IL(f,opt) (NP(opt)+NI(opt)+NT(opt)+(f)) /* receiver h/w bias */
#define IB(s,f,opt) (NR(opt)+MAXSAT*(f)+(s)-1) /* phase bias (s:satno,f:freq) */
#ifdef EXTGSI
extern int resamb_WLNL(rtk_t *rtk, const obsd_t *obs, const int *sat,
const int *iu, const int *ir, int ns, const nav_t *nav,
const double *azel);
extern int resamb_TCAR(rtk_t *rtk, const obsd_t *obs, const int *sat,
const int *iu, const int *ir, int ns, const nav_t *nav,
const double *azel);
#else
extern int resamb_WLNL(rtk_t *rtk, const obsd_t *obs, const int *sat,
const int *iu, const int *ir, int ns, const nav_t *nav,
const double *azel) {return 0;}
extern int resamb_TCAR(rtk_t *rtk, const obsd_t *obs, const int *sat,
const int *iu, const int *ir, int ns, const nav_t *nav,
const double *azel) {return 0;}
#endif
/* global variables ----------------------------------------------------------*/
static int statlevel=0; /* rtk status output level (0:off) */
static FILE *fp_stat=NULL; /* rtk status file pointer */
static char file_stat[1024]=""; /* rtk status file original path */
static gtime_t time_stat={0}; /* rtk status file time */
/* open solution status file ---------------------------------------------------
* open solution status file and set output level
* args : char *file I rtk status file
* int level I rtk status level (0: off)
* return : status (1:ok,0:error)
* notes : file can constain time keywords (%Y,%y,%m...) defined in reppath().
* The time to replace keywords is based on UTC of CPU time.
* output : solution status file record format
*
* $POS,week,tow,stat,posx,posy,posz,posxf,posyf,poszf
* week/tow : gps week no/time of week (s)
* stat : solution status
* posx/posy/posz : position x/y/z ecef (m) float
* posxf/posyf/poszf : position x/y/z ecef (m) fixed
*
* $VELACC,week,tow,stat,vele,veln,velu,acce,accn,accu,velef,velnf,veluf,accef,accnf,accuf
* week/tow : gps week no/time of week (s)
* stat : solution status
* vele/veln/velu : velocity e/n/u (m/s) float
* acce/accn/accu : acceleration e/n/u (m/s^2) float
* velef/velnf/veluf : velocity e/n/u (m/s) fixed
* accef/accnf/accuf : acceleration e/n/u (m/s^2) fixed
*
* $CLK,week,tow,stat,clk1,clk2,clk3,clk4
* week/tow : gps week no/time of week (s)
* stat : solution status
* clk1 : receiver clock bias GPS (ns)
* clk2 : receiver clock bias GLO-GPS (ns)
* clk3 : receiver clock bias GAL-GPS (ns)
* clk4 : receiver clock bias BDS-GPS (ns)
*
* $ION,week,tow,stat,sat,az,el,ion,ion-fixed
* week/tow : gps week no/time of week (s)
* stat : solution status
* sat : satellite id
* az/el : azimuth/elevation angle(deg)
* ion : vertical ionospheric delay L1 (m) float
* ion-fixed: vertical ionospheric delay L1 (m) fixed
*
* $TROP,week,tow,stat,rcv,ztd,ztdf
* week/tow : gps week no/time of week (s)
* stat : solution status
* rcv : receiver (1:rover,2:base station)
* ztd : zenith total delay (m) float
* ztdf : zenith total delay (m) fixed
*
* $HWBIAS,week,tow,stat,frq,bias,biasf
* week/tow : gps week no/time of week (s)
* stat : solution status
* frq : frequency (1:L1,2:L2,...)
* bias : h/w bias coefficient (m/MHz) float
* biasf : h/w bias coefficient (m/MHz) fixed
*
* $SAT,week,tow,sat,frq,az,el,resp,resc,vsat,snr,fix,slip,lock,outc,slipc,rejc
* week/tow : gps week no/time of week (s)
* sat/frq : satellite id/frequency (1:L1,2:L2,...)
* az/el : azimuth/elevation angle (deg)
* resp : pseudorange residual (m)
* resc : carrier-phase residual (m)
* vsat : valid data flag (0:invalid,1:valid)
* snr : signal strength (dbHz)
* fix : ambiguity flag (0:no data,1:float,2:fixed,3:hold,4:ppp)
* slip : cycle-slip flag (bit1:slip,bit2:parity unknown)
* lock : carrier-lock count
* outc : data outage count
* slipc : cycle-slip count
* rejc : data reject (outlier) count
*
*-----------------------------------------------------------------------------*/
extern int rtkopenstat(const char *file, int level)
{
gtime_t time=utc2gpst(timeget());
char path[1024];
trace(3,"rtkopenstat: file=%s level=%d\n",file,level);
if (level<=0) return 0;
reppath(file,path,time,"","");
if (!(fp_stat=fopen(path,"w"))) {
trace(1,"rtkopenstat: file open error path=%s\n",path);
return 0;
}
strcpy(file_stat,file);
time_stat=time;
statlevel=level;
return 1;
}
/* close solution status file --------------------------------------------------
* close solution status file
* args : none
* return : none
*-----------------------------------------------------------------------------*/
extern void rtkclosestat(void)
{
trace(3,"rtkclosestat:\n");
if (fp_stat) fclose(fp_stat);
fp_stat=NULL;
file_stat[0]='\0';
statlevel=0;
}
/* swap solution status file -------------------------------------------------*/
static void swapsolstat(void)
{
gtime_t time=utc2gpst(timeget());
char path[1024];
if ((int)(time2gpst(time ,NULL)/INT_SWAP_STAT)==
(int)(time2gpst(time_stat,NULL)/INT_SWAP_STAT)) {
return;
}
time_stat=time;
if (!reppath(file_stat,path,time,"","")) {
return;
}
if (fp_stat) fclose(fp_stat);
if (!(fp_stat=fopen(path,"w"))) {
trace(2,"swapsolstat: file open error path=%s\n",path);
return;
}
trace(3,"swapsolstat: path=%s\n",path);
}
/* output solution status ----------------------------------------------------*/
static void outsolstat(rtk_t *rtk)
{
ssat_t *ssat;
double tow,pos[3],vel[3],acc[3],vela[3]={0},acca[3]={0},xa[3];
int i,j,week,est,nfreq,nf=NF(&rtk->opt);
char id[32];
if (statlevel<=0||!fp_stat) return;
trace(3,"outsolstat:\n");
/* swap solution status file */
swapsolstat();
est=rtk->opt.mode>=PMODE_DGPS;
nfreq=est?nf:1;
tow=time2gpst(rtk->sol.time,&week);
/* receiver position */
if (est) {
for (i=0;i<3;i++) xa[i]=i<rtk->na?rtk->xa[i]:0.0;
fprintf(fp_stat,"$POS,%d,%.3f,%d,%.4f,%.4f,%.4f,%.4f,%.4f,%.4f\n",week,tow,
rtk->sol.stat,rtk->x[0],rtk->x[1],rtk->x[2],xa[0],xa[1],xa[2]);
}
else {
fprintf(fp_stat,"$POS,%d,%.3f,%d,%.4f,%.4f,%.4f,%.4f,%.4f,%.4f\n",week,tow,
rtk->sol.stat,rtk->sol.rr[0],rtk->sol.rr[1],rtk->sol.rr[2],
0.0,0.0,0.0);
}
/* receiver velocity and acceleration */
if (est&&rtk->opt.dynamics) {
ecef2pos(rtk->sol.rr,pos);
ecef2enu(pos,rtk->x+3,vel);
ecef2enu(pos,rtk->x+6,acc);
if (rtk->na>=6) ecef2enu(pos,rtk->xa+3,vela);
if (rtk->na>=9) ecef2enu(pos,rtk->xa+6,acca);
fprintf(fp_stat,"$VELACC,%d,%.3f,%d,%.4f,%.4f,%.4f,%.5f,%.5f,%.5f,%.4f,%.4f,%.4f,%.5f,%.5f,%.5f\n",
week,tow,rtk->sol.stat,vel[0],vel[1],vel[2],acc[0],acc[1],acc[2],
vela[0],vela[1],vela[2],acca[0],acca[1],acca[2]);
}
else {
ecef2pos(rtk->sol.rr,pos);
ecef2enu(pos,rtk->sol.rr+3,vel);
fprintf(fp_stat,"$VELACC,%d,%.3f,%d,%.4f,%.4f,%.4f,%.5f,%.5f,%.5f,%.4f,%.4f,%.4f,%.5f,%.5f,%.5f\n",
week,tow,rtk->sol.stat,vel[0],vel[1],vel[2],
0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0,0.0);
}
/* receiver clocks */
fprintf(fp_stat,"$CLK,%d,%.3f,%d,%d,%.3f,%.3f,%.3f,%.3f\n",
week,tow,rtk->sol.stat,1,rtk->sol.dtr[0]*1E9,rtk->sol.dtr[1]*1E9,
rtk->sol.dtr[2]*1E9,rtk->sol.dtr[3]*1E9);
/* ionospheric parameters */
if (est&&rtk->opt.ionoopt==IONOOPT_EST) {
for (i=0;i<MAXSAT;i++) {
ssat=rtk->ssat+i;
if (!ssat->vs) continue;
satno2id(i+1,id);
j=II(i+1,&rtk->opt);
xa[0]=j<rtk->na?rtk->xa[j]:0.0;
fprintf(fp_stat,"$ION,%d,%.3f,%d,%s,%.1f,%.1f,%.4f,%.4f\n",week,tow,rtk->sol.stat,
id,ssat->azel[0]*R2D,ssat->azel[1]*R2D,rtk->x[j],xa[0]);
}
}
/* tropospheric parameters */
if (est&&(rtk->opt.tropopt==TROPOPT_EST||rtk->opt.tropopt==TROPOPT_ESTG)) {
for (i=0;i<2;i++) {
j=IT(i,&rtk->opt);
xa[0]=j<rtk->na?rtk->xa[j]:0.0;
fprintf(fp_stat,"$TROP,%d,%.3f,%d,%d,%.4f,%.4f\n",week,tow,rtk->sol.stat,
i+1,rtk->x[j],xa[0]);
}
}
/* receiver h/w bias */
if (est&&rtk->opt.glomodear==2) {
for (i=0;i<nfreq;i++) {
j=IL(i,&rtk->opt);
xa[0]=j<rtk->na?rtk->xa[j]:0.0;
fprintf(fp_stat,"$HWBIAS,%d,%.3f,%d,%d,%.4f,%.4f\n",week,tow,rtk->sol.stat,
i+1,rtk->x[j],xa[0]);
}
}
if (rtk->sol.stat==SOLQ_NONE||statlevel<=1) return;
/* residuals and status */
for (i=0;i<MAXSAT;i++) {
ssat=rtk->ssat+i;
if (!ssat->vs) continue;
satno2id(i+1,id);
for (j=0;j<nfreq;j++) {
fprintf(fp_stat,"$SAT,%d,%.3f,%s,%d,%.1f,%.1f,%.4f,%.4f,%d,%.0f,%d,%d,%d,%d,%d,%d\n",
week,tow,id,j+1,ssat->azel[0]*R2D,ssat->azel[1]*R2D,
ssat->resp [j],ssat->resc[j], ssat->vsat[j],ssat->snr[j]*0.25,
ssat->fix [j],ssat->slip[j]&3,ssat->lock[j],ssat->outc[j],
ssat->slipc[j],ssat->rejc[j]);
}
}
}
/* save error message --------------------------------------------------------*/
static void errmsg(rtk_t *rtk, const char *format, ...)
{
char buff[256],tstr[32];
int n;
va_list ap;
time2str(rtk->sol.time,tstr,2);
n=sprintf(buff,"%s: ",tstr+11);
va_start(ap,format);
n+=vsprintf(buff+n,format,ap);
va_end(ap);
n=n<MAXERRMSG-rtk->neb?n:MAXERRMSG-rtk->neb;
memcpy(rtk->errbuf+rtk->neb,buff,n);
rtk->neb+=n;
trace(2,"%s",buff);
}
/* single-differenced observable ---------------------------------------------*/
static double sdobs(const obsd_t *obs, int i, int j, int f)
{
double pi=f<NFREQ?obs[i].L[f]:obs[i].P[f-NFREQ];
double pj=f<NFREQ?obs[j].L[f]:obs[j].P[f-NFREQ];
return pi==0.0||pj==0.0?0.0:pi-pj;
}
/* single-differenced geometry-free linear combination of phase --------------*/
static double gfobs_L1L2(const obsd_t *obs, int i, int j, const double *lam)
{
double pi=sdobs(obs,i,j,0)*lam[0],pj=sdobs(obs,i,j,1)*lam[1];
return pi==0.0||pj==0.0?0.0:pi-pj;
}
static double gfobs_L1L5(const obsd_t *obs, int i, int j, const double *lam)
{
double pi=sdobs(obs,i,j,0)*lam[0],pj=sdobs(obs,i,j,2)*lam[2];
return pi==0.0||pj==0.0?0.0:pi-pj;
}
/* single-differenced measurement error variance -----------------------------*/
static double varerr(int sat, int sys, double el, double bl, double dt, int f,
const prcopt_t *opt)
{
double a,b,c=opt->err[3]*bl/1E4,d=CLIGHT*opt->sclkstab*dt,fact=1.0;
double sinel=sin(el);
int i=sys==SYS_GLO?1:(sys==SYS_GAL?2:0),nf=NF(opt);
/* extended error model */
if (f>=nf&&opt->exterr.ena[0]) { /* code */
a=opt->exterr.cerr[i][ (f-nf)*2];
b=opt->exterr.cerr[i][1+(f-nf)*2];
if (sys==SYS_SBS) {a*=EFACT_SBS; b*=EFACT_SBS;}
}
else if (f<nf&&opt->exterr.ena[1]) { /* phase */
a=opt->exterr.perr[i][ f*2];
b=opt->exterr.perr[i][1+f*2];
if (sys==SYS_SBS) {a*=EFACT_SBS; b*=EFACT_SBS;}
}
else { /* normal error model */
if (f>=nf) fact=opt->eratio[f-nf];
if (fact<=0.0) fact=opt->eratio[0];
fact*=sys==SYS_GLO?EFACT_GLO:(sys==SYS_SBS?EFACT_SBS:EFACT_GPS);
a=fact*opt->err[1];
b=fact*opt->err[2];
}
return 2.0*(opt->ionoopt==IONOOPT_IFLC?3.0:1.0)*(a*a+b*b/sinel/sinel+c*c)+d*d;
}
/* baseline length -----------------------------------------------------------*/
static double baseline(const double *ru, const double *rb, double *dr)
{
int i;
for (i=0;i<3;i++) dr[i]=ru[i]-rb[i];
return norm(dr,3);
}
/* initialize state and covariance -------------------------------------------*/
static void initx(rtk_t *rtk, double xi, double var, int i)
{
int j;
rtk->x[i]=xi;
for (j=0;j<rtk->nx;j++) {
rtk->P[i+j*rtk->nx]=rtk->P[j+i*rtk->nx]=i==j?var:0.0;
}
}
/* select common satellites between rover and reference station --------------*/
static int selsat(const obsd_t *obs, double *azel, int nu, int nr,
const prcopt_t *opt, int *sat, int *iu, int *ir)
{
int i,j,k=0;
trace(3,"selsat : nu=%d nr=%d\n",nu,nr);
for (i=0,j=nu;i<nu&&j<nu+nr;i++,j++) {
if (obs[i].sat<obs[j].sat) j--;
else if (obs[i].sat>obs[j].sat) i--;
else if (azel[1+j*2]>=opt->elmin) { /* elevation at base station */
sat[k]=obs[i].sat; iu[k]=i; ir[k++]=j;
trace(4,"(%2d) sat=%3d iu=%2d ir=%2d\n",k-1,obs[i].sat,i,j);
}
}
return k;
}
/* temporal update of position/velocity/acceleration -------------------------*/
static void udpos(rtk_t *rtk, double tt)
{
double *F,*FP,*xp,pos[3],Q[9]={0},Qv[9],var=0.0;
int i,j;
trace(3,"udpos : tt=%.3f\n",tt);
/* fixed mode */
if (rtk->opt.mode==PMODE_FIXED) {
for (i=0;i<3;i++) initx(rtk,rtk->opt.ru[i],1E-8,i);
return;
}
/* initialize position for first epoch */
if (norm(rtk->x,3)<=0.0) {
for (i=0;i<3;i++) initx(rtk,rtk->sol.rr[i],VAR_POS,i);
if (rtk->opt.dynamics) {
for (i=3;i<6;i++) initx(rtk,rtk->sol.rr[i],VAR_VEL,i);
for (i=6;i<9;i++) initx(rtk,1E-6,VAR_ACC,i);
}
}
/* static mode */
if (rtk->opt.mode==PMODE_STATIC) return;
/* kinmatic mode without dynamics */
if (!rtk->opt.dynamics) {
for (i=0;i<3;i++) initx(rtk,rtk->sol.rr[i],VAR_POS,i);
return;
}
/* check variance of estimated postion */
for (i=0;i<3;i++) var+=rtk->P[i+i*rtk->nx]; var/=3.0;
if (var>VAR_POS) {
/* reset position with large variance */
for (i=0;i<3;i++) initx(rtk,rtk->sol.rr[i],VAR_POS,i);
for (i=3;i<6;i++) initx(rtk,rtk->sol.rr[i],VAR_VEL,i);
for (i=6;i<9;i++) initx(rtk,1E-6,VAR_ACC,i);
trace(2,"reset rtk position due to large variance: var=%.3f\n",var);
return;
}
/* state transition of position/velocity/acceleration */
F=eye(rtk->nx); FP=mat(rtk->nx,rtk->nx); xp=mat(rtk->nx,1);
for (i=0;i<6;i++) {
F[i+(i+3)*rtk->nx]=tt;
}
/* x=F*x, P=F*P*F+Q */
matmul("NN",rtk->nx,1,rtk->nx,1.0,F,rtk->x,0.0,xp);
matcpy(rtk->x,xp,rtk->nx,1);
matmul("NN",rtk->nx,rtk->nx,rtk->nx,1.0,F,rtk->P,0.0,FP);
matmul("NT",rtk->nx,rtk->nx,rtk->nx,1.0,FP,F,0.0,rtk->P);
/* process noise added to only acceleration */
Q[0]=Q[4]=SQR(rtk->opt.prn[3]); Q[8]=SQR(rtk->opt.prn[4]);
ecef2pos(rtk->x,pos);
covecef(pos,Q,Qv);
for (i=0;i<3;i++) for (j=0;j<3;j++) {
rtk->P[i+6+(j+6)*rtk->nx]+=Qv[i+j*3];
}
free(F); free(FP); free(xp);
}
/* temporal update of ionospheric parameters ---------------------------------*/
static void udion(rtk_t *rtk, double tt, double bl, const int *sat, int ns)
{
double el,fact;
int i,j;
trace(3,"udion : tt=%.1f bl=%.0f ns=%d\n",tt,bl,ns);
for (i=1;i<=MAXSAT;i++) {
j=II(i,&rtk->opt);
if (rtk->x[j]!=0.0&&
rtk->ssat[i-1].outc[0]>GAP_RESION&&rtk->ssat[i-1].outc[1]>GAP_RESION)
rtk->x[j]=0.0;
}
for (i=0;i<ns;i++) {
j=II(sat[i],&rtk->opt);
if (rtk->x[j]==0.0) {
initx(rtk,1E-6,SQR(rtk->opt.std[1]*bl/1E4),j);
}
else {
/* elevation dependent factor of process noise */
el=rtk->ssat[sat[i]-1].azel[1];
fact=cos(el);
rtk->P[j+j*rtk->nx]+=SQR(rtk->opt.prn[1]*bl/1E4*fact)*tt;
}
}
}
/* temporal update of tropospheric parameters --------------------------------*/
static void udtrop(rtk_t *rtk, double tt, double bl)
{
int i,j,k;
trace(3,"udtrop : tt=%.1f\n",tt);
for (i=0;i<2;i++) {
j=IT(i,&rtk->opt);
if (rtk->x[j]==0.0) {
initx(rtk,INIT_ZWD,SQR(rtk->opt.std[2]),j); /* initial zwd */
if (rtk->opt.tropopt>=TROPOPT_ESTG) {
for (k=0;k<2;k++) initx(rtk,1E-6,VAR_GRA,++j);
}
}
else {
rtk->P[j+j*rtk->nx]+=SQR(rtk->opt.prn[2])*tt;
if (rtk->opt.tropopt>=TROPOPT_ESTG) {
for (k=0;k<2;k++) {
rtk->P[++j*(1+rtk->nx)]+=SQR(rtk->opt.prn[2]*0.3)*fabs(rtk->tt);
}
}
}
}
}
/* temporal update of receiver h/w biases ------------------------------------*/
static void udrcvbias(rtk_t *rtk, double tt)
{
int i,j;
trace(3,"udrcvbias: tt=%.1f\n",tt);
for (i=0;i<NFREQGLO;i++) {
j=IL(i,&rtk->opt);
if (rtk->x[j]==0.0) {
initx(rtk,1E-6,VAR_HWBIAS,j);
}
/* hold to fixed solution */
else if (rtk->nfix>=rtk->opt.minfix&&rtk->sol.ratio>rtk->opt.thresar[0]) {
initx(rtk,rtk->xa[j],rtk->Pa[j+j*rtk->na],j);
}
else {
rtk->P[j+j*rtk->nx]+=SQR(PRN_HWBIAS)*tt;
}
}
}
/* detect cycle slip by LLI --------------------------------------------------*/
static void detslp_ll(rtk_t *rtk, const obsd_t *obs, int i, int rcv)
{
unsigned char slip,LLI1,LLI2,LLI;
int f,sat=obs[i].sat;
trace(3,"detslp_ll: i=%d rcv=%d\n",i,rcv);
for (f=0;f<rtk->opt.nf;f++) {
if (obs[i].L[f]==0.0) continue;
/* restore previous LLI */
LLI1=(rtk->ssat[sat-1].slip[f]>>6)&3;
LLI2=(rtk->ssat[sat-1].slip[f]>>4)&3;
LLI=rcv==1?LLI1:LLI2;
/* detect slip by cycle slip flag */
slip=(rtk->ssat[sat-1].slip[f]|obs[i].LLI[f])&3;
if (obs[i].LLI[f]&1) {
errmsg(rtk,"slip detected (sat=%2d rcv=%d LLI%d=%x)\n",
sat,rcv,f+1,obs[i].LLI[f]);
}
/* detect slip by parity unknown flag transition */
if (((LLI&2)&&!(obs[i].LLI[f]&2))||(!(LLI&2)&&(obs[i].LLI[f]&2))) {
errmsg(rtk,"slip detected (sat=%2d rcv=%d LLI%d=%x->%x)\n",
sat,rcv,f+1,LLI,obs[i].LLI[f]);
slip|=1;
}
/* save current LLI and slip flag */
if (rcv==1) rtk->ssat[sat-1].slip[f]=(obs[i].LLI[f]<<6)|(LLI2<<4)|slip;
else rtk->ssat[sat-1].slip[f]=(obs[i].LLI[f]<<4)|(LLI1<<6)|slip;
}
}
/* detect cycle slip by L1-L2 geometry free phase jump -----------------------*/
static void detslp_gf_L1L2(rtk_t *rtk, const obsd_t *obs, int i, int j,
const nav_t *nav)
{
int sat=obs[i].sat;
double g0,g1;
trace(3,"detslp_gf_L1L2: i=%d j=%d\n",i,j);
if (rtk->opt.nf<=1||(g1=gfobs_L1L2(obs,i,j,nav->lam[sat-1]))==0.0) return;
g0=rtk->ssat[sat-1].gf; rtk->ssat[sat-1].gf=g1;
if (g0!=0.0&&fabs(g1-g0)>rtk->opt.thresslip) {
rtk->ssat[sat-1].slip[0]|=1;
rtk->ssat[sat-1].slip[1]|=1;
errmsg(rtk,"slip detected (sat=%2d GF_L1_L2=%.3f %.3f)\n",sat,g0,g1);
}
}
/* detect cycle slip by L1-L5 geometry free phase jump -----------------------*/
static void detslp_gf_L1L5(rtk_t *rtk, const obsd_t *obs, int i, int j,
const nav_t *nav)
{
int sat=obs[i].sat;
double g0,g1;
trace(3,"detslp_gf_L1L5: i=%d j=%d\n",i,j);
if (rtk->opt.nf<=2||(g1=gfobs_L1L5(obs,i,j,nav->lam[sat-1]))==0.0) return;
g0=rtk->ssat[sat-1].gf2; rtk->ssat[sat-1].gf2=g1;
if (g0!=0.0&&fabs(g1-g0)>rtk->opt.thresslip) {
rtk->ssat[sat-1].slip[0]|=1;
rtk->ssat[sat-1].slip[2]|=1;
errmsg(rtk,"slip detected (sat=%2d GF_L1_L5=%.3f %.3f)\n",sat,g0,g1);
}
}
/* detect cycle slip by doppler and phase difference -------------------------*/
static void detslp_dop(rtk_t *rtk, const obsd_t *obs, int i, int rcv,
const nav_t *nav)
{
/* detection with doppler disabled because of clock-jump issue (v.2.3.0) */
#if 0
int f,sat=obs[i].sat;
double tt,dph,dpt,lam,thres;
trace(3,"detslp_dop: i=%d rcv=%d\n",i,rcv);
for (f=0;f<rtk->opt.nf;f++) {
if (obs[i].L[f]==0.0||obs[i].D[f]==0.0||rtk->ph[rcv-1][sat-1][f]==0.0) {
continue;
}
if (fabs(tt=timediff(obs[i].time,rtk->pt[rcv-1][sat-1][f]))<DTTOL) continue;
if ((lam=nav->lam[sat-1][f])<=0.0) continue;
/* cycle slip threshold (cycle) */
thres=MAXACC*tt*tt/2.0/lam+rtk->opt.err[4]*fabs(tt)*4.0;
/* phase difference and doppler x time (cycle) */
dph=obs[i].L[f]-rtk->ph[rcv-1][sat-1][f];
dpt=-obs[i].D[f]*tt;
if (fabs(dph-dpt)<=thres) continue;
rtk->slip[sat-1][f]|=1;
errmsg(rtk,"slip detected (sat=%2d rcv=%d L%d=%.3f %.3f thres=%.3f)\n",
sat,rcv,f+1,dph,dpt,thres);
}
#endif
}
/* temporal update of phase biases -------------------------------------------*/
static void udbias(rtk_t *rtk, double tt, const obsd_t *obs, const int *sat,
const int *iu, const int *ir, int ns, const nav_t *nav)
{
double cp,pr,cp1,cp2,pr1,pr2,*bias,offset,lami,lam1,lam2,C1,C2;
int i,j,f,slip,reset,nf=NF(&rtk->opt);
trace(3,"udbias : tt=%.1f ns=%d\n",tt,ns);
for (i=0;i<ns;i++) {
/* detect cycle slip by LLI */
for (f=0;f<rtk->opt.nf;f++) rtk->ssat[sat[i]-1].slip[f]&=0xFC;
detslp_ll(rtk,obs,iu[i],1);
detslp_ll(rtk,obs,ir[i],2);
/* detect cycle slip by geometry-free phase jump */
detslp_gf_L1L2(rtk,obs,iu[i],ir[i],nav);
detslp_gf_L1L5(rtk,obs,iu[i],ir[i],nav);
/* detect cycle slip by doppler and phase difference */
detslp_dop(rtk,obs,iu[i],1,nav);
detslp_dop(rtk,obs,ir[i],2,nav);
}
for (f=0;f<nf;f++) {
/* reset phase-bias if instantaneous AR or expire obs outage counter */
for (i=1;i<=MAXSAT;i++) {
reset=++rtk->ssat[i-1].outc[f]>(unsigned int)rtk->opt.maxout;
if (rtk->opt.modear==ARMODE_INST&&rtk->x[IB(i,f,&rtk->opt)]!=0.0) {
initx(rtk,0.0,0.0,IB(i,f,&rtk->opt));
}
else if (reset&&rtk->x[IB(i,f,&rtk->opt)]!=0.0) {
initx(rtk,0.0,0.0,IB(i,f,&rtk->opt));
trace(3,"udbias : obs outage counter overflow (sat=%3d L%d n=%d)\n",
i,f+1,rtk->ssat[i-1].outc[f]);
}
if (rtk->opt.modear!=ARMODE_INST&&reset) {
rtk->ssat[i-1].lock[f]=-rtk->opt.minlock;
}
}
/* reset phase-bias if detecting cycle slip */
for (i=0;i<ns;i++) {
j=IB(sat[i],f,&rtk->opt);
rtk->P[j+j*rtk->nx]+=rtk->opt.prn[0]*rtk->opt.prn[0]*tt;
slip=rtk->ssat[sat[i]-1].slip[f];
if (rtk->opt.ionoopt==IONOOPT_IFLC) slip|=rtk->ssat[sat[i]-1].slip[1];
if (rtk->opt.modear==ARMODE_INST||!(slip&1)) continue;
rtk->x[j]=0.0;
rtk->ssat[sat[i]-1].lock[f]=-rtk->opt.minlock;
}
bias=zeros(ns,1);
/* estimate approximate phase-bias by phase - code */
for (i=j=0,offset=0.0;i<ns;i++) {
if (rtk->opt.ionoopt!=IONOOPT_IFLC) {
cp=sdobs(obs,iu[i],ir[i],f); /* cycle */
pr=sdobs(obs,iu[i],ir[i],f+NFREQ);
lami=nav->lam[sat[i]-1][f];
if (cp==0.0||pr==0.0||lami<=0.0) continue;
bias[i]=cp-pr/lami;
}
else {
cp1=sdobs(obs,iu[i],ir[i],0);
cp2=sdobs(obs,iu[i],ir[i],1);
pr1=sdobs(obs,iu[i],ir[i],NFREQ);
pr2=sdobs(obs,iu[i],ir[i],NFREQ+1);
lam1=nav->lam[sat[i]-1][0];
lam2=nav->lam[sat[i]-1][1];
if (cp1==0.0||cp2==0.0||pr1==0.0||pr2==0.0||lam1<=0.0||lam2<=0.0) continue;
C1= SQR(lam2)/(SQR(lam2)-SQR(lam1));
C2=-SQR(lam1)/(SQR(lam2)-SQR(lam1));
bias[i]=(C1*lam1*cp1+C2*lam2*cp2)-(C1*pr1+C2*pr2);
}
if (rtk->x[IB(sat[i],f,&rtk->opt)]!=0.0) {
offset+=bias[i]-rtk->x[IB(sat[i],f,&rtk->opt)];
j++;
}
}
/* correct phase-bias offset to enssure phase-code coherency */
if (j>0) {
for (i=1;i<=MAXSAT;i++) {
if (rtk->x[IB(i,f,&rtk->opt)]!=0.0) rtk->x[IB(i,f,&rtk->opt)]+=offset/j;
}
}
/* set initial states of phase-bias */
for (i=0;i<ns;i++) {
if (bias[i]==0.0||rtk->x[IB(sat[i],f,&rtk->opt)]!=0.0) continue;
initx(rtk,bias[i],SQR(rtk->opt.std[0]),IB(sat[i],f,&rtk->opt));
}
free(bias);
}
}
/* temporal update of states --------------------------------------------------*/
static void udstate(rtk_t *rtk, const obsd_t *obs, const int *sat,
const int *iu, const int *ir, int ns, const nav_t *nav)
{
double tt=fabs(rtk->tt),bl,dr[3];
trace(3,"udstate : ns=%d\n",ns);
/* temporal update of position/velocity/acceleration */
udpos(rtk,tt);
/* temporal update of ionospheric parameters */
if (rtk->opt.ionoopt>=IONOOPT_EST) {
bl=baseline(rtk->x,rtk->rb,dr);
udion(rtk,tt,bl,sat,ns);
}
/* temporal update of tropospheric parameters */
if (rtk->opt.tropopt>=TROPOPT_EST) {
udtrop(rtk,tt,bl);
}
/* temporal update of eceiver h/w bias */
if (rtk->opt.glomodear==2&&(rtk->opt.navsys&SYS_GLO)) {
udrcvbias(rtk,tt);
}
/* temporal update of phase-bias */
if (rtk->opt.mode>PMODE_DGPS) {
udbias(rtk,tt,obs,sat,iu,ir,ns,nav);
}
}
/* undifferenced phase/code residual for satellite ---------------------------*/
static void zdres_sat(int base, double r, const obsd_t *obs, const nav_t *nav,
const double *azel, const double *dant,
const prcopt_t *opt, double *y)
{
const double *lam=nav->lam[obs->sat-1];
double f1,f2,C1,C2,dant_if;
int i,nf=NF(opt);
if (opt->ionoopt==IONOOPT_IFLC) { /* iono-free linear combination */
if (lam[0]==0.0||lam[1]==0.0) return;
if (testsnr(base,0,azel[1],obs->SNR[0]*0.25,&opt->snrmask)||
testsnr(base,1,azel[1],obs->SNR[1]*0.25,&opt->snrmask)) return;
f1=CLIGHT/lam[0];
f2=CLIGHT/lam[1];
C1= SQR(f1)/(SQR(f1)-SQR(f2));
C2=-SQR(f2)/(SQR(f1)-SQR(f2));
dant_if=C1*dant[0]+C2*dant[1];
if (obs->L[0]!=0.0&&obs->L[1]!=0.0) {
y[0]=C1*obs->L[0]*lam[0]+C2*obs->L[1]*lam[1]-r-dant_if;
}
if (obs->P[0]!=0.0&&obs->P[1]!=0.0) {
y[1]=C1*obs->P[0]+C2*obs->P[1]-r-dant_if;
}
}
else {
for (i=0;i<nf;i++) {
if (lam[i]==0.0) continue;
/* check snr mask */
if (testsnr(base,i,azel[1],obs->SNR[i]*0.25,&opt->snrmask)) {
continue;
}
/* residuals = observable - pseudorange */
if (obs->L[i]!=0.0) y[i ]=obs->L[i]*lam[i]-r-dant[i];
if (obs->P[i]!=0.0) y[i+nf]=obs->P[i] -r-dant[i];
}
}
}
/* undifferenced phase/code residuals ----------------------------------------*/
static int zdres(int base, const obsd_t *obs, int n, const double *rs,
const double *dts, const int *svh, const nav_t *nav,
const double *rr, const prcopt_t *opt, int index, double *y,
double *e, double *azel)
{
double r,rr_[3],pos[3],dant[NFREQ]={0},disp[3];
double zhd,zazel[]={0.0,90.0*D2R};
int i,nf=NF(opt);
trace(3,"zdres : n=%d\n",n);
for (i=0;i<n*nf*2;i++) y[i]=0.0;
if (norm(rr,3)<=0.0) return 0; /* no receiver position */
for (i=0;i<3;i++) rr_[i]=rr[i];
/* earth tide correction */
if (opt->tidecorr) {
tidedisp(gpst2utc(obs[0].time),rr_,opt->tidecorr,&nav->erp,
opt->odisp[base],disp);
for (i=0;i<3;i++) rr_[i]+=disp[i];
}
ecef2pos(rr_,pos);
for (i=0;i<n;i++) {
/* compute geometric-range and azimuth/elevation angle */
if ((r=geodist(rs+i*6,rr_,e+i*3))<=0.0) continue;
if (satazel(pos,e+i*3,azel+i*2)<opt->elmin) continue;
/* excluded satellite? */
if (satexclude(obs[i].sat,svh[i],opt)) continue;
/* satellite clock-bias */
r+=-CLIGHT*dts[i*2];
/* troposphere delay model (hydrostatic) */
zhd=tropmodel(obs[0].time,pos,zazel,0.0);
r+=tropmapf(obs[i].time,pos,azel+i*2,NULL)*zhd;
/* receiver antenna phase center correction */
antmodel(opt->pcvr+index,opt->antdel[index],azel+i*2,opt->posopt[1],
dant);
/* undifferenced phase/code residual for satellite */
zdres_sat(base,r,obs+i,nav,azel+i*2,dant,opt,y+i*nf*2);
}
trace(4,"rr_=%.3f %.3f %.3f\n",rr_[0],rr_[1],rr_[2]);
trace(4,"pos=%.9f %.9f %.3f\n",pos[0]*R2D,pos[1]*R2D,pos[2]);
for (i=0;i<n;i++) {
trace(4,"sat=%2d %13.3f %13.3f %13.3f %13.10f %6.1f %5.1f\n",
obs[i].sat,rs[i*6],rs[1+i*6],rs[2+i*6],dts[i*2],azel[i*2]*R2D,
azel[1+i*2]*R2D);
}
trace(4,"y=\n"); tracemat(4,y,nf*2,n,13,3);
return 1;
}
/* test valid observation data -----------------------------------------------*/
static int validobs(int i, int j, int f, int nf, double *y)
{
/* if no phase observable, psudorange is also unusable */
return y[f+i*nf*2]!=0.0&&y[f+j*nf*2]!=0.0&&
(f<nf||(y[f-nf+i*nf*2]!=0.0&&y[f-nf+j*nf*2]!=0.0));
}
/* double-differenced measurement error covariance ---------------------------*/
static void ddcov(const int *nb, int n, const double *Ri, const double *Rj,
int nv, double *R)
{
int i,j,k=0,b;
trace(3,"ddcov : n=%d\n",n);
for (i=0;i<nv*nv;i++) R[i]=0.0;
for (b=0;b<n;k+=nb[b++]) {
for (i=0;i<nb[b];i++) for (j=0;j<nb[b];j++) {
R[k+i+(k+j)*nv]=Ri[k+i]+(i==j?Rj[k+i]:0.0);
}
}
trace(5,"R=\n"); tracemat(5,R,nv,nv,8,6);
}
/* baseline length constraint ------------------------------------------------*/
static int constbl(rtk_t *rtk, const double *x, const double *P, double *v,
double *H, double *Ri, double *Rj, int index)
{
const double thres=0.1; /* threshold for nonliearity (v.2.3.0) */
double xb[3],b[3],bb,var=0.0;
int i;
trace(3,"constbl : \n");
/* no constraint */
if (rtk->opt.baseline[0]<=0.0) return 0;
/* time-adjusted baseline vector and length */
for (i=0;i<3;i++) {
xb[i]=rtk->rb[i]+rtk->rb[i+3]*rtk->sol.age;
b[i]=x[i]-xb[i];
}
bb=norm(b,3);
/* approximate variance of solution */
if (P) {
for (i=0;i<3;i++) var+=P[i+i*rtk->nx];
var/=3.0;
}
/* check nonlinearity */
if (var>thres*thres*bb*bb) {
trace(3,"constbl : equation nonlinear (bb=%.3f var=%.3f)\n",bb,var);
return 0;
}
/* constraint to baseline length */
v[index]=rtk->opt.baseline[0]-bb;
if (H) {
for (i=0;i<3;i++) H[i+index*rtk->nx]=b[i]/bb;
}
Ri[index]=0.0;
Rj[index]=SQR(rtk->opt.baseline[1]);
trace(4,"baseline len v=%13.3f R=%8.6f %8.6f\n",v[index],Ri[index],Rj[index]);
return 1;
}
/* precise tropspheric model -------------------------------------------------*/
static double prectrop(gtime_t time, const double *pos, int r,
const double *azel, const prcopt_t *opt, const double *x,
double *dtdx)
{
double m_w=0.0,cotz,grad_n,grad_e;
int i=IT(r,opt);
/* wet mapping function */
tropmapf(time,pos,azel,&m_w);
if (opt->tropopt>=TROPOPT_ESTG&&azel[1]>0.0) {
/* m_w=m_0+m_0*cot(el)*(Gn*cos(az)+Ge*sin(az)): ref [6] */
cotz=1.0/tan(azel[1]);
grad_n=m_w*cotz*cos(azel[0]);
grad_e=m_w*cotz*sin(azel[0]);
m_w+=grad_n*x[i+1]+grad_e*x[i+2];
dtdx[1]=grad_n*x[i];
dtdx[2]=grad_e*x[i];
}
else dtdx[1]=dtdx[2]=0.0;
dtdx[0]=m_w;
return m_w*x[i];
}
/* glonass inter-channel bias correction -------------------------------------*/
static double gloicbcorr(int sat1, int sat2, const prcopt_t *opt, double lam1,
double lam2, int f)
{
double dfreq;