ViewGraphComparison.py
Go to the documentation of this file.
1 """
2  Compare several methods for optimizing the view-graph.
3  We measure the distance from the ground truth in terms of the norm of
4  local coordinates (geodesic distance) on the F-manifold.
5  We also plot the final error of the optimization.
6 
7  Author: Frank Dellaert (with heavy assist from ChatGPT)
8  Date: October 2024
9 """
10 
11 import argparse
12 
13 import matplotlib.pyplot as plt
14 import numpy as np
15 from gtsam.examples import SFMdata
16 
17 import gtsam
18 from gtsam import (
19  Cal3f,
20  EdgeKey,
21  EssentialMatrix,
22  FundamentalMatrix,
23  LevenbergMarquardtOptimizer,
24  LevenbergMarquardtParams,
25  NonlinearFactorGraph,
26  PinholeCameraCal3f,
27  SimpleFundamentalMatrix,
28  Values,
29 )
30 
31 # For symbol shorthand (e.g., K(0), K(1))
32 K = gtsam.symbol_shorthand.K
33 
34 # Methods to compare
35 methods = ["SimpleF", "Fundamental", "Essential+Ks", "Calibrated"]
36 
37 
38 # Formatter function for printing keys
39 def formatter(key):
40  sym = gtsam.Symbol(key)
41  if sym.chr() == ord("k"):
42  return f"k{sym.index()}"
43  else:
44  edge = EdgeKey(key)
45  return f"({edge.i()},{edge.j()})"
46 
47 
48 def simulate_geometry(num_cameras, rng, num_random_points=12):
49  """simulate geometry (points and poses)"""
50  # Define the camera calibration parameters
51  cal = Cal3f(50.0, 50.0, 50.0)
52 
53  # Create the set of 8 ground-truth landmarks
54  points = SFMdata.createPoints()
55 
56  # Create extra random points in the -10,10 cube around the origin
57  extra_points = rng.uniform(-10, 10, (num_random_points, 3))
58  points.extend([gtsam.Point3(p) for p in extra_points])
59 
60  # Create the set of ground-truth poses
61  poses = SFMdata.posesOnCircle(num_cameras, 30)
62 
63  return points, poses, cal
64 
65 
66 def simulate_data(points, poses, cal, rng, noise_std):
67  """Simulate measurements from each camera pose"""
68  measurements = [[None for _ in points] for _ in poses]
69  for i, pose in enumerate(poses):
70  camera = PinholeCameraCal3f(pose, cal)
71  for j, point in enumerate(points):
72  projection = camera.project(point)
73  noise = rng.normal(0, noise_std, size=2)
74  measurements[i][j] = projection + noise
75 
76  return measurements
77 
78 
79 # Function to compute ground truth matrices
80 def compute_ground_truth(method, poses, cal):
81  E1 = EssentialMatrix.FromPose3(poses[0].between(poses[1]))
82  E2 = EssentialMatrix.FromPose3(poses[0].between(poses[2]))
83  F1 = FundamentalMatrix(cal.K(), E1, cal.K())
84  F2 = FundamentalMatrix(cal.K(), E2, cal.K())
85  if method == "Fundamental":
86  return F1, F2
87  elif method == "SimpleF":
88  f = cal.fx()
89  c = cal.principalPoint()
90  SF1 = SimpleFundamentalMatrix(E1, f, f, c, c)
91  SF2 = SimpleFundamentalMatrix(E2, f, f, c, c)
92  return SF1, SF2
93  elif method == "Essential+Ks" or method == "Calibrated":
94  return E1, E2
95  else:
96  raise ValueError(f"Unknown method {method}")
97 
98 
99 def build_factor_graph(method, num_cameras, measurements, cal):
100  """build the factor graph"""
101  graph = NonlinearFactorGraph()
102 
103  if method == "Fundamental":
104  FactorClass = gtsam.TransferFactorFundamentalMatrix
105  elif method == "SimpleF":
106  FactorClass = gtsam.TransferFactorSimpleFundamentalMatrix
107  elif method == "Essential+Ks":
108  FactorClass = gtsam.EssentialTransferFactorKCal3f
109  # add priors on all calibrations:
110  for i in range(num_cameras):
111  model = gtsam.noiseModel.Isotropic.Sigma(1, 10.0)
112  graph.addPriorCal3f(K(i), cal, model)
113  elif method == "Calibrated":
114  FactorClass = gtsam.EssentialTransferFactorCal3f
115  # No priors on calibration needed
116  else:
117  raise ValueError(f"Unknown method {method}")
118 
119  z = measurements # shorthand
120 
121  for a in range(num_cameras):
122  b = (a + 1) % num_cameras # Next camera
123  c = (a + 2) % num_cameras # Camera after next
124 
125  # Vectors to collect tuples for each factor
126  tuples1 = []
127  tuples2 = []
128  tuples3 = []
129 
130  # Collect data for the three factors
131  for j in range(len(measurements[0])):
132  tuples1.append((z[a][j], z[b][j], z[c][j]))
133  tuples2.append((z[a][j], z[c][j], z[b][j]))
134  tuples3.append((z[c][j], z[b][j], z[a][j]))
135 
136  # Add transfer factors between views a, b, and c.
137  if method in ["Calibrated"]:
138  graph.add(FactorClass(EdgeKey(a, c), EdgeKey(b, c), tuples1, cal))
139  graph.add(FactorClass(EdgeKey(a, b), EdgeKey(b, c), tuples2, cal))
140  graph.add(FactorClass(EdgeKey(a, c), EdgeKey(a, b), tuples3, cal))
141  else:
142  graph.add(FactorClass(EdgeKey(a, c), EdgeKey(b, c), tuples1))
143  graph.add(FactorClass(EdgeKey(a, b), EdgeKey(b, c), tuples2))
144  graph.add(FactorClass(EdgeKey(a, c), EdgeKey(a, b), tuples3))
145 
146  return graph
147 
148 
149 def get_initial_estimate(method, num_cameras, ground_truth, cal):
150  """get initial estimate for method"""
151  initialEstimate = Values()
152  total_dimension = 0
153 
154  if method in ["Fundamental", "SimpleF"]:
155  F1, F2 = ground_truth
156  for a in range(num_cameras):
157  b = (a + 1) % num_cameras # Next camera
158  c = (a + 2) % num_cameras # Camera after next
159  initialEstimate.insert(EdgeKey(a, b).key(), F1)
160  initialEstimate.insert(EdgeKey(a, c).key(), F2)
161  total_dimension += F1.dim() + F2.dim()
162  elif method in ["Essential+Ks", "Calibrated"]:
163  E1, E2 = ground_truth
164  for a in range(num_cameras):
165  b = (a + 1) % num_cameras # Next camera
166  c = (a + 2) % num_cameras # Camera after next
167  initialEstimate.insert(EdgeKey(a, b).key(), E1)
168  initialEstimate.insert(EdgeKey(a, c).key(), E2)
169  total_dimension += E1.dim() + E2.dim()
170  else:
171  raise ValueError(f"Unknown method {method}")
172 
173  if method == "Essential+Ks":
174  # Insert initial calibrations
175  for i in range(num_cameras):
176  initialEstimate.insert(K(i), cal)
177  total_dimension += cal.dim()
178 
179  print(f"Total dimension of the problem: {total_dimension}")
180  return initialEstimate
181 
182 
183 def optimize(graph, initialEstimate, method):
184  """optimize the graph"""
185  params = LevenbergMarquardtParams()
186  params.setlambdaInitial(1e10) # Initialize lambda to a high value
187  params.setlambdaUpperBound(1e10)
188  # params.setAbsoluteErrorTol(0.1)
189  params.setVerbosityLM("SUMMARY")
190  optimizer = LevenbergMarquardtOptimizer(graph, initialEstimate, params)
191  result = optimizer.optimize()
192  iterations = optimizer.iterations()
193  return result, iterations
194 
195 
196 def compute_distances(method, result, ground_truth, num_cameras, cal):
197  """Compute geodesic distances from ground truth"""
198  distances = []
199 
200  F1, F2 = ground_truth["Fundamental"]
201 
202  for a in range(num_cameras):
203  b = (a + 1) % num_cameras
204  c = (a + 2) % num_cameras
205  key_ab = EdgeKey(a, b).key()
206  key_ac = EdgeKey(a, c).key()
207 
208  if method in ["Essential+Ks", "Calibrated"]:
209  E_est_ab = result.atEssentialMatrix(key_ab)
210  E_est_ac = result.atEssentialMatrix(key_ac)
211 
212  # Compute estimated FundamentalMatrices
213  if method == "Fundamental":
214  F_est_ab = result.atFundamentalMatrix(key_ab)
215  F_est_ac = result.atFundamentalMatrix(key_ac)
216  elif method == "SimpleF":
217  SF_est_ab = result.atSimpleFundamentalMatrix(key_ab).matrix()
218  SF_est_ac = result.atSimpleFundamentalMatrix(key_ac).matrix()
219  F_est_ab = FundamentalMatrix(SF_est_ab)
220  F_est_ac = FundamentalMatrix(SF_est_ac)
221  elif method == "Essential+Ks":
222  # Retrieve calibrations from result:
223  cal_a = result.atCal3f(K(a))
224  cal_b = result.atCal3f(K(b))
225  cal_c = result.atCal3f(K(c))
226 
227  # Convert estimated EssentialMatrices to FundamentalMatrices
228  F_est_ab = FundamentalMatrix(cal_a.K(), E_est_ab, cal_b.K())
229  F_est_ac = FundamentalMatrix(cal_a.K(), E_est_ac, cal_c.K())
230  elif method == "Calibrated":
231  # Use ground truth calibration
232  F_est_ab = FundamentalMatrix(cal.K(), E_est_ab, cal.K())
233  F_est_ac = FundamentalMatrix(cal.K(), E_est_ac, cal.K())
234  else:
235  raise ValueError(f"Unknown method {method}")
236 
237  # Compute local coordinates (geodesic distance on the F-manifold)
238  dist_ab = np.linalg.norm(F1.localCoordinates(F_est_ab))
239  dist_ac = np.linalg.norm(F2.localCoordinates(F_est_ac))
240  distances.append(dist_ab)
241  distances.append(dist_ac)
242 
243  return distances
244 
245 
246 def plot_results(results):
247  """plot results"""
248  methods = list(results.keys())
249  final_errors = [results[method]["final_error"] for method in methods]
250  distances = [results[method]["distances"] for method in methods]
251  iterations = [results[method]["iterations"] for method in methods]
252 
253  fig, ax1 = plt.subplots()
254 
255  color = "tab:red"
256  ax1.set_xlabel("Method")
257  ax1.set_ylabel("Median Error (log scale)", color=color)
258  ax1.set_yscale("log")
259  ax1.bar(methods, final_errors, color=color, alpha=0.6)
260  ax1.tick_params(axis="y", labelcolor=color)
261 
262  ax2 = ax1.twinx()
263  color = "tab:blue"
264  ax2.set_ylabel("Median Geodesic Distance", color=color)
265  ax2.plot(methods, distances, color=color, marker="o", linestyle="-")
266  ax2.tick_params(axis="y", labelcolor=color)
267 
268  # Annotate the blue data points with the average number of iterations
269  for i, method in enumerate(methods):
270  ax2.annotate(
271  f"{iterations[i]:.1f}",
272  (i, distances[i]),
273  textcoords="offset points",
274  xytext=(0, 10),
275  ha="center",
276  color=color,
277  )
278 
279  plt.title("Comparison of Methods (Labels show avg iterations)")
280  fig.tight_layout()
281  plt.show()
282 
283 
284 # Main function
285 def main():
286  # Parse command line arguments
287  parser = argparse.ArgumentParser(description="Compare Fundamental and Essential Matrix Methods")
288  parser.add_argument("--num_cameras", type=int, default=4, help="Number of cameras (default: 4)")
289  parser.add_argument("--num_extra_points", type=int, default=12, help="Number of extra random points (default: 12)")
290  parser.add_argument("--num_trials", type=int, default=5, help="Number of trials (default: 5)")
291  parser.add_argument("--seed", type=int, default=42, help="Random seed (default: 42)")
292  parser.add_argument("--noise_std", type=float, default=0.5, help="Standard deviation of noise (default: 0.5)")
293  args = parser.parse_args()
294 
295  # Initialize the random number generator
296  rng = np.random.default_rng(seed=args.seed)
297 
298  # Initialize results dictionary
299  results = {method: {"distances": [], "final_error": [], "iterations": []} for method in methods}
300 
301  # Simulate geometry
302  points, poses, cal = simulate_geometry(args.num_cameras, rng, args.num_extra_points)
303 
304  # Compute ground truth matrices
305  ground_truth = {method: compute_ground_truth(method, poses, cal) for method in methods}
306 
307  # Get initial estimates
308  initial_estimate: dict[Values] = {
309  method: get_initial_estimate(method, args.num_cameras, ground_truth[method], cal) for method in methods
310  }
311  simple_f_result: Values = Values()
312 
313  for trial in range(args.num_trials):
314  print(f"\nTrial {trial + 1}/{args.num_trials}")
315 
316  # Simulate data
317  measurements = simulate_data(points, poses, cal, rng, args.noise_std)
318 
319  for method in methods:
320  print(f"\nRunning method: {method}")
321 
322  # Build the factor graph
323  graph = build_factor_graph(method, args.num_cameras, measurements, cal)
324 
325  # For F, initialize from SimpleF:
326  if method == "Fundamental":
327  initial_estimate[method] = simple_f_result
328 
329  # Optimize the graph
330  result, iterations = optimize(graph, initial_estimate[method], method)
331 
332  # Store SimpleF result as a set of FundamentalMatrices
333  if method == "SimpleF":
334  simple_f_result = Values()
335  for a in range(args.num_cameras):
336  b = (a + 1) % args.num_cameras # Next camera
337  c = (a + 2) % args.num_cameras # Camera after next
338  key_ab = EdgeKey(a, b).key()
339  key_ac = EdgeKey(a, c).key()
340  F1 = result.atSimpleFundamentalMatrix(key_ab).matrix()
341  F2 = result.atSimpleFundamentalMatrix(key_ac).matrix()
342  simple_f_result.insert(key_ab, FundamentalMatrix(F1))
343  simple_f_result.insert(key_ac, FundamentalMatrix(F2))
344 
345  # Compute distances from ground truth
346  distances = compute_distances(method, result, ground_truth, args.num_cameras, cal)
347 
348  # Compute final error
349  final_error = graph.error(result)
350 
351  # Store results
352  results[method]["distances"].extend(distances)
353  results[method]["final_error"].append(final_error)
354  results[method]["iterations"].append(iterations)
355 
356  print(f"Method: {method}")
357  print(f"Final Error: {final_error:.3f}")
358  print(f"Mean Geodesic Distance: {np.mean(distances):.3f}")
359  print(f"Number of Iterations: {iterations}\n")
360 
361  # Average results over trials
362  for method in methods:
363  results[method]["final_error"] = np.median(results[method]["final_error"])
364  results[method]["distances"] = np.median(results[method]["distances"])
365  results[method]["iterations"] = np.median(results[method]["iterations"])
366 
367  # Plot results
368  plot_results(results)
369 
370 
371 if __name__ == "__main__":
372  main()
gtsam.examples.ViewGraphComparison.compute_distances
def compute_distances(method, result, ground_truth, num_cameras, cal)
Definition: ViewGraphComparison.py:196
gtsam.examples.ViewGraphComparison.simulate_geometry
def simulate_geometry(num_cameras, rng, num_random_points=12)
Definition: ViewGraphComparison.py:48
gtsam.examples.ViewGraphComparison.optimize
def optimize(graph, initialEstimate, method)
Definition: ViewGraphComparison.py:183
gtsam::between
Expression< T > between(const Expression< T > &t1, const Expression< T > &t2)
Definition: nonlinear/expressions.h:17
list
Definition: pytypes.h:2166
gtsam::Cal3f
Calibration model with a single focal length and zero skew.
Definition: Cal3f.h:35
gtsam.examples.ViewGraphComparison.build_factor_graph
def build_factor_graph(method, num_cameras, measurements, cal)
Definition: ViewGraphComparison.py:99
gtsam.examples
Definition: python/gtsam/examples/__init__.py:1
gtsam.examples.ViewGraphComparison.K
K
Definition: ViewGraphComparison.py:32
gtsam::range
Double_ range(const Point2_ &p, const Point2_ &q)
Definition: slam/expressions.h:30
gtsam::print
void print(const Matrix &A, const string &s, ostream &stream)
Definition: Matrix.cpp:155
gtsam.examples.ViewGraphComparison.main
def main()
Definition: ViewGraphComparison.py:285
gtsam.examples.ViewGraphComparison.compute_ground_truth
def compute_ground_truth(method, poses, cal)
Definition: ViewGraphComparison.py:80
gtsam::NonlinearFactorGraph
Definition: NonlinearFactorGraph.h:55
gtsam.examples.ViewGraphComparison.get_initial_estimate
def get_initial_estimate(method, num_cameras, ground_truth, cal)
Definition: ViewGraphComparison.py:149
gtsam.examples.ViewGraphComparison.plot_results
def plot_results(results)
Definition: ViewGraphComparison.py:246
gtsam::LevenbergMarquardtOptimizer
Definition: LevenbergMarquardtOptimizer.h:35
matrix
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: gtsam/3rdparty/Eigen/blas/common.h:110
key
const gtsam::Symbol key('X', 0)
gtsam::EdgeKey
Definition: EdgeKey.h:25
gtsam::SimpleFundamentalMatrix
Class for representing a simple fundamental matrix.
Definition: FundamentalMatrix.h:128
gtsam::Values
Definition: Values.h:65
gtsam::LevenbergMarquardtParams
Definition: LevenbergMarquardtParams.h:35
gtsam::Point3
Vector3 Point3
Definition: Point3.h:38
gtsam::noiseModel::Isotropic::Sigma
static shared_ptr Sigma(size_t dim, double sigma, bool smart=true)
Definition: NoiseModel.cpp:624
len
size_t len(handle h)
Get the length of a Python object.
Definition: pytypes.h:2446
gtsam.examples.ViewGraphComparison.formatter
def formatter(key)
Definition: ViewGraphComparison.py:39
gtsam::FundamentalMatrix
Represents a fundamental matrix in computer vision, which encodes the epipolar geometry between two v...
Definition: FundamentalMatrix.h:28
gtsam.examples.ViewGraphComparison.simulate_data
def simulate_data(points, poses, cal, rng, noise_std)
Definition: ViewGraphComparison.py:66
gtsam::Symbol
Definition: inference/Symbol.h:37


gtsam
Author(s):
autogenerated on Sat Nov 16 2024 04:09:42