-
Notifications
You must be signed in to change notification settings - Fork 0
/
bfgs_funcs.cpp
279 lines (240 loc) · 11 KB
/
bfgs_funcs.cpp
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
/*************************************************************
Generalized-ICP Copyright (c) 2009 Aleksandr Segal.
All rights reserved.
Redistribution and use in source and binary forms, with
or without modification, are permitted provided that the
following conditions are met:
* Redistributions of source code must retain the above
copyright notice, this list of conditions and the
following disclaimer.
* Redistributions in binary form must reproduce the above
copyright notice, this list of conditions and the
following disclaimer in the documentation and/or other
materials provided with the distribution.
* The names of the contributors may not be used to endorse
or promote products derived from this software
without specific prior written permission.
THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND
CONTRIBUTORS "AS IS" AND ANY EXPRESS OR IMPLIED
WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A
PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT,
INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE
OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH
DAMAGE.
*************************************************************/
#include "optimize.h"
namespace dgc {
namespace gicp {
// GICP cost function
double GICPOptimizer::f(const gsl_vector *x, void *params) {
GICPOptData *opt_data = (GICPOptData *)params;
double pt1[3];
double pt2[3];
double res[3]; // residual
double temp[3];
gsl_vector_view gsl_pt1 = gsl_vector_view_array(pt1, 3);
gsl_vector_view gsl_pt2 = gsl_vector_view_array(pt2, 3);
gsl_vector_view gsl_res = gsl_vector_view_array(res, 3);
gsl_vector_view gsl_temp = gsl_vector_view_array(temp, 3);
gsl_matrix_view gsl_M;
dgc_transform_t t;
// initialize the temp variable; if it happens to be NaN at start, bad things will happen in blas routines below
temp[0] = 0;
temp[1] = 0;
temp[2] = 0;
// take the base transformation
dgc_transform_copy(t, opt_data->base_t);
// apply the current state
apply_state(t, x);
double f = 0;
double temp_double = 0;
int N = opt_data->p1->Size();
for(int i = 0; i < N; i++) {
int j = opt_data->nn_indecies[i];
if(j != -1) {
// get point 1
pt1[0] = (*opt_data->p1)[i].x;
pt1[1] = (*opt_data->p1)[i].y;
pt1[2] = (*opt_data->p1)[i].z;
// get point 2
pt2[0] = (*opt_data->p2)[j].x;
pt2[1] = (*opt_data->p2)[j].y;
pt2[2] = (*opt_data->p2)[j].z;
//get M-matrix
gsl_M = gsl_matrix_view_array(&opt_data->M[i][0][0], 3, 3);
//transform point 1
dgc_transform_point(&pt1[0], &pt1[1], &pt1[2], t);
res[0] = pt1[0] - pt2[0];
res[1] = pt1[1] - pt2[1];
res[2] = pt1[2] - pt2[2];
//cout << "res: (" << res[0] << ", " <<res[1] <<", " << res[2] << ")" << endl;
// temp := M*res
gsl_blas_dsymv(CblasLower, 1., &gsl_M.matrix, &gsl_res.vector, 0., &gsl_temp.vector);
// temp_double := res'*temp = temp'*M*temp
gsl_blas_ddot(&gsl_res.vector, &gsl_temp.vector, &temp_double);
// increment total error
f += temp_double/(double)opt_data->num_matches;
//cout << "temp: " << temp_double << endl;
//cout << "f: " << f << "\t (" << opt_data->num_matches << ")" << endl;
//print_gsl_matrix(&gsl_M.matrix, "M");
}
}
return f;
}
void GICPOptimizer::df(const gsl_vector *x, void *params, gsl_vector *g) {
GICPOptData *opt_data = (GICPOptData *)params;
double pt1[3];
double pt2[3];
double res[3]; // residual
double temp[3]; // temp local vector
double temp_mat[9]; // temp matrix used for accumulating the rotation gradient
gsl_vector_view gsl_pt1 = gsl_vector_view_array(pt1, 3);
gsl_vector_view gsl_pt2 = gsl_vector_view_array(pt2, 3);
gsl_vector_view gsl_res = gsl_vector_view_array(res, 3);
gsl_vector_view gsl_temp = gsl_vector_view_array(temp, 3);
gsl_vector_view gsl_gradient_t = gsl_vector_subvector(g, 0, 3); // translation comp. of gradient
gsl_matrix_view gsl_temp_mat_r = gsl_matrix_view_array(temp_mat, 3, 3);
gsl_matrix_view gsl_M;
dgc_transform_t t;
double temp_double;
// take the base transformation
dgc_transform_copy(t, opt_data->base_t);
// apply the current state
apply_state(t, x);
// zero all accumulator variables
gsl_vector_set_zero(g);
gsl_vector_set_zero(&gsl_temp.vector);
gsl_matrix_set_zero(&gsl_temp_mat_r.matrix);
for(int i = 0; i < opt_data->p1->Size(); i++) {
int j = opt_data->nn_indecies[i];
if(j != -1) {
// get point 1
pt1[0] = (*opt_data->p1)[i].x;
pt1[1] = (*opt_data->p1)[i].y;
pt1[2] = (*opt_data->p1)[i].z;
// get point 2
pt2[0] = (*opt_data->p2)[j].x;
pt2[1] = (*opt_data->p2)[j].y;
pt2[2] = (*opt_data->p2)[j].z;
//get M-matrix
gsl_M = gsl_matrix_view_array(&opt_data->M[i][0][0], 3, 3);
//transform point 1
dgc_transform_point(&pt1[0], &pt1[1], &pt1[2], t);
res[0] = pt1[0] - pt2[0];
res[1] = pt1[1] - pt2[1];
res[2] = pt1[2] - pt2[2];
// temp := M*res
gsl_blas_dsymv(CblasLower, 1., &gsl_M.matrix, &gsl_res.vector, 0., &gsl_temp.vector);
// temp_double := res'*temp = temp'*M*res
gsl_blas_ddot(&gsl_res.vector, &gsl_temp.vector, &temp_double);
// increment total translation gradient:
// gsl_gradient_t += 2*M*res/num_matches
gsl_blas_dsymv(CblasLower, 2./(double)opt_data->num_matches, &gsl_M.matrix, &gsl_res.vector, 1., &gsl_gradient_t.vector);
if(opt_data->solve_rotation) {
// compute rotation gradient here
// get back the original untransformed point to compute the rotation gradient
pt1[0] = (*opt_data->p1)[i].x;
pt1[1] = (*opt_data->p1)[i].y;
pt1[2] = (*opt_data->p1)[i].z;
dgc_transform_point(&pt1[0], &pt1[1], &pt1[2], opt_data->base_t);
gsl_blas_dger(2./(double)opt_data->num_matches, &gsl_pt1.vector, &gsl_temp.vector, &gsl_temp_mat_r.matrix);
}
}
}
// the above loop sets up the gradient with respect to the translation, and the matrix derivative w.r.t. the rotation matrix
// this code sets up the matrix derivatives dR/dPhi, dR/dPsi, dR/dTheta. i.e. the derivatives of the whole rotation matrix with respect to the euler angles
// note that this code assumes the XYZ order of euler angles, with the Z rotation corresponding to bearing. This means the Z angle is negative of what it would be
// in the regular XYZ euler-angle convention.
// now use the d/dR matrix to compute the derivative with respect to euler angles and put it directly into g[3], g[4], g[5];
if(opt_data->solve_rotation) {
compute_dr(x, &gsl_temp_mat_r.matrix, g);
}
}
void GICPOptimizer::fdf(const gsl_vector *x, void *params, double * f, gsl_vector *g) {
GICPOptData *opt_data = (GICPOptData *)params;
double pt1[3];
double pt2[3];
double res[3]; // residual
double temp[3]; // temp local vector
double temp_mat[9]; // temp matrix used for accumulating the rotation gradient
gsl_vector_view gsl_pt1 = gsl_vector_view_array(pt1, 3);
gsl_vector_view gsl_pt2 = gsl_vector_view_array(pt2, 3);
gsl_vector_view gsl_res = gsl_vector_view_array(res, 3);
gsl_vector_view gsl_temp = gsl_vector_view_array(temp, 3);
gsl_vector_view gsl_gradient_t = gsl_vector_subvector(g, 0, 3); // translation comp. of gradient
gsl_vector_view gsl_gradient_r = gsl_vector_subvector(g, 3, 3); // rotation comp. of gradient
gsl_matrix_view gsl_temp_mat_r = gsl_matrix_view_array(temp_mat, 3, 3);
gsl_matrix_view gsl_M;
dgc_transform_t t;
double temp_double;
// take the base transformation
dgc_transform_copy(t, opt_data->base_t);
// apply the current state
apply_state(t, x);
// zero all accumulator variables
*f = 0;
gsl_vector_set_zero(g);
gsl_vector_set_zero(&gsl_temp.vector);
gsl_matrix_set_zero(&gsl_temp_mat_r.matrix);
for(int i = 0; i < opt_data->p1->Size(); i++) {
int j = opt_data->nn_indecies[i];
if(j != -1) {
// get point 1
pt1[0] = (*opt_data->p1)[i].x;
pt1[1] = (*opt_data->p1)[i].y;
pt1[2] = (*opt_data->p1)[i].z;
// get point 2
pt2[0] = (*opt_data->p2)[j].x;
pt2[1] = (*opt_data->p2)[j].y;
pt2[2] = (*opt_data->p2)[j].z;
//cout << "accessing " << i << " of " << opt_data->p1->Size() << ", " << opt_data->p2->Size() << endl;
//get M-matrix
gsl_M = gsl_matrix_view_array(&opt_data->M[i][0][0], 3, 3);
//transform point 1
dgc_transform_point(&pt1[0], &pt1[1], &pt1[2], t);
res[0] = pt1[0] - pt2[0];
res[1] = pt1[1] - pt2[1];
res[2] = pt1[2] - pt2[2];
// compute the transformed residual
// temp := M*res
//print_gsl_matrix(&gsl_M.matrix, "gsl_m");
gsl_blas_dsymv(CblasLower, 1., &gsl_M.matrix, &gsl_res.vector, 0., &gsl_temp.vector);
// compute M-norm of the residual
// temp_double := res'*temp = temp'*M*res
gsl_blas_ddot(&gsl_res.vector, &gsl_temp.vector, &temp_double);
// accumulate total error: f += res'*M*res
*f += temp_double/(double)opt_data->num_matches;
// accumulate translation gradient:
// gsl_gradient_t += 2*M*res
gsl_blas_dsymv(CblasLower, 2./(double)opt_data->num_matches, &gsl_M.matrix, &gsl_res.vector, 1., &gsl_gradient_t.vector);
if(opt_data->solve_rotation) {
// accumulate the rotation gradient matrix
// get back the original untransformed point to compute the rotation gradient
pt1[0] = (*opt_data->p1)[i].x;
pt1[1] = (*opt_data->p1)[i].y;
pt1[2] = (*opt_data->p1)[i].z;
dgc_transform_point(&pt1[0], &pt1[1], &pt1[2], opt_data->base_t);
// gsl_temp_mat_r += 2*(gsl_temp).(gsl_pt1)' [ = (2*M*residual).(gsl_pt1)' ]
gsl_blas_dger(2./(double)opt_data->num_matches, &gsl_pt1.vector, &gsl_temp.vector, &gsl_temp_mat_r.matrix);
}
}
}
// the above loop sets up the gradient with respect to the translation, and the matrix derivative w.r.t. the rotation matrix
// this code sets up the matrix derivatives dR/dPhi, dR/dPsi, dR/dTheta. i.e. the derivatives of the whole rotation matrix with respect to the euler angles
// note that this code assumes the XYZ order of euler angles, with the Z rotation corresponding to bearing. This means the Z angle is negative of what it would be
// in the regular XYZ euler-angle convention.
if(opt_data->solve_rotation) {
// now use the d/dR matrix to compute the derivative with respect to euler angles and put it directly into g[3], g[4], g[5];
compute_dr(x, &gsl_temp_mat_r.matrix, g);
}
}
}
}