[37] | 1 | ///////////////////////////////////////////////////////////////////////////////// |
---|
| 2 | //// |
---|
| 3 | //// Simple drivers for sparse bundle adjustment based on the |
---|
| 4 | //// Levenberg - Marquardt minimization algorithm |
---|
| 5 | //// This file provides simple wrappers to the functions defined in sba_levmar.c |
---|
| 6 | //// Copyright (C) 2004 Manolis Lourakis (lourakis@ics.forth.gr) |
---|
| 7 | //// Institute of Computer Science, Foundation for Research & Technology - Hellas |
---|
| 8 | //// Heraklion, Crete, Greece. |
---|
| 9 | //// |
---|
| 10 | //// This program is free software; you can redistribute it and/or modify |
---|
| 11 | //// it under the terms of the GNU General Public License as published by |
---|
| 12 | //// the Free Software Foundation; either version 2 of the License, or |
---|
| 13 | //// (at your option) any later version. |
---|
| 14 | //// |
---|
| 15 | //// This program is distributed in the hope that it will be useful, |
---|
| 16 | //// but WITHOUT ANY WARRANTY; without even the implied warranty of |
---|
| 17 | //// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
---|
| 18 | //// GNU General Public License for more details. |
---|
| 19 | //// |
---|
| 20 | /////////////////////////////////////////////////////////////////////////////////// |
---|
| 21 | |
---|
| 22 | #include <stdio.h> |
---|
| 23 | #include <stdlib.h> |
---|
| 24 | #include <string.h> |
---|
| 25 | #include <math.h> |
---|
| 26 | #include <float.h> |
---|
| 27 | |
---|
| 28 | #include "sba.h" |
---|
| 29 | |
---|
| 30 | |
---|
| 31 | #define FABS(x) (((x)>=0)? (x) : -(x)) |
---|
| 32 | |
---|
| 33 | struct wrap_motstr_data_ { |
---|
| 34 | void (*proj)(int j, int i, double *aj, double *bi, double *xij, void *adata); // Q |
---|
| 35 | void (*projac)(int j, int i, double *aj, double *bi, double *Aij, double *Bij, void *adata); // dQ/da, dQ/db |
---|
| 36 | int cnp, pnp, mnp; /* parameter numbers */ |
---|
| 37 | void *adata; |
---|
| 38 | }; |
---|
| 39 | |
---|
| 40 | struct wrap_mot_data_ { |
---|
| 41 | void (*proj)(int j, int i, double *aj, double *xij, void *adata); // Q |
---|
| 42 | void (*projac)(int j, int i, double *aj, double *Aij, void *adata); // dQ/da |
---|
| 43 | int cnp, mnp; /* parameter numbers */ |
---|
| 44 | void *adata; |
---|
| 45 | }; |
---|
| 46 | |
---|
| 47 | struct wrap_str_data_ { |
---|
| 48 | void (*proj)(int j, int i, double *bi, double *xij, void *adata); // Q |
---|
| 49 | void (*projac)(int j, int i, double *bi, double *Bij, void *adata); // dQ/db |
---|
| 50 | int pnp, mnp; /* parameter numbers */ |
---|
| 51 | void *adata; |
---|
| 52 | }; |
---|
| 53 | |
---|
| 54 | /* Routines to estimate the estimated measurement vector (i.e. "func") and |
---|
| 55 | * its sparse jacobian (i.e. "fjac") needed by BA expert drivers. Code below |
---|
| 56 | * makes use of user-supplied functions computing "Q", "dQ/da", d"Q/db", |
---|
| 57 | * i.e. predicted projection and associated jacobians for a SINGLE image measurement. |
---|
| 58 | * Notice also that what follows is two pairs of "func" and corresponding "fjac" routines. |
---|
| 59 | * The first is to be used in full (i.e. motion + structure) BA, the second in |
---|
| 60 | * motion only BA. |
---|
| 61 | */ |
---|
| 62 | |
---|
| 63 | /* FULL BUNDLE ADJUSTMENT */ |
---|
| 64 | |
---|
| 65 | /* Given a parameter vector p made up of the 3D coordinates of n points and the parameters of m cameras, compute in |
---|
| 66 | * hx the prediction of the measurements, i.e. the projections of 3D points in the m images. The measurements |
---|
| 67 | * are returned in the order (hx_11^T, .. hx_1m^T, ..., hx_n1^T, .. hx_nm^T)^T, where hx_ij is the predicted |
---|
| 68 | * projection of the i-th point on the j-th camera. |
---|
| 69 | * Caller supplies rcidxs and rcsubs which can be used as working memory. |
---|
| 70 | * Notice that depending on idxij, some of the hx_ij might be missing |
---|
| 71 | * |
---|
| 72 | */ |
---|
| 73 | static void sba_motstr_Qs(double *p, struct sba_crsm *idxij, int *rcidxs, int *rcsubs, double *hx, void *adata) |
---|
| 74 | { |
---|
| 75 | register int i, j; |
---|
| 76 | int cnp, pnp, mnp; |
---|
| 77 | double *pa, *pb, *paj, *pbi, *pxij; |
---|
| 78 | int n, m, nnz; |
---|
| 79 | struct wrap_motstr_data_ *wdata; |
---|
| 80 | void (*proj)(int j, int i, double *aj, double *bi, double *xij, void *proj_adata); |
---|
| 81 | void *proj_adata; |
---|
| 82 | |
---|
| 83 | wdata=(struct wrap_motstr_data_ *)adata; |
---|
| 84 | cnp=wdata->cnp; pnp=wdata->pnp; mnp=wdata->mnp; |
---|
| 85 | proj=wdata->proj; |
---|
| 86 | proj_adata=wdata->adata; |
---|
| 87 | |
---|
| 88 | n=idxij->nr; m=idxij->nc; |
---|
| 89 | pa=p; pb=p+m*cnp; |
---|
| 90 | |
---|
| 91 | for(j=0; j<m; ++j){ |
---|
| 92 | /* j-th camera parameters */ |
---|
| 93 | paj=pa+j*cnp; |
---|
| 94 | |
---|
| 95 | nnz=sba_crsm_col_elmidxs(idxij, j, rcidxs, rcsubs); /* find nonzero hx_ij, i=0...n-1 */ |
---|
| 96 | |
---|
| 97 | for(i=0; i<nnz; ++i){ |
---|
| 98 | pbi=pb + rcsubs[i]*pnp; |
---|
| 99 | pxij=hx + idxij->val[rcidxs[i]]*mnp; // set pxij to point to hx_ij |
---|
| 100 | |
---|
| 101 | (*proj)(j, rcsubs[i], paj, pbi, pxij, proj_adata); // evaluate Q in pxij |
---|
| 102 | } |
---|
| 103 | } |
---|
| 104 | } |
---|
| 105 | |
---|
| 106 | /* Given a parameter vector p made up of the 3D coordinates of n points and the parameters of m cameras, compute in |
---|
| 107 | * jac the jacobian of the predicted measurements, i.e. the jacobian of the projections of 3D points in the m images. |
---|
| 108 | * The jacobian is returned in the order (A_11, ..., A_1m, ..., A_n1, ..., A_nm, B_11, ..., B_1m, ..., B_n1, ..., B_nm), |
---|
| 109 | * where A_ij=dx_ij/db_j and B_ij=dx_ij/db_i (see HZ). |
---|
| 110 | * Caller supplies rcidxs and rcsubs which can be used as working memory. |
---|
| 111 | * Notice that depending on idxij, some of the A_ij, B_ij might be missing |
---|
| 112 | * |
---|
| 113 | */ |
---|
| 114 | static void sba_motstr_Qs_jac(double *p, struct sba_crsm *idxij, int *rcidxs, int *rcsubs, double *jac, void *adata) |
---|
| 115 | { |
---|
| 116 | register int i, j; |
---|
| 117 | int cnp, pnp, mnp; |
---|
| 118 | double *pa, *pb, *paj, *pbi, *jaca, *jacb, *pAij, *pBij; |
---|
| 119 | int n, m, nnz, Asz, Bsz, idx; |
---|
| 120 | struct wrap_motstr_data_ *wdata; |
---|
| 121 | void (*projac)(int j, int i, double *aj, double *bi, double *Aij, double *Bij, void *projac_adata); |
---|
| 122 | void *projac_adata; |
---|
| 123 | |
---|
| 124 | |
---|
| 125 | wdata=(struct wrap_motstr_data_ *)adata; |
---|
| 126 | cnp=wdata->cnp; pnp=wdata->pnp; mnp=wdata->mnp; |
---|
| 127 | projac=wdata->projac; |
---|
| 128 | projac_adata=wdata->adata; |
---|
| 129 | |
---|
| 130 | n=idxij->nr; m=idxij->nc; |
---|
| 131 | pa=p; pb=p+m*cnp; |
---|
| 132 | Asz=mnp*cnp; Bsz=mnp*pnp; |
---|
| 133 | jaca=jac; jacb=jac+idxij->nnz*Asz; |
---|
| 134 | |
---|
| 135 | for(j=0; j<m; ++j){ |
---|
| 136 | /* j-th camera parameters */ |
---|
| 137 | paj=pa+j*cnp; |
---|
| 138 | |
---|
| 139 | nnz=sba_crsm_col_elmidxs(idxij, j, rcidxs, rcsubs); /* find nonzero hx_ij, i=0...n-1 */ |
---|
| 140 | |
---|
| 141 | for(i=0; i<nnz; ++i){ |
---|
| 142 | pbi=pb + rcsubs[i]*pnp; |
---|
| 143 | idx=idxij->val[rcidxs[i]]; |
---|
| 144 | pAij=jaca + idx*Asz; // set pAij to point to A_ij |
---|
| 145 | pBij=jacb + idx*Bsz; // set pBij to point to B_ij |
---|
| 146 | |
---|
| 147 | (*projac)(j, rcsubs[i], paj, pbi, pAij, pBij, projac_adata); // evaluate dQ/da, dQ/db in pAij, pBij |
---|
| 148 | } |
---|
| 149 | } |
---|
| 150 | } |
---|
| 151 | |
---|
| 152 | /* Given a parameter vector p made up of the 3D coordinates of n points and the parameters of m cameras, compute in |
---|
| 153 | * jac the jacobian of the predicted measurements, i.e. the jacobian of the projections of 3D points in the m images. |
---|
| 154 | * The jacobian is approximated with the aid of finite differences and is returned in the order |
---|
| 155 | * (A_11, ..., A_1m, ..., A_n1, ..., A_nm, B_11, ..., B_1m, ..., B_n1, ..., B_nm), |
---|
| 156 | * where A_ij=dx_ij/da_j and B_ij=dx_ij/db_i (see HZ). |
---|
| 157 | * Notice that depending on idxij, some of the A_ij, B_ij might be missing |
---|
| 158 | * |
---|
| 159 | * Problem-specific information is assumed to be stored in a structure pointed to by "dat". |
---|
| 160 | * |
---|
| 161 | * NOTE: This function is provided mainly for illustration purposes; in case that execution time is a concern, |
---|
| 162 | * the jacobian should be computed analytically |
---|
| 163 | */ |
---|
| 164 | static void sba_motstr_Qs_fdjac( |
---|
| 165 | double *p, /* I: current parameter estimate, (m*cnp+n*pnp)x1 */ |
---|
| 166 | struct sba_crsm *idxij, /* I: sparse matrix containing the location of x_ij in hx */ |
---|
| 167 | int *rcidxs, /* work array for the indexes of nonzero elements of a single sparse matrix row/column */ |
---|
| 168 | int *rcsubs, /* work array for the subscripts of nonzero elements in a single sparse matrix row/column */ |
---|
| 169 | double *jac, /* O: array for storing the approximated jacobian */ |
---|
| 170 | void *dat) /* I: points to a "wrap_motstr_data_" structure */ |
---|
| 171 | { |
---|
| 172 | register int i, j, ii, jj; |
---|
| 173 | double *pa, *pb, *paj, *pbi, *jaca, *jacb; |
---|
| 174 | register double *pAB; |
---|
| 175 | int n, m, nnz, Asz, Bsz; |
---|
| 176 | |
---|
| 177 | double tmp; |
---|
| 178 | register double d, d1; |
---|
| 179 | |
---|
| 180 | struct wrap_motstr_data_ *fdjd; |
---|
| 181 | void (*proj)(int j, int i, double *aj, double *bi, double *xij, void *adata); |
---|
| 182 | double *hxij, *hxxij; |
---|
| 183 | int cnp, pnp, mnp; |
---|
| 184 | void *adata; |
---|
| 185 | |
---|
| 186 | /* retrieve problem-specific information passed in *dat */ |
---|
| 187 | fdjd=(struct wrap_motstr_data_ *)dat; |
---|
| 188 | proj=fdjd->proj; |
---|
| 189 | cnp=fdjd->cnp; pnp=fdjd->pnp; mnp=fdjd->mnp; |
---|
| 190 | adata=fdjd->adata; |
---|
| 191 | |
---|
| 192 | n=idxij->nr; m=idxij->nc; |
---|
| 193 | pa=p; pb=p+m*cnp; |
---|
| 194 | Asz=mnp*cnp; Bsz=mnp*pnp; |
---|
| 195 | jaca=jac; jacb=jac+idxij->nnz*Asz; |
---|
| 196 | |
---|
| 197 | /* allocate memory for hxij, hxxij */ |
---|
| 198 | if((hxij=malloc(2*mnp*sizeof(double)))==NULL){ |
---|
| 199 | fprintf(stderr, "memory allocation request failed in sba_motstr_Qs_fdjac()!\n"); |
---|
| 200 | exit(1); |
---|
| 201 | } |
---|
| 202 | hxxij=hxij+mnp; |
---|
| 203 | |
---|
| 204 | if(cnp){ // is motion varying? |
---|
| 205 | /* compute A_ij */ |
---|
| 206 | for(j=0; j<m; ++j){ |
---|
| 207 | paj=pa+j*cnp; // j-th camera parameters |
---|
| 208 | |
---|
| 209 | nnz=sba_crsm_col_elmidxs(idxij, j, rcidxs, rcsubs); /* find nonzero A_ij, i=0...n-1 */ |
---|
| 210 | for(jj=0; jj<cnp; ++jj){ |
---|
| 211 | /* determine d=max(SBA_DELTA_SCALE*|paj[jj]|, SBA_MIN_DELTA), see HZ */ |
---|
| 212 | d=(double)(SBA_DELTA_SCALE)*paj[jj]; // force evaluation |
---|
| 213 | d=FABS(d); |
---|
| 214 | if(d<SBA_MIN_DELTA) d=SBA_MIN_DELTA; |
---|
| 215 | d1=1.0/d; /* invert so that divisions can be carried out faster as multiplications */ |
---|
| 216 | |
---|
| 217 | for(i=0; i<nnz; ++i){ |
---|
| 218 | pbi=pb + rcsubs[i]*pnp; // i-th point parameters |
---|
| 219 | (*proj)(j, rcsubs[i], paj, pbi, hxij, adata); // evaluate supplied function on current solution |
---|
| 220 | |
---|
| 221 | tmp=paj[jj]; |
---|
| 222 | paj[jj]+=d; |
---|
| 223 | (*proj)(j, rcsubs[i], paj, pbi, hxxij, adata); |
---|
| 224 | paj[jj]=tmp; /* restore */ |
---|
| 225 | |
---|
| 226 | pAB=jaca + idxij->val[rcidxs[i]]*Asz; // set pAB to point to A_ij |
---|
| 227 | for(ii=0; ii<mnp; ++ii) |
---|
| 228 | pAB[ii*cnp+jj]=(hxxij[ii]-hxij[ii])*d1; |
---|
| 229 | } |
---|
| 230 | } |
---|
| 231 | } |
---|
| 232 | } |
---|
| 233 | |
---|
| 234 | if(pnp){ // is structure varying? |
---|
| 235 | /* compute B_ij */ |
---|
| 236 | for(i=0; i<n; ++i){ |
---|
| 237 | pbi=pb+i*pnp; // i-th point parameters |
---|
| 238 | |
---|
| 239 | nnz=sba_crsm_row_elmidxs(idxij, i, rcidxs, rcsubs); /* find nonzero B_ij, j=0...m-1 */ |
---|
| 240 | for(jj=0; jj<pnp; ++jj){ |
---|
| 241 | /* determine d=max(SBA_DELTA_SCALE*|pbi[jj]|, SBA_MIN_DELTA), see HZ */ |
---|
| 242 | d=(double)(SBA_DELTA_SCALE)*pbi[jj]; // force evaluation |
---|
| 243 | d=FABS(d); |
---|
| 244 | if(d<SBA_MIN_DELTA) d=SBA_MIN_DELTA; |
---|
| 245 | d1=1.0/d; /* invert so that divisions can be carried out faster as multiplications */ |
---|
| 246 | |
---|
| 247 | for(j=0; j<nnz; ++j){ |
---|
| 248 | paj=pa + rcsubs[j]*cnp; // j-th camera parameters |
---|
| 249 | (*proj)(rcsubs[j], i, paj, pbi, hxij, adata); // evaluate supplied function on current solution |
---|
| 250 | |
---|
| 251 | tmp=pbi[jj]; |
---|
| 252 | pbi[jj]+=d; |
---|
| 253 | (*proj)(rcsubs[j], i, paj, pbi, hxxij, adata); |
---|
| 254 | pbi[jj]=tmp; /* restore */ |
---|
| 255 | |
---|
| 256 | pAB=jacb + idxij->val[rcidxs[j]]*Bsz; // set pAB to point to B_ij |
---|
| 257 | for(ii=0; ii<mnp; ++ii) |
---|
| 258 | pAB[ii*pnp+jj]=(hxxij[ii]-hxij[ii])*d1; |
---|
| 259 | } |
---|
| 260 | } |
---|
| 261 | } |
---|
| 262 | } |
---|
| 263 | |
---|
| 264 | free(hxij); |
---|
| 265 | } |
---|
| 266 | |
---|
| 267 | /* BUNDLE ADJUSTMENT FOR CAMERA PARAMETERS ONLY */ |
---|
| 268 | |
---|
| 269 | /* Given a parameter vector p made up of the parameters of m cameras, compute in |
---|
| 270 | * hx the prediction of the measurements, i.e. the projections of 3D points in the m images. |
---|
| 271 | * The measurements are returned in the order (hx_11^T, .. hx_1m^T, ..., hx_n1^T, .. hx_nm^T)^T, |
---|
| 272 | * where hx_ij is the predicted projection of the i-th point on the j-th camera. |
---|
| 273 | * Caller supplies rcidxs and rcsubs which can be used as working memory. |
---|
| 274 | * Notice that depending on idxij, some of the hx_ij might be missing |
---|
| 275 | * |
---|
| 276 | */ |
---|
| 277 | static void sba_mot_Qs(double *p, struct sba_crsm *idxij, int *rcidxs, int *rcsubs, double *hx, void *adata) |
---|
| 278 | { |
---|
| 279 | register int i, j; |
---|
| 280 | int cnp, mnp; |
---|
| 281 | double *paj, *pxij; |
---|
| 282 | //int n; |
---|
| 283 | int m, nnz; |
---|
| 284 | struct wrap_mot_data_ *wdata; |
---|
| 285 | void (*proj)(int j, int i, double *aj, double *xij, void *proj_adata); |
---|
| 286 | void *proj_adata; |
---|
| 287 | |
---|
| 288 | wdata=(struct wrap_mot_data_ *)adata; |
---|
| 289 | cnp=wdata->cnp; mnp=wdata->mnp; |
---|
| 290 | proj=wdata->proj; |
---|
| 291 | proj_adata=wdata->adata; |
---|
| 292 | |
---|
| 293 | //n=idxij->nr; |
---|
| 294 | m=idxij->nc; |
---|
| 295 | |
---|
| 296 | for(j=0; j<m; ++j){ |
---|
| 297 | /* j-th camera parameters */ |
---|
| 298 | paj=p+j*cnp; |
---|
| 299 | |
---|
| 300 | nnz=sba_crsm_col_elmidxs(idxij, j, rcidxs, rcsubs); /* find nonzero hx_ij, i=0...n-1 */ |
---|
| 301 | |
---|
| 302 | for(i=0; i<nnz; ++i){ |
---|
| 303 | pxij=hx + idxij->val[rcidxs[i]]*mnp; // set pxij to point to hx_ij |
---|
| 304 | |
---|
| 305 | (*proj)(j, rcsubs[i], paj, pxij, proj_adata); // evaluate Q in pxij |
---|
| 306 | } |
---|
| 307 | } |
---|
| 308 | } |
---|
| 309 | |
---|
| 310 | /* Given a parameter vector p made up of the parameters of m cameras, compute in jac |
---|
| 311 | * the jacobian of the predicted measurements, i.e. the jacobian of the projections of 3D points in the m images. |
---|
| 312 | * The jacobian is returned in the order (A_11, ..., A_1m, ..., A_n1, ..., A_nm), |
---|
| 313 | * where A_ij=dx_ij/db_j (see HZ). |
---|
| 314 | * Caller supplies rcidxs and rcsubs which can be used as working memory. |
---|
| 315 | * Notice that depending on idxij, some of the A_ij might be missing |
---|
| 316 | * |
---|
| 317 | */ |
---|
| 318 | static void sba_mot_Qs_jac(double *p, struct sba_crsm *idxij, int *rcidxs, int *rcsubs, double *jac, void *adata) |
---|
| 319 | { |
---|
| 320 | register int i, j; |
---|
| 321 | int cnp, mnp; |
---|
| 322 | double *paj, *jaca, *pAij; |
---|
| 323 | //int n; |
---|
| 324 | int m, nnz, Asz, idx; |
---|
| 325 | struct wrap_mot_data_ *wdata; |
---|
| 326 | void (*projac)(int j, int i, double *aj, double *Aij, void *projac_adata); |
---|
| 327 | void *projac_adata; |
---|
| 328 | |
---|
| 329 | wdata=(struct wrap_mot_data_ *)adata; |
---|
| 330 | cnp=wdata->cnp; mnp=wdata->mnp; |
---|
| 331 | projac=wdata->projac; |
---|
| 332 | projac_adata=wdata->adata; |
---|
| 333 | |
---|
| 334 | //n=idxij->nr; |
---|
| 335 | m=idxij->nc; |
---|
| 336 | Asz=mnp*cnp; |
---|
| 337 | jaca=jac; |
---|
| 338 | |
---|
| 339 | for(j=0; j<m; ++j){ |
---|
| 340 | /* j-th camera parameters */ |
---|
| 341 | paj=p+j*cnp; |
---|
| 342 | |
---|
| 343 | nnz=sba_crsm_col_elmidxs(idxij, j, rcidxs, rcsubs); /* find nonzero hx_ij, i=0...n-1 */ |
---|
| 344 | |
---|
| 345 | for(i=0; i<nnz; ++i){ |
---|
| 346 | idx=idxij->val[rcidxs[i]]; |
---|
| 347 | pAij=jaca + idx*Asz; // set pAij to point to A_ij |
---|
| 348 | |
---|
| 349 | (*projac)(j, rcsubs[i], paj, pAij, projac_adata); // evaluate dQ/da in pAij |
---|
| 350 | } |
---|
| 351 | } |
---|
| 352 | } |
---|
| 353 | |
---|
| 354 | /* Given a parameter vector p made up of the parameters of m cameras, compute in jac the jacobian |
---|
| 355 | * of the predicted measurements, i.e. the jacobian of the projections of 3D points in the m images. |
---|
| 356 | * The jacobian is approximated with the aid of finite differences and is returned in the order |
---|
| 357 | * (A_11, ..., A_1m, ..., A_n1, ..., A_nm), where A_ij=dx_ij/da_j (see HZ). |
---|
| 358 | * Notice that depending on idxij, some of the A_ij might be missing |
---|
| 359 | * |
---|
| 360 | * Problem-specific information is assumed to be stored in a structure pointed to by "dat". |
---|
| 361 | * |
---|
| 362 | * NOTE: This function is provided mainly for illustration purposes; in case that execution time is a concern, |
---|
| 363 | * the jacobian should be computed analytically |
---|
| 364 | */ |
---|
| 365 | static void sba_mot_Qs_fdjac( |
---|
| 366 | double *p, /* I: current parameter estimate, (m*cnp)x1 */ |
---|
| 367 | struct sba_crsm *idxij, /* I: sparse matrix containing the location of x_ij in hx */ |
---|
| 368 | int *rcidxs, /* work array for the indexes of nonzero elements of a single sparse matrix row/column */ |
---|
| 369 | int *rcsubs, /* work array for the subscripts of nonzero elements in a single sparse matrix row/column */ |
---|
| 370 | double *jac, /* O: array for storing the approximated jacobian */ |
---|
| 371 | void *dat) /* I: points to a "wrap_mot_data_" structure */ |
---|
| 372 | { |
---|
| 373 | register int i, j, ii, jj; |
---|
| 374 | double *paj, *jaca; |
---|
| 375 | register double *pA; |
---|
| 376 | //int n; |
---|
| 377 | int m, nnz, Asz; |
---|
| 378 | |
---|
| 379 | double tmp; |
---|
| 380 | register double d, d1; |
---|
| 381 | |
---|
| 382 | struct wrap_mot_data_ *fdjd; |
---|
| 383 | void (*proj)(int j, int i, double *aj, double *xij, void *adata); |
---|
| 384 | double *hxij, *hxxij; |
---|
| 385 | int cnp, mnp; |
---|
| 386 | void *adata; |
---|
| 387 | |
---|
| 388 | /* retrieve problem-specific information passed in *dat */ |
---|
| 389 | fdjd=(struct wrap_mot_data_ *)dat; |
---|
| 390 | proj=fdjd->proj; |
---|
| 391 | cnp=fdjd->cnp; mnp=fdjd->mnp; |
---|
| 392 | adata=fdjd->adata; |
---|
| 393 | |
---|
| 394 | //n=idxij->nr; |
---|
| 395 | m=idxij->nc; |
---|
| 396 | Asz=mnp*cnp; |
---|
| 397 | jaca=jac; |
---|
| 398 | |
---|
| 399 | /* allocate memory for hxij, hxxij */ |
---|
| 400 | if((hxij=malloc(2*mnp*sizeof(double)))==NULL){ |
---|
| 401 | fprintf(stderr, "memory allocation request failed in sba_mot_Qs_fdjac()!\n"); |
---|
| 402 | exit(1); |
---|
| 403 | } |
---|
| 404 | hxxij=hxij+mnp; |
---|
| 405 | |
---|
| 406 | /* compute A_ij */ |
---|
| 407 | for(j=0; j<m; ++j){ |
---|
| 408 | paj=p+j*cnp; // j-th camera parameters |
---|
| 409 | |
---|
| 410 | nnz=sba_crsm_col_elmidxs(idxij, j, rcidxs, rcsubs); /* find nonzero A_ij, i=0...n-1 */ |
---|
| 411 | for(jj=0; jj<cnp; ++jj){ |
---|
| 412 | /* determine d=max(SBA_DELTA_SCALE*|paj[jj]|, SBA_MIN_DELTA), see HZ */ |
---|
| 413 | d=(double)(SBA_DELTA_SCALE)*paj[jj]; // force evaluation |
---|
| 414 | d=FABS(d); |
---|
| 415 | if(d<SBA_MIN_DELTA) d=SBA_MIN_DELTA; |
---|
| 416 | d1=1.0/d; /* invert so that divisions can be carried out faster as multiplications */ |
---|
| 417 | |
---|
| 418 | for(i=0; i<nnz; ++i){ |
---|
| 419 | (*proj)(j, rcsubs[i], paj, hxij, adata); // evaluate supplied function on current solution |
---|
| 420 | |
---|
| 421 | tmp=paj[jj]; |
---|
| 422 | paj[jj]+=d; |
---|
| 423 | (*proj)(j, rcsubs[i], paj, hxxij, adata); |
---|
| 424 | paj[jj]=tmp; /* restore */ |
---|
| 425 | |
---|
| 426 | pA=jaca + idxij->val[rcidxs[i]]*Asz; // set pA to point to A_ij |
---|
| 427 | for(ii=0; ii<mnp; ++ii) |
---|
| 428 | pA[ii*cnp+jj]=(hxxij[ii]-hxij[ii])*d1; |
---|
| 429 | } |
---|
| 430 | } |
---|
| 431 | } |
---|
| 432 | |
---|
| 433 | free(hxij); |
---|
| 434 | } |
---|
| 435 | |
---|
| 436 | /* BUNDLE ADJUSTMENT FOR STRUCTURE PARAMETERS ONLY */ |
---|
| 437 | |
---|
| 438 | /* Given a parameter vector p made up of the 3D coordinates of n points, compute in |
---|
| 439 | * hx the prediction of the measurements, i.e. the projections of 3D points in the m images. The measurements |
---|
| 440 | * are returned in the order (hx_11^T, .. hx_1m^T, ..., hx_n1^T, .. hx_nm^T)^T, where hx_ij is the predicted |
---|
| 441 | * projection of the i-th point on the j-th camera. |
---|
| 442 | * Caller supplies rcidxs and rcsubs which can be used as working memory. |
---|
| 443 | * Notice that depending on idxij, some of the hx_ij might be missing |
---|
| 444 | * |
---|
| 445 | */ |
---|
| 446 | static void sba_str_Qs(double *p, struct sba_crsm *idxij, int *rcidxs, int *rcsubs, double *hx, void *adata) |
---|
| 447 | { |
---|
| 448 | register int i, j; |
---|
| 449 | int pnp, mnp; |
---|
| 450 | double *pbi, *pxij; |
---|
| 451 | //int n; |
---|
| 452 | int m, nnz; |
---|
| 453 | struct wrap_str_data_ *wdata; |
---|
| 454 | void (*proj)(int j, int i, double *bi, double *xij, void *proj_adata); |
---|
| 455 | void *proj_adata; |
---|
| 456 | |
---|
| 457 | wdata=(struct wrap_str_data_ *)adata; |
---|
| 458 | pnp=wdata->pnp; mnp=wdata->mnp; |
---|
| 459 | proj=wdata->proj; |
---|
| 460 | proj_adata=wdata->adata; |
---|
| 461 | |
---|
| 462 | //n=idxij->nr; |
---|
| 463 | m=idxij->nc; |
---|
| 464 | |
---|
| 465 | for(j=0; j<m; ++j){ |
---|
| 466 | nnz=sba_crsm_col_elmidxs(idxij, j, rcidxs, rcsubs); /* find nonzero hx_ij, i=0...n-1 */ |
---|
| 467 | |
---|
| 468 | for(i=0; i<nnz; ++i){ |
---|
| 469 | pbi=p + rcsubs[i]*pnp; |
---|
| 470 | pxij=hx + idxij->val[rcidxs[i]]*mnp; // set pxij to point to hx_ij |
---|
| 471 | |
---|
| 472 | (*proj)(j, rcsubs[i], pbi, pxij, proj_adata); // evaluate Q in pxij |
---|
| 473 | } |
---|
| 474 | } |
---|
| 475 | } |
---|
| 476 | |
---|
| 477 | /* Given a parameter vector p made up of the 3D coordinates of n points, compute in |
---|
| 478 | * jac the jacobian of the predicted measurements, i.e. the jacobian of the projections of 3D points in the m images. |
---|
| 479 | * The jacobian is returned in the order (B_11, ..., B_1m, ..., B_n1, ..., B_nm), where B_ij=dx_ij/db_i (see HZ). |
---|
| 480 | * Caller supplies rcidxs and rcsubs which can be used as working memory. |
---|
| 481 | * Notice that depending on idxij, some of the B_ij might be missing |
---|
| 482 | * |
---|
| 483 | */ |
---|
| 484 | static void sba_str_Qs_jac(double *p, struct sba_crsm *idxij, int *rcidxs, int *rcsubs, double *jac, void *adata) |
---|
| 485 | { |
---|
| 486 | register int i, j; |
---|
| 487 | int pnp, mnp; |
---|
| 488 | double *pbi, *pBij; |
---|
| 489 | //int n; |
---|
| 490 | int m, nnz, Bsz, idx; |
---|
| 491 | struct wrap_str_data_ *wdata; |
---|
| 492 | void (*projac)(int j, int i, double *bi, double *Bij, void *projac_adata); |
---|
| 493 | void *projac_adata; |
---|
| 494 | |
---|
| 495 | |
---|
| 496 | wdata=(struct wrap_str_data_ *)adata; |
---|
| 497 | pnp=wdata->pnp; mnp=wdata->mnp; |
---|
| 498 | projac=wdata->projac; |
---|
| 499 | projac_adata=wdata->adata; |
---|
| 500 | |
---|
| 501 | //n=idxij->nr; |
---|
| 502 | m=idxij->nc; |
---|
| 503 | Bsz=mnp*pnp; |
---|
| 504 | |
---|
| 505 | for(j=0; j<m; ++j){ |
---|
| 506 | |
---|
| 507 | nnz=sba_crsm_col_elmidxs(idxij, j, rcidxs, rcsubs); /* find nonzero hx_ij, i=0...n-1 */ |
---|
| 508 | |
---|
| 509 | for(i=0; i<nnz; ++i){ |
---|
| 510 | pbi=p + rcsubs[i]*pnp; |
---|
| 511 | idx=idxij->val[rcidxs[i]]; |
---|
| 512 | pBij=jac + idx*Bsz; // set pBij to point to B_ij |
---|
| 513 | |
---|
| 514 | (*projac)(j, rcsubs[i], pbi, pBij, projac_adata); // evaluate dQ/db in pBij |
---|
| 515 | } |
---|
| 516 | } |
---|
| 517 | } |
---|
| 518 | |
---|
| 519 | /* Given a parameter vector p made up of the 3D coordinates of n points, compute in |
---|
| 520 | * jac the jacobian of the predicted measurements, i.e. the jacobian of the projections of 3D points in the m images. |
---|
| 521 | * The jacobian is approximated with the aid of finite differences and is returned in the order |
---|
| 522 | * (B_11, ..., B_1m, ..., B_n1, ..., B_nm), where B_ij=dx_ij/db_i (see HZ). |
---|
| 523 | * Notice that depending on idxij, some of the B_ij might be missing |
---|
| 524 | * |
---|
| 525 | * Problem-specific information is assumed to be stored in a structure pointed to by "dat". |
---|
| 526 | * |
---|
| 527 | * NOTE: This function is provided mainly for illustration purposes; in case that execution time is a concern, |
---|
| 528 | * the jacobian should be computed analytically |
---|
| 529 | */ |
---|
| 530 | static void sba_str_Qs_fdjac( |
---|
| 531 | double *p, /* I: current parameter estimate, (n*pnp)x1 */ |
---|
| 532 | struct sba_crsm *idxij, /* I: sparse matrix containing the location of x_ij in hx */ |
---|
| 533 | int *rcidxs, /* work array for the indexes of nonzero elements of a single sparse matrix row/column */ |
---|
| 534 | int *rcsubs, /* work array for the subscripts of nonzero elements in a single sparse matrix row/column */ |
---|
| 535 | double *jac, /* O: array for storing the approximated jacobian */ |
---|
| 536 | void *dat) /* I: points to a "wrap_str_data_" structure */ |
---|
| 537 | { |
---|
| 538 | register int i, j, ii, jj; |
---|
| 539 | double *pbi; |
---|
| 540 | register double *pB; |
---|
| 541 | //int m; |
---|
| 542 | int n, nnz, Bsz; |
---|
| 543 | |
---|
| 544 | double tmp; |
---|
| 545 | register double d, d1; |
---|
| 546 | |
---|
| 547 | struct wrap_str_data_ *fdjd; |
---|
| 548 | void (*proj)(int j, int i, double *bi, double *xij, void *adata); |
---|
| 549 | double *hxij, *hxxij; |
---|
| 550 | int pnp, mnp; |
---|
| 551 | void *adata; |
---|
| 552 | |
---|
| 553 | /* retrieve problem-specific information passed in *dat */ |
---|
| 554 | fdjd=(struct wrap_str_data_ *)dat; |
---|
| 555 | proj=fdjd->proj; |
---|
| 556 | pnp=fdjd->pnp; mnp=fdjd->mnp; |
---|
| 557 | adata=fdjd->adata; |
---|
| 558 | |
---|
| 559 | n=idxij->nr; |
---|
| 560 | //m=idxij->nc; |
---|
| 561 | Bsz=mnp*pnp; |
---|
| 562 | |
---|
| 563 | /* allocate memory for hxij, hxxij */ |
---|
| 564 | if((hxij=malloc(2*mnp*sizeof(double)))==NULL){ |
---|
| 565 | fprintf(stderr, "memory allocation request failed in sba_str_Qs_fdjac()!\n"); |
---|
| 566 | exit(1); |
---|
| 567 | } |
---|
| 568 | hxxij=hxij+mnp; |
---|
| 569 | |
---|
| 570 | /* compute B_ij */ |
---|
| 571 | for(i=0; i<n; ++i){ |
---|
| 572 | pbi=p+i*pnp; // i-th point parameters |
---|
| 573 | |
---|
| 574 | nnz=sba_crsm_row_elmidxs(idxij, i, rcidxs, rcsubs); /* find nonzero B_ij, j=0...m-1 */ |
---|
| 575 | for(jj=0; jj<pnp; ++jj){ |
---|
| 576 | /* determine d=max(SBA_DELTA_SCALE*|pbi[jj]|, SBA_MIN_DELTA), see HZ */ |
---|
| 577 | d=(double)(SBA_DELTA_SCALE)*pbi[jj]; // force evaluation |
---|
| 578 | d=FABS(d); |
---|
| 579 | if(d<SBA_MIN_DELTA) d=SBA_MIN_DELTA; |
---|
| 580 | d1=1.0/d; /* invert so that divisions can be carried out faster as multiplications */ |
---|
| 581 | |
---|
| 582 | for(j=0; j<nnz; ++j){ |
---|
| 583 | (*proj)(rcsubs[j], i, pbi, hxij, adata); // evaluate supplied function on current solution |
---|
| 584 | |
---|
| 585 | tmp=pbi[jj]; |
---|
| 586 | pbi[jj]+=d; |
---|
| 587 | (*proj)(rcsubs[j], i, pbi, hxxij, adata); |
---|
| 588 | pbi[jj]=tmp; /* restore */ |
---|
| 589 | |
---|
| 590 | pB=jac + idxij->val[rcidxs[j]]*Bsz; // set pB to point to B_ij |
---|
| 591 | for(ii=0; ii<mnp; ++ii) |
---|
| 592 | pB[ii*pnp+jj]=(hxxij[ii]-hxij[ii])*d1; |
---|
| 593 | } |
---|
| 594 | } |
---|
| 595 | } |
---|
| 596 | |
---|
| 597 | free(hxij); |
---|
| 598 | } |
---|
| 599 | |
---|
| 600 | |
---|
| 601 | /* |
---|
| 602 | * Simple driver to sba_motstr_levmar_x for bundle adjustment on camera and structure parameters. |
---|
| 603 | */ |
---|
| 604 | |
---|
| 605 | int sba_motstr_levmar( |
---|
| 606 | const int n, /* number of points */ |
---|
| 607 | const int m, /* number of images */ |
---|
| 608 | const int mcon,/* number of images (starting from the 1st) whose parameters should not be modified. |
---|
| 609 | * All A_ij (see below) with j<mcon are assumed to be zero |
---|
| 610 | */ |
---|
| 611 | char *vmask, /* visibility mask: vmask[i][j]=1 if point i visible in image j, 0 otherwise. nxm */ |
---|
| 612 | double *p, /* initial parameter vector p0: (a1, ..., am, b1, ..., bn). |
---|
| 613 | * aj are the image j parameters, bi are the i-th point parameters, |
---|
| 614 | * size m*cnp + n*pnp |
---|
| 615 | */ |
---|
| 616 | const int cnp,/* number of parameters for ONE camera; e.g. 6 for Euclidean cameras */ |
---|
| 617 | const int pnp,/* number of parameters for ONE point; e.g. 3 for Euclidean points */ |
---|
| 618 | double *x, /* measurements vector: (x_11^T, .. x_1m^T, ..., x_n1^T, .. x_nm^T)^T where |
---|
| 619 | * x_ij is the projection of the i-th point on the j-th image. |
---|
| 620 | * NOTE: some of the x_ij might be missing, if point i is not visible in image j; |
---|
| 621 | * see vmask[i][j], max. size n*m*mnp |
---|
| 622 | */ |
---|
| 623 | const int mnp,/* number of parameters for EACH measurement; usually 2 */ |
---|
| 624 | void (*proj)(int j, int i, double *aj, double *bi, double *xij, void *adata), |
---|
| 625 | /* functional relation computing a SINGLE image measurement. Assuming that |
---|
| 626 | * the parameters of point i are bi and the parameters of camera j aj, |
---|
| 627 | * computes a prediction of \hat{x}_{ij}. aj is cnp x 1, bi is pnp x 1 and |
---|
| 628 | * xij is mnp x 1. This function is called only if point i is visible in |
---|
| 629 | * image j (i.e. vmask[i][j]==1) |
---|
| 630 | */ |
---|
| 631 | void (*projac)(int j, int i, double *aj, double *bi, double *Aij, double *Bij, void *adata), |
---|
| 632 | /* functional relation to evaluate d x_ij / d a_j and |
---|
| 633 | * d x_ij / d b_i in Aij and Bij resp. |
---|
| 634 | * This function is called only if point i is visible in * image j |
---|
| 635 | * (i.e. vmask[i][j]==1). Also, A_ij and B_ij are mnp x cnp and mnp x pnp |
---|
| 636 | * matrices resp. and they should be stored in row-major order. |
---|
| 637 | * |
---|
| 638 | * If NULL, the jacobians are approximated by repetitive proj calls |
---|
| 639 | * and finite differences. |
---|
| 640 | */ |
---|
| 641 | void *adata, /* pointer to possibly additional data, passed uninterpreted to proj, projac */ |
---|
| 642 | |
---|
| 643 | int itmax, /* I: maximum number of iterations. itmax==0 signals jacobian verification followed by immediate return */ |
---|
| 644 | int verbose, /* I: verbosity */ |
---|
| 645 | double opts[SBA_OPTSSZ], |
---|
| 646 | /* I: minim. options [\mu, \epsilon1, \epsilon2]. Respectively the scale factor for initial \mu, |
---|
| 647 | * stoping thresholds for ||J^T e||_inf, ||dp||_2 and ||e||_2 |
---|
| 648 | */ |
---|
| 649 | double info[SBA_INFOSZ] |
---|
| 650 | /* O: information regarding the minimization. Set to NULL if don't care |
---|
| 651 | * info[0]=||e||_2 at initial p. |
---|
| 652 | * info[1-4]=[ ||e||_2, ||J^T e||_inf, ||dp||_2, mu/max[J^T J]_ii ], all computed at estimated p. |
---|
| 653 | * info[5]= # iterations, |
---|
| 654 | * info[6]=reason for terminating: 1 - stopped by small gradient J^T e |
---|
| 655 | * 2 - stopped by small dp |
---|
| 656 | * 3 - stopped by itmax |
---|
| 657 | * 4 - singular matrix. Restart from current p with increased mu |
---|
| 658 | * 5 - too many attempts to increase damping. Restart with increased mu |
---|
| 659 | * 6 - stopped by small ||e||_2 |
---|
| 660 | * info[7]= # function evaluations |
---|
| 661 | * info[8]= # jacobian evaluations |
---|
| 662 | * info[9]= # number of linear systems solved, i.e. number of attempts for reducing error |
---|
| 663 | */ |
---|
| 664 | ) |
---|
| 665 | { |
---|
| 666 | int retval; |
---|
| 667 | struct wrap_motstr_data_ wdata; |
---|
| 668 | static void (*fjac)(double *p, struct sba_crsm *idxij, int *rcidxs, int *rcsubs, double *jac, void *adata); |
---|
| 669 | |
---|
| 670 | wdata.proj=proj; |
---|
| 671 | wdata.projac=projac; |
---|
| 672 | wdata.cnp=cnp; |
---|
| 673 | wdata.pnp=pnp; |
---|
| 674 | wdata.mnp=mnp; |
---|
| 675 | wdata.adata=adata; |
---|
| 676 | |
---|
| 677 | fjac=(projac)? sba_motstr_Qs_jac : sba_motstr_Qs_fdjac; |
---|
| 678 | retval=sba_motstr_levmar_x(n, m, mcon, vmask, p, cnp, pnp, x, mnp, sba_motstr_Qs, fjac, &wdata, itmax, verbose, opts, info); |
---|
| 679 | |
---|
| 680 | if(info){ |
---|
| 681 | register int i; |
---|
| 682 | int nvis; |
---|
| 683 | |
---|
| 684 | /* count visible image points */ |
---|
| 685 | for(i=nvis=0; i<n*m; ++i) |
---|
| 686 | nvis+=vmask[i]; |
---|
| 687 | |
---|
| 688 | /* each "func" & "fjac" evaluation requires nvis "proj" & "projac" evaluations */ |
---|
| 689 | info[7]*=nvis; |
---|
| 690 | info[8]*=nvis; |
---|
| 691 | } |
---|
| 692 | |
---|
| 693 | return retval; |
---|
| 694 | } |
---|
| 695 | |
---|
| 696 | |
---|
| 697 | /* |
---|
| 698 | * Simple driver to sba_mot_levmar_x for bundle adjustment on camera parameters. |
---|
| 699 | */ |
---|
| 700 | |
---|
| 701 | int sba_mot_levmar( |
---|
| 702 | const int n, /* number of points */ |
---|
| 703 | const int m, /* number of images */ |
---|
| 704 | const int mcon,/* number of images (starting from the 1st) whose parameters should not be modified. |
---|
| 705 | * All A_ij (see below) with j<mcon are assumed to be zero |
---|
| 706 | */ |
---|
| 707 | char *vmask, /* visibility mask: vmask[i][j]=1 if point i visible in image j, 0 otherwise. nxm */ |
---|
| 708 | double *p, /* initial parameter vector p0: (a1, ..., am). |
---|
| 709 | * aj are the image j parameters, size m*cnp */ |
---|
| 710 | const int cnp,/* number of parameters for ONE camera; e.g. 6 for Euclidean cameras */ |
---|
| 711 | double *x, /* measurements vector: (x_11^T, .. x_1m^T, ..., x_n1^T, .. x_nm^T)^T where |
---|
| 712 | * x_ij is the projection of the i-th point on the j-th image. |
---|
| 713 | * NOTE: some of the x_ij might be missing, if point i is not visible in image j; |
---|
| 714 | * see vmask[i][j], max. size n*m*mnp |
---|
| 715 | */ |
---|
| 716 | const int mnp,/* number of parameters for EACH measurement; usually 2 */ |
---|
| 717 | void (*proj)(int j, int i, double *aj, double *xij, void *adata), |
---|
| 718 | /* functional relation computing a SINGLE image measurement. Assuming that |
---|
| 719 | * the parameters of camera j are aj, computes a prediction of \hat{x}_{ij} |
---|
| 720 | * for point i. aj is cnp x 1 and xij is mnp x 1. |
---|
| 721 | * This function is called only if point i is visible in image j (i.e. vmask[i][j]==1) |
---|
| 722 | */ |
---|
| 723 | void (*projac)(int j, int i, double *aj, double *Aij, void *adata), |
---|
| 724 | /* functional relation to evaluate d x_ij / d a_j in Aij |
---|
| 725 | * This function is called only if point i is visible in image j |
---|
| 726 | * (i.e. vmask[i][j]==1). Also, A_ij are a mnp x cnp matrices |
---|
| 727 | * and should be stored in row-major order. |
---|
| 728 | * |
---|
| 729 | * If NULL, the jacobian is approximated by repetitive proj calls |
---|
| 730 | * and finite differences. |
---|
| 731 | */ |
---|
| 732 | void *adata, /* pointer to possibly additional data, passed uninterpreted to proj, projac */ |
---|
| 733 | |
---|
| 734 | int itmax, /* I: maximum number of iterations. itmax==0 signals jacobian verification followed by immediate return */ |
---|
| 735 | int verbose, /* I: verbosity */ |
---|
| 736 | double opts[SBA_OPTSSZ], |
---|
| 737 | /* I: minim. options [\mu, \epsilon1, \epsilon2]. Respectively the scale factor for initial \mu, |
---|
| 738 | * stoping thresholds for ||J^T e||_inf, ||dp||_2 and ||e||_2 |
---|
| 739 | */ |
---|
| 740 | double info[SBA_INFOSZ] |
---|
| 741 | /* O: information regarding the minimization. Set to NULL if don't care |
---|
| 742 | * info[0]=||e||_2 at initial p. |
---|
| 743 | * info[1-4]=[ ||e||_2, ||J^T e||_inf, ||dp||_2, mu/max[J^T J]_ii ], all computed at estimated p. |
---|
| 744 | * info[5]= # iterations, |
---|
| 745 | * info[6]=reason for terminating: 1 - stopped by small gradient J^T e |
---|
| 746 | * 2 - stopped by small dp |
---|
| 747 | * 3 - stopped by itmax |
---|
| 748 | * 4 - singular matrix. Restart from current p with increased mu |
---|
| 749 | * 5 - too many attempts to increase damping. Restart with increased mu |
---|
| 750 | * 6 - stopped by small ||e||_2 |
---|
| 751 | * info[7]= # function evaluations |
---|
| 752 | * info[8]= # jacobian evaluations |
---|
| 753 | * info[9]= # number of linear systems solved, i.e. number of attempts for reducing error |
---|
| 754 | */ |
---|
| 755 | ) |
---|
| 756 | { |
---|
| 757 | int retval; |
---|
| 758 | struct wrap_mot_data_ wdata; |
---|
| 759 | void (*fjac)(double *p, struct sba_crsm *idxij, int *rcidxs, int *rcsubs, double *jac, void *adata); |
---|
| 760 | |
---|
| 761 | wdata.proj=proj; |
---|
| 762 | wdata.projac=projac; |
---|
| 763 | wdata.cnp=cnp; |
---|
| 764 | wdata.mnp=mnp; |
---|
| 765 | wdata.adata=adata; |
---|
| 766 | |
---|
| 767 | fjac=(projac)? sba_mot_Qs_jac : sba_mot_Qs_fdjac; |
---|
| 768 | retval=sba_mot_levmar_x(n, m, mcon, vmask, p, cnp, x, mnp, sba_mot_Qs, fjac, &wdata, itmax, verbose, opts, info); |
---|
| 769 | |
---|
| 770 | if(info){ |
---|
| 771 | register int i; |
---|
| 772 | int nvis; |
---|
| 773 | |
---|
| 774 | /* count visible image points */ |
---|
| 775 | for(i=nvis=0; i<n*m; ++i) |
---|
| 776 | nvis+=vmask[i]; |
---|
| 777 | |
---|
| 778 | /* each "func" & "fjac" evaluation requires nvis "proj" & "projac" evaluations */ |
---|
| 779 | info[7]*=nvis; |
---|
| 780 | info[8]*=nvis; |
---|
| 781 | } |
---|
| 782 | |
---|
| 783 | return retval; |
---|
| 784 | } |
---|
| 785 | |
---|
| 786 | /* |
---|
| 787 | * Simple driver to sba_str_levmar_x for bundle adjustment on structure parameters. |
---|
| 788 | */ |
---|
| 789 | |
---|
| 790 | int sba_str_levmar( |
---|
| 791 | const int n, /* number of points */ |
---|
| 792 | const int m, /* number of images */ |
---|
| 793 | char *vmask, /* visibility mask: vmask[i][j]=1 if point i visible in image j, 0 otherwise. nxm */ |
---|
| 794 | double *p, /* initial parameter vector p0: (b1, ..., bn). |
---|
| 795 | * bi are the i-th point parameters, size n*pnp |
---|
| 796 | */ |
---|
| 797 | const int pnp,/* number of parameters for ONE point; e.g. 3 for Euclidean points */ |
---|
| 798 | double *x, /* measurements vector: (x_11^T, .. x_1m^T, ..., x_n1^T, .. x_nm^T)^T where |
---|
| 799 | * x_ij is the projection of the i-th point on the j-th image. |
---|
| 800 | * NOTE: some of the x_ij might be missing, if point i is not visible in image j; |
---|
| 801 | * see vmask[i][j], max. size n*m*mnp |
---|
| 802 | */ |
---|
| 803 | const int mnp,/* number of parameters for EACH measurement; usually 2 */ |
---|
| 804 | void (*proj)(int j, int i, double *bi, double *xij, void *adata), |
---|
| 805 | /* functional relation computing a SINGLE image measurement. Assuming that |
---|
| 806 | * the parameters of point i are bi, computes a prediction of \hat{x}_{ij}. |
---|
| 807 | * bi is pnp x 1 and xij is mnp x 1. This function is called only if point |
---|
| 808 | * i is visible in image j (i.e. vmask[i][j]==1) |
---|
| 809 | */ |
---|
| 810 | void (*projac)(int j, int i, double *bi, double *Bij, void *adata), |
---|
| 811 | /* functional relation to evaluate d x_ij / d b_i in Bij. |
---|
| 812 | * This function is called only if point i is visible in image j |
---|
| 813 | * (i.e. vmask[i][j]==1). Also, B_ij are mnp x pnp matrices |
---|
| 814 | * and they should be stored in row-major order. |
---|
| 815 | * |
---|
| 816 | * If NULL, the jacobians are approximated by repetitive proj calls |
---|
| 817 | * and finite differences. |
---|
| 818 | */ |
---|
| 819 | void *adata, /* pointer to possibly additional data, passed uninterpreted to proj, projac */ |
---|
| 820 | |
---|
| 821 | int itmax, /* I: maximum number of iterations. itmax==0 signals jacobian verification followed by immediate return */ |
---|
| 822 | int verbose, /* I: verbosity */ |
---|
| 823 | double opts[SBA_OPTSSZ], |
---|
| 824 | /* I: minim. options [\mu, \epsilon1, \epsilon2]. Respectively the scale factor for initial \mu, |
---|
| 825 | * stoping thresholds for ||J^T e||_inf, ||dp||_2 and ||e||_2 |
---|
| 826 | */ |
---|
| 827 | double info[SBA_INFOSZ] |
---|
| 828 | /* O: information regarding the minimization. Set to NULL if don't care |
---|
| 829 | * info[0]=||e||_2 at initial p. |
---|
| 830 | * info[1-4]=[ ||e||_2, ||J^T e||_inf, ||dp||_2, mu/max[J^T J]_ii ], all computed at estimated p. |
---|
| 831 | * info[5]= # iterations, |
---|
| 832 | * info[6]=reason for terminating: 1 - stopped by small gradient J^T e |
---|
| 833 | * 2 - stopped by small dp |
---|
| 834 | * 3 - stopped by itmax |
---|
| 835 | * 4 - singular matrix. Restart from current p with increased mu |
---|
| 836 | * 5 - too many attempts to increase damping. Restart with increased mu |
---|
| 837 | * 6 - stopped by small ||e||_2 |
---|
| 838 | * info[7]= # function evaluations |
---|
| 839 | * info[8]= # jacobian evaluations |
---|
| 840 | * info[9]= # number of linear systems solved, i.e. number of attempts for reducing error |
---|
| 841 | */ |
---|
| 842 | ) |
---|
| 843 | { |
---|
| 844 | int retval; |
---|
| 845 | struct wrap_str_data_ wdata; |
---|
| 846 | static void (*fjac)(double *p, struct sba_crsm *idxij, int *rcidxs, int *rcsubs, double *jac, void *adata); |
---|
| 847 | |
---|
| 848 | wdata.proj=proj; |
---|
| 849 | wdata.projac=projac; |
---|
| 850 | wdata.pnp=pnp; |
---|
| 851 | wdata.mnp=mnp; |
---|
| 852 | wdata.adata=adata; |
---|
| 853 | |
---|
| 854 | fjac=(projac)? sba_str_Qs_jac : sba_str_Qs_fdjac; |
---|
| 855 | retval=sba_str_levmar_x(n, m, vmask, p, pnp, x, mnp, sba_str_Qs, fjac, &wdata, itmax, verbose, opts, info); |
---|
| 856 | |
---|
| 857 | if(info){ |
---|
| 858 | register int i; |
---|
| 859 | int nvis; |
---|
| 860 | |
---|
| 861 | /* count visible image points */ |
---|
| 862 | for(i=nvis=0; i<n*m; ++i) |
---|
| 863 | nvis+=vmask[i]; |
---|
| 864 | |
---|
| 865 | /* each "func" & "fjac" evaluation requires nvis "proj" & "projac" evaluations */ |
---|
| 866 | info[7]*=nvis; |
---|
| 867 | info[8]*=nvis; |
---|
| 868 | } |
---|
| 869 | |
---|
| 870 | return retval; |
---|
| 871 | } |
---|