HybridGaussianFactorGraph.cpp
Go to the documentation of this file.
1 /* ----------------------------------------------------------------------------
2 
3  * GTSAM Copyright 2010, Georgia Tech Research Corporation,
4  * Atlanta, Georgia 30332-0415
5  * All Rights Reserved
6  * Authors: Frank Dellaert, et al. (see THANKS for the full author list)
7 
8  * See LICENSE for the license information
9 
10  * -------------------------------------------------------------------------- */
11 
21 #include <gtsam/base/utilities.h>
38 #include <gtsam/inference/Key.h>
45 
46 #include <cstddef>
47 #include <iostream>
48 #include <memory>
49 #include <stdexcept>
50 #include <utility>
51 #include <vector>
52 
53 #define GTSAM_HYBRID_WITH_TABLEFACTOR 1
54 
55 namespace gtsam {
56 
58 template class EliminateableFactorGraph<HybridGaussianFactorGraph>;
59 
60 using std::dynamic_pointer_cast;
62 
64 struct Result {
65  // Gaussian conditional resulting from elimination.
67  double negLogK; // Negative log of the normalization constant K.
68  GaussianFactor::shared_ptr factor; // Leftover factor 𝜏.
69  double scalar; // Scalar value associated with factor 𝜏.
70 
71  bool operator==(const Result &other) const {
72  return conditional == other.conditional && negLogK == other.negLogK &&
73  factor == other.factor && scalar == other.scalar;
74  }
75 };
77 
78 static const VectorValues kEmpty;
79 
80 /* ************************************************************************ */
81 // Throw a runtime exception for method specified in string s, and factor f:
82 static void throwRuntimeError(const std::string &s,
83  const std::shared_ptr<Factor> &f) {
84  auto &fr = *f;
85  throw std::runtime_error(s + " not implemented for factor type " +
86  demangle(typeid(fr).name()) + ".");
87 }
88 
89 /* ************************************************************************ */
91  KeySet discrete_keys = graph.discreteKeySet();
92  const VariableIndex index(graph);
94  index, KeyVector(discrete_keys.begin(), discrete_keys.end()), true);
95 }
96 
97 /* ************************************************************************ */
98 static void printFactor(const std::shared_ptr<Factor> &factor,
99  const DiscreteValues &assignment,
100  const KeyFormatter &keyFormatter) {
101  if (auto hgf = dynamic_pointer_cast<HybridGaussianFactor>(factor)) {
102  if (assignment.empty()) {
103  hgf->print("HybridGaussianFactor:", keyFormatter);
104  } else {
105  hgf->operator()(assignment)
106  .first->print("HybridGaussianFactor, component:", keyFormatter);
107  }
108  } else if (auto gf = dynamic_pointer_cast<GaussianFactor>(factor)) {
109  factor->print("GaussianFactor:\n", keyFormatter);
110 
111  } else if (auto df = dynamic_pointer_cast<DiscreteFactor>(factor)) {
112  factor->print("DiscreteFactor:\n", keyFormatter);
113  } else if (auto hc = dynamic_pointer_cast<HybridConditional>(factor)) {
114  if (hc->isContinuous()) {
115  factor->print("GaussianConditional:\n", keyFormatter);
116  } else if (hc->isDiscrete()) {
117  factor->print("DiscreteConditional:\n", keyFormatter);
118  } else {
119  if (assignment.empty()) {
120  hc->print("HybridConditional:", keyFormatter);
121  } else {
122  hc->asHybrid()
123  ->choose(assignment)
124  ->print("HybridConditional, component:\n", keyFormatter);
125  }
126  }
127  } else {
128  factor->print("Unknown factor type\n", keyFormatter);
129  }
130 }
131 
132 /* ************************************************************************ */
133 void HybridGaussianFactorGraph::print(const std::string &s,
134  const KeyFormatter &keyFormatter) const {
135  std::cout << (s.empty() ? "" : s + " ") << std::endl;
136  std::cout << "size: " << size() << std::endl;
137 
138  for (size_t i = 0; i < factors_.size(); i++) {
139  auto &&factor = factors_[i];
140  if (factor == nullptr) {
141  std::cout << "Factor " << i << ": nullptr\n";
142  continue;
143  }
144  // Print the factor
145  std::cout << "Factor " << i << "\n";
146  printFactor(factor, {}, keyFormatter);
147  std::cout << "\n";
148  }
149  std::cout.flush();
150 }
151 
152 /* ************************************************************************ */
154  const HybridValues &values, const std::string &str,
155  const KeyFormatter &keyFormatter,
156  const std::function<bool(const Factor * /*factor*/,
157  double /*whitenedError*/, size_t /*index*/)>
158  &printCondition) const {
159  std::cout << str << "size: " << size() << std::endl << std::endl;
160 
161  for (size_t i = 0; i < factors_.size(); i++) {
162  auto &&factor = factors_[i];
163  if (factor == nullptr) {
164  std::cout << "Factor " << i << ": nullptr\n";
165  continue;
166  }
167  const double errorValue = factor->error(values);
168  if (!printCondition(factor.get(), errorValue, i))
169  continue; // User-provided filter did not pass
170 
171  // Print the factor
172  std::cout << "Factor " << i << ", error = " << errorValue << "\n";
173  printFactor(factor, values.discrete(), keyFormatter);
174  std::cout << "\n";
175  }
176  std::cout.flush();
177 }
178 
179 /* ************************************************************************ */
181  const {
183 
184  for (auto &f : factors_) {
185  // TODO(dellaert): can we make this cleaner and less error-prone?
186  if (auto orphan = dynamic_pointer_cast<OrphanWrapper>(f)) {
187  continue; // Ignore OrphanWrapper
188  } else if (auto gf = dynamic_pointer_cast<GaussianFactor>(f)) {
189  result += gf;
190  } else if (auto gc = dynamic_pointer_cast<GaussianConditional>(f)) {
191  result += gc;
192  } else if (auto gmf = dynamic_pointer_cast<HybridGaussianFactor>(f)) {
193  result += *gmf;
194  } else if (auto gm = dynamic_pointer_cast<HybridGaussianConditional>(f)) {
195  result += *gm; // handled above already?
196  } else if (auto hc = dynamic_pointer_cast<HybridConditional>(f)) {
197  if (auto gm = hc->asHybrid()) {
198  result += *gm;
199  } else if (auto g = hc->asGaussian()) {
200  result += g;
201  } else {
202  // Has to be discrete.
203  // TODO(dellaert): in C++20, we can use std::visit.
204  continue;
205  }
206  } else if (dynamic_pointer_cast<DiscreteFactor>(f)) {
207  // Don't do anything for discrete-only factors
208  // since we want to eliminate continuous values only.
209  continue;
210  } else {
211  // TODO(dellaert): there was an unattributed comment here: We need to
212  // handle the case where the object is actually an BayesTreeOrphanWrapper!
213  throwRuntimeError("gtsam::collectProductFactor", f);
214  }
215  }
216 
217  return result;
218 }
219 
220 /* ************************************************************************ */
221 static std::pair<HybridConditional::shared_ptr, std::shared_ptr<Factor>>
223  const Ordering &frontalKeys) {
225  for (auto &f : factors) {
226  if (auto gf = dynamic_pointer_cast<GaussianFactor>(f)) {
227  gfg.push_back(gf);
228  } else if (auto orphan = dynamic_pointer_cast<OrphanWrapper>(f)) {
229  // Ignore orphaned clique.
230  // TODO(dellaert): is this correct? If so explain here.
231  } else if (auto hc = dynamic_pointer_cast<HybridConditional>(f)) {
232  auto gc = hc->asGaussian();
233  if (!gc) throwRuntimeError("continuousElimination", gc);
234  gfg.push_back(gc);
235  } else {
236  throwRuntimeError("continuousElimination", f);
237  }
238  }
239 
240  auto result = EliminatePreferCholesky(gfg, frontalKeys);
241  return {std::make_shared<HybridConditional>(result.first), result.second};
242 }
243 
244 /* ************************************************************************ */
253  const DiscreteKeys &discreteKeys,
254  const AlgebraicDecisionTree<Key> &errors) {
255  double min_log = errors.min();
256  AlgebraicDecisionTree<Key> potentials(
257  errors, [&min_log](const double x) { return exp(-(x - min_log)); });
258 #if GTSAM_HYBRID_WITH_TABLEFACTOR
259  return std::make_shared<TableFactor>(discreteKeys, potentials);
260 #else
261  return std::make_shared<DecisionTreeFactor>(discreteKeys, potentials);
262 #endif
263 }
264 
265 /* ************************************************************************ */
269 
270  for (auto &f : factors) {
271  if (auto df = dynamic_pointer_cast<DiscreteFactor>(f)) {
272  dfg.push_back(df);
273 
274  } else if (auto gmf = dynamic_pointer_cast<HybridGaussianFactor>(f)) {
275  // Case where we have a HybridGaussianFactor with no continuous keys.
276  // In this case, compute a discrete factor from the remaining error.
277  auto calculateError = [&](const auto &pair) -> double {
278  auto [factor, scalar] = pair;
279  // If factor is null, it has been pruned, hence return infinite error
280  if (!factor) return std::numeric_limits<double>::infinity();
281  return scalar + factor->error(kEmpty);
282  };
283  AlgebraicDecisionTree<Key> errors(gmf->factors(), calculateError);
284  dfg.push_back(DiscreteFactorFromErrors(gmf->discreteKeys(), errors));
285 
286  } else if (auto orphan = dynamic_pointer_cast<OrphanWrapper>(f)) {
287  // Ignore orphaned clique.
288  // TODO(dellaert): is this correct? If so explain here.
289  } else if (auto hc = dynamic_pointer_cast<HybridConditional>(f)) {
290  auto dc = hc->asDiscrete();
291  if (!dc) throwRuntimeError("discreteElimination", dc);
292 #if GTSAM_HYBRID_TIMING
293  gttic_(ConvertConditionalToTableFactor);
294 #endif
295  if (auto dtc = std::dynamic_pointer_cast<TableDistribution>(dc)) {
297  dfg.push_back(dtc->table());
298  } else {
299 #if GTSAM_HYBRID_WITH_TABLEFACTOR
300  // Convert DiscreteConditional to TableFactor
301  auto tdc = std::make_shared<TableFactor>(*dc);
302  dfg.push_back(tdc);
303 #else
304  dfg.push_back(dc);
305 #endif
306  }
307 #if GTSAM_HYBRID_TIMING
308  gttoc_(ConvertConditionalToTableFactor);
309 #endif
310  } else {
311  throwRuntimeError("discreteElimination", f);
312  }
313  }
314 
315  return dfg;
316 }
317 
318 /* ************************************************************************ */
319 static std::pair<HybridConditional::shared_ptr, std::shared_ptr<Factor>>
321  const Ordering &frontalKeys) {
322 #if GTSAM_HYBRID_TIMING
324 #endif
326 #if GTSAM_HYBRID_TIMING
328 #endif
329 
330 #if GTSAM_HYBRID_TIMING
332 #endif
333 #if GTSAM_HYBRID_WITH_TABLEFACTOR
334  // Check if separator is empty.
335  // This is the same as checking if the number of frontal variables
336  // is the same as the number of variables in the DiscreteFactorGraph.
337  // If the separator is empty, we have a clique of all the discrete variables
338  // so we can use the TableFactor for efficiency.
339  if (frontalKeys.size() == dfg.keys().size()) {
340  // Get product factor
342 
343  // Check type of product, and get as TableFactor for efficiency.
344  // Use object instead of pointer since we need it
345  // for the TableDistribution constructor.
346  TableFactor p;
347  if (auto tf = std::dynamic_pointer_cast<TableFactor>(product)) {
348  p = *tf;
349  } else {
350  p = TableFactor(product->toDecisionTreeFactor());
351  }
352  auto conditional = std::make_shared<TableDistribution>(p);
353 
354  DiscreteFactor::shared_ptr sum = p.sum(frontalKeys);
355 
356  return {std::make_shared<HybridConditional>(conditional), sum};
357 
358  } else {
359 #endif
360  // Perform sum-product.
361  auto result = EliminateDiscrete(dfg, frontalKeys);
362  return {std::make_shared<HybridConditional>(result.first), result.second};
363 #if GTSAM_HYBRID_WITH_TABLEFACTOR
364  }
365 #endif
366 #if GTSAM_HYBRID_TIMING
368 #endif
369 }
370 
371 /* ************************************************************************ */
378 static std::shared_ptr<Factor> createDiscreteFactor(
379  const ResultTree &eliminationResults,
380  const DiscreteKeys &discreteSeparator) {
381  auto calculateError = [&](const Result &result) -> double {
382  if (result.conditional && result.factor) {
383  // `error` has the following contributions:
384  // - the scalar is the sum of all mode-dependent constants
385  // - factor->error(kempty) is the error remaining after elimination
386  // - negLogK is what is given to the conditional to normalize
387  return result.scalar + result.factor->error(kEmpty) - result.negLogK;
388  } else if (!result.conditional && !result.factor) {
389  // If the factor has been pruned, return infinite error
390  return std::numeric_limits<double>::infinity();
391  } else {
392  throw std::runtime_error("createDiscreteFactor has mixed NULLs");
393  }
394  };
395 
396 #if GTSAM_HYBRID_TIMING
397  gttic_(DiscreteBoundaryErrors);
398 #endif
399  AlgebraicDecisionTree<Key> errors(eliminationResults, calculateError);
400 #if GTSAM_HYBRID_TIMING
401  gttoc_(DiscreteBoundaryErrors);
402  gttic_(DiscreteBoundaryResult);
403 #endif
404  auto result = DiscreteFactorFromErrors(discreteSeparator, errors);
405 #if GTSAM_HYBRID_TIMING
406  gttoc_(DiscreteBoundaryResult);
407 #endif
408  return result;
409 }
410 
411 /* *******************************************************************************/
412 // Create HybridGaussianFactor on the separator, taking care to correct
413 // for conditional constants.
414 static std::shared_ptr<Factor> createHybridGaussianFactor(
415  const ResultTree &eliminationResults,
416  const DiscreteKeys &discreteSeparator) {
417  // Correct for the normalization constant used up by the conditional
418  auto correct = [&](const Result &result) -> GaussianFactorValuePair {
419  if (result.conditional && result.factor) {
420  return {result.factor, result.scalar - result.negLogK};
421  } else if (!result.conditional && !result.factor) {
422  return {nullptr, std::numeric_limits<double>::infinity()};
423  } else {
424  throw std::runtime_error("createHybridGaussianFactors has mixed NULLs");
425  }
426  };
427 #if GTSAM_HYBRID_TIMING
428  gttic_(HybridCreateGaussianFactor);
429 #endif
430  DecisionTree<Key, GaussianFactorValuePair> newFactors(eliminationResults,
431  correct);
432 #if GTSAM_HYBRID_TIMING
433  gttoc_(HybridCreateGaussianFactor);
434 #endif
435 
436  return std::make_shared<HybridGaussianFactor>(discreteSeparator, newFactors);
437 }
438 
439 /* *******************************************************************************/
441 static auto GetDiscreteKeys =
442  [](const HybridGaussianFactorGraph &hfg) -> DiscreteKeys {
443  const std::set<DiscreteKey> discreteKeySet = hfg.discreteKeys();
444  return {discreteKeySet.begin(), discreteKeySet.end()};
445 };
446 
447 /* *******************************************************************************/
448 std::pair<HybridConditional::shared_ptr, std::shared_ptr<Factor>>
450  // Since we eliminate all continuous variables first,
451  // the discrete separator will be *all* the discrete keys.
452  DiscreteKeys discreteSeparator = GetDiscreteKeys(*this);
453 
454 #if GTSAM_HYBRID_TIMING
455  gttic_(HybridCollectProductFactor);
456 #endif
457  // Collect all the factors to create a set of Gaussian factor graphs in a
458  // decision tree indexed by all discrete keys involved. Just like any hybrid
459  // factor, every assignment also has a scalar error, in this case the sum of
460  // all errors in the graph. This error is assignment-specific and accounts for
461  // any difference in noise models used.
463 #if GTSAM_HYBRID_TIMING
464  gttoc_(HybridCollectProductFactor);
465 #endif
466 
467  // Check if a factor is null
468  auto isNull = [](const GaussianFactor::shared_ptr &ptr) { return !ptr; };
469 
470  // This is the elimination method on the leaf nodes
471  bool someContinuousLeft = false;
472  auto eliminate =
473  [&](const std::pair<GaussianFactorGraph, double> &pair) -> Result {
474  const auto &[graph, scalar] = pair;
475 
476  // If any product contains a pruned factor, prune it here. Done here as it's
477  // non non-trivial to do within collectProductFactor.
478  if (graph.empty() || std::any_of(graph.begin(), graph.end(), isNull)) {
479  return {nullptr, 0.0, nullptr, 0.0};
480  }
481 
482  // Expensive elimination of product factor.
483  auto [conditional, factor] =
485 
486  // Record whether there any continuous variables left
487  someContinuousLeft |= !factor->empty();
488 
489  // We pass on the scalar unmodified.
490  return {conditional, conditional->negLogConstant(), factor, scalar};
491  };
492 
493 #if GTSAM_HYBRID_TIMING
494  gttic_(HybridEliminate);
495 #endif
496  // Perform elimination!
497  const ResultTree eliminationResults(productFactor, eliminate);
498 #if GTSAM_HYBRID_TIMING
499  gttoc_(HybridEliminate);
500 #endif
501 
502  // If there are no more continuous parents we create a DiscreteFactor with the
503  // error for each discrete choice. Otherwise, create a HybridGaussianFactor
504  // on the separator, taking care to correct for conditional constants.
505  auto newFactor =
506  someContinuousLeft
507  ? createHybridGaussianFactor(eliminationResults, discreteSeparator)
508  : createDiscreteFactor(eliminationResults, discreteSeparator);
509 
510  // Create the HybridGaussianConditional without re-calculating constants:
512  eliminationResults, [](const Result &result) -> GaussianFactorValuePair {
513  return {result.conditional, result.negLogK};
514  });
515  auto hybridGaussian =
516  std::make_shared<HybridGaussianConditional>(discreteSeparator, pairs);
517 
518  return {std::make_shared<HybridConditional>(hybridGaussian), newFactor};
519 }
520 
521 /* ************************************************************************
522  * Function to eliminate variables **under the following assumptions**:
523  * 1. When the ordering is fully continuous, and the graph only contains
524  * continuous and hybrid factors
525  * 2. When the ordering is fully discrete, and the graph only contains discrete
526  * factors
527  *
528  * Any usage outside of this is considered incorrect.
529  *
530  * \warning This function is not meant to be used with arbitrary hybrid factor
531  * graphs. For example, if there exists continuous parents, and one tries to
532  * eliminate a discrete variable (as specified in the ordering), the result will
533  * be INCORRECT and there will be NO error raised.
534  */
535 std::pair<HybridConditional::shared_ptr, std::shared_ptr<Factor>> //
537  const Ordering &keys) {
538  // NOTE: Because we are in the Conditional Gaussian regime there are only
539  // a few cases:
540  // 1. continuous variable, make a hybrid Gaussian conditional if there are
541  // hybrid factors;
542  // 2. continuous variable, we make a Gaussian Factor if there are no hybrid
543  // factors;
544  // 3. discrete variable, no continuous factor is allowed
545  // (escapes Conditional Gaussian regime), if discrete only we do the discrete
546  // elimination.
547 
548  // However it is not that simple. During elimination it is possible that the
549  // multifrontal needs to eliminate an ordering that contains both Gaussian and
550  // hybrid variables, for example x1, c1.
551  // In this scenario, we will have a density P(x1, c1) that is a Conditional
552  // Linear Gaussian P(x1|c1)P(c1) (see Murphy02).
553 
554  // The issue here is that, how can we know which variable is discrete if we
555  // unify Values? Obviously we can tell using the factors, but is that fast?
556 
557  // In the case of multifrontal, we will need to use a constrained ordering
558  // so that the discrete parts will be guaranteed to be eliminated last!
559  // Because of all these reasons, we carefully consider how to
560  // implement the hybrid factors so that we do not get poor performance.
561 
562  // The first thing is how to represent the HybridGaussianConditional.
563  // A very possible scenario is that the incoming factors will have different
564  // levels of discrete keys. For example, imagine we are going to eliminate the
565  // fragment: $\phi(x1,c1,c2)$, $\phi(x1,c2,c3)$, which is perfectly valid.
566  // Now we will need to know how to retrieve the corresponding continuous
567  // densities for the assignment (c1,c2,c3) (OR (c2,c3,c1), note there is NO
568  // defined order!). We also need to consider when there is pruning. Two
569  // hybrid factors could have different pruning patterns - one could have
570  // (c1=0,c2=1) pruned, and another could have (c2=0,c3=1) pruned, and this
571  // creates a big problem in how to identify the intersection of non-pruned
572  // branches.
573 
574  // Our approach is first building the collection of all discrete keys. After
575  // that we enumerate the space of all key combinations *lazily* so that the
576  // exploration branch terminates whenever an assignment yields NULL in any of
577  // the hybrid factors.
578 
579  // When the number of assignments is large we may encounter stack overflows.
580  // However this is also the case with iSAM2, so no pressure :)
581 
582  // Check the factors:
583  // 1. if all factors are discrete, then we can do discrete elimination:
584  // 2. if all factors are continuous, then we can do continuous elimination:
585  // 3. if not, we do hybrid elimination:
586 
587  bool only_discrete = true, only_continuous = true;
588  for (auto &&factor : factors) {
589  if (auto hybrid_factor = std::dynamic_pointer_cast<HybridFactor>(factor)) {
590  if (hybrid_factor->isDiscrete()) {
591  only_continuous = false;
592  } else if (hybrid_factor->isContinuous()) {
593  only_discrete = false;
594  } else if (hybrid_factor->isHybrid()) {
595  only_continuous = false;
596  only_discrete = false;
597  break;
598  }
599  } else if (auto cont_factor =
600  std::dynamic_pointer_cast<GaussianFactor>(factor)) {
601  only_discrete = false;
602  } else if (auto discrete_factor =
603  std::dynamic_pointer_cast<DiscreteFactor>(factor)) {
604  only_continuous = false;
605  }
606  }
607 
608  // NOTE: We should really defer the product here because of pruning
609 
610  if (only_discrete) {
611  // Case 1: we are only dealing with discrete
613  } else if (only_continuous) {
614  // Case 2: we are only dealing with continuous
616  } else {
617  // Case 3: We are now in the hybrid land!
618  return factors.eliminate(keys);
619  }
620 }
621 
622 /* ************************************************************************ */
624  const VectorValues &continuousValues) const {
626  // Iterate over each factor.
627  for (auto &factor : factors_) {
628  if (auto hf = std::dynamic_pointer_cast<HybridFactor>(factor)) {
629  // Add errorTree for hybrid factors, includes HybridGaussianConditionals!
630  result = result + hf->errorTree(continuousValues);
631  } else if (auto df = std::dynamic_pointer_cast<DiscreteFactor>(factor)) {
632  // If discrete, just add its errorTree as well
633  result = result + df->errorTree();
634  } else if (auto gf = dynamic_pointer_cast<GaussianFactor>(factor)) {
635  // For a continuous only factor, just add its error
636  result = result + gf->error(continuousValues);
637  } else {
638  throwRuntimeError("HybridGaussianFactorGraph::errorTree", factor);
639  }
640  }
641  return result;
642 }
643 
644 /* ************************************************************************ */
646  double error = this->error(values);
647  // NOTE: The 0.5 term is handled by each factor
648  return std::exp(-error);
649 }
650 
651 /* ************************************************************************ */
653  const VectorValues &continuousValues) const {
654  AlgebraicDecisionTree<Key> errors = this->errorTree(continuousValues);
655  AlgebraicDecisionTree<Key> p = errors.apply([](double error) {
656  // NOTE: The 0.5 term is handled by each factor
657  return exp(-error);
658  });
659  return p / p.sum();
660 }
661 
662 /* ************************************************************************ */
664  const DiscreteValues &assignment) const {
666  for (auto &&f : *this) {
667  if (auto gf = std::dynamic_pointer_cast<GaussianFactor>(f)) {
668  gfg.push_back(gf);
669  } else if (auto gc = std::dynamic_pointer_cast<GaussianConditional>(f)) {
670  gfg.push_back(gf);
671  } else if (auto hgf = std::dynamic_pointer_cast<HybridGaussianFactor>(f)) {
672  gfg.push_back((*hgf)(assignment).first);
673  } else if (auto hgc = dynamic_pointer_cast<HybridGaussianConditional>(f)) {
674  gfg.push_back((*hgc)(assignment));
675  } else if (auto hc = dynamic_pointer_cast<HybridConditional>(f)) {
676  if (auto gc = hc->asGaussian())
677  gfg.push_back(gc);
678  else if (auto hgc = hc->asHybrid())
679  gfg.push_back((*hgc)(assignment));
680  } else {
681  continue;
682  }
683  }
684  return gfg;
685 }
686 
687 /* ************************************************************************ */
690  for (auto &&f : factors_) {
691  if (auto discreteFactor = std::dynamic_pointer_cast<DiscreteFactor>(f)) {
692  dfg.push_back(discreteFactor);
693  }
694  }
695  return dfg;
696 }
697 
698 } // namespace gtsam
gtsam::BayesTreeOrphanWrapper
Definition: BayesTree.h:283
gtsam::EliminateDiscrete
std::pair< DiscreteConditional::shared_ptr, DiscreteFactor::shared_ptr > EliminateDiscrete(const DiscreteFactorGraph &factors, const Ordering &frontalKeys)
Main elimination function for DiscreteFactorGraph.
Definition: DiscreteFactorGraph.cpp:202
HybridJunctionTree.h
gtsam::Ordering::ColamdConstrainedLast
static Ordering ColamdConstrainedLast(const FACTOR_GRAPH &graph, const KeyVector &constrainLast, bool forceOrder=false)
Definition: inference/Ordering.h:112
TableDistribution.h
gtsam::TableFactor
Definition: TableFactor.h:54
gtsam::Result
Result from elimination.
Definition: HybridGaussianFactorGraph.cpp:64
GaussianFactorGraph.h
Linear Factor Graph where all factors are Gaussians.
gtsam::HybridValues
Definition: HybridValues.h:37
GaussianConditional.h
Conditional Gaussian Base class.
gtsam::FactorGraph< Factor >::factors_
FastVector< sharedFactor > factors_
Definition: FactorGraph.h:92
gtsam::DiscreteFactorGraph
Definition: DiscreteFactorGraph.h:99
gtsam::Result::conditional
GaussianConditional::shared_ptr conditional
Definition: HybridGaussianFactorGraph.cpp:66
HybridGaussianConditional.h
A hybrid conditional in the Conditional Linear Gaussian scheme.
s
RealScalar s
Definition: level1_cplx_impl.h:126
gtsam::FactorGraph< Factor >::error
double error(const HybridValues &values) const
Definition: FactorGraph-inst.h:66
keys
const KeyVector keys
Definition: testRegularImplicitSchurFactor.cpp:40
DiscreteFactorGraph.h
gtsam::DiscreteFactorGraph::keys
KeySet keys() const
Definition: DiscreteFactorGraph.cpp:45
gtsam::DiscreteFactorFromErrors
static DiscreteFactor::shared_ptr DiscreteFactorFromErrors(const DiscreteKeys &discreteKeys, const AlgebraicDecisionTree< Key > &errors)
Take negative log-values, shift them so that the minimum value is 0, and then exponentiate to create ...
Definition: HybridGaussianFactorGraph.cpp:252
simple_graph::factors
const GaussianFactorGraph factors
Definition: testJacobianFactor.cpp:213
x
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
Definition: gnuplot_common_settings.hh:12
EliminateableFactorGraph-inst.h
gtsam::GetDiscreteKeys
static auto GetDiscreteKeys
Get the discrete keys from the HybridGaussianFactorGraph as DiscreteKeys.
Definition: HybridGaussianFactorGraph.cpp:441
gtsam::DecisionTreeFactor::error
double error(const DiscreteValues &values) const override
Calculate error for DiscreteValues x, is -log(probability).
Definition: DecisionTreeFactor.cpp:57
gtsam::HybridGaussianProductFactor
Alias for DecisionTree of GaussianFactorGraphs and their scalar sums.
Definition: HybridGaussianProductFactor.h:34
HessianFactor.h
Contains the HessianFactor class, a general quadratic factor.
different_sigmas::values
HybridValues values
Definition: testHybridBayesNet.cpp:247
gtsam::Factor
Definition: Factor.h:70
gtsam::DiscreteKeys
DiscreteKeys is a set of keys that can be assembled using the & operator.
Definition: DiscreteKey.h:41
gtsam::FastSet
Definition: FastSet.h:53
any_of
negation< all_of< negation< Ts >... > > any_of
Definition: wrap/pybind11/include/pybind11/detail/common.h:761
GaussianJunctionTree.h
exp
const EIGEN_DEVICE_FUNC ExpReturnType exp() const
Definition: ArrayCwiseUnaryOps.h:97
gtsam::HybridGaussianFactorGraph::discreteFactors
DiscreteFactorGraph discreteFactors() const
Helper method to get all the discrete factors as a DiscreteFactorGraph.
Definition: HybridGaussianFactorGraph.cpp:688
gtsam::KeyVector
FastVector< Key > KeyVector
Define collection type once and for all - also used in wrappers.
Definition: Key.h:92
result
Values result
Definition: OdometryOptimize.cpp:8
HybridGaussianFactor.h
A set of GaussianFactors, indexed by a set of discrete keys.
gtsam::HybridGaussianFactorGraph::collectProductFactor
HybridGaussianProductFactor collectProductFactor() const
Create a decision tree of factor graphs out of this hybrid factor graph.
Definition: HybridGaussianFactorGraph.cpp:180
utilities.h
hc
Vector3d hc
Definition: Tridiagonalization_householderCoefficients.cpp:5
name
static char name[]
Definition: rgamma.c:72
gtsam::AlgebraicDecisionTree< Key >
pruning_fixture::factor
DecisionTreeFactor factor(D &C &B &A, "0.0 0.0 0.0 0.60658897 0.61241912 0.61241969 0.61247685 0.61247742 0.0 " "0.0 0.0 0.99995287 1.0 1.0 1.0 1.0")
gtsam::DiscreteFactorGraph::scaledProduct
DiscreteFactor::shared_ptr scaledProduct() const
Return product of all factors as a single factor, which is scaled by the max value to prevent underfl...
Definition: DiscreteFactorGraph.cpp:117
gtsam::EliminatePreferCholesky
std::pair< std::shared_ptr< GaussianConditional >, std::shared_ptr< GaussianFactor > > EliminatePreferCholesky(const GaussianFactorGraph &factors, const Ordering &keys)
Definition: HessianFactor.cpp:539
gtsam::Result::scalar
double scalar
Definition: HybridGaussianFactorGraph.cpp:69
Key.h
gtsam::GaussianFactorGraph
Definition: GaussianFactorGraph.h:73
GaussianEliminationTree.h
HybridFactor.h
gtsam::GaussianFactorValuePair
std::pair< GaussianFactor::shared_ptr, double > GaussianFactorValuePair
Alias for pair of GaussianFactor::shared_pointer and a double value.
Definition: HybridGaussianFactor.h:38
gtsam::VectorValues
Definition: VectorValues.h:74
gttoc_
#define gttoc_(label)
Definition: timing.h:273
gttic_
#define gttic_(label)
Definition: timing.h:268
gtsam::throwRuntimeError
static void throwRuntimeError(const std::string &s, const std::shared_ptr< Factor > &f)
Definition: HybridGaussianFactorGraph.cpp:82
gtsam::createDiscreteFactor
static std::shared_ptr< Factor > createDiscreteFactor(const ResultTree &eliminationResults, const DiscreteKeys &discreteSeparator)
Definition: HybridGaussianFactorGraph.cpp:378
gtsam::KeyFormatter
std::function< std::string(Key)> KeyFormatter
Typedef for a function to format a key, i.e. to convert it to a string.
Definition: Key.h:35
gtsam::CollectDiscreteFactors
static DiscreteFactorGraph CollectDiscreteFactors(const HybridGaussianFactorGraph &factors)
Definition: HybridGaussianFactorGraph.cpp:266
DiscreteJunctionTree.h
scalar
mxArray * scalar(mxClassID classid)
Definition: matlab.h:82
gtsam::GaussianFactor::shared_ptr
std::shared_ptr< This > shared_ptr
shared_ptr to this class
Definition: GaussianFactor.h:42
gtsam::continuousElimination
static std::pair< HybridConditional::shared_ptr, std::shared_ptr< Factor > > continuousElimination(const HybridGaussianFactorGraph &factors, const Ordering &frontalKeys)
Definition: HybridGaussianFactorGraph.cpp:222
gtsam::FactorGraph< Factor >::keys
KeySet keys() const
Definition: FactorGraph-inst.h:85
Assignment.h
An assignment from labels to a discrete value index (size_t)
gtsam::discreteElimination
static std::pair< HybridConditional::shared_ptr, std::shared_ptr< Factor > > discreteElimination(const HybridGaussianFactorGraph &factors, const Ordering &frontalKeys)
Definition: HybridGaussianFactorGraph.cpp:320
gtsam::HybridGaussianFactorGraph::eliminate
std::pair< std::shared_ptr< HybridConditional >, std::shared_ptr< Factor > > eliminate(const Ordering &keys) const
Eliminate the given continuous keys.
Definition: HybridGaussianFactorGraph.cpp:449
gtsam::HybridGaussianFactorGraph
Definition: HybridGaussianFactorGraph.h:106
gtsam::kEmpty
static const VectorValues kEmpty
Definition: HybridGaussianFactorGraph.cpp:78
gtsam::DecisionTreeFactor::print
void print(const std::string &s="DecisionTreeFactor:\n", const KeyFormatter &formatter=DefaultKeyFormatter) const override
print
Definition: DecisionTreeFactor.cpp:118
gtsam::EliminateHybrid
std::pair< HybridConditional::shared_ptr, std::shared_ptr< Factor > > EliminateHybrid(const HybridGaussianFactorGraph &factors, const Ordering &keys)
Main elimination function for HybridGaussianFactorGraph.
Definition: HybridGaussianFactorGraph.cpp:536
gtsam::printFactor
static void printFactor(const std::shared_ptr< Factor > &factor, const DiscreteValues &assignment, const KeyFormatter &keyFormatter)
Definition: HybridGaussianFactorGraph.cpp:98
gtsam::DecisionTree::apply
DecisionTree apply(const Unary &op) const
Definition: DecisionTree-inl.h:1003
DiscreteEliminationTree.h
TableFactor.h
gtsam::VariableIndex
Definition: VariableIndex.h:41
g
void g(const string &key, int i)
Definition: testBTree.cpp:41
DiscreteKey.h
specialized key for discrete variables
str
Definition: pytypes.h:1560
gtsam::Result::operator==
bool operator==(const Result &other) const
Definition: HybridGaussianFactorGraph.cpp:71
JacobianFactor.h
tree::f
Point2(* f)(const Point3 &, OptionalJacobian< 2, 3 >)
Definition: testExpression.cpp:218
gtsam::FactorGraph< Factor >::size
size_t size() const
Definition: FactorGraph.h:297
gtsam::DecisionTree
a decision tree is a function from assignments to values.
Definition: DecisionTree.h:62
HybridEliminationTree.h
different_sigmas::gc
const auto gc
Definition: testHybridBayesNet.cpp:233
HybridGaussianFactorGraph.h
Linearized Hybrid factor graph that uses type erasure.
gtsam
traits
Definition: SFMdata.h:40
gtsam::HybridGaussianFactorGraph::choose
GaussianFactorGraph choose(const DiscreteValues &assignment) const
Get the GaussianFactorGraph at a given discrete assignment. Note this corresponds to the Gaussian pos...
Definition: HybridGaussianFactorGraph.cpp:663
DiscreteValues.h
gtsam::DiscreteFactor::shared_ptr
std::shared_ptr< DiscreteFactor > shared_ptr
shared_ptr to this class
Definition: DiscreteFactor.h:45
gtsam::HybridGaussianFactorGraph::printErrors
void printErrors(const HybridValues &values, const std::string &str="HybridGaussianFactorGraph: ", const KeyFormatter &keyFormatter=DefaultKeyFormatter, const std::function< bool(const Factor *, double, size_t)> &printCondition=[](const Factor *, double, size_t) { return true;}) const
Print the errors of each factor in the hybrid factor graph.
Definition: HybridGaussianFactorGraph.cpp:153
gtsam::HybridGaussianFactorGraph::probPrime
double probPrime(const HybridValues &values) const
Compute the unnormalized posterior probability for a continuous vector values given a specific assign...
Definition: HybridGaussianFactorGraph.cpp:645
gtsam::DiscreteValues
Definition: DiscreteValues.h:34
gtsam::FactorGraph::push_back
IsDerived< DERIVEDFACTOR > push_back(std::shared_ptr< DERIVEDFACTOR > factor)
Add a factor directly using a shared_ptr.
Definition: FactorGraph.h:147
gtsam::createHybridGaussianFactor
static std::shared_ptr< Factor > createHybridGaussianFactor(const ResultTree &eliminationResults, const DiscreteKeys &discreteSeparator)
Definition: HybridGaussianFactorGraph.cpp:414
gtsam::demangle
std::string demangle(const char *name)
Pretty print Value type name.
Definition: types.cpp:37
gtsam::Result::factor
GaussianFactor::shared_ptr factor
Definition: HybridGaussianFactorGraph.cpp:68
p
float * p
Definition: Tutorial_Map_using.cpp:9
gtsam::GaussianConditional::shared_ptr
std::shared_ptr< This > shared_ptr
shared_ptr to this class
Definition: GaussianConditional.h:46
gtsam::HybridGaussianFactorGraph::errorTree
AlgebraicDecisionTree< Key > errorTree(const VectorValues &continuousValues) const
Compute error for each discrete assignment, and return as a tree.
Definition: HybridGaussianFactorGraph.cpp:623
gtsam::HybridGaussianFactorGraph::print
void print(const std::string &s="HybridGaussianFactorGraph", const KeyFormatter &keyFormatter=DefaultKeyFormatter) const override
Print out graph to std::cout, with optional key formatter.
Definition: HybridGaussianFactorGraph.cpp:133
product
void product(const MatrixType &m)
Definition: product.h:20
gtsam::HybridOrdering
const Ordering HybridOrdering(const HybridGaussianFactorGraph &graph)
Return a Colamd constrained ordering where the discrete keys are eliminated after the continuous keys...
Definition: HybridGaussianFactorGraph.cpp:90
gtsam::AlgebraicDecisionTree::min
double min() const
Find the minimum values amongst all leaves.
Definition: AlgebraicDecisionTree.h:227
different_sigmas::hgc
const auto hgc
Definition: testHybridBayesNet.cpp:238
gtsam::HybridGaussianFactorGraph::discretePosterior
AlgebraicDecisionTree< Key > discretePosterior(const VectorValues &continuousValues) const
Computer posterior P(M|X=x) when all continuous values X are given. This is efficient as this simply ...
Definition: HybridGaussianFactorGraph.cpp:652
graph
NonlinearFactorGraph graph
Definition: doc/Code/OdometryExample.cpp:2
gtsam::Ordering
Definition: inference/Ordering.h:33
i
int i
Definition: BiCGSTAB_step_by_step.cpp:9
pybind_wrapper_test_script.other
other
Definition: pybind_wrapper_test_script.py:42
HybridConditional.h
gtsam::Result::negLogK
double negLogK
Definition: HybridGaussianFactorGraph.cpp:67


gtsam
Author(s):
autogenerated on Wed Mar 19 2025 03:01:48