001/******************************************************************************* 002 * This software is provided as a supplement to the authors' textbooks on digital 003 * image processing published by Springer-Verlag in various languages and editions. 004 * Permission to use and distribute this software is granted under the BSD 2-Clause 005 * "Simplified" License (see http://opensource.org/licenses/BSD-2-Clause). 006 * Copyright (c) 2006-2023 Wilhelm Burger, Mark J. Burge. All rights reserved. 007 * Visit https://imagingbook.com for additional details. 008 ******************************************************************************/ 009package imagingbook.common.geometry.fitting.ellipse.algebraic; 010 011import imagingbook.common.geometry.basic.Pnt2d; 012import imagingbook.common.math.Matrix; 013import imagingbook.common.math.eigen.GeneralizedSymmetricEigenDecomposition; 014import org.apache.commons.math3.linear.MatrixUtils; 015import org.apache.commons.math3.linear.RealMatrix; 016 017import static imagingbook.common.math.Arithmetic.sqr; 018 019/** 020 * <p> 021 * Algebraic ellipse fit based on Fitzgibbon's original method [1], as described in [2] (WITHOUT the additional 022 * numerical improvements). Based on generalized symmetric eigenproblem solution. See [3, Sec. 11.2.1] for a detailed 023 * description. Note: This implementation does not use data centering nor accepts a specific reference point. With 024 * exactly 5 input points (generally sufficient for ellipse fitting) the scatter matrix X is singular and thus the 025 * Cholesky decomposition used by the {@link GeneralizedSymmetricEigenDecomposition} cannot be applied. Thus at least 6 026 * distinct input points are required (i.e., no duplicate points are allowed). 027 * </p> 028 * <p> 029 * [1] A. W. Fitzgibbon, M. Pilu, and R. B. Fisher. Direct least- squares fitting of ellipses. IEEE Transactions on 030 * Pattern Analysis and Machine Intelligence 21(5), 476-480 (1999). <br> [2] R. Halíř and J. Flusser. Numerically stable 031 * direct least squares fitting of ellipses. In "Proceedings of the 6th International Conference in Central Europe on 032 * Computer Graphics and Visualization (WSCG’98)", pp. 125-132, Plzeň, CZ (February 1998). <br> [3] W. Burger, M.J. 033 * Burge, <em>Digital Image Processing – An Algorithmic Introduction</em>, 3rd ed, Springer (2022). 034 * </p> 035 * 036 * @author WB 037 * @version 2022/11/17 038 */ 039public class EllipseFitFitzgibbonOriginal implements EllipseFitAlgebraic { 040 041 // constraint matrix 042 // constraint matrix 043 private static final RealMatrix C = Matrix.makeRealMatrix(6, 6, 044 0, 0, 2, 0, 0, 0, 045 0,-1, 0, 0, 0, 0, 046 2, 0, 0, 0, 0, 0, 047 0, 0, 0, 0, 0, 0, 048 0, 0, 0, 0, 0, 0, 049 0, 0, 0, 0, 0, 0 ); 050 051 private final double[] q; // = (A,B,C,D,E,F) algebraic ellipse parameters 052 053 public EllipseFitFitzgibbonOriginal(Pnt2d[] points) { 054 if (points.length < 6) { 055 throw new IllegalArgumentException("fitter requires at least 6 sample points instead of " + points.length); 056 } 057 this.q = fit(points); 058 } 059 060 @Override 061 public double[] getParameters() { 062 return this.q; 063 } 064 065 private double[] fit(Pnt2d[] points) { 066 final int n = points.length; 067 068 // design matrix X: 069 RealMatrix X = MatrixUtils.createRealMatrix(n, 6); 070 071 for (int i = 0; i < n; i++) { 072 double x = points[i].getX(); 073 double y = points[i].getY(); 074 double[] ri = {sqr(x), x * y, sqr(y), x, y, 1}; 075 X.setRow(i, ri); 076 } 077 078 // scatter matrix S: 079 RealMatrix S = X.transpose().multiply(X); 080 081 // solve C*p = lambda*S*p, 082 // equiv. to A*x = lambda*B*x (with A, B symmetric, B positive definite): 083 GeneralizedSymmetricEigenDecomposition eigen = 084 new GeneralizedSymmetricEigenDecomposition(C, S, 1e-15, 1e-15); 085 086 if (eigen.hasComplexEigenvalues()) { 087 throw new RuntimeException("encountered complex eigenvalues"); 088 } 089 090 double[] evals = eigen.getRealEigenvalues(); 091 int k = Matrix.idxMax(evals); // index of the largest eigenvalue 092 return eigen.getEigenvector(k).toArray(); // associated eigenvector 093 } 094 095}