1 /*M///////////////////////////////////////////////////////////////////////////////////////
2 //
3 //  IMPORTANT: READ BEFORE DOWNLOADING, COPYING, INSTALLING OR USING.
4 //
5 //  By downloading, copying, installing or using the software you agree to this license.
6 //  If you do not agree to this license, do not download, install,
7 //  copy or use the software.
8 //
9 //
10 //                           License Agreement
11 //                For Open Source Computer Vision Library
12 //
13 // Copyright (C) 2000-2008, Intel Corporation, all rights reserved.
14 // Copyright (C) 2009-2011, Willow Garage Inc., all rights reserved.
15 // Third party copyrights are property of their respective owners.
16 //
17 // Redistribution and use in source and binary forms, with or without modification,
18 // are permitted provided that the following conditions are met:
19 //
20 //   * Redistribution's of source code must retain the above copyright notice,
21 //     this list of conditions and the following disclaimer.
22 //
23 //   * Redistribution's in binary form must reproduce the above copyright notice,
24 //     this list of conditions and the following disclaimer in the documentation
25 //     and/or other materials provided with the distribution.
26 //
27 //   * The name of the copyright holders may not be used to endorse or promote products
28 //     derived from this software without specific prior written permission.
29 //
30 // This software is provided by the copyright holders and contributors "as is" and
31 // any express or implied warranties, including, but not limited to, the implied
32 // warranties of merchantability and fitness for a particular purpose are disclaimed.
33 // In no event shall the Intel Corporation or contributors be liable for any direct,
34 // indirect, incidental, special, exemplary, or consequential damages
35 // (including, but not limited to, procurement of substitute goods or services;
36 // loss of use, data, or profits; or business interruption) however caused
37 // and on any theory of liability, whether in contract, strict liability,
38 // or tort (including negligence or otherwise) arising in any way out of
39 // the use of this software, even if advised of the possibility of such damage.
40 //
41 //M*/
42 
43 #include "precomp.hpp"
44 #include "fisheye.hpp"
45 #include <limits>
46 
47 namespace cv { namespace
48 {
49     struct JacobianRow
50     {
51         Vec2d df, dc;
52         Vec4d dk;
53         Vec3d dom, dT;
54         double dalpha;
55     };
56 
57     void subMatrix(const Mat& src, Mat& dst, const std::vector<uchar>& cols, const std::vector<uchar>& rows);
58 }}
59 
60 //////////////////////////////////////////////////////////////////////////////////////////////////////////////
61 /// cv::fisheye::projectPoints
62 
projectPoints(InputArray objectPoints,OutputArray imagePoints,const Affine3d & affine,InputArray K,InputArray D,double alpha,OutputArray jacobian)63 void cv::fisheye::projectPoints(InputArray objectPoints, OutputArray imagePoints, const Affine3d& affine,
64     InputArray K, InputArray D, double alpha, OutputArray jacobian)
65 {
66     CV_INSTRUMENT_REGION();
67 
68     projectPoints(objectPoints, imagePoints, affine.rvec(), affine.translation(), K, D, alpha, jacobian);
69 }
70 
projectPoints(InputArray objectPoints,OutputArray imagePoints,InputArray _rvec,InputArray _tvec,InputArray _K,InputArray _D,double alpha,OutputArray jacobian)71 void cv::fisheye::projectPoints(InputArray objectPoints, OutputArray imagePoints, InputArray _rvec,
72         InputArray _tvec, InputArray _K, InputArray _D, double alpha, OutputArray jacobian)
73 {
74     CV_INSTRUMENT_REGION();
75 
76     // will support only 3-channel data now for points
77     CV_Assert(objectPoints.type() == CV_32FC3 || objectPoints.type() == CV_64FC3);
78     imagePoints.create(objectPoints.size(), CV_MAKETYPE(objectPoints.depth(), 2));
79     size_t n = objectPoints.total();
80 
81     CV_Assert(_rvec.total() * _rvec.channels() == 3 && (_rvec.depth() == CV_32F || _rvec.depth() == CV_64F));
82     CV_Assert(_tvec.total() * _tvec.channels() == 3 && (_tvec.depth() == CV_32F || _tvec.depth() == CV_64F));
83     CV_Assert(_tvec.getMat().isContinuous() && _rvec.getMat().isContinuous());
84 
85     Vec3d om = _rvec.depth() == CV_32F ? (Vec3d)*_rvec.getMat().ptr<Vec3f>() : *_rvec.getMat().ptr<Vec3d>();
86     Vec3d T  = _tvec.depth() == CV_32F ? (Vec3d)*_tvec.getMat().ptr<Vec3f>() : *_tvec.getMat().ptr<Vec3d>();
87 
88     CV_Assert(_K.size() == Size(3,3) && (_K.type() == CV_32F || _K.type() == CV_64F) && _D.type() == _K.type() && _D.total() == 4);
89 
90     cv::Vec2d f, c;
91     if (_K.depth() == CV_32F)
92     {
93 
94         Matx33f K = _K.getMat();
95         f = Vec2f(K(0, 0), K(1, 1));
96         c = Vec2f(K(0, 2), K(1, 2));
97     }
98     else
99     {
100         Matx33d K = _K.getMat();
101         f = Vec2d(K(0, 0), K(1, 1));
102         c = Vec2d(K(0, 2), K(1, 2));
103     }
104 
105     Vec4d k = _D.depth() == CV_32F ? (Vec4d)*_D.getMat().ptr<Vec4f>(): *_D.getMat().ptr<Vec4d>();
106 
107     const bool isJacobianNeeded = jacobian.needed();
108     JacobianRow *Jn = 0;
109     if (isJacobianNeeded)
110     {
111         int nvars = 2 + 2 + 1 + 4 + 3 + 3; // f, c, alpha, k, om, T,
112         jacobian.create(2*(int)n, nvars, CV_64F);
113         Jn = jacobian.getMat().ptr<JacobianRow>(0);
114     }
115 
116     Matx33d R;
117     Matx<double, 3, 9> dRdom;
118     Rodrigues(om, R, dRdom);
119     Affine3d aff(om, T);
120 
121     const Vec3f* Xf = objectPoints.getMat().ptr<Vec3f>();
122     const Vec3d* Xd = objectPoints.getMat().ptr<Vec3d>();
123     Vec2f *xpf = imagePoints.getMat().ptr<Vec2f>();
124     Vec2d *xpd = imagePoints.getMat().ptr<Vec2d>();
125 
126     for(size_t i = 0; i < n; ++i)
127     {
128         Vec3d Xi = objectPoints.depth() == CV_32F ? (Vec3d)Xf[i] : Xd[i];
129         Vec3d Y = aff*Xi;
130         if (fabs(Y[2]) < DBL_MIN)
131             Y[2] = 1;
132         Vec2d x(Y[0]/Y[2], Y[1]/Y[2]);
133 
134         double r2 = x.dot(x);
135         double r = std::sqrt(r2);
136 
137         // Angle of the incoming ray:
138         double theta = atan(r);
139 
140         double theta2 = theta*theta, theta3 = theta2*theta, theta4 = theta2*theta2, theta5 = theta4*theta,
141                 theta6 = theta3*theta3, theta7 = theta6*theta, theta8 = theta4*theta4, theta9 = theta8*theta;
142 
143         double theta_d = theta + k[0]*theta3 + k[1]*theta5 + k[2]*theta7 + k[3]*theta9;
144 
145         double inv_r = r > 1e-8 ? 1.0/r : 1;
146         double cdist = r > 1e-8 ? theta_d * inv_r : 1;
147 
148         Vec2d xd1 = x * cdist;
149         Vec2d xd3(xd1[0] + alpha*xd1[1], xd1[1]);
150         Vec2d final_point(xd3[0] * f[0] + c[0], xd3[1] * f[1] + c[1]);
151 
152         if (objectPoints.depth() == CV_32F)
153             xpf[i] = final_point;
154         else
155             xpd[i] = final_point;
156 
157         if (isJacobianNeeded)
158         {
159             //Vec3d Xi = pdepth == CV_32F ? (Vec3d)Xf[i] : Xd[i];
160             //Vec3d Y = aff*Xi;
161             double dYdR[] = { Xi[0], Xi[1], Xi[2], 0, 0, 0, 0, 0, 0,
162                               0, 0, 0, Xi[0], Xi[1], Xi[2], 0, 0, 0,
163                               0, 0, 0, 0, 0, 0, Xi[0], Xi[1], Xi[2] };
164 
165             Matx33d dYdom_data = Matx<double, 3, 9>(dYdR) * dRdom.t();
166             const Vec3d *dYdom = (Vec3d*)dYdom_data.val;
167 
168             Matx33d dYdT_data = Matx33d::eye();
169             const Vec3d *dYdT = (Vec3d*)dYdT_data.val;
170 
171             //Vec2d x(Y[0]/Y[2], Y[1]/Y[2]);
172             Vec3d dxdom[2];
173             dxdom[0] = (1.0/Y[2]) * dYdom[0] - x[0]/Y[2] * dYdom[2];
174             dxdom[1] = (1.0/Y[2]) * dYdom[1] - x[1]/Y[2] * dYdom[2];
175 
176             Vec3d dxdT[2];
177             dxdT[0]  = (1.0/Y[2]) * dYdT[0] - x[0]/Y[2] * dYdT[2];
178             dxdT[1]  = (1.0/Y[2]) * dYdT[1] - x[1]/Y[2] * dYdT[2];
179 
180             //double r2 = x.dot(x);
181             Vec3d dr2dom = 2 * x[0] * dxdom[0] + 2 * x[1] * dxdom[1];
182             Vec3d dr2dT  = 2 * x[0] *  dxdT[0] + 2 * x[1] *  dxdT[1];
183 
184             //double r = std::sqrt(r2);
185             double drdr2 = r > 1e-8 ? 1.0/(2*r) : 1;
186             Vec3d drdom = drdr2 * dr2dom;
187             Vec3d drdT  = drdr2 * dr2dT;
188 
189             // Angle of the incoming ray:
190             //double theta = atan(r);
191             double dthetadr = 1.0/(1+r2);
192             Vec3d dthetadom = dthetadr * drdom;
193             Vec3d dthetadT  = dthetadr *  drdT;
194 
195             //double theta_d = theta + k[0]*theta3 + k[1]*theta5 + k[2]*theta7 + k[3]*theta9;
196             double dtheta_ddtheta = 1 + 3*k[0]*theta2 + 5*k[1]*theta4 + 7*k[2]*theta6 + 9*k[3]*theta8;
197             Vec3d dtheta_ddom = dtheta_ddtheta * dthetadom;
198             Vec3d dtheta_ddT  = dtheta_ddtheta * dthetadT;
199             Vec4d dtheta_ddk  = Vec4d(theta3, theta5, theta7, theta9);
200 
201             //double inv_r = r > 1e-8 ? 1.0/r : 1;
202             //double cdist = r > 1e-8 ? theta_d / r : 1;
203             Vec3d dcdistdom = inv_r * (dtheta_ddom - cdist*drdom);
204             Vec3d dcdistdT  = inv_r * (dtheta_ddT  - cdist*drdT);
205             Vec4d dcdistdk  = inv_r *  dtheta_ddk;
206 
207             //Vec2d xd1 = x * cdist;
208             Vec4d dxd1dk[2];
209             Vec3d dxd1dom[2], dxd1dT[2];
210             dxd1dom[0] = x[0] * dcdistdom + cdist * dxdom[0];
211             dxd1dom[1] = x[1] * dcdistdom + cdist * dxdom[1];
212             dxd1dT[0]  = x[0] * dcdistdT  + cdist * dxdT[0];
213             dxd1dT[1]  = x[1] * dcdistdT  + cdist * dxdT[1];
214             dxd1dk[0]  = x[0] * dcdistdk;
215             dxd1dk[1]  = x[1] * dcdistdk;
216 
217             //Vec2d xd3(xd1[0] + alpha*xd1[1], xd1[1]);
218             Vec4d dxd3dk[2];
219             Vec3d dxd3dom[2], dxd3dT[2];
220             dxd3dom[0] = dxd1dom[0] + alpha * dxd1dom[1];
221             dxd3dom[1] = dxd1dom[1];
222             dxd3dT[0]  = dxd1dT[0]  + alpha * dxd1dT[1];
223             dxd3dT[1]  = dxd1dT[1];
224             dxd3dk[0]  = dxd1dk[0]  + alpha * dxd1dk[1];
225             dxd3dk[1]  = dxd1dk[1];
226 
227             Vec2d dxd3dalpha(xd1[1], 0);
228 
229             //final jacobian
230             Jn[0].dom = f[0] * dxd3dom[0];
231             Jn[1].dom = f[1] * dxd3dom[1];
232 
233             Jn[0].dT = f[0] * dxd3dT[0];
234             Jn[1].dT = f[1] * dxd3dT[1];
235 
236             Jn[0].dk = f[0] * dxd3dk[0];
237             Jn[1].dk = f[1] * dxd3dk[1];
238 
239             Jn[0].dalpha = f[0] * dxd3dalpha[0];
240             Jn[1].dalpha = 0; //f[1] * dxd3dalpha[1];
241 
242             Jn[0].df = Vec2d(xd3[0], 0);
243             Jn[1].df = Vec2d(0, xd3[1]);
244 
245             Jn[0].dc = Vec2d(1, 0);
246             Jn[1].dc = Vec2d(0, 1);
247 
248             //step to jacobian rows for next point
249             Jn += 2;
250         }
251     }
252 }
253 
254 //////////////////////////////////////////////////////////////////////////////////////////////////////////////
255 /// cv::fisheye::distortPoints
256 
distortPoints(InputArray undistorted,OutputArray distorted,InputArray K,InputArray D,double alpha)257 void cv::fisheye::distortPoints(InputArray undistorted, OutputArray distorted, InputArray K, InputArray D, double alpha)
258 {
259     CV_INSTRUMENT_REGION();
260 
261     // will support only 2-channel data now for points
262     CV_Assert(undistorted.type() == CV_32FC2 || undistorted.type() == CV_64FC2);
263     distorted.create(undistorted.size(), undistorted.type());
264     size_t n = undistorted.total();
265 
266     CV_Assert(K.size() == Size(3,3) && (K.type() == CV_32F || K.type() == CV_64F) && D.total() == 4);
267 
268     cv::Vec2d f, c;
269     if (K.depth() == CV_32F)
270     {
271         Matx33f camMat = K.getMat();
272         f = Vec2f(camMat(0, 0), camMat(1, 1));
273         c = Vec2f(camMat(0, 2), camMat(1, 2));
274     }
275     else
276     {
277         Matx33d camMat = K.getMat();
278         f = Vec2d(camMat(0, 0), camMat(1, 1));
279         c = Vec2d(camMat(0 ,2), camMat(1, 2));
280     }
281 
282     Vec4d k = D.depth() == CV_32F ? (Vec4d)*D.getMat().ptr<Vec4f>(): *D.getMat().ptr<Vec4d>();
283 
284     const Vec2f* Xf = undistorted.getMat().ptr<Vec2f>();
285     const Vec2d* Xd = undistorted.getMat().ptr<Vec2d>();
286     Vec2f *xpf = distorted.getMat().ptr<Vec2f>();
287     Vec2d *xpd = distorted.getMat().ptr<Vec2d>();
288 
289     for(size_t i = 0; i < n; ++i)
290     {
291         Vec2d x = undistorted.depth() == CV_32F ? (Vec2d)Xf[i] : Xd[i];
292 
293         double r2 = x.dot(x);
294         double r = std::sqrt(r2);
295 
296         // Angle of the incoming ray:
297         double theta = atan(r);
298 
299         double theta2 = theta*theta, theta3 = theta2*theta, theta4 = theta2*theta2, theta5 = theta4*theta,
300                 theta6 = theta3*theta3, theta7 = theta6*theta, theta8 = theta4*theta4, theta9 = theta8*theta;
301 
302         double theta_d = theta + k[0]*theta3 + k[1]*theta5 + k[2]*theta7 + k[3]*theta9;
303 
304         double inv_r = r > 1e-8 ? 1.0/r : 1;
305         double cdist = r > 1e-8 ? theta_d * inv_r : 1;
306 
307         Vec2d xd1 = x * cdist;
308         Vec2d xd3(xd1[0] + alpha*xd1[1], xd1[1]);
309         Vec2d final_point(xd3[0] * f[0] + c[0], xd3[1] * f[1] + c[1]);
310 
311         if (undistorted.depth() == CV_32F)
312             xpf[i] = final_point;
313         else
314             xpd[i] = final_point;
315     }
316 }
317 
318 //////////////////////////////////////////////////////////////////////////////////////////////////////////////
319 /// cv::fisheye::undistortPoints
320 
undistortPoints(InputArray distorted,OutputArray undistorted,InputArray K,InputArray D,InputArray R,InputArray P)321 void cv::fisheye::undistortPoints( InputArray distorted, OutputArray undistorted, InputArray K, InputArray D, InputArray R, InputArray P)
322 {
323     CV_INSTRUMENT_REGION();
324 
325     // will support only 2-channel data now for points
326     CV_Assert(distorted.type() == CV_32FC2 || distorted.type() == CV_64FC2);
327     undistorted.create(distorted.size(), distorted.type());
328 
329     CV_Assert(P.empty() || P.size() == Size(3, 3) || P.size() == Size(4, 3));
330     CV_Assert(R.empty() || R.size() == Size(3, 3) || R.total() * R.channels() == 3);
331     CV_Assert(D.total() == 4 && K.size() == Size(3, 3) && (K.depth() == CV_32F || K.depth() == CV_64F));
332 
333     cv::Vec2d f, c;
334     if (K.depth() == CV_32F)
335     {
336         Matx33f camMat = K.getMat();
337         f = Vec2f(camMat(0, 0), camMat(1, 1));
338         c = Vec2f(camMat(0, 2), camMat(1, 2));
339     }
340     else
341     {
342         Matx33d camMat = K.getMat();
343         f = Vec2d(camMat(0, 0), camMat(1, 1));
344         c = Vec2d(camMat(0, 2), camMat(1, 2));
345     }
346 
347     Vec4d k = D.depth() == CV_32F ? (Vec4d)*D.getMat().ptr<Vec4f>(): *D.getMat().ptr<Vec4d>();
348 
349     cv::Matx33d RR = cv::Matx33d::eye();
350     if (!R.empty() && R.total() * R.channels() == 3)
351     {
352         cv::Vec3d rvec;
353         R.getMat().convertTo(rvec, CV_64F);
354         RR = cv::Affine3d(rvec).rotation();
355     }
356     else if (!R.empty() && R.size() == Size(3, 3))
357         R.getMat().convertTo(RR, CV_64F);
358 
359     if(!P.empty())
360     {
361         cv::Matx33d PP;
362         P.getMat().colRange(0, 3).convertTo(PP, CV_64F);
363         RR = PP * RR;
364     }
365 
366     // start undistorting
367     const cv::Vec2f* srcf = distorted.getMat().ptr<cv::Vec2f>();
368     const cv::Vec2d* srcd = distorted.getMat().ptr<cv::Vec2d>();
369     cv::Vec2f* dstf = undistorted.getMat().ptr<cv::Vec2f>();
370     cv::Vec2d* dstd = undistorted.getMat().ptr<cv::Vec2d>();
371 
372     size_t n = distorted.total();
373     int sdepth = distorted.depth();
374 
375     for(size_t i = 0; i < n; i++ )
376     {
377         Vec2d pi = sdepth == CV_32F ? (Vec2d)srcf[i] : srcd[i];  // image point
378         Vec2d pw((pi[0] - c[0])/f[0], (pi[1] - c[1])/f[1]);      // world point
379 
380         double theta_d = sqrt(pw[0]*pw[0] + pw[1]*pw[1]);
381 
382         // the current camera model is only valid up to 180 FOV
383         // for larger FOV the loop below does not converge
384         // clip values so we still get plausible results for super fisheye images > 180 grad
385         theta_d = min(max(-CV_PI/2., theta_d), CV_PI/2.);
386 
387         bool converged = false;
388         double theta = theta_d;
389 
390         double scale = 0.0;
391 
392         if (fabs(theta_d) > 1e-8)
393         {
394             // compensate distortion iteratively
395 
396             const double EPS = 1e-8; // or std::numeric_limits<double>::epsilon();
397 
398             for (int j = 0; j < 10; j++)
399             {
400                 double theta2 = theta*theta, theta4 = theta2*theta2, theta6 = theta4*theta2, theta8 = theta6*theta2;
401                 double k0_theta2 = k[0] * theta2, k1_theta4 = k[1] * theta4, k2_theta6 = k[2] * theta6, k3_theta8 = k[3] * theta8;
402                 /* new_theta = theta - theta_fix, theta_fix = f0(theta) / f0'(theta) */
403                 double theta_fix = (theta * (1 + k0_theta2 + k1_theta4 + k2_theta6 + k3_theta8) - theta_d) /
404                                    (1 + 3*k0_theta2 + 5*k1_theta4 + 7*k2_theta6 + 9*k3_theta8);
405                 theta = theta - theta_fix;
406                 if (fabs(theta_fix) < EPS)
407                 {
408                     converged = true;
409                     break;
410                 }
411             }
412 
413             scale = std::tan(theta) / theta_d;
414         }
415         else
416         {
417             converged = true;
418         }
419 
420         // theta is monotonously increasing or decreasing depending on the sign of theta
421         // if theta has flipped, it might converge due to symmetry but on the opposite of the camera center
422         // so we can check whether theta has changed the sign during the optimization
423         bool theta_flipped = ((theta_d < 0 && theta > 0) || (theta_d > 0 && theta < 0));
424 
425         if (converged && !theta_flipped)
426         {
427             Vec2d pu = pw * scale; //undistorted point
428 
429             // reproject
430             Vec3d pr = RR * Vec3d(pu[0], pu[1], 1.0); // rotated point optionally multiplied by new camera matrix
431             Vec2d fi(pr[0]/pr[2], pr[1]/pr[2]);       // final
432 
433             if( sdepth == CV_32F )
434                 dstf[i] = fi;
435             else
436                 dstd[i] = fi;
437         }
438         else
439         {
440             // Vec2d fi(std::numeric_limits<double>::quiet_NaN(), std::numeric_limits<double>::quiet_NaN());
441             Vec2d fi(-1000000.0, -1000000.0);
442 
443             if( sdepth == CV_32F )
444                 dstf[i] = fi;
445             else
446                 dstd[i] = fi;
447         }
448     }
449 }
450 
451 //////////////////////////////////////////////////////////////////////////////////////////////////////////////
452 /// cv::fisheye::undistortPoints
453 
initUndistortRectifyMap(InputArray K,InputArray D,InputArray R,InputArray P,const cv::Size & size,int m1type,OutputArray map1,OutputArray map2)454 void cv::fisheye::initUndistortRectifyMap( InputArray K, InputArray D, InputArray R, InputArray P,
455     const cv::Size& size, int m1type, OutputArray map1, OutputArray map2 )
456 {
457     CV_INSTRUMENT_REGION();
458 
459     CV_Assert( m1type == CV_16SC2 || m1type == CV_32F || m1type <=0 );
460     map1.create( size, m1type <= 0 ? CV_16SC2 : m1type );
461     map2.create( size, map1.type() == CV_16SC2 ? CV_16UC1 : CV_32F );
462 
463     CV_Assert((K.depth() == CV_32F || K.depth() == CV_64F) && (D.depth() == CV_32F || D.depth() == CV_64F));
464     CV_Assert((P.empty() || P.depth() == CV_32F || P.depth() == CV_64F) && (R.empty() || R.depth() == CV_32F || R.depth() == CV_64F));
465     CV_Assert(K.size() == Size(3, 3) && (D.empty() || D.total() == 4));
466     CV_Assert(R.empty() || R.size() == Size(3, 3) || R.total() * R.channels() == 3);
467     CV_Assert(P.empty() || P.size() == Size(3, 3) || P.size() == Size(4, 3));
468 
469     cv::Vec2d f, c;
470     if (K.depth() == CV_32F)
471     {
472         Matx33f camMat = K.getMat();
473         f = Vec2f(camMat(0, 0), camMat(1, 1));
474         c = Vec2f(camMat(0, 2), camMat(1, 2));
475     }
476     else
477     {
478         Matx33d camMat = K.getMat();
479         f = Vec2d(camMat(0, 0), camMat(1, 1));
480         c = Vec2d(camMat(0, 2), camMat(1, 2));
481     }
482 
483     Vec4d k = Vec4d::all(0);
484     if (!D.empty())
485         k = D.depth() == CV_32F ? (Vec4d)*D.getMat().ptr<Vec4f>(): *D.getMat().ptr<Vec4d>();
486 
487     cv::Matx33d RR  = cv::Matx33d::eye();
488     if (!R.empty() && R.total() * R.channels() == 3)
489     {
490         cv::Vec3d rvec;
491         R.getMat().convertTo(rvec, CV_64F);
492         RR = Affine3d(rvec).rotation();
493     }
494     else if (!R.empty() && R.size() == Size(3, 3))
495         R.getMat().convertTo(RR, CV_64F);
496 
497     cv::Matx33d PP = cv::Matx33d::eye();
498     if (!P.empty())
499         P.getMat().colRange(0, 3).convertTo(PP, CV_64F);
500 
501     cv::Matx33d iR = (PP * RR).inv(cv::DECOMP_SVD);
502 
503     for( int i = 0; i < size.height; ++i)
504     {
505         float* m1f = map1.getMat().ptr<float>(i);
506         float* m2f = map2.getMat().ptr<float>(i);
507         short*  m1 = (short*)m1f;
508         ushort* m2 = (ushort*)m2f;
509 
510         double _x = i*iR(0, 1) + iR(0, 2),
511                _y = i*iR(1, 1) + iR(1, 2),
512                _w = i*iR(2, 1) + iR(2, 2);
513 
514         for( int j = 0; j < size.width; ++j)
515         {
516             double u, v;
517             if( _w <= 0)
518             {
519                 u = (_x > 0) ? -std::numeric_limits<double>::infinity() : std::numeric_limits<double>::infinity();
520                 v = (_y > 0) ? -std::numeric_limits<double>::infinity() : std::numeric_limits<double>::infinity();
521             }
522             else
523             {
524                 double x = _x/_w, y = _y/_w;
525 
526                 double r = sqrt(x*x + y*y);
527                 double theta = atan(r);
528 
529                 double theta2 = theta*theta, theta4 = theta2*theta2, theta6 = theta4*theta2, theta8 = theta4*theta4;
530                 double theta_d = theta * (1 + k[0]*theta2 + k[1]*theta4 + k[2]*theta6 + k[3]*theta8);
531 
532                 double scale = (r == 0) ? 1.0 : theta_d / r;
533                 u = f[0]*x*scale + c[0];
534                 v = f[1]*y*scale + c[1];
535             }
536 
537             if( m1type == CV_16SC2 )
538             {
539                 int iu = cv::saturate_cast<int>(u*cv::INTER_TAB_SIZE);
540                 int iv = cv::saturate_cast<int>(v*cv::INTER_TAB_SIZE);
541                 m1[j*2+0] = (short)(iu >> cv::INTER_BITS);
542                 m1[j*2+1] = (short)(iv >> cv::INTER_BITS);
543                 m2[j] = (ushort)((iv & (cv::INTER_TAB_SIZE-1))*cv::INTER_TAB_SIZE + (iu & (cv::INTER_TAB_SIZE-1)));
544             }
545             else if( m1type == CV_32FC1 )
546             {
547                 m1f[j] = (float)u;
548                 m2f[j] = (float)v;
549             }
550 
551             _x += iR(0, 0);
552             _y += iR(1, 0);
553             _w += iR(2, 0);
554         }
555     }
556 }
557 
558 //////////////////////////////////////////////////////////////////////////////////////////////////////////////
559 /// cv::fisheye::undistortImage
560 
undistortImage(InputArray distorted,OutputArray undistorted,InputArray K,InputArray D,InputArray Knew,const Size & new_size)561 void cv::fisheye::undistortImage(InputArray distorted, OutputArray undistorted,
562         InputArray K, InputArray D, InputArray Knew, const Size& new_size)
563 {
564     CV_INSTRUMENT_REGION();
565 
566     Size size = !new_size.empty() ? new_size : distorted.size();
567 
568     cv::Mat map1, map2;
569     fisheye::initUndistortRectifyMap(K, D, cv::Matx33d::eye(), Knew, size, CV_16SC2, map1, map2 );
570     cv::remap(distorted, undistorted, map1, map2, INTER_LINEAR, BORDER_CONSTANT);
571 }
572 
573 
574 //////////////////////////////////////////////////////////////////////////////////////////////////////////////
575 /// cv::fisheye::estimateNewCameraMatrixForUndistortRectify
576 
estimateNewCameraMatrixForUndistortRectify(InputArray K,InputArray D,const Size & image_size,InputArray R,OutputArray P,double balance,const Size & new_size,double fov_scale)577 void cv::fisheye::estimateNewCameraMatrixForUndistortRectify(InputArray K, InputArray D, const Size &image_size, InputArray R,
578     OutputArray P, double balance, const Size& new_size, double fov_scale)
579 {
580     CV_INSTRUMENT_REGION();
581 
582     CV_Assert( K.size() == Size(3, 3)       && (K.depth() == CV_32F || K.depth() == CV_64F));
583     CV_Assert(D.empty() || ((D.total() == 4) && (D.depth() == CV_32F || D.depth() == CV_64F)));
584 
585     int w = image_size.width, h = image_size.height;
586     balance = std::min(std::max(balance, 0.0), 1.0);
587 
588     cv::Mat points(1, 4, CV_64FC2);
589     Vec2d* pptr = points.ptr<Vec2d>();
590     pptr[0] = Vec2d(w/2, 0);
591     pptr[1] = Vec2d(w, h/2);
592     pptr[2] = Vec2d(w/2, h);
593     pptr[3] = Vec2d(0, h/2);
594 
595     fisheye::undistortPoints(points, points, K, D, R);
596     cv::Scalar center_mass = mean(points);
597     cv::Vec2d cn(center_mass.val);
598 
599     double aspect_ratio = (K.depth() == CV_32F) ? K.getMat().at<float >(0,0)/K.getMat().at<float> (1,1)
600                                                 : K.getMat().at<double>(0,0)/K.getMat().at<double>(1,1);
601 
602     // convert to identity ratio
603     cn[0] *= aspect_ratio;
604     for(size_t i = 0; i < points.total(); ++i)
605         pptr[i][1] *= aspect_ratio;
606 
607     double minx = DBL_MAX, miny = DBL_MAX, maxx = -DBL_MAX, maxy = -DBL_MAX;
608     for(size_t i = 0; i < points.total(); ++i)
609     {
610         miny = std::min(miny, pptr[i][1]);
611         maxy = std::max(maxy, pptr[i][1]);
612         minx = std::min(minx, pptr[i][0]);
613         maxx = std::max(maxx, pptr[i][0]);
614     }
615 
616     double f1 = w * 0.5/(cn[0] - minx);
617     double f2 = w * 0.5/(maxx - cn[0]);
618     double f3 = h * 0.5 * aspect_ratio/(cn[1] - miny);
619     double f4 = h * 0.5 * aspect_ratio/(maxy - cn[1]);
620 
621     double fmin = std::min(f1, std::min(f2, std::min(f3, f4)));
622     double fmax = std::max(f1, std::max(f2, std::max(f3, f4)));
623 
624     double f = balance * fmin + (1.0 - balance) * fmax;
625     f *= fov_scale > 0 ? 1.0/fov_scale : 1.0;
626 
627     cv::Vec2d new_f(f, f), new_c = -cn * f + Vec2d(w, h * aspect_ratio) * 0.5;
628 
629     // restore aspect ratio
630     new_f[1] /= aspect_ratio;
631     new_c[1] /= aspect_ratio;
632 
633     if (!new_size.empty())
634     {
635         double rx = new_size.width /(double)image_size.width;
636         double ry = new_size.height/(double)image_size.height;
637 
638         new_f[0] *= rx;  new_f[1] *= ry;
639         new_c[0] *= rx;  new_c[1] *= ry;
640     }
641 
642     Mat(Matx33d(new_f[0], 0, new_c[0],
643                 0, new_f[1], new_c[1],
644                 0,        0,       1)).convertTo(P, P.empty() ? K.type() : P.type());
645 }
646 
647 
648 //////////////////////////////////////////////////////////////////////////////////////////////////////////////
649 /// cv::fisheye::stereoRectify
650 
stereoRectify(InputArray K1,InputArray D1,InputArray K2,InputArray D2,const Size & imageSize,InputArray _R,InputArray _tvec,OutputArray R1,OutputArray R2,OutputArray P1,OutputArray P2,OutputArray Q,int flags,const Size & newImageSize,double balance,double fov_scale)651 void cv::fisheye::stereoRectify( InputArray K1, InputArray D1, InputArray K2, InputArray D2, const Size& imageSize,
652         InputArray _R, InputArray _tvec, OutputArray R1, OutputArray R2, OutputArray P1, OutputArray P2,
653         OutputArray Q, int flags, const Size& newImageSize, double balance, double fov_scale)
654 {
655     CV_INSTRUMENT_REGION();
656 
657     CV_Assert((_R.size() == Size(3, 3) || _R.total() * _R.channels() == 3) && (_R.depth() == CV_32F || _R.depth() == CV_64F));
658     CV_Assert(_tvec.total() * _tvec.channels() == 3 && (_tvec.depth() == CV_32F || _tvec.depth() == CV_64F));
659 
660 
661     cv::Mat aaa = _tvec.getMat().reshape(3, 1);
662 
663     Vec3d rvec; // Rodrigues vector
664     if (_R.size() == Size(3, 3))
665     {
666         cv::Matx33d rmat;
667         _R.getMat().convertTo(rmat, CV_64F);
668         rvec = Affine3d(rmat).rvec();
669     }
670     else if (_R.total() * _R.channels() == 3)
671         _R.getMat().convertTo(rvec, CV_64F);
672 
673     Vec3d tvec;
674     _tvec.getMat().convertTo(tvec, CV_64F);
675 
676     // rectification algorithm
677     rvec *= -0.5;              // get average rotation
678 
679     Matx33d r_r;
680     Rodrigues(rvec, r_r);  // rotate cameras to same orientation by averaging
681 
682     Vec3d t = r_r * tvec;
683     Vec3d uu(t[0] > 0 ? 1 : -1, 0, 0);
684 
685     // calculate global Z rotation
686     Vec3d ww = t.cross(uu);
687     double nw = norm(ww);
688     if (nw > 0.0)
689         ww *= acos(fabs(t[0])/cv::norm(t))/nw;
690 
691     Matx33d wr;
692     Rodrigues(ww, wr);
693 
694     // apply to both views
695     Matx33d ri1 = wr * r_r.t();
696     Mat(ri1, false).convertTo(R1, R1.empty() ? CV_64F : R1.type());
697     Matx33d ri2 = wr * r_r;
698     Mat(ri2, false).convertTo(R2, R2.empty() ? CV_64F : R2.type());
699     Vec3d tnew = ri2 * tvec;
700 
701     // calculate projection/camera matrices. these contain the relevant rectified image internal params (fx, fy=fx, cx, cy)
702     Matx33d newK1, newK2;
703     estimateNewCameraMatrixForUndistortRectify(K1, D1, imageSize, R1, newK1, balance, newImageSize, fov_scale);
704     estimateNewCameraMatrixForUndistortRectify(K2, D2, imageSize, R2, newK2, balance, newImageSize, fov_scale);
705 
706     double fc_new = std::min(newK1(1,1), newK2(1,1));
707     Point2d cc_new[2] = { Vec2d(newK1(0, 2), newK1(1, 2)), Vec2d(newK2(0, 2), newK2(1, 2)) };
708 
709     // Vertical focal length must be the same for both images to keep the epipolar constraint use fy for fx also.
710     // For simplicity, set the principal points for both cameras to be the average
711     // of the two principal points (either one of or both x- and y- coordinates)
712     if( flags & cv::CALIB_ZERO_DISPARITY )
713         cc_new[0] = cc_new[1] = (cc_new[0] + cc_new[1]) * 0.5;
714     else
715         cc_new[0].y = cc_new[1].y = (cc_new[0].y + cc_new[1].y)*0.5;
716 
717     Mat(Matx34d(fc_new, 0, cc_new[0].x, 0,
718                 0, fc_new, cc_new[0].y, 0,
719                 0,      0,           1, 0), false).convertTo(P1, P1.empty() ? CV_64F : P1.type());
720 
721     Mat(Matx34d(fc_new, 0, cc_new[1].x, tnew[0]*fc_new, // baseline * focal length;,
722                 0, fc_new, cc_new[1].y,              0,
723                 0,      0,           1,              0), false).convertTo(P2, P2.empty() ? CV_64F : P2.type());
724 
725     if (Q.needed())
726         Mat(Matx44d(1, 0, 0,           -cc_new[0].x,
727                     0, 1, 0,           -cc_new[0].y,
728                     0, 0, 0,            fc_new,
729                     0, 0, -1./tnew[0], (cc_new[0].x - cc_new[1].x)/tnew[0]), false).convertTo(Q, Q.empty() ? CV_64F : Q.depth());
730 }
731 
732 //////////////////////////////////////////////////////////////////////////////////////////////////////////////
733 /// cv::fisheye::calibrate
734 
calibrate(InputArrayOfArrays objectPoints,InputArrayOfArrays imagePoints,const Size & image_size,InputOutputArray K,InputOutputArray D,OutputArrayOfArrays rvecs,OutputArrayOfArrays tvecs,int flags,cv::TermCriteria criteria)735 double cv::fisheye::calibrate(InputArrayOfArrays objectPoints, InputArrayOfArrays imagePoints, const Size& image_size,
736                                     InputOutputArray K, InputOutputArray D, OutputArrayOfArrays rvecs, OutputArrayOfArrays tvecs,
737                                     int flags , cv::TermCriteria criteria)
738 {
739     CV_INSTRUMENT_REGION();
740 
741     CV_Assert(!objectPoints.empty() && !imagePoints.empty() && objectPoints.total() == imagePoints.total());
742     CV_Assert(objectPoints.type() == CV_32FC3 || objectPoints.type() == CV_64FC3);
743     CV_Assert(imagePoints.type() == CV_32FC2 || imagePoints.type() == CV_64FC2);
744     CV_Assert(K.empty() || (K.size() == Size(3,3)));
745     CV_Assert(D.empty() || (D.total() == 4));
746     CV_Assert(rvecs.empty() || (rvecs.channels() == 3));
747     CV_Assert(tvecs.empty() || (tvecs.channels() == 3));
748 
749     CV_Assert((!K.empty() && !D.empty()) || !(flags & CALIB_USE_INTRINSIC_GUESS));
750 
751     using namespace cv::internal;
752     //-------------------------------Initialization
753     IntrinsicParams finalParam;
754     IntrinsicParams currentParam;
755     IntrinsicParams errors;
756 
757     finalParam.isEstimate[0] = flags & CALIB_FIX_FOCAL_LENGTH ? 0 : 1;
758     finalParam.isEstimate[1] = flags & CALIB_FIX_FOCAL_LENGTH ? 0 : 1;
759     finalParam.isEstimate[2] = flags & CALIB_FIX_PRINCIPAL_POINT ? 0 : 1;
760     finalParam.isEstimate[3] = flags & CALIB_FIX_PRINCIPAL_POINT ? 0 : 1;
761     finalParam.isEstimate[4] = flags & CALIB_FIX_SKEW ? 0 : 1;
762     finalParam.isEstimate[5] = flags & CALIB_FIX_K1 ? 0 : 1;
763     finalParam.isEstimate[6] = flags & CALIB_FIX_K2 ? 0 : 1;
764     finalParam.isEstimate[7] = flags & CALIB_FIX_K3 ? 0 : 1;
765     finalParam.isEstimate[8] = flags & CALIB_FIX_K4 ? 0 : 1;
766 
767     const int recompute_extrinsic = flags & CALIB_RECOMPUTE_EXTRINSIC ? 1: 0;
768     const int check_cond = flags & CALIB_CHECK_COND ? 1 : 0;
769 
770     const double alpha_smooth = 0.4;
771     const double thresh_cond = 1e6;
772     double change = 1;
773     Vec2d err_std;
774 
775     Matx33d _K;
776     Vec4d _D;
777     if (flags & CALIB_USE_INTRINSIC_GUESS)
778     {
779         K.getMat().convertTo(_K, CV_64FC1);
780         D.getMat().convertTo(_D, CV_64FC1);
781         finalParam.Init(Vec2d(_K(0,0), _K(1, 1)),
782                         Vec2d(_K(0,2), _K(1, 2)),
783                         Vec4d(flags & CALIB_FIX_K1 ? 0 : _D[0],
784                               flags & CALIB_FIX_K2 ? 0 : _D[1],
785                               flags & CALIB_FIX_K3 ? 0 : _D[2],
786                               flags & CALIB_FIX_K4 ? 0 : _D[3]),
787                         _K(0, 1) / _K(0, 0));
788     }
789     else
790     {
791         finalParam.Init(Vec2d(max(image_size.width, image_size.height) / CV_PI, max(image_size.width, image_size.height) / CV_PI),
792                         Vec2d(image_size.width  / 2.0 - 0.5, image_size.height / 2.0 - 0.5));
793     }
794 
795     errors.isEstimate = finalParam.isEstimate;
796 
797     std::vector<Vec3d> omc(objectPoints.total()), Tc(objectPoints.total());
798 
799     CalibrateExtrinsics(objectPoints, imagePoints, finalParam, check_cond, thresh_cond, omc, Tc);
800 
801 
802     //-------------------------------Optimization
803     for(int iter = 0; iter < std::numeric_limits<int>::max(); ++iter)
804     {
805         if ((criteria.type == 1 && iter >= criteria.maxCount)  ||
806             (criteria.type == 2 && change <= criteria.epsilon) ||
807             (criteria.type == 3 && (change <= criteria.epsilon || iter >= criteria.maxCount)))
808             break;
809 
810         double alpha_smooth2 = 1 - std::pow(1 - alpha_smooth, iter + 1.0);
811 
812         Mat JJ2, ex3;
813         ComputeJacobians(objectPoints, imagePoints, finalParam, omc, Tc, check_cond,thresh_cond, JJ2, ex3);
814 
815         Mat G;
816         solve(JJ2, ex3, G);
817         currentParam = finalParam + alpha_smooth2*G;
818 
819         change = norm(Vec4d(currentParam.f[0], currentParam.f[1], currentParam.c[0], currentParam.c[1]) -
820                 Vec4d(finalParam.f[0], finalParam.f[1], finalParam.c[0], finalParam.c[1]))
821                 / norm(Vec4d(currentParam.f[0], currentParam.f[1], currentParam.c[0], currentParam.c[1]));
822 
823         finalParam = currentParam;
824 
825         if (recompute_extrinsic)
826         {
827             CalibrateExtrinsics(objectPoints,  imagePoints, finalParam, check_cond,
828                                     thresh_cond, omc, Tc);
829         }
830     }
831 
832     //-------------------------------Validation
833     double rms;
834     EstimateUncertainties(objectPoints, imagePoints, finalParam,  omc, Tc, errors, err_std, thresh_cond,
835                               check_cond, rms);
836 
837     //-------------------------------
838     _K = Matx33d(finalParam.f[0], finalParam.f[0] * finalParam.alpha, finalParam.c[0],
839             0,                    finalParam.f[1], finalParam.c[1],
840             0,                                  0,               1);
841 
842     if (K.needed()) cv::Mat(_K).convertTo(K, K.empty() ? CV_64FC1 : K.type());
843     if (D.needed()) cv::Mat(finalParam.k).convertTo(D, D.empty() ? CV_64FC1 : D.type());
844     if (rvecs.isMatVector())
845     {
846         int N = (int)objectPoints.total();
847 
848         if(rvecs.empty())
849             rvecs.create(N, 1, CV_64FC3);
850 
851         if(tvecs.empty())
852             tvecs.create(N, 1, CV_64FC3);
853 
854         for(int i = 0; i < N; i++ )
855         {
856             rvecs.create(3, 1, CV_64F, i, true);
857             tvecs.create(3, 1, CV_64F, i, true);
858             memcpy(rvecs.getMat(i).ptr(), omc[i].val, sizeof(Vec3d));
859             memcpy(tvecs.getMat(i).ptr(), Tc[i].val, sizeof(Vec3d));
860         }
861     }
862     else
863     {
864         if (rvecs.needed()) cv::Mat(omc).convertTo(rvecs, rvecs.empty() ? CV_64FC3 : rvecs.type());
865         if (tvecs.needed()) cv::Mat(Tc).convertTo(tvecs, tvecs.empty() ? CV_64FC3 : tvecs.type());
866     }
867 
868     return rms;
869 }
870 
871 //////////////////////////////////////////////////////////////////////////////////////////////////////////////
872 /// cv::fisheye::stereoCalibrate
873 
stereoCalibrate(InputArrayOfArrays objectPoints,InputArrayOfArrays imagePoints1,InputArrayOfArrays imagePoints2,InputOutputArray K1,InputOutputArray D1,InputOutputArray K2,InputOutputArray D2,Size imageSize,OutputArray R,OutputArray T,int flags,TermCriteria criteria)874 double cv::fisheye::stereoCalibrate(InputArrayOfArrays objectPoints, InputArrayOfArrays imagePoints1, InputArrayOfArrays imagePoints2,
875                                     InputOutputArray K1, InputOutputArray D1, InputOutputArray K2, InputOutputArray D2, Size imageSize,
876                                     OutputArray R, OutputArray T, int flags, TermCriteria criteria)
877 {
878     CV_INSTRUMENT_REGION();
879 
880     CV_Assert(!objectPoints.empty() && !imagePoints1.empty() && !imagePoints2.empty());
881     CV_Assert(objectPoints.total() == imagePoints1.total() || imagePoints1.total() == imagePoints2.total());
882     CV_Assert(objectPoints.type() == CV_32FC3 || objectPoints.type() == CV_64FC3);
883     CV_Assert(imagePoints1.type() == CV_32FC2 || imagePoints1.type() == CV_64FC2);
884     CV_Assert(imagePoints2.type() == CV_32FC2 || imagePoints2.type() == CV_64FC2);
885 
886     CV_Assert(K1.empty() || (K1.size() == Size(3,3)));
887     CV_Assert(D1.empty() || (D1.total() == 4));
888     CV_Assert(K2.empty() || (K2.size() == Size(3,3)));
889     CV_Assert(D2.empty() || (D2.total() == 4));
890 
891     CV_Assert((!K1.empty() && !K2.empty() && !D1.empty() && !D2.empty()) || !(flags & CALIB_FIX_INTRINSIC));
892 
893     //-------------------------------Initialization
894 
895     const int threshold = 50;
896     const double thresh_cond = 1e6;
897     const int check_cond = 1;
898 
899     int n_points = (int)objectPoints.getMat(0).total();
900     int n_images = (int)objectPoints.total();
901 
902     double change = 1;
903 
904     cv::internal::IntrinsicParams intrinsicLeft;
905     cv::internal::IntrinsicParams intrinsicRight;
906 
907     cv::internal::IntrinsicParams intrinsicLeft_errors;
908     cv::internal::IntrinsicParams intrinsicRight_errors;
909 
910     Matx33d _K1, _K2;
911     Vec4d _D1, _D2;
912     if (!K1.empty()) K1.getMat().convertTo(_K1, CV_64FC1);
913     if (!D1.empty()) D1.getMat().convertTo(_D1, CV_64FC1);
914     if (!K2.empty()) K2.getMat().convertTo(_K2, CV_64FC1);
915     if (!D2.empty()) D2.getMat().convertTo(_D2, CV_64FC1);
916 
917     std::vector<Vec3d> rvecs1(n_images), tvecs1(n_images), rvecs2(n_images), tvecs2(n_images);
918 
919     if (!(flags & CALIB_FIX_INTRINSIC))
920     {
921         calibrate(objectPoints, imagePoints1, imageSize, _K1, _D1, rvecs1, tvecs1, flags, TermCriteria(3, 20, 1e-6));
922         calibrate(objectPoints, imagePoints2, imageSize, _K2, _D2, rvecs2, tvecs2, flags, TermCriteria(3, 20, 1e-6));
923     }
924 
925     intrinsicLeft.Init(Vec2d(_K1(0,0), _K1(1, 1)), Vec2d(_K1(0,2), _K1(1, 2)),
926                        Vec4d(_D1[0], _D1[1], _D1[2], _D1[3]), _K1(0, 1) / _K1(0, 0));
927 
928     intrinsicRight.Init(Vec2d(_K2(0,0), _K2(1, 1)), Vec2d(_K2(0,2), _K2(1, 2)),
929                         Vec4d(_D2[0], _D2[1], _D2[2], _D2[3]), _K2(0, 1) / _K2(0, 0));
930 
931     if ((flags & CALIB_FIX_INTRINSIC))
932     {
933         cv::internal::CalibrateExtrinsics(objectPoints,  imagePoints1, intrinsicLeft, check_cond, thresh_cond, rvecs1, tvecs1);
934         cv::internal::CalibrateExtrinsics(objectPoints,  imagePoints2, intrinsicRight, check_cond, thresh_cond, rvecs2, tvecs2);
935     }
936 
937     intrinsicLeft.isEstimate[0] = flags & CALIB_FIX_INTRINSIC ? 0 : 1;
938     intrinsicLeft.isEstimate[1] = flags & CALIB_FIX_INTRINSIC ? 0 : 1;
939     intrinsicLeft.isEstimate[2] = flags & CALIB_FIX_INTRINSIC ? 0 : 1;
940     intrinsicLeft.isEstimate[3] = flags & CALIB_FIX_INTRINSIC ? 0 : 1;
941     intrinsicLeft.isEstimate[4] = flags & (CALIB_FIX_SKEW | CALIB_FIX_INTRINSIC) ? 0 : 1;
942     intrinsicLeft.isEstimate[5] = flags & (CALIB_FIX_K1 | CALIB_FIX_INTRINSIC) ? 0 : 1;
943     intrinsicLeft.isEstimate[6] = flags & (CALIB_FIX_K2 | CALIB_FIX_INTRINSIC) ? 0 : 1;
944     intrinsicLeft.isEstimate[7] = flags & (CALIB_FIX_K3 | CALIB_FIX_INTRINSIC) ? 0 : 1;
945     intrinsicLeft.isEstimate[8] = flags & (CALIB_FIX_K4 | CALIB_FIX_INTRINSIC) ? 0 : 1;
946 
947     intrinsicRight.isEstimate[0] = flags & CALIB_FIX_INTRINSIC ? 0 : 1;
948     intrinsicRight.isEstimate[1] = flags & CALIB_FIX_INTRINSIC ? 0 : 1;
949     intrinsicRight.isEstimate[2] = flags & CALIB_FIX_INTRINSIC ? 0 : 1;
950     intrinsicRight.isEstimate[3] = flags & CALIB_FIX_INTRINSIC ? 0 : 1;
951     intrinsicRight.isEstimate[4] = flags & (CALIB_FIX_SKEW | CALIB_FIX_INTRINSIC) ? 0 : 1;
952     intrinsicRight.isEstimate[5] = flags & (CALIB_FIX_K1 | CALIB_FIX_INTRINSIC) ? 0 : 1;
953     intrinsicRight.isEstimate[6] = flags & (CALIB_FIX_K2 | CALIB_FIX_INTRINSIC) ? 0 : 1;
954     intrinsicRight.isEstimate[7] = flags & (CALIB_FIX_K3 | CALIB_FIX_INTRINSIC) ? 0 : 1;
955     intrinsicRight.isEstimate[8] = flags & (CALIB_FIX_K4 | CALIB_FIX_INTRINSIC) ? 0 : 1;
956 
957     intrinsicLeft_errors.isEstimate = intrinsicLeft.isEstimate;
958     intrinsicRight_errors.isEstimate = intrinsicRight.isEstimate;
959 
960     std::vector<uchar> selectedParams;
961     std::vector<uchar> tmp(6 * (n_images + 1), 1);
962     selectedParams.insert(selectedParams.end(), intrinsicLeft.isEstimate.begin(), intrinsicLeft.isEstimate.end());
963     selectedParams.insert(selectedParams.end(), intrinsicRight.isEstimate.begin(), intrinsicRight.isEstimate.end());
964     selectedParams.insert(selectedParams.end(), tmp.begin(), tmp.end());
965 
966     //Init values for rotation and translation between two views
967     cv::Mat om_list(1, n_images, CV_64FC3), T_list(1, n_images, CV_64FC3);
968     cv::Mat om_ref, R_ref, T_ref, R1, R2;
969     for (int image_idx = 0; image_idx < n_images; ++image_idx)
970     {
971         cv::Rodrigues(rvecs1[image_idx], R1);
972         cv::Rodrigues(rvecs2[image_idx], R2);
973         R_ref = R2 * R1.t();
974         T_ref = cv::Mat(tvecs2[image_idx]) - R_ref * cv::Mat(tvecs1[image_idx]);
975         cv::Rodrigues(R_ref, om_ref);
976         om_ref.reshape(3, 1).copyTo(om_list.col(image_idx));
977         T_ref.reshape(3, 1).copyTo(T_list.col(image_idx));
978     }
979     cv::Vec3d omcur = cv::internal::median3d(om_list);
980     cv::Vec3d Tcur  = cv::internal::median3d(T_list);
981 
982     cv::Mat J = cv::Mat::zeros(4 * n_points * n_images, 18 + 6 * (n_images + 1), CV_64FC1),
983             e = cv::Mat::zeros(4 * n_points * n_images, 1, CV_64FC1), Jkk, ekk;
984 
985     for(int iter = 0; ; ++iter)
986     {
987         if ((criteria.type == 1 && iter >= criteria.maxCount)  ||
988             (criteria.type == 2 && change <= criteria.epsilon) ||
989             (criteria.type == 3 && (change <= criteria.epsilon || iter >= criteria.maxCount)))
990             break;
991 
992         J.create(4 * n_points * n_images, 18 + 6 * (n_images + 1), CV_64FC1);
993         e.create(4 * n_points * n_images, 1, CV_64FC1);
994         Jkk.create(4 * n_points, 18 + 6 * (n_images + 1), CV_64FC1);
995         ekk.create(4 * n_points, 1, CV_64FC1);
996 
997         cv::Mat omr, Tr, domrdomckk, domrdTckk, domrdom, domrdT, dTrdomckk, dTrdTckk, dTrdom, dTrdT;
998 
999         for (int image_idx = 0; image_idx < n_images; ++image_idx)
1000         {
1001             Jkk = cv::Mat::zeros(4 * n_points, 18 + 6 * (n_images + 1), CV_64FC1);
1002 
1003             cv::Mat object  = objectPoints.getMat(image_idx).clone();
1004             cv::Mat imageLeft  = imagePoints1.getMat(image_idx).clone();
1005             cv::Mat imageRight  = imagePoints2.getMat(image_idx).clone();
1006             cv::Mat jacobians, projected;
1007 
1008             //left camera jacobian
1009             cv::Mat rvec = cv::Mat(rvecs1[image_idx]);
1010             cv::Mat tvec  = cv::Mat(tvecs1[image_idx]);
1011             cv::internal::projectPoints(object, projected, rvec, tvec, intrinsicLeft, jacobians);
1012             cv::Mat(cv::Mat((imageLeft - projected).t()).reshape(1, 1).t()).copyTo(ekk.rowRange(0, 2 * n_points));
1013             jacobians.colRange(8, 11).copyTo(Jkk.colRange(24 + image_idx * 6, 27 + image_idx * 6).rowRange(0, 2 * n_points));
1014             jacobians.colRange(11, 14).copyTo(Jkk.colRange(27 + image_idx * 6, 30 + image_idx * 6).rowRange(0, 2 * n_points));
1015             jacobians.colRange(0, 2).copyTo(Jkk.colRange(0, 2).rowRange(0, 2 * n_points));
1016             jacobians.colRange(2, 4).copyTo(Jkk.colRange(2, 4).rowRange(0, 2 * n_points));
1017             jacobians.colRange(4, 8).copyTo(Jkk.colRange(5, 9).rowRange(0, 2 * n_points));
1018             jacobians.col(14).copyTo(Jkk.col(4).rowRange(0, 2 * n_points));
1019 
1020             //right camera jacobian
1021             cv::internal::compose_motion(rvec, tvec, omcur, Tcur, omr, Tr, domrdomckk, domrdTckk, domrdom, domrdT, dTrdomckk, dTrdTckk, dTrdom, dTrdT);
1022             rvec = cv::Mat(rvecs2[image_idx]);
1023             tvec  = cv::Mat(tvecs2[image_idx]);
1024 
1025             cv::internal::projectPoints(object, projected, omr, Tr, intrinsicRight, jacobians);
1026             cv::Mat(cv::Mat((imageRight - projected).t()).reshape(1, 1).t()).copyTo(ekk.rowRange(2 * n_points, 4 * n_points));
1027             cv::Mat dxrdom = jacobians.colRange(8, 11) * domrdom + jacobians.colRange(11, 14) * dTrdom;
1028             cv::Mat dxrdT = jacobians.colRange(8, 11) * domrdT + jacobians.colRange(11, 14)* dTrdT;
1029             cv::Mat dxrdomckk = jacobians.colRange(8, 11) * domrdomckk + jacobians.colRange(11, 14) * dTrdomckk;
1030             cv::Mat dxrdTckk = jacobians.colRange(8, 11) * domrdTckk + jacobians.colRange(11, 14) * dTrdTckk;
1031 
1032             dxrdom.copyTo(Jkk.colRange(18, 21).rowRange(2 * n_points, 4 * n_points));
1033             dxrdT.copyTo(Jkk.colRange(21, 24).rowRange(2 * n_points, 4 * n_points));
1034             dxrdomckk.copyTo(Jkk.colRange(24 + image_idx * 6, 27 + image_idx * 6).rowRange(2 * n_points, 4 * n_points));
1035             dxrdTckk.copyTo(Jkk.colRange(27 + image_idx * 6, 30 + image_idx * 6).rowRange(2 * n_points, 4 * n_points));
1036             jacobians.colRange(0, 2).copyTo(Jkk.colRange(9 + 0, 9 + 2).rowRange(2 * n_points, 4 * n_points));
1037             jacobians.colRange(2, 4).copyTo(Jkk.colRange(9 + 2, 9 + 4).rowRange(2 * n_points, 4 * n_points));
1038             jacobians.colRange(4, 8).copyTo(Jkk.colRange(9 + 5, 9 + 9).rowRange(2 * n_points, 4 * n_points));
1039             jacobians.col(14).copyTo(Jkk.col(9 + 4).rowRange(2 * n_points, 4 * n_points));
1040 
1041             //check goodness of sterepair
1042             double abs_max  = 0;
1043             for (int i = 0; i < 4 * n_points; i++)
1044             {
1045                 if (fabs(ekk.at<double>(i)) > abs_max)
1046                 {
1047                     abs_max = fabs(ekk.at<double>(i));
1048                 }
1049             }
1050 
1051             CV_Assert(abs_max < threshold); // bad stereo pair
1052 
1053             Jkk.copyTo(J.rowRange(image_idx * 4 * n_points, (image_idx + 1) * 4 * n_points));
1054             ekk.copyTo(e.rowRange(image_idx * 4 * n_points, (image_idx + 1) * 4 * n_points));
1055         }
1056 
1057         cv::Vec6d oldTom(Tcur[0], Tcur[1], Tcur[2], omcur[0], omcur[1], omcur[2]);
1058 
1059         //update all parameters
1060         cv::subMatrix(J, J, selectedParams, std::vector<uchar>(J.rows, 1));
1061         int a = cv::countNonZero(intrinsicLeft.isEstimate);
1062         int b = cv::countNonZero(intrinsicRight.isEstimate);
1063         cv::Mat deltas;
1064         solve(J.t() * J, J.t()*e, deltas);
1065         if (a > 0)
1066             intrinsicLeft = intrinsicLeft + deltas.rowRange(0, a);
1067         if (b > 0)
1068             intrinsicRight = intrinsicRight + deltas.rowRange(a, a + b);
1069         omcur = omcur + cv::Vec3d(deltas.rowRange(a + b, a + b + 3));
1070         Tcur = Tcur + cv::Vec3d(deltas.rowRange(a + b + 3, a + b + 6));
1071         for (int image_idx = 0; image_idx < n_images; ++image_idx)
1072         {
1073             rvecs1[image_idx] = cv::Mat(cv::Mat(rvecs1[image_idx]) + deltas.rowRange(a + b + 6 + image_idx * 6, a + b + 9 + image_idx * 6));
1074             tvecs1[image_idx] = cv::Mat(cv::Mat(tvecs1[image_idx]) + deltas.rowRange(a + b + 9 + image_idx * 6, a + b + 12 + image_idx * 6));
1075         }
1076 
1077         cv::Vec6d newTom(Tcur[0], Tcur[1], Tcur[2], omcur[0], omcur[1], omcur[2]);
1078         change = cv::norm(newTom - oldTom) / cv::norm(newTom);
1079     }
1080 
1081     double rms = 0;
1082     const Vec2d* ptr_e = e.ptr<Vec2d>();
1083     for (size_t i = 0; i < e.total() / 2; i++)
1084     {
1085         rms += ptr_e[i][0] * ptr_e[i][0] + ptr_e[i][1] * ptr_e[i][1];
1086     }
1087 
1088     rms /= ((double)e.total() / 2.0);
1089     rms = sqrt(rms);
1090 
1091     _K1 = Matx33d(intrinsicLeft.f[0], intrinsicLeft.f[0] * intrinsicLeft.alpha, intrinsicLeft.c[0],
1092                                        0,                       intrinsicLeft.f[1], intrinsicLeft.c[1],
1093                                        0,                                        0,                 1);
1094 
1095     _K2 = Matx33d(intrinsicRight.f[0], intrinsicRight.f[0] * intrinsicRight.alpha, intrinsicRight.c[0],
1096                                         0,                        intrinsicRight.f[1], intrinsicRight.c[1],
1097                                         0,                                          0,                  1);
1098 
1099     Mat _R;
1100     Rodrigues(omcur, _R);
1101 
1102     if (K1.needed()) cv::Mat(_K1).convertTo(K1, K1.empty() ? CV_64FC1 : K1.type());
1103     if (K2.needed()) cv::Mat(_K2).convertTo(K2, K2.empty() ? CV_64FC1 : K2.type());
1104     if (D1.needed()) cv::Mat(intrinsicLeft.k).convertTo(D1, D1.empty() ? CV_64FC1 : D1.type());
1105     if (D2.needed()) cv::Mat(intrinsicRight.k).convertTo(D2, D2.empty() ? CV_64FC1 : D2.type());
1106     if (R.needed()) _R.convertTo(R, R.empty() ? CV_64FC1 : R.type());
1107     if (T.needed()) cv::Mat(Tcur).convertTo(T, T.empty() ? CV_64FC1 : T.type());
1108 
1109     return rms;
1110 }
1111 
1112 namespace cv{ namespace {
subMatrix(const Mat & src,Mat & dst,const std::vector<uchar> & cols,const std::vector<uchar> & rows)1113 void subMatrix(const Mat& src, Mat& dst, const std::vector<uchar>& cols, const std::vector<uchar>& rows)
1114 {
1115     CV_Assert(src.channels() == 1);
1116 
1117     int nonzeros_cols = cv::countNonZero(cols);
1118     Mat tmp(src.rows, nonzeros_cols, CV_64F);
1119 
1120     for (int i = 0, j = 0; i < (int)cols.size(); i++)
1121     {
1122         if (cols[i])
1123         {
1124             src.col(i).copyTo(tmp.col(j++));
1125         }
1126     }
1127 
1128     int nonzeros_rows  = cv::countNonZero(rows);
1129     dst.create(nonzeros_rows, nonzeros_cols, CV_64F);
1130     for (int i = 0, j = 0; i < (int)rows.size(); i++)
1131     {
1132         if (rows[i])
1133         {
1134             tmp.row(i).copyTo(dst.row(j++));
1135         }
1136     }
1137 }
1138 
1139 }}
1140 
IntrinsicParams()1141 cv::internal::IntrinsicParams::IntrinsicParams():
1142     f(Vec2d::all(0)), c(Vec2d::all(0)), k(Vec4d::all(0)), alpha(0), isEstimate(9,0)
1143 {
1144 }
1145 
IntrinsicParams(Vec2d _f,Vec2d _c,Vec4d _k,double _alpha)1146 cv::internal::IntrinsicParams::IntrinsicParams(Vec2d _f, Vec2d _c, Vec4d _k, double _alpha):
1147     f(_f), c(_c), k(_k), alpha(_alpha), isEstimate(9,0)
1148 {
1149 }
1150 
operator +(const Mat & a)1151 cv::internal::IntrinsicParams cv::internal::IntrinsicParams::operator+(const Mat& a)
1152 {
1153     CV_Assert(a.type() == CV_64FC1);
1154     IntrinsicParams tmp;
1155     const double* ptr = a.ptr<double>();
1156 
1157     int j = 0;
1158     tmp.f[0]    = this->f[0]    + (isEstimate[0] ? ptr[j++] : 0);
1159     tmp.f[1]    = this->f[1]    + (isEstimate[1] ? ptr[j++] : 0);
1160     tmp.c[0]    = this->c[0]    + (isEstimate[2] ? ptr[j++] : 0);
1161     tmp.c[1]    = this->c[1]    + (isEstimate[3] ? ptr[j++] : 0);
1162     tmp.alpha   = this->alpha   + (isEstimate[4] ? ptr[j++] : 0);
1163     tmp.k[0]    = this->k[0]    + (isEstimate[5] ? ptr[j++] : 0);
1164     tmp.k[1]    = this->k[1]    + (isEstimate[6] ? ptr[j++] : 0);
1165     tmp.k[2]    = this->k[2]    + (isEstimate[7] ? ptr[j++] : 0);
1166     tmp.k[3]    = this->k[3]    + (isEstimate[8] ? ptr[j++] : 0);
1167 
1168     tmp.isEstimate = isEstimate;
1169     return tmp;
1170 }
1171 
operator =(const Mat & a)1172 cv::internal::IntrinsicParams& cv::internal::IntrinsicParams::operator =(const Mat& a)
1173 {
1174     CV_Assert(a.type() == CV_64FC1);
1175     const double* ptr = a.ptr<double>();
1176 
1177     int j = 0;
1178 
1179     this->f[0]  = isEstimate[0] ? ptr[j++] : 0;
1180     this->f[1]  = isEstimate[1] ? ptr[j++] : 0;
1181     this->c[0]  = isEstimate[2] ? ptr[j++] : 0;
1182     this->c[1]  = isEstimate[3] ? ptr[j++] : 0;
1183     this->alpha = isEstimate[4] ? ptr[j++] : 0;
1184     this->k[0]  = isEstimate[5] ? ptr[j++] : 0;
1185     this->k[1]  = isEstimate[6] ? ptr[j++] : 0;
1186     this->k[2]  = isEstimate[7] ? ptr[j++] : 0;
1187     this->k[3]  = isEstimate[8] ? ptr[j++] : 0;
1188 
1189     return *this;
1190 }
1191 
Init(const cv::Vec2d & _f,const cv::Vec2d & _c,const cv::Vec4d & _k,const double & _alpha)1192 void cv::internal::IntrinsicParams::Init(const cv::Vec2d& _f, const cv::Vec2d& _c, const cv::Vec4d& _k, const double& _alpha)
1193 {
1194     this->c = _c;
1195     this->f = _f;
1196     this->k = _k;
1197     this->alpha = _alpha;
1198 }
1199 
projectPoints(cv::InputArray objectPoints,cv::OutputArray imagePoints,cv::InputArray _rvec,cv::InputArray _tvec,const IntrinsicParams & param,cv::OutputArray jacobian)1200 void cv::internal::projectPoints(cv::InputArray objectPoints, cv::OutputArray imagePoints,
1201                    cv::InputArray _rvec,cv::InputArray _tvec,
1202                    const IntrinsicParams& param, cv::OutputArray jacobian)
1203 {
1204     CV_INSTRUMENT_REGION();
1205 
1206     CV_Assert(!objectPoints.empty() && (objectPoints.type() == CV_32FC3 || objectPoints.type() == CV_64FC3));
1207     Matx33d K(param.f[0], param.f[0] * param.alpha, param.c[0],
1208                        0,               param.f[1], param.c[1],
1209                        0,                        0,         1);
1210     fisheye::projectPoints(objectPoints, imagePoints, _rvec, _tvec, K, param.k, param.alpha, jacobian);
1211 }
1212 
ComputeExtrinsicRefine(const Mat & imagePoints,const Mat & objectPoints,Mat & rvec,Mat & tvec,Mat & J,const int MaxIter,const IntrinsicParams & param,const double thresh_cond)1213 void cv::internal::ComputeExtrinsicRefine(const Mat& imagePoints, const Mat& objectPoints, Mat& rvec,
1214                             Mat&  tvec, Mat& J, const int MaxIter,
1215                             const IntrinsicParams& param, const double thresh_cond)
1216 {
1217     CV_Assert(!objectPoints.empty() && objectPoints.type() == CV_64FC3);
1218     CV_Assert(!imagePoints.empty() && imagePoints.type() == CV_64FC2);
1219     CV_Assert(rvec.total() > 2 && tvec.total() > 2);
1220     Vec6d extrinsics(rvec.at<double>(0), rvec.at<double>(1), rvec.at<double>(2),
1221                     tvec.at<double>(0), tvec.at<double>(1), tvec.at<double>(2));
1222     double change = 1;
1223     int iter = 0;
1224 
1225     while (change > 1e-10 && iter < MaxIter)
1226     {
1227         std::vector<Point2d> x;
1228         Mat jacobians;
1229         projectPoints(objectPoints, x, rvec, tvec, param, jacobians);
1230 
1231         Mat ex = imagePoints - Mat(x).t();
1232         ex = ex.reshape(1, 2);
1233 
1234         J = jacobians.colRange(8, 14).clone();
1235 
1236         SVD svd(J, SVD::NO_UV);
1237         double condJJ = svd.w.at<double>(0)/svd.w.at<double>(5);
1238 
1239         if (condJJ > thresh_cond)
1240             change = 0;
1241         else
1242         {
1243             Vec6d param_innov;
1244             solve(J, ex.reshape(1, (int)ex.total()), param_innov, DECOMP_SVD + DECOMP_NORMAL);
1245 
1246             Vec6d param_up = extrinsics + param_innov;
1247             change = norm(param_innov)/norm(param_up);
1248             extrinsics = param_up;
1249             iter = iter + 1;
1250 
1251             rvec = Mat(Vec3d(extrinsics.val));
1252             tvec = Mat(Vec3d(extrinsics.val+3));
1253         }
1254     }
1255 }
1256 
ComputeHomography(Mat m,Mat M)1257 cv::Mat cv::internal::ComputeHomography(Mat m, Mat M)
1258 {
1259     CV_INSTRUMENT_REGION();
1260 
1261     int Np = m.cols;
1262 
1263     if (m.rows < 3)
1264     {
1265         vconcat(m, Mat::ones(1, Np, CV_64FC1), m);
1266     }
1267     if (M.rows < 3)
1268     {
1269         vconcat(M, Mat::ones(1, Np, CV_64FC1), M);
1270     }
1271 
1272     divide(m, Mat::ones(3, 1, CV_64FC1) * m.row(2), m);
1273     divide(M, Mat::ones(3, 1, CV_64FC1) * M.row(2), M);
1274 
1275     Mat ax = m.row(0).clone();
1276     Mat ay = m.row(1).clone();
1277 
1278     double mxx = mean(ax)[0];
1279     double myy = mean(ay)[0];
1280 
1281     ax = ax - mxx;
1282     ay = ay - myy;
1283 
1284     double scxx = mean(abs(ax))[0];
1285     double scyy = mean(abs(ay))[0];
1286 
1287     Mat Hnorm (Matx33d( 1/scxx,        0.0,     -mxx/scxx,
1288                          0.0,     1/scyy,     -myy/scyy,
1289                          0.0,        0.0,           1.0 ));
1290 
1291     Mat inv_Hnorm (Matx33d( scxx,     0,   mxx,
1292                                     0,  scyy,   myy,
1293                                     0,     0,     1 ));
1294     Mat mn =  Hnorm * m;
1295 
1296     Mat L = Mat::zeros(2*Np, 9, CV_64FC1);
1297 
1298     for (int i = 0; i < Np; ++i)
1299     {
1300         for (int j = 0; j < 3; j++)
1301         {
1302             L.at<double>(2 * i, j) = M.at<double>(j, i);
1303             L.at<double>(2 * i + 1, j + 3) = M.at<double>(j, i);
1304             L.at<double>(2 * i, j + 6) = -mn.at<double>(0,i) * M.at<double>(j, i);
1305             L.at<double>(2 * i + 1, j + 6) = -mn.at<double>(1,i) * M.at<double>(j, i);
1306         }
1307     }
1308 
1309     if (Np > 4) L = L.t() * L;
1310     SVD svd(L);
1311     Mat hh = svd.vt.row(8) / svd.vt.row(8).at<double>(8);
1312     Mat Hrem = hh.reshape(1, 3);
1313     Mat H = inv_Hnorm * Hrem;
1314 
1315     if (Np > 4)
1316     {
1317         Mat hhv = H.reshape(1, 9)(Rect(0, 0, 1, 8)).clone();
1318         for (int iter = 0; iter < 10; iter++)
1319         {
1320             Mat mrep = H * M;
1321             Mat J = Mat::zeros(2 * Np, 8, CV_64FC1);
1322             Mat MMM;
1323             divide(M, Mat::ones(3, 1, CV_64FC1) * mrep(Rect(0, 2, mrep.cols, 1)), MMM);
1324             divide(mrep, Mat::ones(3, 1, CV_64FC1) * mrep(Rect(0, 2, mrep.cols, 1)), mrep);
1325             Mat m_err = m(Rect(0,0, m.cols, 2)) - mrep(Rect(0,0, mrep.cols, 2));
1326             m_err = Mat(m_err.t()).reshape(1, m_err.cols * m_err.rows);
1327             Mat MMM2, MMM3;
1328             multiply(Mat::ones(3, 1, CV_64FC1) * mrep(Rect(0, 0, mrep.cols, 1)), MMM, MMM2);
1329             multiply(Mat::ones(3, 1, CV_64FC1) * mrep(Rect(0, 1, mrep.cols, 1)), MMM, MMM3);
1330 
1331             for (int i = 0; i < Np; ++i)
1332             {
1333                 for (int j = 0; j < 3; ++j)
1334                 {
1335                     J.at<double>(2 * i, j)         = -MMM.at<double>(j, i);
1336                     J.at<double>(2 * i + 1, j + 3) = -MMM.at<double>(j, i);
1337                 }
1338 
1339                 for (int j = 0; j < 2; ++j)
1340                 {
1341                     J.at<double>(2 * i, j + 6)     = MMM2.at<double>(j, i);
1342                     J.at<double>(2 * i + 1, j + 6) = MMM3.at<double>(j, i);
1343                 }
1344             }
1345             divide(M, Mat::ones(3, 1, CV_64FC1) * mrep(Rect(0,2,mrep.cols,1)), MMM);
1346             Mat hh_innov = (J.t() * J).inv() * (J.t()) * m_err;
1347             Mat hhv_up = hhv - hh_innov;
1348             Mat tmp;
1349             vconcat(hhv_up, Mat::ones(1,1,CV_64FC1), tmp);
1350             Mat H_up = tmp.reshape(1,3);
1351             hhv = hhv_up;
1352             H = H_up;
1353         }
1354     }
1355     return H;
1356 }
1357 
NormalizePixels(const Mat & imagePoints,const IntrinsicParams & param)1358 cv::Mat cv::internal::NormalizePixels(const Mat& imagePoints, const IntrinsicParams& param)
1359 {
1360     CV_INSTRUMENT_REGION();
1361 
1362     CV_Assert(!imagePoints.empty() && imagePoints.type() == CV_64FC2);
1363 
1364     Mat distorted((int)imagePoints.total(), 1, CV_64FC2), undistorted;
1365     const Vec2d* ptr   = imagePoints.ptr<Vec2d>();
1366     Vec2d* ptr_d = distorted.ptr<Vec2d>();
1367     for (size_t i = 0; i < imagePoints.total(); ++i)
1368     {
1369         ptr_d[i] = (ptr[i] - param.c).mul(Vec2d(1.0 / param.f[0], 1.0 / param.f[1]));
1370         ptr_d[i][0] -= param.alpha * ptr_d[i][1];
1371     }
1372     cv::fisheye::undistortPoints(distorted, undistorted, Matx33d::eye(), param.k);
1373     return undistorted;
1374 }
1375 
InitExtrinsics(const Mat & _imagePoints,const Mat & _objectPoints,const IntrinsicParams & param,Mat & omckk,Mat & Tckk)1376 void cv::internal::InitExtrinsics(const Mat& _imagePoints, const Mat& _objectPoints, const IntrinsicParams& param, Mat& omckk, Mat& Tckk)
1377 {
1378     CV_Assert(!_objectPoints.empty() && _objectPoints.type() == CV_64FC3);
1379     CV_Assert(!_imagePoints.empty() && _imagePoints.type() == CV_64FC2);
1380 
1381     Mat imagePointsNormalized = NormalizePixels(_imagePoints, param).reshape(1).t();
1382     Mat objectPoints = _objectPoints.reshape(1).t();
1383     Mat objectPointsMean, covObjectPoints;
1384     Mat Rckk;
1385     int Np = imagePointsNormalized.cols;
1386     calcCovarMatrix(objectPoints, covObjectPoints, objectPointsMean, COVAR_NORMAL | COVAR_COLS);
1387     SVD svd(covObjectPoints);
1388     Mat R(svd.vt);
1389     if (norm(R(Rect(2, 0, 1, 2))) < 1e-6)
1390         R = Mat::eye(3,3, CV_64FC1);
1391     if (determinant(R) < 0)
1392         R = -R;
1393     Mat T = -R * objectPointsMean;
1394     Mat X_new = R * objectPoints + T * Mat::ones(1, Np, CV_64FC1);
1395     Mat H = ComputeHomography(imagePointsNormalized, X_new(Rect(0,0,X_new.cols,2)));
1396     double sc = .5 * (norm(H.col(0)) + norm(H.col(1)));
1397     H = H / sc;
1398     Mat u1 = H.col(0).clone();
1399     double norm_u1 = norm(u1);
1400     CV_Assert(fabs(norm_u1) > 0);
1401     u1  = u1 / norm_u1;
1402     Mat u2 = H.col(1).clone() - u1.dot(H.col(1).clone()) * u1;
1403     double norm_u2 = norm(u2);
1404     CV_Assert(fabs(norm_u2) > 0);
1405     u2 = u2 / norm_u2;
1406     Mat u3 = u1.cross(u2);
1407     Mat RRR;
1408     hconcat(u1, u2, RRR);
1409     hconcat(RRR, u3, RRR);
1410     Rodrigues(RRR, omckk);
1411     Rodrigues(omckk, Rckk);
1412     Tckk = H.col(2).clone();
1413     Tckk = Tckk + Rckk * T;
1414     Rckk = Rckk * R;
1415     Rodrigues(Rckk, omckk);
1416 }
1417 
CalibrateExtrinsics(InputArrayOfArrays objectPoints,InputArrayOfArrays imagePoints,const IntrinsicParams & param,const int check_cond,const double thresh_cond,InputOutputArray omc,InputOutputArray Tc)1418 void cv::internal::CalibrateExtrinsics(InputArrayOfArrays objectPoints, InputArrayOfArrays imagePoints,
1419                          const IntrinsicParams& param, const int check_cond,
1420                          const double thresh_cond, InputOutputArray omc, InputOutputArray Tc)
1421 {
1422     CV_Assert(!objectPoints.empty() && (objectPoints.type() == CV_32FC3 || objectPoints.type() == CV_64FC3));
1423     CV_Assert(!imagePoints.empty() && (imagePoints.type() == CV_32FC2 || imagePoints.type() == CV_64FC2));
1424     CV_Assert(omc.type() == CV_64FC3 || Tc.type() == CV_64FC3);
1425 
1426     if (omc.empty()) omc.create(1, (int)objectPoints.total(), CV_64FC3);
1427     if (Tc.empty()) Tc.create(1, (int)objectPoints.total(), CV_64FC3);
1428 
1429     const int maxIter = 20;
1430 
1431     for(int image_idx = 0; image_idx < (int)imagePoints.total(); ++image_idx)
1432     {
1433         Mat omckk, Tckk, JJ_kk;
1434         Mat image, object;
1435 
1436         objectPoints.getMat(image_idx).convertTo(object,  CV_64FC3);
1437         imagePoints.getMat (image_idx).convertTo(image, CV_64FC2);
1438 
1439         bool imT = image.rows < image.cols;
1440         bool obT = object.rows < object.cols;
1441 
1442         InitExtrinsics(imT ? image.t() : image, obT ? object.t() : object, param, omckk, Tckk);
1443 
1444         ComputeExtrinsicRefine(!imT ? image.t() : image, !obT ? object.t() : object, omckk, Tckk, JJ_kk, maxIter, param, thresh_cond);
1445         if (check_cond)
1446         {
1447             SVD svd(JJ_kk, SVD::NO_UV);
1448             if(svd.w.at<double>(0) / svd.w.at<double>((int)svd.w.total() - 1) > thresh_cond )
1449                 CV_Error( cv::Error::StsInternal, format("CALIB_CHECK_COND - Ill-conditioned matrix for input array %d",image_idx));
1450         }
1451         omckk.reshape(3,1).copyTo(omc.getMat().col(image_idx));
1452         Tckk.reshape(3,1).copyTo(Tc.getMat().col(image_idx));
1453     }
1454 }
1455 
ComputeJacobians(InputArrayOfArrays objectPoints,InputArrayOfArrays imagePoints,const IntrinsicParams & param,InputArray omc,InputArray Tc,const int & check_cond,const double & thresh_cond,Mat & JJ2,Mat & ex3)1456 void cv::internal::ComputeJacobians(InputArrayOfArrays objectPoints, InputArrayOfArrays imagePoints,
1457                       const IntrinsicParams& param,  InputArray omc, InputArray Tc,
1458                       const int& check_cond, const double& thresh_cond, Mat& JJ2, Mat& ex3)
1459 {
1460     CV_Assert(!objectPoints.empty() && (objectPoints.type() == CV_32FC3 || objectPoints.type() == CV_64FC3));
1461     CV_Assert(!imagePoints.empty() && (imagePoints.type() == CV_32FC2 || imagePoints.type() == CV_64FC2));
1462 
1463     CV_Assert(!omc.empty() && omc.type() == CV_64FC3);
1464     CV_Assert(!Tc.empty() && Tc.type() == CV_64FC3);
1465 
1466     int n = (int)objectPoints.total();
1467 
1468     JJ2 = Mat::zeros(9 + 6 * n, 9 + 6 * n, CV_64FC1);
1469     ex3 = Mat::zeros(9 + 6 * n, 1, CV_64FC1 );
1470 
1471     for (int image_idx = 0; image_idx < n; ++image_idx)
1472     {
1473         Mat image, object;
1474         objectPoints.getMat(image_idx).convertTo(object, CV_64FC3);
1475         imagePoints.getMat (image_idx).convertTo(image, CV_64FC2);
1476 
1477         bool imT = image.rows < image.cols;
1478         Mat om(omc.getMat().col(image_idx)), T(Tc.getMat().col(image_idx));
1479 
1480         std::vector<Point2d> x;
1481         Mat jacobians;
1482         projectPoints(object, x, om, T, param, jacobians);
1483         Mat exkk = (imT ? image.t() : image) - Mat(x);
1484 
1485         Mat A(jacobians.rows, 9, CV_64FC1);
1486         jacobians.colRange(0, 4).copyTo(A.colRange(0, 4));
1487         jacobians.col(14).copyTo(A.col(4));
1488         jacobians.colRange(4, 8).copyTo(A.colRange(5, 9));
1489 
1490         A = A.t();
1491 
1492         Mat B = jacobians.colRange(8, 14).clone();
1493         B = B.t();
1494 
1495         JJ2(Rect(0, 0, 9, 9)) += A * A.t();
1496         JJ2(Rect(9 + 6 * image_idx, 9 + 6 * image_idx, 6, 6)) = B * B.t();
1497 
1498         JJ2(Rect(9 + 6 * image_idx, 0, 6, 9)) = A * B.t();
1499         JJ2(Rect(0, 9 + 6 * image_idx, 9, 6)) = JJ2(Rect(9 + 6 * image_idx, 0, 6, 9)).t();
1500 
1501         ex3.rowRange(0, 9) += A * exkk.reshape(1, 2 * exkk.rows);
1502         ex3.rowRange(9 + 6 * image_idx, 9 + 6 * (image_idx + 1)) = B * exkk.reshape(1, 2 * exkk.rows);
1503 
1504         if (check_cond)
1505         {
1506             Mat JJ_kk = B.t();
1507             SVD svd(JJ_kk, SVD::NO_UV);
1508             CV_Assert(svd.w.at<double>(0) / svd.w.at<double>(svd.w.rows - 1) < thresh_cond);
1509         }
1510     }
1511 
1512     std::vector<uchar> idxs(param.isEstimate);
1513     idxs.insert(idxs.end(), 6 * n, 1);
1514 
1515     subMatrix(JJ2, JJ2, idxs, idxs);
1516     subMatrix(ex3, ex3, std::vector<uchar>(1, 1), idxs);
1517 }
1518 
EstimateUncertainties(InputArrayOfArrays objectPoints,InputArrayOfArrays imagePoints,const IntrinsicParams & params,InputArray omc,InputArray Tc,IntrinsicParams & errors,Vec2d & std_err,double thresh_cond,int check_cond,double & rms)1519 void cv::internal::EstimateUncertainties(InputArrayOfArrays objectPoints, InputArrayOfArrays imagePoints,
1520                            const IntrinsicParams& params, InputArray omc, InputArray Tc,
1521                            IntrinsicParams& errors, Vec2d& std_err, double thresh_cond, int check_cond, double& rms)
1522 {
1523     CV_Assert(!objectPoints.empty() && (objectPoints.type() == CV_32FC3 || objectPoints.type() == CV_64FC3));
1524     CV_Assert(!imagePoints.empty() && (imagePoints.type() == CV_32FC2 || imagePoints.type() == CV_64FC2));
1525 
1526     CV_Assert(!omc.empty() && omc.type() == CV_64FC3);
1527     CV_Assert(!Tc.empty() && Tc.type() == CV_64FC3);
1528 
1529     int total_ex = 0;
1530     for (int image_idx = 0; image_idx < (int)objectPoints.total(); ++image_idx)
1531     {
1532         total_ex += (int)objectPoints.getMat(image_idx).total();
1533     }
1534     Mat ex(total_ex, 1, CV_64FC2);
1535     int insert_idx = 0;
1536     for (int image_idx = 0; image_idx < (int)objectPoints.total(); ++image_idx)
1537     {
1538         Mat image, object;
1539         objectPoints.getMat(image_idx).convertTo(object, CV_64FC3);
1540         imagePoints.getMat (image_idx).convertTo(image, CV_64FC2);
1541 
1542         bool imT = image.rows < image.cols;
1543 
1544         Mat om(omc.getMat().col(image_idx)), T(Tc.getMat().col(image_idx));
1545 
1546         std::vector<Point2d> x;
1547         projectPoints(object, x, om, T, params, noArray());
1548         Mat ex_ = (imT ? image.t() : image) - Mat(x);
1549         ex_.copyTo(ex.rowRange(insert_idx, insert_idx + ex_.rows));
1550         insert_idx += ex_.rows;
1551     }
1552 
1553     meanStdDev(ex, noArray(), std_err);
1554     std_err *= sqrt((double)ex.total()/((double)ex.total() - 1.0));
1555 
1556     Vec<double, 1> sigma_x;
1557     meanStdDev(ex.reshape(1, 1), noArray(), sigma_x);
1558     sigma_x  *= sqrt(2.0 * (double)ex.total()/(2.0 * (double)ex.total() - 1.0));
1559 
1560     Mat JJ2, ex3;
1561     ComputeJacobians(objectPoints, imagePoints, params, omc, Tc, check_cond, thresh_cond, JJ2, ex3);
1562 
1563     sqrt(JJ2.inv(), JJ2);
1564 
1565     errors = 3 * sigma_x(0) * JJ2.diag();
1566     rms = sqrt(norm(ex, NORM_L2SQR)/ex.total());
1567 }
1568 
dAB(InputArray A,InputArray B,OutputArray dABdA,OutputArray dABdB)1569 void cv::internal::dAB(InputArray A, InputArray B, OutputArray dABdA, OutputArray dABdB)
1570 {
1571     CV_Assert(A.getMat().cols == B.getMat().rows);
1572     CV_Assert(A.type() == CV_64FC1 && B.type() == CV_64FC1);
1573 
1574     int p = A.getMat().rows;
1575     int n = A.getMat().cols;
1576     int q = B.getMat().cols;
1577 
1578     dABdA.create(p * q, p * n, CV_64FC1);
1579     dABdB.create(p * q, q * n, CV_64FC1);
1580 
1581     dABdA.getMat() = Mat::zeros(p * q, p * n, CV_64FC1);
1582     dABdB.getMat() = Mat::zeros(p * q, q * n, CV_64FC1);
1583 
1584     for (int i = 0; i < q; ++i)
1585     {
1586         for (int j = 0; j < p; ++j)
1587         {
1588             int ij = j + i * p;
1589             for (int k = 0; k < n; ++k)
1590             {
1591                 int kj = j + k * p;
1592                 dABdA.getMat().at<double>(ij, kj) = B.getMat().at<double>(k, i);
1593             }
1594         }
1595     }
1596 
1597     for (int i = 0; i < q; ++i)
1598     {
1599         A.getMat().copyTo(dABdB.getMat().rowRange(i * p, i * p + p).colRange(i * n, i * n + n));
1600     }
1601 }
1602 
JRodriguesMatlab(const Mat & src,Mat & dst)1603 void cv::internal::JRodriguesMatlab(const Mat& src, Mat& dst)
1604 {
1605     Mat tmp(src.cols, src.rows, src.type());
1606     if (src.rows == 9)
1607     {
1608         Mat(src.row(0).t()).copyTo(tmp.col(0));
1609         Mat(src.row(1).t()).copyTo(tmp.col(3));
1610         Mat(src.row(2).t()).copyTo(tmp.col(6));
1611         Mat(src.row(3).t()).copyTo(tmp.col(1));
1612         Mat(src.row(4).t()).copyTo(tmp.col(4));
1613         Mat(src.row(5).t()).copyTo(tmp.col(7));
1614         Mat(src.row(6).t()).copyTo(tmp.col(2));
1615         Mat(src.row(7).t()).copyTo(tmp.col(5));
1616         Mat(src.row(8).t()).copyTo(tmp.col(8));
1617     }
1618     else
1619     {
1620         Mat(src.col(0).t()).copyTo(tmp.row(0));
1621         Mat(src.col(1).t()).copyTo(tmp.row(3));
1622         Mat(src.col(2).t()).copyTo(tmp.row(6));
1623         Mat(src.col(3).t()).copyTo(tmp.row(1));
1624         Mat(src.col(4).t()).copyTo(tmp.row(4));
1625         Mat(src.col(5).t()).copyTo(tmp.row(7));
1626         Mat(src.col(6).t()).copyTo(tmp.row(2));
1627         Mat(src.col(7).t()).copyTo(tmp.row(5));
1628         Mat(src.col(8).t()).copyTo(tmp.row(8));
1629     }
1630     dst = tmp.clone();
1631 }
1632 
compose_motion(InputArray _om1,InputArray _T1,InputArray _om2,InputArray _T2,Mat & om3,Mat & T3,Mat & dom3dom1,Mat & dom3dT1,Mat & dom3dom2,Mat & dom3dT2,Mat & dT3dom1,Mat & dT3dT1,Mat & dT3dom2,Mat & dT3dT2)1633 void cv::internal::compose_motion(InputArray _om1, InputArray _T1, InputArray _om2, InputArray _T2,
1634                     Mat& om3, Mat& T3, Mat& dom3dom1, Mat& dom3dT1, Mat& dom3dom2,
1635                     Mat& dom3dT2, Mat& dT3dom1, Mat& dT3dT1, Mat& dT3dom2, Mat& dT3dT2)
1636 {
1637     Mat om1 = _om1.getMat();
1638     Mat om2 = _om2.getMat();
1639     Mat T1 = _T1.getMat().reshape(1, 3);
1640     Mat T2 = _T2.getMat().reshape(1, 3);
1641 
1642     //% Rotations:
1643     Mat R1, R2, R3, dR1dom1(9, 3, CV_64FC1), dR2dom2;
1644     Rodrigues(om1, R1, dR1dom1);
1645     Rodrigues(om2, R2, dR2dom2);
1646     JRodriguesMatlab(dR1dom1, dR1dom1);
1647     JRodriguesMatlab(dR2dom2, dR2dom2);
1648     R3 = R2 * R1;
1649     Mat dR3dR2, dR3dR1;
1650     dAB(R2, R1, dR3dR2, dR3dR1);
1651     Mat dom3dR3;
1652     Rodrigues(R3, om3, dom3dR3);
1653     JRodriguesMatlab(dom3dR3, dom3dR3);
1654     dom3dom1 = dom3dR3 * dR3dR1 * dR1dom1;
1655     dom3dom2 = dom3dR3 * dR3dR2 * dR2dom2;
1656     dom3dT1 = Mat::zeros(3, 3, CV_64FC1);
1657     dom3dT2 = Mat::zeros(3, 3, CV_64FC1);
1658 
1659     //% Translations:
1660     Mat T3t = R2 * T1;
1661     Mat dT3tdR2, dT3tdT1;
1662     dAB(R2, T1, dT3tdR2, dT3tdT1);
1663     Mat dT3tdom2 = dT3tdR2 * dR2dom2;
1664     T3 = T3t + T2;
1665     dT3dT1 = dT3tdT1;
1666     dT3dT2 = Mat::eye(3, 3, CV_64FC1);
1667     dT3dom2 = dT3tdom2;
1668     dT3dom1 = Mat::zeros(3, 3, CV_64FC1);
1669 }
1670 
median(const Mat & row)1671 double cv::internal::median(const Mat& row)
1672 {
1673     CV_Assert(row.type() == CV_64FC1);
1674     CV_Assert(!row.empty() && row.rows == 1);
1675     Mat tmp = row.clone();
1676     sort(tmp, tmp, 0);
1677     if ((int)tmp.total() % 2) return tmp.at<double>((int)tmp.total() / 2);
1678     else return 0.5 *(tmp.at<double>((int)tmp.total() / 2) + tmp.at<double>((int)tmp.total() / 2 - 1));
1679 }
1680 
median3d(InputArray m)1681 cv::Vec3d cv::internal::median3d(InputArray m)
1682 {
1683     CV_Assert(m.depth() == CV_64F && m.getMat().rows == 1);
1684     Mat M = Mat(m.getMat().t()).reshape(1).t();
1685     return Vec3d(median(M.row(0)), median(M.row(1)), median(M.row(2)));
1686 }
1687