40FastList<Key> createKeyList(
const Vector& I) {
42 for (
int i = 0; i < I.size(); i++)
48FastList<Key> createKeyList(std::string s,
const Vector& I) {
51 for (
int i = 0; i < I.size(); i++)
52 set.push_back(Symbol(c, I[i]));
57KeyVector createKeyVector(
const Vector& I) {
59 for (
int i = 0; i < I.size(); i++)
65KeyVector createKeyVector(std::string s,
const Vector& I) {
68 for (
int i = 0; i < I.size(); i++)
69 set.push_back(Symbol(c, I[i]));
74KeySet createKeySet(
const Vector& I) {
76 for (
int i = 0; i < I.size(); i++)
82KeySet createKeySet(std::string s,
const Vector& I) {
85 for (
int i = 0; i < I.size(); i++)
86 set.insert(
symbol(c, I[i]));
96 Matrix result(points.
size() + points2.size(), 2);
99 for (
const auto& key_value : points) {
100 result.row(j++) = key_value.value;
102 for (
const auto& key_value : points2) {
103 if (key_value.value.rows() == 2) {
104 result.row(j++) = key_value.value;
116 Matrix result(points.
size() + points2.size(), 3);
119 for (
const auto& key_value : points) {
120 result.row(j++) = key_value.value;
122 for (
const auto& key_value : points2) {
123 if (key_value.value.rows() == 3) {
124 result.row(j++) = key_value.value;
138 Matrix result(poses.
size(), 3);
140 for(
const auto& key_value: poses)
141 result.row(j++) << key_value.value.x(), key_value.value.y(), key_value.value.theta();
153 Matrix result(poses.
size(), 12);
155 for(
const auto& key_value: poses) {
156 result.row(j).segment(0, 3) << key_value.value.rotation().matrix().row(0);
157 result.row(j).segment(3, 3) << key_value.value.rotation().matrix().row(1);
158 result.row(j).segment(6, 3) << key_value.value.rotation().matrix().row(2);
159 result.row(j).tail(3) = key_value.value.translation();
167 noiseModel::Isotropic::shared_ptr model =
174 for (
const auto& key_value : values.
filter<gtsam::Vector>()) {
175 if (key_value.value.rows() == 2) {
176 values.
update<gtsam::Vector>(key_value.key,
186 Vector3(sigmaT, sigmaT, sigmaR));
188 for(
const auto& key_value: values.
filter<
Pose2>()) {
195 noiseModel::Isotropic::shared_ptr model =
202 for (
const auto& key_value : values.
filter<gtsam::Vector>()) {
203 if (key_value.value.rows() == 3) {
204 values.
update<gtsam::Vector>(key_value.key,
220 const Vector& J,
const Matrix& Z,
double depth) {
222 throw std::invalid_argument(
"insertBackProjections: Z must be 2*J");
223 if (Z.cols() != J.size())
224 throw std::invalid_argument(
225 "insertBackProjections: J and Z must have same number of entries");
226 for (
int k = 0; k < Z.cols(); k++) {
227 Point2 p(Z(0, k), Z(1, k));
246 const Cal3_S2::shared_ptr K,
const Pose3& body_P_sensor =
Pose3()) {
248 throw std::invalid_argument(
"addMeasurements: Z must be 2*K");
249 if (Z.cols() != J.size())
250 throw std::invalid_argument(
251 "addMeasurements: J and Z must have same number of entries");
252 for (
int k = 0; k < Z.cols(); k++) {
255 Point2(Z(0, k), Z(1, k)), model, i,
Key(J(k)), K, body_P_sensor));
264 for(
const NonlinearFactor::shared_ptr& f: graph)
270 for(
const NonlinearFactor::shared_ptr& f: graph) {
271 boost::shared_ptr<const GenericProjectionFactor<Pose3, Point3> > p =
272 boost::dynamic_pointer_cast<const GenericProjectionFactor<Pose3, Point3> >(
275 errors.col(k++) = p->unwhitenedError(values);
296 world.
insert(key, base.compose(pose));
297 }
catch (
const std::exception& e1) {
302 }
catch (
const std::exception& e2) {
305 std::cerr <<
"Values[key] is neither Pose2 nor Point2, so skip" << std::endl;
The most common 5DOF 3D->2D calibration.
Base class for all pinhole cameras.
sampling from a NoiseModel
A non-templated config holding any types of Manifold-group elements.
Factor Graph consisting of non-linear factors.
Non-linear factor base classes.
Reprojection of a LANDMARK to a 2D point.
Global functions in a separate testing namespace.
Definition: chartTesting.h:28
FastVector< Key > KeyVector
Define collection type once and for all - also used in wrappers.
Definition: Key.h:86
Key symbol(unsigned char c, std::uint64_t j)
Create a symbol key from a character and index, i.e.
Definition: Symbol.h:135
Vector2 Point2
As of GTSAM 4, in order to make GTSAM more lean, it is now possible to just typedef Point2 to Vector2...
Definition: Point2.h:27
Vector3 Point3
As of GTSAM 4, in order to make GTSAM more lean, it is now possible to just typedef Point3 to Vector3...
Definition: Point3.h:35
noiseModel::Base::shared_ptr SharedNoiseModel
Note, deliberately not in noiseModel namespace.
Definition: NoiseModel.h:736
gtsam::enable_if_t< needs_eigen_aligned_allocator< T >::value, boost::shared_ptr< T > > make_shared(Args &&... args)
Add our own make_shared as a layer of wrapping on boost::make_shared This solves the problem with the...
Definition: make_shared.h:57
std::uint64_t Key
Integer nonlinear key type.
Definition: types.h:69
Definition: PinholeCamera.h:33
Point3 backproject(const Point2 &p, double depth, OptionalJacobian< 3, 6 > Dresult_dpose=boost::none, OptionalJacobian< 3, 2 > Dresult_dp=boost::none, OptionalJacobian< 3, 1 > Dresult_ddepth=boost::none, OptionalJacobian< 3, DimK > Dresult_dcal=boost::none) const
backproject a 2-dimensional point to a 3-dimensional point at given depth
Definition: PinholePose.h:132
GTSAM_EXPORT Point2 transformFrom(const Point2 &point, OptionalJacobian< 2, 3 > Dpose=boost::none, OptionalJacobian< 2, 2 > Dpoint=boost::none) const
Return point coordinates in global frame.
Definition: Pose2.cpp:218
static shared_ptr Sigmas(const Vector &sigmas, bool smart=true)
A diagonal noise model created by specifying a Vector of sigmas, i.e.
Definition: NoiseModel.cpp:270
static shared_ptr Sigma(size_t dim, double sigma, bool smart=true)
An isotropic noise model created by specifying a standard devation sigma.
Definition: NoiseModel.cpp:567
Sampling structure that keeps internal random number generators for diagonal distributions specified ...
Definition: Sampler.h:31
Vector sample() const
sample from distribution
Definition: Sampler.cpp:51
A non-linear factor graph is a graph of non-Gaussian, i.e.
Definition: NonlinearFactorGraph.h:78
A filtered view of a const Values, returned from Values::filter.
Definition: Values-inl.h:169
size_t size() const
Returns the number of values in this view.
Definition: Values-inl.h:189
A non-templated config holding any types of Manifold-group elements.
Definition: Values.h:63
const ValueType at(Key j) const
Retrieve a variable by key j.
Definition: Values-inl.h:346
KeyVector keys() const
Returns a set of keys in the config Note: by construction, the list is ordered.
Definition: Values.cpp:183
void insert(Key j, const Value &val)
Add a variable with the given j, throws KeyAlreadyExists<J> if j is already present.
Definition: Values.cpp:132
void update(Key j, const Value &val)
single element change of existing element
Definition: Values.cpp:153
Filtered< Value > filter(const std::function< bool(Key)> &filterFcn)
Return a filtered view of this Values class, without copying any data.
Definition: Values-inl.h:239
Definition: ProjectionFactor.h:40
Matrix extractPoint2(const Values &values)
Extract all Point2 values into a single matrix [x y].
Definition: utilities.h:91
Values localToWorld(const Values &local, const Pose2 &base, const KeyVector user_keys=KeyVector())
Convert from local to world coordinates.
Definition: utilities.h:281
Matrix extractPoint3(const Values &values)
Extract all Point3 values into a single matrix [x y z].
Definition: utilities.h:111
Matrix reprojectionErrors(const NonlinearFactorGraph &graph, const Values &values)
Calculate the errors of all projection factors in a graph.
Definition: utilities.h:260
void perturbPoint2(Values &values, double sigma, int32_t seed=42u)
Perturb all Point2 values using normally distributed noise.
Definition: utilities.h:166
void perturbPose2(Values &values, double sigmaT, double sigmaR, int32_t seed=42u)
Perturb all Pose2 values using normally distributed noise.
Definition: utilities.h:183
Values allPose2s(const Values &values)
Extract all Pose3 values.
Definition: utilities.h:131
void insertProjectionFactors(NonlinearFactorGraph &graph, Key i, const Vector &J, const Matrix &Z, const SharedNoiseModel &model, const Cal3_S2::shared_ptr K, const Pose3 &body_P_sensor=Pose3())
Insert multiple projection factors for a single pose key.
Definition: utilities.h:244
Values allPose3s(const Values &values)
Extract all Pose3 values.
Definition: utilities.h:146
void insertBackprojections(Values &values, const PinholeCamera< Cal3_S2 > &camera, const Vector &J, const Matrix &Z, double depth)
Insert a number of initial point values by backprojecting.
Definition: utilities.h:219
Matrix extractPose3(const Values &values)
Extract all Pose3 values into a single matrix [r11 r12 r13 r21 r22 r23 r31 r32 r33 x y z].
Definition: utilities.h:151
Matrix extractPose2(const Values &values)
Extract all Pose2 values into a single matrix [x y theta].
Definition: utilities.h:136
void perturbPoint3(Values &values, double sigma, int32_t seed=42u)
Perturb all Point3 values using normally distributed noise.
Definition: utilities.h:194