g2o
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | Protected Attributes | List of all members
g2o::SparseOptimizerOnline Class Reference

#include <graph_optimizer_sparse_online.h>

Inheritance diagram for g2o::SparseOptimizerOnline:
Inheritance graph
[legend]
Collaboration diagram for g2o::SparseOptimizerOnline:
Collaboration graph
[legend]

Public Member Functions

 SparseOptimizerOnline (bool pcg=false)
 
virtual ~SparseOptimizerOnline ()
 
int optimize (int iterations, bool online=false)
 
virtual bool updateInitialization (HyperGraph::VertexSet &vset, HyperGraph::EdgeSet &eset)
 
void update (double *update)
 
virtual bool initSolver (int dimension, int batchEveryN)
 
virtual void gnuplotVisualization ()
 
- Public Member Functions inherited from g2o::SparseOptimizer
 SparseOptimizer ()
 
virtual ~SparseOptimizer ()
 
virtual bool initializeOptimization (HyperGraph::EdgeSet &eset)
 
virtual bool initializeOptimization (HyperGraph::VertexSet &vset, int level=0)
 
virtual bool initializeOptimization (int level=0)
 
virtual void computeInitialGuess ()
 
virtual void computeInitialGuess (EstimatePropagatorCost &propagator)
 
virtual void setToOrigin ()
 
bool computeMarginals (SparseBlockMatrix< MatrixX > &spinv, const std::vector< std::pair< int, int > > &blockIndices)
 
bool computeMarginals (SparseBlockMatrix< MatrixX > &spinv, const Vertex *vertex)
 
bool computeMarginals (SparseBlockMatrix< MatrixX > &spinv, const VertexContainer &vertices)
 
virtual VertexfindGauge ()
 finds a gauge in the graph to remove the undefined dof.
 
bool gaugeFreedom ()
 
double activeChi2 () const
 
double activeRobustChi2 () const
 
bool verbose () const
 verbose information during optimization
 
void setVerbose (bool verbose)
 
void setForceStopFlag (bool *flag)
 
bool * forceStopFlag () const
 
bool terminate ()
 if external stop flag is given, return its state. False otherwise
 
const VertexContainerindexMapping () const
 the index mapping of the vertices
 
const VertexContaineractiveVertices () const
 the vertices active in the current optimization
 
const EdgeContaineractiveEdges () const
 the edges active in the current optimization
 
virtual bool removeVertex (HyperGraph::Vertex *v, bool detach=false)
 
VertexContainer::const_iterator findActiveVertex (const OptimizableGraph::Vertex *v) const
 
EdgeContainer::const_iterator findActiveEdge (const OptimizableGraph::Edge *e) const
 
const OptimizationAlgorithmalgorithm () const
 the solver used by the optimizer
 
OptimizationAlgorithmsolver ()
 
void setAlgorithm (OptimizationAlgorithm *algorithm)
 
void push (SparseOptimizer::VertexContainer &vlist)
 push the estimate of a subset of the variables onto a stack
 
void push (HyperGraph::VertexSet &vlist)
 push the estimate of a subset of the variables onto a stack
 
void push ()
 push all the active vertices onto a stack
 
void pop (SparseOptimizer::VertexContainer &vlist)
 pop (restore) the estimate a subset of the variables from the stack
 
void pop (HyperGraph::VertexSet &vlist)
 pop (restore) the estimate a subset of the variables from the stack
 
void pop ()
 pop (restore) the estimate of the active vertices from the stack
 
void discardTop (SparseOptimizer::VertexContainer &vlist)
 
void discardTop ()
 same as above, but for the active vertices
 
virtual void clear ()
 
void computeActiveErrors ()
 
 G2O_ATTRIBUTE_DEPRECATED (void linearizeSystem())
 
void update (const double *update)
 
const BatchStatisticsContainerbatchStatistics () const
 
BatchStatisticsContainerbatchStatistics ()
 
void setComputeBatchStatistics (bool computeBatchStatistics)
 
bool computeBatchStatistics () const
 
bool addComputeErrorAction (HyperGraphAction *action)
 add an action to be executed before the error vectors are computed
 
bool removeComputeErrorAction (HyperGraphAction *action)
 
virtual void discardTop ()
 
virtual void discardTop (HyperGraph::VertexSet &vset)
 
- Public Member Functions inherited from g2o::OptimizableGraph
Vertexvertex (int id)
 returns the vertex number id appropriately casted
 
const Vertexvertex (int id) const
 returns the vertex number id appropriately casted
 
 OptimizableGraph ()
 empty constructor
 
virtual ~OptimizableGraph ()
 
virtual bool addVertex (HyperGraph::Vertex *v, Data *userData)
 
virtual bool addVertex (HyperGraph::Vertex *v)
 
bool addVertex (OptimizableGraph::Vertex *v, Data *userData)
 
bool addVertex (OptimizableGraph::Vertex *v)
 
virtual bool addEdge (HyperGraph::Edge *e)
 
bool addEdge (OptimizableGraph::Edge *e)
 
virtual bool setEdgeVertex (HyperGraph::Edge *e, int pos, HyperGraph::Vertex *v)
 
double chi2 () const
 returns the chi2 of the current configuration
 
int maxDimension () const
 return the maximum dimension of all vertices in the graph
 
void recomputeJacobianWorkspaceSize ()
 
std::set< int > dimensions () const
 
virtual void preIteration (int)
 
virtual void postIteration (int)
 
bool addPreIterationAction (HyperGraphAction *action)
 add an action to be executed before each iteration
 
bool addPostIterationAction (HyperGraphAction *action)
 add an action to be executed after each iteration
 
bool removePreIterationAction (HyperGraphAction *action)
 remove an action that should no longer be executed before each iteration
 
bool removePostIterationAction (HyperGraphAction *action)
 remove an action that should no longer be executed after each iteration
 
virtual bool load (std::istream &is)
 
bool load (const char *filename)
 
virtual bool save (std::ostream &os, int level=0) const
 save the graph to a stream. Again uses the Factory system.
 
bool save (const char *filename, int level=0) const
 function provided for convenience, see save() above
 
bool saveSubset (std::ostream &os, HyperGraph::VertexSet &vset, int level=0)
 save a subgraph to a stream. Again uses the Factory system.
 
bool saveSubset (std::ostream &os, HyperGraph::EdgeSet &eset)
 save a subgraph to a stream. Again uses the Factory system.
 
virtual void setFixed (HyperGraph::VertexSet &vset, bool fixed)
 fixes/releases a set of vertices
 
void setRenamedTypesFromString (const std::string &types)
 
bool isSolverSuitable (const OptimizationAlgorithmProperty &solverProperty, const std::set< int > &vertDims=std::set< int >()) const
 
virtual void clearParameters ()
 
bool addParameter (Parameter *p)
 
Parameterparameter (int id)
 
bool verifyInformationMatrices (bool verbose=false) const
 
bool saveVertex (std::ostream &os, Vertex *v) const
 
bool saveParameter (std::ostream &os, Parameter *v) const
 
bool saveEdge (std::ostream &os, Edge *e) const
 
bool saveUserData (std::ostream &os, HyperGraph::Data *v) const
 
JacobianWorkspacejacobianWorkspace ()
 the workspace for storing the Jacobians of the graph
 
const JacobianWorkspacejacobianWorkspace () const
 
ParameterContainerparameters ()
 
const ParameterContainerparameters () const
 
void forEachVertex (std::function< void(OptimizableGraph::Vertex *)> fn)
 apply a unary function to all vertices
 
void forEachVertex (HyperGraph::VertexSet &vset, std::function< void(OptimizableGraph::Vertex *)> fn)
 apply a unary function to the vertices in vset
 
- Public Member Functions inherited from g2o::HyperGraph
 HyperGraph ()
 constructs an empty hyper graph
 
virtual ~HyperGraph ()
 destroys the hyper-graph and all the vertices of the graph
 
Vertexvertex (int id)
 
const Vertexvertex (int id) const
 
virtual bool removeEdge (Edge *e)
 
const VertexIDMapvertices () const
 
VertexIDMapvertices ()
 
const EdgeSetedges () const
 
EdgeSetedges ()
 
virtual bool mergeVertices (Vertex *vBig, Vertex *vSmall, bool erase)
 
virtual bool detachVertex (Vertex *v)
 
virtual bool changeId (Vertex *v, int newId)
 

Public Attributes

int slamDimension
 
HyperGraph::EdgeSetnewEdges
 
bool batchStep
 
bool vizWithGnuplot
 
- Public Attributes inherited from g2o::OptimizableGraph
class G2O_CORE_API Vertex
 
class G2O_CORE_API Edge
 
- Public Attributes inherited from g2o::HyperGraph
class G2O_CORE_API Data
 
class G2O_CORE_API DataContainer
 
class G2O_CORE_API Vertex
 
class G2O_CORE_API Edge
 

Protected Attributes

FILE * _gnuplot
 
bool _usePcg
 
Solver_underlyingSolver
 
- Protected Attributes inherited from g2o::SparseOptimizer
bool * _forceStopFlag
 
bool _verbose
 
VertexContainer _ivMap
 
VertexContainer _activeVertices
 sorted according to VertexIDCompare
 
EdgeContainer _activeEdges
 sorted according to EdgeIDCompare
 
OptimizationAlgorithm_algorithm
 
BatchStatisticsContainer _batchStatistics
 
bool _computeBatchStatistics
 
- Protected Attributes inherited from g2o::OptimizableGraph
std::map< std::string, std::string > _renamedTypesLookup
 
long long _nextEdgeId
 
std::vector< HyperGraphActionSet_graphActions
 
ParameterContainer _parameters
 
JacobianWorkspace _jacobianWorkspace
 
- Protected Attributes inherited from g2o::HyperGraph
VertexIDMap _vertices
 
EdgeSet _edges
 

Additional Inherited Members

- Public Types inherited from g2o::SparseOptimizer
enum  { AT_COMPUTEACTIVERROR = OptimizableGraph::AT_NUM_ELEMENTS , AT_NUM_ELEMENTS }
 
- Public Types inherited from g2o::OptimizableGraph
enum  ActionType { AT_PREITERATION , AT_POSTITERATION , AT_NUM_ELEMENTS }
 
typedef std::set< HyperGraphAction * > HyperGraphActionSet
 
typedef std::vector< OptimizableGraph::Vertex * > VertexContainer
 vector container for vertices
 
typedef std::vector< OptimizableGraph::Edge * > EdgeContainer
 vector container for edges
 
- Public Types inherited from g2o::HyperGraph
typedef std::bitset< HyperGraph::HGET_NUM_ELEMS > GraphElemBitset
 
typedef std::set< Edge * > EdgeSet
 
typedef std::set< Vertex * > VertexSet
 
typedef std::unordered_map< int, Vertex * > VertexIDMap
 
typedef std::vector< Vertex * > VertexContainer
 
- Static Public Member Functions inherited from g2o::OptimizableGraph
static bool initMultiThreading ()
 
- Static Public Attributes inherited from g2o::HyperGraph
static const int UnassignedId = -1
 
static const int InvalidId = -2
 
- Protected Member Functions inherited from g2o::SparseOptimizer
void sortVectorContainers ()
 
bool buildIndexMapping (SparseOptimizer::VertexContainer &vlist)
 
void clearIndexMapping ()
 
- Protected Member Functions inherited from g2o::OptimizableGraph
void performActions (int iter, HyperGraphActionSet &actions)
 

Detailed Description

Definition at line 37 of file graph_optimizer_sparse_online.h.

Constructor & Destructor Documentation

◆ SparseOptimizerOnline()

g2o::SparseOptimizerOnline::SparseOptimizerOnline ( bool  pcg = false)
explicit

◆ ~SparseOptimizerOnline()

g2o::SparseOptimizerOnline::~SparseOptimizerOnline ( )
virtual

Definition at line 73 of file graph_optimizer_sparse_online.cpp.

73 {
74 if (_gnuplot) {
75#ifdef WINDOWS
76 _pclose(_gnuplot);
77#else
78 pclose(_gnuplot);
79#endif
80 }
81}

References _gnuplot.

Member Function Documentation

◆ gnuplotVisualization()

void g2o::SparseOptimizerOnline::gnuplotVisualization ( )
virtual

Definition at line 226 of file graph_optimizer_sparse_online.cpp.

226 {
227 if (slamDimension == 3) {
228 if (!_gnuplot) {
229#ifdef WINDOWS
230 _gnuplot = _popen("gnuplot -persistent", "w");
231#else
232 _gnuplot = popen("gnuplot -persistent", "w");
233#endif
234 if (_gnuplot == 0) return;
235 fprintf(_gnuplot, "set terminal X11 noraise\n");
236 fprintf(_gnuplot, "set size ratio -1\n");
237 }
238 fprintf(_gnuplot, "plot \"-\" w l\n");
239 for (EdgeSet::iterator it = edges().begin(); it != edges().end(); ++it) {
240 OnlineEdgeSE2* e = static_cast<OnlineEdgeSE2*>(*it);
241 OnlineVertexSE2* v1 = static_cast<OnlineVertexSE2*>(e->vertices()[0]);
242 OnlineVertexSE2* v2 = static_cast<OnlineVertexSE2*>(e->vertices()[1]);
243 fprintf(_gnuplot, "%f %f\n", v1->updatedEstimate.translation().x(),
244 v1->updatedEstimate.translation().y());
245 fprintf(_gnuplot, "%f %f\n\n", v2->updatedEstimate.translation().x(),
246 v2->updatedEstimate.translation().y());
247 }
248 fprintf(_gnuplot, "e\n");
249 }
250 if (slamDimension == 6) {
251 if (!_gnuplot) {
252#ifdef WINDOWS
253 _gnuplot = _popen("gnuplot -persistent", "w");
254#else
255 _gnuplot = popen("gnuplot -persistent", "w");
256#endif
257 if (_gnuplot == 0) return;
258 fprintf(_gnuplot, "set terminal X11 noraise\n");
259 }
260 fprintf(_gnuplot, "splot \"-\" w l\n");
261 for (EdgeSet::iterator it = edges().begin(); it != edges().end(); ++it) {
262 OnlineEdgeSE3* e = (OnlineEdgeSE3*)*it;
263 OnlineVertexSE3* v1 = static_cast<OnlineVertexSE3*>(e->vertices()[0]);
264 OnlineVertexSE3* v2 = static_cast<OnlineVertexSE3*>(e->vertices()[1]);
265 fprintf(_gnuplot, "%f %f %f\n", v1->updatedEstimate.translation().x(),
266 v1->updatedEstimate.translation().y(),
267 v1->updatedEstimate.translation().z());
268 fprintf(_gnuplot, "%f %f %f \n\n\n",
269 v2->updatedEstimate.translation().x(),
270 v2->updatedEstimate.translation().y(),
271 v2->updatedEstimate.translation().z());
272 }
273 fprintf(_gnuplot, "e\n");
274 }
275}
const EdgeSet & edges() const

References _gnuplot, g2o::HyperGraph::edges(), slamDimension, g2o::SE2::translation(), g2o::OnlineVertexSE2::updatedEstimate, g2o::OnlineVertexSE3::updatedEstimate, and g2o::HyperGraph::Edge::vertices().

Referenced by g2o::SparseOptimizerIncremental::optimize(), and optimize().

◆ initSolver()

bool g2o::SparseOptimizerOnline::initSolver ( int  dimension,
int  batchEveryN 
)
virtual

Reimplemented in g2o::SparseOptimizerIncremental.

Definition at line 188 of file graph_optimizer_sparse_online.cpp.

188 {
189 slamDimension = dimension;
190 OptimizationAlgorithmFactory* solverFactory =
192 OptimizationAlgorithmProperty solverProperty;
193 if (_usePcg) {
194 std::unique_ptr<Solver> s;
195 if (dimension == 3) {
196 s = AllocatePCGSolver<3, 2>();
197 } else {
198 s = AllocatePCGSolver<6, 3>();
199 }
200
201 OptimizationAlgorithmGaussNewton* gaussNewton =
202 new OptimizationAlgorithmGaussNewton(std::move(s));
203 setAlgorithm(gaussNewton);
204 } else {
205 if (dimension == 3) {
207 solverFactory->construct("gn_fix3_2_cholmod", solverProperty));
208 } else {
210 solverFactory->construct("gn_fix6_3_cholmod", solverProperty));
211 }
212 }
213
214 OptimizationAlgorithmGaussNewton* gaussNewton =
215 dynamic_cast<OptimizationAlgorithmGaussNewton*>(solver());
216 _underlyingSolver = &gaussNewton->solver();
217
218 if (!solver()) {
219 cerr << "Error allocating solver. Allocating CHOLMOD solver failed!"
220 << endl;
221 return false;
222 }
223 return true;
224}
static OptimizationAlgorithmFactory * instance()
return the instance
void setAlgorithm(OptimizationAlgorithm *algorithm)
OptimizationAlgorithm * solver()

References _underlyingSolver, _usePcg, g2o::OptimizationAlgorithmFactory::construct(), g2o::OptimizationAlgorithmFactory::instance(), g2o::SparseOptimizer::setAlgorithm(), slamDimension, g2o::OptimizationAlgorithmWithHessian::solver(), and g2o::SparseOptimizer::solver().

Referenced by g2o::G2oSlamInterface::addNode().

◆ optimize()

int g2o::SparseOptimizerOnline::optimize ( int  iterations,
bool  online = false 
)
virtual

starts one optimization run given the current configuration of the graph, and the current settings stored in the class instance. It can be called only after initializeOptimization

Reimplemented from g2o::SparseOptimizer.

Definition at line 83 of file graph_optimizer_sparse_online.cpp.

83 {
84 // return SparseOptimizer::optimize(iterations, online);
85
86 (void)iterations; // we only do one iteration anyhow
88
89 bool ok = true;
90
91 solver->init(online);
92 if (!online) {
94 if (!ok) {
95 cerr << __PRETTY_FUNCTION__ << ": Failure while building CCS structure"
96 << endl;
97 return 0;
98 }
99 }
100
101 if (_usePcg) batchStep = true;
102
103 if (!online || batchStep) {
104 // cerr << "BATCH" << endl;
105 //_underlyingSolver->buildStructure();
106 // copy over the updated estimate as new linearization point
107 if (slamDimension == 3) {
108 for (size_t i = 0; i < indexMapping().size(); ++i) {
109 OnlineVertexSE2* v = static_cast<OnlineVertexSE2*>(indexMapping()[i]);
110 v->setEstimate(v->updatedEstimate);
111 }
112 } else if (slamDimension == 6) {
113 for (size_t i = 0; i < indexMapping().size(); ++i) {
114 OnlineVertexSE3* v = static_cast<OnlineVertexSE3*>(indexMapping()[i]);
115 v->setEstimate(v->updatedEstimate);
116 }
117 }
118
120 // SparseOptimizer::linearizeSystem();
122 } else {
123 // cerr << "UPDATE" << endl;
124 // compute the active errors for the required edges
125 for (HyperGraph::EdgeSet::iterator it = newEdges->begin();
126 it != newEdges->end(); ++it) {
127 OptimizableGraph::Edge* e = static_cast<OptimizableGraph::Edge*>(*it);
128 e->computeError();
129 }
130 // linearize the constraints and update the Hessian
131 for (HyperGraph::EdgeSet::iterator it = newEdges->begin();
132 it != newEdges->end(); ++it) {
133 OptimizableGraph::Edge* e = static_cast<OptimizableGraph::Edge*>(*it);
134 e->linearizeOplus(jacobianWorkspace());
135 e->constructQuadraticForm();
136 }
137 // update the b vector
138 for (int i = 0; i < static_cast<int>(indexMapping().size()); ++i) {
139 OptimizableGraph::Vertex* v = indexMapping()[i];
140 int iBase = v->colInHessian();
141 v->copyB(_underlyingSolver->b() + iBase);
142 }
143 }
144 ok = _underlyingSolver->solve();
146
147 if (verbose()) {
149 cerr << "nodes = " << vertices().size()
150 << "\t edges= " << _activeEdges.size()
151 << "\t chi2= " << FIXED(activeChi2()) << endl;
152 }
153
155
156 if (!ok) return 0;
157 return 1;
158}
const VertexIDMap & vertices() const
virtual bool init(bool online=false)=0
double * b()
return b, the right hand side of the system
Definition solver.h:101
virtual bool buildStructure(bool zeroBlocks=false)=0
double * x()
return x, the solution vector
Definition solver.h:98
virtual bool solve()=0
virtual bool buildSystem()=0
EdgeContainer _activeEdges
sorted according to EdgeIDCompare
OptimizationAlgorithm * _algorithm
bool verbose() const
verbose information during optimization
const VertexContainer & indexMapping() const
the index mapping of the vertices
#define __PRETTY_FUNCTION__
Definition macros.h:90
class G2O_CORE_API OptimizationAlgorithm
JacobianWorkspace & jacobianWorkspace()
the workspace for storing the Jacobians of the graph

References __PRETTY_FUNCTION__, g2o::SparseOptimizer::_activeEdges, g2o::SparseOptimizer::_algorithm, _underlyingSolver, _usePcg, g2o::SparseOptimizer::activeChi2(), g2o::Solver::b(), batchStep, g2o::Solver::buildStructure(), g2o::Solver::buildSystem(), g2o::OptimizableGraph::Vertex::colInHessian(), g2o::SparseOptimizer::computeActiveErrors(), g2o::OptimizableGraph::Edge::computeError(), g2o::OptimizableGraph::Edge::constructQuadraticForm(), g2o::OptimizableGraph::Vertex::copyB(), gnuplotVisualization(), g2o::SparseOptimizer::indexMapping(), g2o::OptimizationAlgorithm::init(), g2o::OptimizableGraph::jacobianWorkspace(), g2o::OptimizableGraph::Edge::linearizeOplus(), newEdges, g2o::BaseVertex< D, T >::setEstimate(), slamDimension, g2o::Solver::solve(), g2o::SparseOptimizer::solver(), update(), g2o::OnlineVertexSE2::updatedEstimate, g2o::OnlineVertexSE3::updatedEstimate, g2o::SparseOptimizer::verbose(), g2o::HyperGraph::vertices(), vizWithGnuplot, and g2o::Solver::x().

Referenced by g2o::G2oSlamInterface::solve().

◆ update()

void g2o::SparseOptimizerOnline::update ( double *  update)

Definition at line 160 of file graph_optimizer_sparse_online.cpp.

160 {
161 if (slamDimension == 3) {
162 for (size_t i = 0; i < _ivMap.size(); ++i) {
163 OnlineVertexSE2* v = static_cast<OnlineVertexSE2*>(_ivMap[i]);
164 v->oplusUpdatedEstimate(update);
165 update += 3;
166 }
167 } else if (slamDimension == 6) {
168 for (size_t i = 0; i < _ivMap.size(); ++i) {
169 OnlineVertexSE3* v = static_cast<OnlineVertexSE3*>(_ivMap[i]);
170 v->oplusUpdatedEstimate(update);
171 update += 6;
172 }
173 }
174}

References g2o::SparseOptimizer::_ivMap, g2o::OnlineVertexSE2::oplusUpdatedEstimate(), g2o::OnlineVertexSE3::oplusUpdatedEstimate(), slamDimension, and update().

Referenced by g2o::SparseOptimizerIncremental::optimize(), optimize(), and update().

◆ updateInitialization()

bool g2o::SparseOptimizerOnline::updateInitialization ( HyperGraph::VertexSet vset,
HyperGraph::EdgeSet eset 
)
virtual

HACK updating the internal structures for online processing

Reimplemented from g2o::SparseOptimizer.

Reimplemented in g2o::SparseOptimizerIncremental.

Definition at line 176 of file graph_optimizer_sparse_online.cpp.

177 {
178 newEdges = &eset;
179 bool result = SparseOptimizer::updateInitialization(vset, eset);
180 for (HyperGraph::VertexSet::iterator it = vset.begin(); it != vset.end();
181 ++it) {
182 OptimizableGraph::Vertex* v = static_cast<OptimizableGraph::Vertex*>(*it);
183 v->clearQuadraticForm(); // be sure that b is zero for this vertex
184 }
185 return result;
186}
virtual bool updateInitialization(HyperGraph::VertexSet &vset, HyperGraph::EdgeSet &eset)

References g2o::OptimizableGraph::Vertex::clearQuadraticForm(), newEdges, and g2o::SparseOptimizer::updateInitialization().

Referenced by g2o::G2oSlamInterface::solve(), and g2o::SparseOptimizerIncremental::updateInitialization().

Member Data Documentation

◆ _gnuplot

FILE* g2o::SparseOptimizerOnline::_gnuplot
protected

Definition at line 62 of file graph_optimizer_sparse_online.h.

Referenced by gnuplotVisualization(), and ~SparseOptimizerOnline().

◆ _underlyingSolver

Solver* g2o::SparseOptimizerOnline::_underlyingSolver
protected

◆ _usePcg

bool g2o::SparseOptimizerOnline::_usePcg
protected

Definition at line 63 of file graph_optimizer_sparse_online.h.

Referenced by initSolver(), and optimize().

◆ batchStep

bool g2o::SparseOptimizerOnline::batchStep

◆ newEdges

HyperGraph::EdgeSet* g2o::SparseOptimizerOnline::newEdges

Definition at line 54 of file graph_optimizer_sparse_online.h.

Referenced by optimize(), and updateInitialization().

◆ slamDimension

int g2o::SparseOptimizerOnline::slamDimension

◆ vizWithGnuplot

bool g2o::SparseOptimizerOnline::vizWithGnuplot

The documentation for this class was generated from the following files: