forked from ndhuan/GPSRTK
-
Notifications
You must be signed in to change notification settings - Fork 0
/
lambda.c
190 lines (175 loc) · 6.36 KB
/
lambda.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
/*------------------------------------------------------------------------------
* lambda.c : integer ambiguity resolution
*
* Copyright (C) 2007-2008 by T.TAKASU, All rights reserved.
*
* reference :
* [1] P.J.G.Teunissen, The least-square ambiguity decorrelation adjustment:
* a method for fast GPS ambiguity estimation, J.Geodesy, Vol.70, 65-82,
* 1995
* [2] X.-W.Chang, X.Yang, T.Zhou, MLAMBDA: A modified LAMBDA method for
* integer least-squares estimation, J.Geodesy, Vol.79, 552-565, 2005
*
* version : $Revision: 1.1 $ $Date: 2008/07/17 21:48:06 $
* history : 2007/01/13 1.0 new
*-----------------------------------------------------------------------------*/
#include "rtk.h"
/* constants/macros ----------------------------------------------------------*/
#define LOOPMAX 10000 /* maximum count of search loop */
#define SGN(x) ((x)<=0.0?-1.0:1.0)
#define ROUND(x) (floor((x)+0.5))
#define SWAP(x,y) do {double tmp_; tmp_=x; x=y; y=tmp_;} while (0)
/* LD factorization (Q=L'*diag(D)*L) -----------------------------------------*/
static int LD(int n, const double *Q, double *L, double *D)
{
int i,j,k,info=0;
double a,*A=mat(n,n);
memcpy(A,Q,sizeof(double)*n*n);
for (i=n-1;i>=0;i--) {
if ((D[i]=A[i+i*n])<=0.0) {info=-1; break;}
a=sqrt(D[i]);
for (j=0;j<=i;j++) L[i+j*n]=A[i+j*n]/a;
for (j=0;j<=i-1;j++) for (k=0;k<=j;k++) A[j+k*n]-=L[i+k*n]*L[i+j*n];
for (j=0;j<=i;j++) L[i+j*n]/=L[i+i*n];
}
free(A);
//if (info) fprintf(stderr,"%s : LD factorization error\n",__FILE__);
return info;
}
/* integer gauss transformation ----------------------------------------------*/
static void gauss(int n, double *L, double *Z, int i, int j)
{
int k,mu;
if ((mu=(int)ROUND(L[i+j*n]))!=0) {
for (k=i;k<n;k++) L[k+n*j]-=(double)mu*L[k+i*n];
for (k=0;k<n;k++) Z[k+n*j]-=(double)mu*Z[k+i*n];
}
}
/* permutations --------------------------------------------------------------*/
static void perm(int n, double *L, double *D, int j, double del, double *Z)
{
int k;
double eta,lam,a0,a1;
eta=D[j]/del;
lam=D[j+1]*L[j+1+j*n]/del;
D[j]=eta*D[j+1]; D[j+1]=del;
for (k=0;k<=j-1;k++) {
a0=L[j+k*n]; a1=L[j+1+k*n];
L[j+k*n]=-L[j+1+j*n]*a0+a1;
L[j+1+k*n]=eta*a0+lam*a1;
}
L[j+1+j*n]=lam;
for (k=j+2;k<n;k++) SWAP(L[k+j*n],L[k+(j+1)*n]);
for (k=0;k<n;k++) SWAP(Z[k+j*n],Z[k+(j+1)*n]);
}
/* lambda reduction (z=Z'*a, Qz=Z'*Q*Z=L'*diag(D)*L) (ref.[1]) ---------------*/
static void reduction(int n, double *L, double *D, double *Z)
{
int i,j,k;
double del;
j=n-2; k=n-2;
while (j>=0) {
if (j<=k) for (i=j+1;i<n;i++) gauss(n,L,Z,i,j);
del=D[j]+L[j+1+j*n]*L[j+1+j*n]*D[j+1];
if (del+1E-6<D[j+1]) { /* compared considering numerical error */
perm(n,L,D,j,del,Z);
k=j; j=n-2;
}
else j--;
}
}
/* modified lambda (mlambda) search (ref. [2]) -------------------------------*/
static int search(int n, int m, const double *L, const double *D,
const double *zs, double *zn, double *s)
{
int i,j,k,c,nn=0,imax=0;
double newdist,maxdist=1E99,y;
double *S=zeros(n,n),*dist=mat(n,1),*zb=mat(n,1),*z=mat(n,1),*step=mat(n,1);
k=n-1; dist[k]=0.0;
zb[k]=zs[k];
z[k]=ROUND(zb[k]); y=zb[k]-z[k]; step[k]=SGN(y);
for (c=0;c<LOOPMAX;c++) {
newdist=dist[k]+y*y/D[k];
if (newdist<maxdist) {
if (k!=0) {
dist[--k]=newdist;
for (i=0;i<=k;i++)
S[k+i*n]=S[k+1+i*n]+(z[k+1]-zb[k+1])*L[k+1+i*n];
zb[k]=zs[k]+S[k+k*n];
z[k]=ROUND(zb[k]); y=zb[k]-z[k]; step[k]=SGN(y);
}
else {
if (nn<m) {
if (nn==0||newdist>s[imax]) imax=nn;
for (i=0;i<n;i++) zn[i+nn*n]=z[i];
s[nn++]=newdist;
}
else {
if (newdist<s[imax]) {
for (i=0;i<n;i++) zn[i+imax*n]=z[i];
s[imax]=newdist;
for (i=imax=0;i<m;i++) if (s[imax]<s[i]) imax=i;
}
maxdist=s[imax];
}
z[0]+=step[0]; y=zb[0]-z[0]; step[0]=-step[0]-SGN(step[0]);
}
}
else {
if (k==n-1) break;
else {
k++;
z[k]+=step[k]; y=zb[k]-z[k]; step[k]=-step[k]-SGN(step[k]);
}
}
}
for (i=0;i<m-1;i++) { /* sort by s */
for (j=i+1;j<m;j++) {
if (s[i]<s[j]) continue;
SWAP(s[i],s[j]);
for (k=0;k<n;k++) SWAP(zn[k+i*n],zn[k+j*n]);
}
}
free(S); free(dist); free(zb); free(z); free(step);
if (c>=LOOPMAX) {
//fprintf(stderr,"%s : search loop count overflow\n",__FILE__);
return -1;
}
return 0;
}
/* lambda/mlambda integer least-square estimation ------------------------------
* integer least-square estimation. reduction is performed by lambda (ref.[1]),
* and search by mlambda (ref.[2]).
* args : int n I number of float parameters
* int m I number of fixed solutions
* double *a I float parameters (n x 1)
* double *Q I covariance matrix of float parameters (n x n)
* double *F O fixed solutions (n x m)
* double *s O sum of squared residulas of fixed solutions (1 x m)
* return : status (0:ok,other:error)
* notes : matrix stored by column-major order (fortran convension)
*-----------------------------------------------------------------------------*/
extern int lambda(int n, int m, const double *a, const double *Q, double *f,
double *s)
{
int info;
double *L,*D,*Z,*z,*E;
if (n<=0||m<=0) return -1;
L=zeros(n,n); D=mat(n,1); Z=eye(n); z=mat(n,1),E=mat(n,m);
/* LD factorization */
if (!(LD(n,Q,L,D))) {
/* lambda reduction */
reduction(n,L,D,Z);
matmul("TN",n,1,n,1.0,Z,a,0.0,z); /* z=Z'*a */
/* mlambda search */
if (!(search(n,m,L,D,z,E,s))) {
info=solve("T",Z,E,n,m,f); /* F=Z'\E */
}
else
info = -2;
}
else
info=-3;
free(L); free(D); free(Z); free(z); free(E);
return info;
}