MRPT
2.0.1
|
A class for storing an occupancy grid map.
COccupancyGridMap2D is a class for storing a metric map representation in the form of a probabilistic occupancy grid map: value of 0 means certainly occupied, 1 means a certainly empty cell. Initially 0.5 means uncertainty.
The cells keep the log-odd representation of probabilities instead of the probabilities themselves. More details can be found at https://www.mrpt.org/Occupancy_Grids
The algorithm for updating the grid from a laser scanner can optionally take into account the progressive widening of the beams, as described in this page
Some implemented methods are:
Definition at line 53 of file COccupancyGridMap2D.h.
#include <mrpt/maps/COccupancyGridMap2D.h>
Classes | |
struct | TCriticalPointsList |
The structure used to store the set of Voronoi diagram critical points. More... | |
struct | TEntropyInfo |
Used for returning entropy related information. More... | |
class | TInsertionOptions |
With this struct options are provided to the observation insertion process. More... | |
struct | TLaserSimulUncertaintyParams |
Input params for laserScanSimulatorWithUncertainty() More... | |
struct | TLaserSimulUncertaintyResult |
Output params for laserScanSimulatorWithUncertainty() More... | |
class | TLikelihoodOptions |
With this struct options are provided to the observation likelihood computation process. More... | |
class | TLikelihoodOutput |
Some members of this struct will contain intermediate or output data after calling "computeObservationLikelihood" for some likelihood functions. More... | |
struct | TMapDefinition |
struct | TMapDefinitionBase |
struct | TUpdateCellsInfoChangeOnly |
An internal structure for storing data related to counting the new information apported by some observation. More... | |
Public Types | |
enum | TLikelihoodMethod { lmMeanInformation = 0, lmRayTracing, lmConsensus, lmCellsDifference, lmLikelihoodField_Thrun, lmLikelihoodField_II, lmConsensusOWA } |
The type for selecting a likelihood computation method. More... | |
using | cellType = OccGridCellTraits::cellType |
The type of the map cells: More... | |
using | cellTypeUnsigned = OccGridCellTraits::cellTypeUnsigned |
using | TPairLikelihoodIndex = std::pair< double, mrpt::math::TPoint2D > |
Auxiliary private class. More... | |
using | cell_t = OccGridCellTraits::cellType |
The type of cells. More... | |
using | traits_t = detail::logoddscell_traits< OccGridCellTraits::cellType > |
Public Member Functions | |
const std::vector< cellType > & | getRawMap () const |
Read-only access to the raw cell contents (cells are in log-odd units) More... | |
void | updateCell (int x, int y, float v) |
Performs the Bayesian fusion of a new observation of a cell. More... | |
void | fill (float default_value=0.5f) |
Fills all the cells with a default value. More... | |
COccupancyGridMap2D (float min_x=-20.0f, float max_x=20.0f, float min_y=-20.0f, float max_y=20.0f, float resolution=0.05f) | |
Constructor. More... | |
void | setSize (float x_min, float x_max, float y_min, float y_max, float resolution, float default_value=0.5f) |
Change the size of gridmap, erasing all its previous contents. More... | |
void | resizeGrid (float new_x_min, float new_x_max, float new_y_min, float new_y_max, float new_cells_default_value=0.5f, bool additionalMargin=true) noexcept |
Change the size of gridmap, maintaining previous contents. More... | |
double | getArea () const |
Returns the area of the gridmap, in square meters. More... | |
unsigned int | getSizeX () const |
Returns the horizontal size of grid map in cells count. More... | |
unsigned int | getSizeY () const |
Returns the vertical size of grid map in cells count. More... | |
float | getXMin () const |
Returns the "x" coordinate of left side of grid map. More... | |
float | getXMax () const |
Returns the "x" coordinate of right side of grid map. More... | |
float | getYMin () const |
Returns the "y" coordinate of top side of grid map. More... | |
float | getYMax () const |
Returns the "y" coordinate of bottom side of grid map. More... | |
float | getResolution () const |
Returns the resolution of the grid map. More... | |
int | x2idx (float x) const |
Transform a coordinate value into a cell index. More... | |
int | y2idx (float y) const |
int | x2idx (double x) const |
int | y2idx (double y) const |
float | idx2x (const size_t cx) const |
Transform a cell index into a coordinate value. More... | |
float | idx2y (const size_t cy) const |
int | x2idx (float x, float xmin) const |
Transform a coordinate value into a cell index, using a diferent "x_min" value. More... | |
int | y2idx (float y, float ymin) const |
void | setCell (int x, int y, float value) |
Change the contents [0,1] of a cell, given its index. More... | |
float | getCell (int x, int y) const |
Read the real valued [0,1] contents of a cell, given its index. More... | |
cellType * | getRow (int cy) |
Access to a "row": mainly used for drawing grid as a bitmap efficiently, do not use it normally. More... | |
const cellType * | getRow (int cy) const |
Access to a "row": mainly used for drawing grid as a bitmap efficiently, do not use it normally. More... | |
void | setPos (float x, float y, float value) |
Change the contents [0,1] of a cell, given its coordinates. More... | |
float | getPos (float x, float y) const |
Read the real valued [0,1] contents of a cell, given its coordinates. More... | |
bool | isStaticPos (float x, float y, float threshold=0.7f) const |
Returns "true" if cell is "static", i.e.if its occupancy is below a given threshold. More... | |
bool | isStaticCell (int cx, int cy, float threshold=0.7f) const |
void | setBasisCell (int x, int y, uint8_t value) |
Change a cell in the "basis" maps.Used for Voronoi calculation. More... | |
unsigned char | getBasisCell (int x, int y) const |
Reads a cell in the "basis" maps.Used for Voronoi calculation. More... | |
void | copyMapContentFrom (const COccupancyGridMap2D &otherMap) |
copy the gridmap contents, but not all the options, from another map instance More... | |
void | subSample (int downRatio) |
Performs a downsampling of the gridmap, by a given factor: resolution/=ratio. More... | |
void | computeEntropy (TEntropyInfo &info) const |
Computes the entropy and related values of this grid map. More... | |
int | computeClearance (int cx, int cy, int *basis_x, int *basis_y, int *nBasis, bool GetContourPoint=false) const |
Compute the clearance of a given cell, and returns its two first basis (closest obstacle) points.Used to build Voronoi and critical points. More... | |
float | computeClearance (float x, float y, float maxSearchDistance) const |
An alternative method for computing the clearance of a given location (in meters). More... | |
float | computePathCost (float x1, float y1, float x2, float y2) const |
Compute the 'cost' of traversing a segment of the map according to the occupancy of traversed cells. More... | |
double | computeLikelihoodField_Thrun (const CPointsMap *pm, const mrpt::poses::CPose2D *relativePose=nullptr) |
Computes the likelihood [0,1] of a set of points, given the current grid map as reference. More... | |
double | computeLikelihoodField_II (const CPointsMap *pm, const mrpt::poses::CPose2D *relativePose=nullptr) |
Computes the likelihood [0,1] of a set of points, given the current grid map as reference. More... | |
bool | saveAsBitmapFile (const std::string &file) const |
Saves the gridmap as a graphical file (BMP,PNG,...). More... | |
template<class CLANDMARKSMAP > | |
bool | saveAsBitmapFileWithLandmarks (const std::string &file, const CLANDMARKSMAP *landmarks, bool addTextLabels=false, const mrpt::img::TColor &marks_color=mrpt::img::TColor(0, 0, 255)) const |
Saves the gridmap as a graphical bitmap file, 8 bit gray scale, 1 pixel is 1 cell, and with an overlay of landmarks. More... | |
void | getAsImage (mrpt::img::CImage &img, bool verticalFlip=false, bool forceRGB=false, bool tricolor=false) const |
Returns the grid as a 8-bit graylevel image, where each pixel is a cell (output image is RGB only if forceRGB is true) If "tricolor" is true, only three gray levels will appear in the image: gray for unobserved cells, and black/white for occupied/empty cells respectively. More... | |
void | getAsImageFiltered (img::CImage &img, bool verticalFlip=false, bool forceRGB=false) const |
Returns the grid as a 8-bit graylevel image, where each pixel is a cell (output image is RGB only if forceRGB is true) - This method filters the image for easy feature detection If "tricolor" is true, only three gray levels will appear in the image: gray for unobserved cells, and black/white for occupied/empty cells respectively. More... | |
void | getAs3DObject (mrpt::opengl::CSetOfObjects::Ptr &outObj) const override |
Returns a 3D plane with its texture being the occupancy grid and transparency proportional to "uncertainty" (i.e. More... | |
void | getAsPointCloud (mrpt::maps::CSimplePointsMap &pm, const float occup_threshold=0.5f) const |
Get a point cloud with all (border) occupied cells as points. More... | |
bool | isEmpty () const override |
Returns true upon map construction or after calling clear(), the return changes to false upon successful insertObservation() or any other method to load data in the map. More... | |
bool | loadFromBitmapFile (const std::string &file, float resolution, const mrpt::math::TPoint2D &origin=mrpt::math::TPoint2D(std::numeric_limits< double >::max(), std::numeric_limits< double >::max())) |
Load the gridmap from a image in a file (the format can be any supported by CImage::loadFromFile). More... | |
bool | loadFromBitmap (const mrpt::img::CImage &img, float resolution, const mrpt::math::TPoint2D &origin=mrpt::math::TPoint2D(std::numeric_limits< double >::max(), std::numeric_limits< double >::max())) |
Load the gridmap from a image in a file (the format can be any supported by CImage::loadFromFile). More... | |
void | determineMatching2D (const mrpt::maps::CMetricMap *otherMap, const mrpt::poses::CPose2D &otherMapPose, mrpt::tfest::TMatchingPairList &correspondences, const TMatchingParams ¶ms, TMatchingExtraResults &extraResults) const override |
See the base class for more details: In this class it is implemented as correspondences of the passed points map to occupied cells. More... | |
float | compute3DMatchingRatio (const mrpt::maps::CMetricMap *otherMap, const mrpt::poses::CPose3D &otherMapPose, const TMatchingRatioParams ¶ms) const override |
See docs in base class: in this class this always returns 0. More... | |
void | saveMetricMapRepresentationToFile (const std::string &filNamePrefix) const override |
This virtual method saves the map to a file "filNamePrefix"+< some_file_extension >, as an image or in any other applicable way (Notice that other methods to save the map may be implemented in classes implementing this virtual interface). More... | |
void | clear () |
Erase all the contents of the map. More... | |
void | loadFromProbabilisticPosesAndObservations (const mrpt::maps::CSimpleMap &Map) |
Load the map contents from a CSimpleMap object, erasing all previous content of the map. More... | |
void | loadFromSimpleMap (const mrpt::maps::CSimpleMap &Map) |
! More... | |
bool | insertObservation (const mrpt::obs::CObservation &obs, const mrpt::poses::CPose3D *robotPose=nullptr) |
Insert the observation information into this map. More... | |
bool | insertObservationPtr (const mrpt::obs::CObservation::Ptr &obs, const mrpt::poses::CPose3D *robotPose=nullptr) |
A wrapper for smart pointers, just calls the non-smart pointer version. More... | |
double | computeObservationLikelihood (const mrpt::obs::CObservation &obs, const mrpt::poses::CPose3D &takenFrom) |
Computes the log-likelihood of a given observation given an arbitrary robot 3D pose. More... | |
double | computeObservationLikelihood (const mrpt::obs::CObservation &obs, const mrpt::poses::CPose2D &takenFrom) |
virtual bool | canComputeObservationLikelihood (const mrpt::obs::CObservation &obs) const |
Returns true if this map is able to compute a sensible likelihood function for this observation (i.e. More... | |
double | computeObservationsLikelihood (const mrpt::obs::CSensoryFrame &sf, const mrpt::poses::CPose2D &takenFrom) |
Returns the sum of the log-likelihoods of each individual observation within a mrpt::obs::CSensoryFrame. More... | |
bool | canComputeObservationsLikelihood (const mrpt::obs::CSensoryFrame &sf) const |
Returns true if this map is able to compute a sensible likelihood function for this observation (i.e. More... | |
virtual void | determineMatching3D (const mrpt::maps::CMetricMap *otherMap, const mrpt::poses::CPose3D &otherMapPose, mrpt::tfest::TMatchingPairList &correspondences, const TMatchingParams ¶ms, TMatchingExtraResults &extraResults) const |
Computes the matchings between this and another 3D points map - method used in 3D-ICP. More... | |
virtual void | auxParticleFilterCleanUp () |
This method is called at the end of each "prediction-update-map
insertion" cycle within "mrpt::slam::CMetricMapBuilderRBPF::processActionObservation". More... | |
virtual float | squareDistanceToClosestCorrespondence (float x0, float y0) const |
Returns the square distance from the 2D point (x0,y0) to the closest correspondence in the map. More... | |
virtual const mrpt::maps::CSimplePointsMap * | getAsSimplePointsMap () const |
If the map is a simple points map or it's a multi-metric map that contains EXACTLY one simple points map, return it. More... | |
mrpt::maps::CSimplePointsMap * | getAsSimplePointsMap () |
virtual mxArray * | writeToMatlab () const |
Introduces a pure virtual method responsible for writing to a mxArray Matlab object, typically a MATLAB struct whose contents are documented in each derived class. More... | |
RTTI classes and functions for polymorphic hierarchies | |
mrpt::rtti::CObject::Ptr | duplicateGetSmartPtr () const |
Makes a deep copy of the object and returns a smart pointer to it. More... | |
Static Public Member Functions | |
static float | l2p (const cellType l) |
Scales an integer representation of the log-odd into a real valued probability in [0,1], using p=exp(l)/(1+exp(l)) More... | |
static uint8_t | l2p_255 (const cellType l) |
Scales an integer representation of the log-odd into a linear scale [0,255], using p=exp(l)/(1+exp(l)) More... | |
static cellType | p2l (const float p) |
Scales a real valued probability in [0,1] to an integer representation of: log(p)-log(1-p) in the valid range of cellType. More... | |
static bool | saveAsBitmapTwoMapsWithCorrespondences (const std::string &fileName, const COccupancyGridMap2D *m1, const COccupancyGridMap2D *m2, const mrpt::tfest::TMatchingPairList &corrs) |
Saves a composite image with two gridmaps and lines representing a set of correspondences between them. More... | |
static bool | saveAsEMFTwoMapsWithCorrespondences (const std::string &fileName, const COccupancyGridMap2D *m1, const COccupancyGridMap2D *m2, const mrpt::tfest::TMatchingPairList &corrs) |
Saves a composite image with two gridmaps and numbers for the correspondences between them. More... | |
static void | updateCell_fast_occupied (const unsigned x, const unsigned y, const cell_t logodd_obs, const cell_t thres, cell_t *mapArray, const unsigned _size_x) |
Performs the Bayesian fusion of a new observation of a cell, without checking for grid limits nor updateInfoChangeOnly. More... | |
static void | updateCell_fast_occupied (cell_t *theCell, const cell_t logodd_obs, const cell_t thres) |
Performs the Bayesian fusion of a new observation of a cell, without checking for grid limits nor updateInfoChangeOnly. More... | |
static void | updateCell_fast_free (const unsigned x, const unsigned y, const cell_t logodd_obs, const cell_t thres, cell_t *mapArray, const unsigned _size_x) |
Performs the Bayesian fusion of a new observation of a cell, without checking for grid limits nor updateInfoChangeOnly. More... | |
static void | updateCell_fast_free (cell_t *theCell, const cell_t logodd_obs, const cell_t thres) |
Performs the Bayesian fusion of a new observation of a cell, without checking for grid limits nor updateInfoChangeOnly. More... | |
Public Attributes | |
struct mrpt::maps::COccupancyGridMap2D::TUpdateCellsInfoChangeOnly | updateInfoChangeOnly |
TInsertionOptions | insertionOptions |
With this struct options are provided to the observation insertion process. More... | |
mrpt::maps::COccupancyGridMap2D::TLikelihoodOptions | likelihoodOptions |
class mrpt::maps::COccupancyGridMap2D::TLikelihoodOutput | likelihoodOutputs |
struct mrpt::maps::COccupancyGridMap2D::TCriticalPointsList | CriticalPointsList |
TMapGenericParams | genericMapParams |
Common params to all maps. More... | |
Static Public Attributes | |
static constexpr cellType | OCCGRID_CELLTYPE_MIN |
Discrete to float conversion factors: The min/max values of the integer cell type, eg. More... | |
static constexpr cellType | OCCGRID_CELLTYPE_MAX |
static double | RAYTRACE_STEP_SIZE_IN_CELL_UNITS = 0.8 |
(Default:1.0) Can be set to <1 if a more fine raytracing is needed in sonarSimulator() and laserScanSimulator(), or >1 to speed it up. More... | |
Protected Member Functions | |
void | freeMap () |
Frees the dynamic memory buffers of map. More... | |
void | OnPostSuccesfulInsertObs (const mrpt::obs::CObservation &) override |
See base class. More... | |
void | setCell_nocheck (int x, int y, float value) |
Change the contents [0,1] of a cell, given its index. More... | |
float | getCell_nocheck (int x, int y) const |
Read the real valued [0,1] contents of a cell, given its index. More... | |
void | setRawCell (unsigned int cellIndex, cellType b) |
Changes a cell by its absolute index (Do not use it normally) More... | |
double | computeObservationLikelihood_Consensus (const mrpt::obs::CObservation &obs, const mrpt::poses::CPose2D &takenFrom) |
One of the methods that can be selected for implementing "computeObservationLikelihood" (This method is the Range-Scan Likelihood Consensus for gridmaps, see the ICRA2007 paper by Blanco et al.) More... | |
double | computeObservationLikelihood_ConsensusOWA (const mrpt::obs::CObservation &obs, const mrpt::poses::CPose2D &takenFrom) |
One of the methods that can be selected for implementing "computeObservationLikelihood". More... | |
double | computeObservationLikelihood_CellsDifference (const mrpt::obs::CObservation &obs, const mrpt::poses::CPose2D &takenFrom) |
One of the methods that can be selected for implementing "computeObservationLikelihood". More... | |
double | computeObservationLikelihood_MI (const mrpt::obs::CObservation &obs, const mrpt::poses::CPose2D &takenFrom) |
One of the methods that can be selected for implementing "computeObservationLikelihood". More... | |
double | computeObservationLikelihood_rayTracing (const mrpt::obs::CObservation &obs, const mrpt::poses::CPose2D &takenFrom) |
One of the methods that can be selected for implementing "computeObservationLikelihood". More... | |
double | computeObservationLikelihood_likelihoodField_Thrun (const mrpt::obs::CObservation &obs, const mrpt::poses::CPose2D &takenFrom) |
One of the methods that can be selected for implementing "computeObservationLikelihood". More... | |
double | computeObservationLikelihood_likelihoodField_II (const mrpt::obs::CObservation &obs, const mrpt::poses::CPose2D &takenFrom) |
One of the methods that can be selected for implementing "computeObservationLikelihood". More... | |
void | internal_clear () override |
Clear the map: It set all cells to their default occupancy value (0.5), without changing the resolution (the grid extension is reset to the default values). More... | |
bool | internal_insertObservation (const mrpt::obs::CObservation &obs, const mrpt::poses::CPose3D *robotPose=nullptr) override |
Insert the observation information into this map. More... | |
void | publishEvent (const mrptEvent &e) const |
Called when you want this object to emit an event to all the observers currently subscribed to this object. More... | |
bool | hasSubscribers () const |
Can be called by a derived class before preparing an event for publishing with publishEvent to determine if there is no one subscribed, so it can save the wasted time preparing an event that will be not read. More... | |
CSerializable virtual methods | |
uint8_t | serializeGetVersion () const override |
Must return the current versioning number of the object. More... | |
void | serializeTo (mrpt::serialization::CArchive &out) const override |
Pure virtual method for writing (serializing) to an abstract archive. More... | |
void | serializeFrom (mrpt::serialization::CArchive &in, uint8_t serial_version) override |
Pure virtual method for reading (deserializing) from an abstract archive. More... | |
CSerializable virtual methods | |
virtual void | serializeTo (CSchemeArchiveBase &out) const |
Virtual method for writing (serializing) to an abstract schema based archive. More... | |
virtual void | serializeFrom (CSchemeArchiveBase &in) |
Virtual method for reading (deserializing) from an abstract schema based archive. More... | |
Static Protected Member Functions | |
static CLogOddsGridMapLUT< cellType > & | get_logodd_lut () |
Lookup tables for log-odds. More... | |
template<typename T > | |
static T | H (const T p) |
Entropy computation internal function: More... | |
Protected Attributes | |
std::vector< cellType > | map |
Store of cell occupancy values. More... | |
uint32_t | size_x {0} |
The size of the grid in cells. More... | |
uint32_t | size_y {0} |
float | x_min {} |
The limits of the grid in "units" (meters) More... | |
float | x_max {} |
float | y_min {} |
float | y_max {} |
float | resolution {} |
Cell size, i.e. More... | |
std::vector< double > | precomputedLikelihood |
Auxiliary variables to speed up the computation of observation likelihood values for LF method among others, at a high cost in memory (see TLikelihoodOptions::enableLikelihoodCache). More... | |
bool | m_likelihoodCacheOutDated {true} |
mrpt::containers::CDynamicGrid< uint8_t > | m_basis_map |
Used for Voronoi calculation.Same struct as "map", but contains a "0" if not a basis point. More... | |
mrpt::containers::CDynamicGrid< uint16_t > | m_voronoi_diagram |
Used to store the Voronoi diagram. More... | |
bool | m_is_empty {true} |
True upon construction; used by isEmpty() More... | |
float | voroni_free_threshold {} |
The free-cells threshold used to compute the Voronoi diagram. More... | |
Static Protected Attributes | |
static std::vector< float > | entropyTable |
Internally used to speed-up entropy calculation. More... | |
Private Member Functions | |
double | internal_computeObservationLikelihood (const mrpt::obs::CObservation &obs, const mrpt::poses::CPose3D &takenFrom) override |
Internal method called by computeObservationLikelihood() More... | |
bool | internal_canComputeObservationLikelihood (const mrpt::obs::CObservation &obs) const override |
Returns true if this map is able to compute a sensible likelihood function for this observation (i.e. More... | |
unsigned char | GetNeighborhood (int cx, int cy) const |
Returns a byte with the occupancy of the 8 sorrounding cells. More... | |
int | direction2idx (int dx, int dy) |
Returns the index [0,7] of the given movement, or -1 if invalid one. More... | |
Private Attributes | |
int | direccion_vecino_x [8] |
Used to store the 8 possible movements from a cell to the sorrounding ones.Filled in the constructor. More... | |
int | direccion_vecino_y [8] |
Friends | |
class | CMultiMetricMap |
class | CMultiMetricMapPDF |
RTTI stuff | |
using | Ptr = std::shared_ptr< mrpt::maps ::COccupancyGridMap2D > |
using | ConstPtr = std::shared_ptr< const mrpt::maps ::COccupancyGridMap2D > |
using | UniquePtr = std::unique_ptr< mrpt::maps ::COccupancyGridMap2D > |
using | ConstUniquePtr = std::unique_ptr< const mrpt::maps ::COccupancyGridMap2D > |
static const mrpt::rtti::TRuntimeClassId | runtimeClassId |
static constexpr const char * | className = "mrpt::maps" "::" "COccupancyGridMap2D" |
static const mrpt::rtti::TRuntimeClassId * | _GetBaseClass () |
static constexpr auto | getClassName () |
static const mrpt::rtti::TRuntimeClassId & | GetRuntimeClassIdStatic () |
static std::shared_ptr< CObject > | CreateObject () |
template<typename... Args> | |
static Ptr | Create (Args &&... args) |
template<typename Alloc , typename... Args> | |
static Ptr | CreateAlloc (const Alloc &alloc, Args &&... args) |
template<typename... Args> | |
static UniquePtr | CreateUnique (Args &&... args) |
virtual const mrpt::rtti::TRuntimeClassId * | GetRuntimeClass () const override |
Returns information about the class of an object in runtime. More... | |
virtual mrpt::rtti::CObject * | clone () const override |
Returns a deep copy (clone) of the object, indepently of its class. More... | |
Voronoi methods | |
void | buildVoronoiDiagram (float threshold, float robot_size, int x1=0, int x2=0, int y1=0, int y2=0) |
Build the Voronoi diagram of the grid map. More... | |
uint16_t | getVoroniClearance (int cx, int cy) const |
Reads a the clearance of a cell (in centimeters), after building the Voronoi diagram with buildVoronoiDiagram. More... | |
const mrpt::containers::CDynamicGrid< uint8_t > & | getBasisMap () const |
Return the auxiliary "basis" map built while building the Voronoi diagram. More... | |
const mrpt::containers::CDynamicGrid< uint16_t > & | getVoronoiDiagram () const |
Return the Voronoi diagram; each cell contains the distance to its closer obstacle, or 0 if not part of the Voronoi diagram. More... | |
void | findCriticalPoints (float filter_distance) |
Builds a list with the critical points from Voronoi diagram, which must must be built before calling this method. More... | |
void | setVoroniClearance (int cx, int cy, uint16_t dist) |
Used to set the clearance of a cell, while building the Voronoi diagram. More... | |
Sensor simulators | |
enum | TLaserSimulUncertaintyMethod { sumUnscented = 0, sumMonteCarlo } |
Methods for TLaserSimulUncertaintyParams in laserScanSimulatorWithUncertainty() More... | |
void | laserScanSimulator (mrpt::obs::CObservation2DRangeScan &inout_Scan, const mrpt::poses::CPose2D &robotPose, float threshold=0.6f, size_t N=361, float noiseStd=0, unsigned int decimation=1, float angleNoiseStd=mrpt::DEG2RAD(.0)) const |
Simulates a laser range scan into the current grid map. More... | |
void | sonarSimulator (mrpt::obs::CObservationRange &inout_observation, const mrpt::poses::CPose2D &robotPose, float threshold=0.5f, float rangeNoiseStd=0.f, float angleNoiseStd=mrpt::DEG2RAD(0.f)) const |
Simulates the observations of a sonar rig into the current grid map. More... | |
void | simulateScanRay (const double x, const double y, const double angle_direction, float &out_range, bool &out_valid, const double max_range_meters, const float threshold_free=0.4f, const double noiseStd=.0, const double angleNoiseStd=.0) const |
Simulate just one "ray" in the grid map. More... | |
void | laserScanSimulatorWithUncertainty (const TLaserSimulUncertaintyParams &in_params, TLaserSimulUncertaintyResult &out_results) const |
Like laserScanSimulatorWithUncertainty() (see it for a discussion of most parameters) but taking into account the robot pose uncertainty and generating a vector of predicted variances for each ray. More... | |
Map Definition Interface stuff (see * mrpt::maps::TMetricMapInitializer) @{ | |
static const size_t | m_private_map_register_id = mrpt::maps::internal::TMetricMapTypesRegistry::Instance().doRegister( "mrpt::maps::COccupancyGridMap2D,occupancyGrid" , & mrpt::maps::COccupancyGridMap2D ::MapDefinition, & mrpt::maps::COccupancyGridMap2D ::internal_CreateFromMapDefinition) |
ID used to initialize class registration (just ignore it) More... | |
static mrpt::maps::TMetricMapInitializer * | MapDefinition () |
Returns default map definition initializer. More... | |
static COccupancyGridMap2D * | CreateFromMapDefinition (const mrpt::maps::TMetricMapInitializer &def) |
Constructor from a map definition structure: initializes the map and * its parameters accordingly. More... | |
static mrpt::maps::CMetricMap * | internal_CreateFromMapDefinition (const mrpt::maps::TMetricMapInitializer &def) |
|
inherited |
The type of cells.
Definition at line 29 of file CLogOddsGridMap2D.h.
The type of the map cells:
Definition at line 60 of file COccupancyGridMap2D.h.
Definition at line 61 of file COccupancyGridMap2D.h.
using mrpt::maps::COccupancyGridMap2D::ConstPtr = std::shared_ptr<const mrpt::maps :: COccupancyGridMap2D > |
Definition at line 57 of file COccupancyGridMap2D.h.
using mrpt::maps::COccupancyGridMap2D::ConstUniquePtr = std::unique_ptr<const mrpt::maps :: COccupancyGridMap2D > |
Definition at line 57 of file COccupancyGridMap2D.h.
using mrpt::maps::COccupancyGridMap2D::Ptr = std::shared_ptr< mrpt::maps :: COccupancyGridMap2D > |
A type for the associated smart pointer
Definition at line 57 of file COccupancyGridMap2D.h.
using mrpt::maps::COccupancyGridMap2D::TPairLikelihoodIndex = std::pair<double, mrpt::math::TPoint2D> |
Auxiliary private class.
Definition at line 628 of file COccupancyGridMap2D.h.
|
inherited |
Definition at line 30 of file CLogOddsGridMap2D.h.
using mrpt::maps::COccupancyGridMap2D::UniquePtr = std::unique_ptr< mrpt::maps :: COccupancyGridMap2D > |
Definition at line 57 of file COccupancyGridMap2D.h.
Methods for TLaserSimulUncertaintyParams in laserScanSimulatorWithUncertainty()
Enumerator | |
---|---|
sumUnscented | Performs an unscented transform. |
sumMonteCarlo | Montecarlo-based estimation. |
Definition at line 831 of file COccupancyGridMap2D.h.
The type for selecting a likelihood computation method.
Enumerator | |
---|---|
lmMeanInformation | |
lmRayTracing | |
lmConsensus | |
lmCellsDifference | |
lmLikelihoodField_Thrun | |
lmLikelihoodField_II | |
lmConsensusOWA |
Definition at line 533 of file COccupancyGridMap2D.h.
COccupancyGridMap2D::COccupancyGridMap2D | ( | float | min_x = -20.0f , |
float | max_x = 20.0f , |
||
float | min_y = -20.0f , |
||
float | max_y = 20.0f , |
||
float | resolution = 0.05f |
||
) |
Constructor.
Definition at line 101 of file COccupancyGridMap2D_common.cpp.
References MRPT_END, MRPT_START, and setSize().
|
staticprotected |
|
inlinevirtualinherited |
This method is called at the end of each "prediction-update-map insertion" cycle within "mrpt::slam::CMetricMapBuilderRBPF::processActionObservation".
This method should normally do nothing, but in some cases can be used to free auxiliary cached variables.
Reimplemented in mrpt::maps::CLandmarksMap, and mrpt::maps::CMultiMetricMap.
Definition at line 282 of file CMetricMap.h.
void COccupancyGridMap2D::buildVoronoiDiagram | ( | float | threshold, |
float | robot_size, | ||
int | x1 = 0 , |
||
int | x2 = 0 , |
||
int | y1 = 0 , |
||
int | y2 = 0 |
||
) |
Build the Voronoi diagram of the grid map.
threshold | The threshold for binarizing the map. |
robot_size | Size in "units" (meters) of robot, approx. |
x1 | Left coordinate of area to be computed. Default, entire map. |
x2 | Right coordinate of area to be computed. Default, entire map. |
y1 | Top coordinate of area to be computed. Default, entire map. |
y2 | Bottom coordinate of area to be computed. Default, entire map. |
Definition at line 24 of file COccupancyGridMap2D_voronoi.cpp.
References ASSERT_EQUAL_, and mrpt::round().
|
virtualinherited |
Returns true if this map is able to compute a sensible likelihood function for this observation (i.e.
an occupancy grid map cannot with an image). See: Maps and observations compatibility matrix
obs | The observation. |
Definition at line 161 of file CMetricMap.cpp.
|
inherited |
Returns true if this map is able to compute a sensible likelihood function for this observation (i.e.
an occupancy grid map cannot with an image). See: Maps and observations compatibility matrix
sf | The observations. |
Definition at line 85 of file CMetricMap.cpp.
References mrpt::obs::CSensoryFrame::begin(), and mrpt::obs::CSensoryFrame::end().
|
inherited |
Erase all the contents of the map.
Definition at line 30 of file CMetricMap.cpp.
Referenced by mrpt::maps::CHeightGridMap2D_MRF::CHeightGridMap2D_MRF(), mrpt::maps::CReflectivityGridMap2D::clear(), mrpt::maps::CHeightGridMap2D::clear(), mrpt::maps::CRandomFieldGridMap2D::clear(), mrpt::apps::MonteCarloLocalization_Base::do_pf_localization(), mrpt::maps::CPointsMap::extractCylinder(), mrpt::maps::CPointsMap::extractPoints(), mrpt::ros1bridge::fromROS(), mrpt::opengl::CAngularObservationMesh::generatePointCloud(), getAsPointCloud(), mrpt::nav::CReactiveNavigationSystem3D::loggingGetWSObstaclesAndShape(), mrpt::hmtslam::CLSLAM_RBPF_2DLASER::prediction_and_update_pfOptimalProposal(), mrpt::vision::projectMatchedFeatures(), run_rnav_test(), mrpt::maps::CPointsMap::setFromPCLPointCloud(), mrpt::maps::CColouredPointsMap::setFromPCLPointCloudRGB(), mrpt::maps::CPointsMapXYZI::setFromPCLPointCloudXYZI(), mrpt::maps::CRandomFieldGridMap2D::setSize(), and mrpt::nav::PlannerTPS_VirtualBase::transformPointcloudWithSquareClipping().
|
overridevirtual |
Returns a deep copy (clone) of the object, indepently of its class.
Implements mrpt::rtti::CObject.
|
overridevirtual |
See docs in base class: in this class this always returns 0.
Reimplemented from mrpt::maps::CMetricMap.
Definition at line 762 of file COccupancyGridMap2D_common.cpp.
int COccupancyGridMap2D::computeClearance | ( | int | cx, |
int | cy, | ||
int * | basis_x, | ||
int * | basis_y, | ||
int * | nBasis, | ||
bool | GetContourPoint = false |
||
) | const |
Compute the clearance of a given cell, and returns its two first basis (closest obstacle) points.Used to build Voronoi and critical points.
cx | The cell index |
cy | The cell index |
basis_x | Target buffer for coordinates of basis, having a size of two "ints". |
basis_y | Target buffer for coordinates of basis, having a size of two "ints". |
nBasis | The number of found basis: Can be 0,1 or 2. |
GetContourPoint | If "true" the basis are not returned, but the closest free cells.Default at false. |
Definition at line 255 of file COccupancyGridMap2D_voronoi.cpp.
References dir, M_2PI, M_PI, M_PIf, N_CIRCULOS, and mrpt::round().
float COccupancyGridMap2D::computeClearance | ( | float | x, |
float | y, | ||
float | maxSearchDistance | ||
) | const |
An alternative method for computing the clearance of a given location (in meters).
Definition at line 561 of file COccupancyGridMap2D_voronoi.cpp.
References mrpt::square().
void COccupancyGridMap2D::computeEntropy | ( | TEntropyInfo & | info | ) | const |
Computes the entropy and related values of this grid map.
The entropy is computed as the summed entropy of each cell, taking them as discrete random variables following a Bernoulli distribution:
info | The output information is returned here |
Definition at line 352 of file COccupancyGridMap2D_common.cpp.
References mrpt::d2f(), mrpt::maps::COccupancyGridMap2D::TEntropyInfo::effectiveMappedArea, mrpt::maps::COccupancyGridMap2D::TEntropyInfo::effectiveMappedCells, entropyTable, H(), mrpt::maps::COccupancyGridMap2D::TEntropyInfo::H, mrpt::maps::COccupancyGridMap2D::TEntropyInfo::I, l2p(), map, MAX_H, mrpt::maps::COccupancyGridMap2D::TEntropyInfo::mean_H, mrpt::maps::COccupancyGridMap2D::TEntropyInfo::mean_I, and resolution.
Referenced by mrpt::maps::CMultiMetricMapPDF::getCurrentJointEntropy().
double COccupancyGridMap2D::computeLikelihoodField_II | ( | const CPointsMap * | pm, |
const mrpt::poses::CPose2D * | relativePose = nullptr |
||
) |
Computes the likelihood [0,1] of a set of points, given the current grid map as reference.
pm | The points map |
relativePose | The relative pose of the points map in this map's coordinates, or nullptr for (0,0,0). See "likelihoodOptions" for configuration parameters. |
Definition at line 712 of file COccupancyGridMap2D_likelihood.cpp.
References mrpt::maps::CPointsMap::getPoint(), MRPT_END, MRPT_START, mrpt::maps::CPointsMap::size(), mrpt::square(), mrpt::math::TPoint2D_data< T >::x, and mrpt::math::TPoint2D_data< T >::y.
double COccupancyGridMap2D::computeLikelihoodField_Thrun | ( | const CPointsMap * | pm, |
const mrpt::poses::CPose2D * | relativePose = nullptr |
||
) |
Computes the likelihood [0,1] of a set of points, given the current grid map as reference.
pm | The points map |
relativePose | The relative pose of the points map in this map's coordinates, or nullptr for (0,0,0). See "likelihoodOptions" for configuration parameters. |
Definition at line 517 of file COccupancyGridMap2D_likelihood.cpp.
References mrpt::maps::CPointsMap::getPoint(), mrpt::keep_min(), LIK_LF_CACHE_INVALID, MRPT_END, MRPT_START, mrpt::poses::CPose2D::phi(), mrpt::round(), mrpt::maps::CPointsMap::size(), mrpt::square(), mrpt::math::TPoint2D_data< T >::x, mrpt::poses::CPoseOrPoint< DERIVEDCLASS, DIM >::x(), mrpt::math::TPoint2D_data< T >::y, and mrpt::poses::CPoseOrPoint< DERIVEDCLASS, DIM >::y().
|
inherited |
Computes the log-likelihood of a given observation given an arbitrary robot 3D pose.
See: Maps and observations compatibility matrix
takenFrom | The robot's pose the observation is supposed to be taken from. |
obs | The observation. |
Definition at line 170 of file CMetricMap.cpp.
Referenced by mrpt::maps::CMultiMetricMapPDF::PF_SLAM_computeObservationLikelihoodForParticle().
|
inherited |
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
Definition at line 76 of file CMetricMap.cpp.
|
protected |
One of the methods that can be selected for implementing "computeObservationLikelihood".
Definition at line 241 of file COccupancyGridMap2D_likelihood.cpp.
References IS_CLASS, mrpt::round(), mrpt::poses::CPoseOrPoint< DERIVEDCLASS, DIM >::x(), and mrpt::poses::CPoseOrPoint< DERIVEDCLASS, DIM >::y().
|
protected |
One of the methods that can be selected for implementing "computeObservationLikelihood" (This method is the Range-Scan Likelihood Consensus for gridmaps, see the ICRA2007 paper by Blanco et al.)
Definition at line 85 of file COccupancyGridMap2D_likelihood.cpp.
References mrpt::obs::CObservation2DRangeScan::buildAuxPointsMap(), mrpt::poses::CPose2D::composePoint(), IS_CLASS, mrpt::math::internal::ProvideStaticResize< Derived >::size(), mrpt::math::TPoint2D_data< T >::x, and mrpt::math::TPoint2D_data< T >::y.
|
protected |
One of the methods that can be selected for implementing "computeObservationLikelihood".
TODO: This method is described in....
Definition at line 144 of file COccupancyGridMap2D_likelihood.cpp.
References ASSERT_, mrpt::poses::CPose2D::composePoint(), IS_CLASS, mrpt::maps::CPointsMap::TInsertionOptions::minDistBetweenLaserPoints, mrpt::math::internal::ProvideStaticResize< Derived >::size(), mrpt::math::TPoint2D_data< T >::x, and mrpt::math::TPoint2D_data< T >::y.
|
protected |
One of the methods that can be selected for implementing "computeObservationLikelihood".
Definition at line 480 of file COccupancyGridMap2D_likelihood.cpp.
References IS_CLASS, MRPT_END, and MRPT_START.
|
protected |
One of the methods that can be selected for implementing "computeObservationLikelihood".
Definition at line 426 of file COccupancyGridMap2D_likelihood.cpp.
References mrpt::maps::CPointsMap::TInsertionOptions::horizontalTolerance, mrpt::maps::CPointsMap::insertionOptions, mrpt::maps::CMetricMap::insertObservation(), IS_CLASS, mrpt::maps::CPointsMap::TInsertionOptions::isPlanarMap, mrpt::maps::CPointsMap::TInsertionOptions::minDistBetweenLaserPoints, MRPT_END, and MRPT_START.
|
protected |
One of the methods that can be selected for implementing "computeObservationLikelihood".
Definition at line 308 of file COccupancyGridMap2D_likelihood.cpp.
References MRPT_END, and MRPT_START.
|
protected |
One of the methods that can be selected for implementing "computeObservationLikelihood".
Definition at line 352 of file COccupancyGridMap2D_likelihood.cpp.
References mrpt::obs::CObservation2DRangeScan::aperture, mrpt::obs::CObservation2DRangeScan::getScanRange(), mrpt::obs::CObservation2DRangeScan::getScanSize(), IS_CLASS, mrpt::obs::CObservation2DRangeScan::maxRange, mrpt::obs::CObservation2DRangeScan::rightToLeft, mrpt::obs::CObservation2DRangeScan::sensorPose, and mrpt::square().
|
inherited |
Returns the sum of the log-likelihoods of each individual observation within a mrpt::obs::CSensoryFrame.
See: Maps and observations compatibility matrix
takenFrom | The robot's pose the observation is supposed to be taken from. |
sf | The set of observations in a CSensoryFrame. |
Definition at line 66 of file CMetricMap.cpp.
Referenced by mrpt::hmtslam::CLSLAM_RBPF_2DLASER::auxiliarComputeObservationLikelihood().
float COccupancyGridMap2D::computePathCost | ( | float | x1, |
float | y1, | ||
float | x2, | ||
float | y2 | ||
) | const |
Compute the 'cost' of traversing a segment of the map according to the occupancy of traversed cells.
Definition at line 741 of file COccupancyGridMap2D_common.cpp.
References mrpt::d2f(), getPos(), resolution, mrpt::round(), and mrpt::square().
void COccupancyGridMap2D::copyMapContentFrom | ( | const COccupancyGridMap2D & | otherMap | ) |
copy the gridmap contents, but not all the options, from another map instance
Definition at line 121 of file COccupancyGridMap2D_common.cpp.
References mrpt::containers::CDynamicGrid< T >::clear(), freeMap(), m_basis_map, m_is_empty, m_likelihoodCacheOutDated, m_voronoi_diagram, map, resolution, size_x, size_y, x_max, x_min, y_max, and y_min.
|
inlinestatic |
Definition at line 57 of file COccupancyGridMap2D.h.
Referenced by mrpt::graphslam::CGraphSlamEngine< GRAPH_T >::getMap(), mrpt::graphslam::CGraphSlamEngine< GRAPH_T >::initClass(), and CGraphSlamHandler< GRAPH_T >::saveMap().
|
inlinestatic |
Definition at line 57 of file COccupancyGridMap2D.h.
|
static |
Constructor from a map definition structure: initializes the map and * its parameters accordingly.
Definition at line 32 of file COccupancyGridMap2D_common.cpp.
|
static |
|
inlinestatic |
Definition at line 57 of file COccupancyGridMap2D.h.
|
overridevirtual |
See the base class for more details: In this class it is implemented as correspondences of the passed points map to occupied cells.
NOTICE: That the "z" dimension is ignored in the points. Clip the points as appropiated if needed before calling this method.
Reimplemented from mrpt::maps::CMetricMap.
Definition at line 518 of file COccupancyGridMap2D_common.cpp.
References ASSERT_, ASSERT_ABOVE_, ASSERT_BELOW_, mrpt::poses::CPose2D::asTPose(), CLASS_ID, mrpt::d2f(), mrpt::rtti::TRuntimeClassId::derivedFrom(), mrpt::maps::CMetricMap::GetRuntimeClass(), idx2x(), idx2y(), map, MRPT_END, MRPT_START, mrpt::tfest::TMatchingPair::other_idx, mrpt::tfest::TMatchingPair::other_x, mrpt::tfest::TMatchingPair::other_y, mrpt::tfest::TMatchingPair::other_z, p2l(), params, mrpt::math::TPose2D::phi, resolution, mrpt::round(), mrpt::math::internal::ProvideStaticResize< Derived >::size(), size_x, size_y, mrpt::square(), mrpt::tfest::TMatchingPair::this_idx, mrpt::tfest::TMatchingPair::this_x, mrpt::tfest::TMatchingPair::this_y, mrpt::tfest::TMatchingPair::this_z, mrpt::math::TPose2D::x, x2idx(), x_max, mrpt::math::TPose2D::y, y2idx(), y_max, and y_min.
|
virtualinherited |
Computes the matchings between this and another 3D points map - method used in 3D-ICP.
This method finds the set of point pairs in each map.
The method is the most time critical one into ICP-like algorithms.
The algorithm is:
otherMap | [IN] The other map to compute the matching with. |
otherMapPose | [IN] The pose of the other map as seen from "this". |
params | [IN] Parameters for the determination of pairings. |
correspondences | [OUT] The detected matchings pairs. |
extraResults | [OUT] Other results. |
Reimplemented in mrpt::maps::CPointsMap.
Definition at line 131 of file CMetricMap.cpp.
References MRPT_END, MRPT_START, and THROW_EXCEPTION.
Referenced by mrpt::slam::CICP::ICP3D_Method_Classic().
|
private |
Returns the index [0,7] of the given movement, or -1 if invalid one.
Definition at line 515 of file COccupancyGridMap2D_voronoi.cpp.
|
inlineinherited |
Makes a deep copy of the object and returns a smart pointer to it.
Definition at line 204 of file CObject.h.
References mrpt::rtti::CObject::clone().
Referenced by mrpt::obs::CRawlog::insert().
void COccupancyGridMap2D::fill | ( | float | default_value = 0.5f | ) |
Fills all the cells with a default value.
Definition at line 429 of file COccupancyGridMap2D_common.cpp.
References m_likelihoodCacheOutDated, map, and p2l().
Referenced by run_rnav_test().
void COccupancyGridMap2D::findCriticalPoints | ( | float | filter_distance | ) |
Builds a list with the critical points from Voronoi diagram, which must must be built before calling this method.
filter_distance | The minimum distance between two critical points. |
Definition at line 97 of file COccupancyGridMap2D_voronoi.cpp.
References ASSERT_EQUAL_, mrpt::containers::clear(), and mrpt::round().
|
protected |
Frees the dynamic memory buffers of map.
Definition at line 318 of file COccupancyGridMap2D_common.cpp.
References mrpt::containers::CDynamicGrid< T >::clear(), m_basis_map, m_is_empty, m_likelihoodCacheOutDated, m_voronoi_diagram, map, MRPT_END, MRPT_START, size_x, and size_y.
Referenced by copyMapContentFrom(), and setSize().
|
staticprotected |
Lookup tables for log-odds.
Definition at line 93 of file COccupancyGridMap2D_common.cpp.
References logodd_lut.
Referenced by l2p(), l2p_255(), and p2l().
|
inline |
Returns the area of the gridmap, in square meters.
Definition at line 268 of file COccupancyGridMap2D.h.
References resolution, size_x, size_y, and mrpt::square().
Referenced by mrpt::slam::CGridMapAligner::AlignPDF_robustMatch().
|
overridevirtual |
Returns a 3D plane with its texture being the occupancy grid and transparency proportional to "uncertainty" (i.e.
a value of 0.5 is fully transparent)
Implements mrpt::maps::CMetricMap.
Definition at line 148 of file COccupancyGridMap2D_getAs.cpp.
References mrpt::img::CH_GRAY, mrpt::opengl::CTexturedPlane::Create(), MRPT_END, and MRPT_START.
void COccupancyGridMap2D::getAsImage | ( | mrpt::img::CImage & | img, |
bool | verticalFlip = false , |
||
bool | forceRGB = false , |
||
bool | tricolor = false |
||
) | const |
Returns the grid as a 8-bit graylevel image, where each pixel is a cell (output image is RGB only if forceRGB is true) If "tricolor" is true, only three gray levels will appear in the image: gray for unobserved cells, and black/white for occupied/empty cells respectively.
Definition at line 29 of file COccupancyGridMap2D_getAs.cpp.
References mrpt::img::CH_GRAY, mrpt::img::CH_RGB, and mrpt::img::CImage::resize().
Referenced by mrpt::slam::CGridMapAligner::AlignPDF_correlation(), saveAsBitmapTwoMapsWithCorrespondences(), saveAsEMFTwoMapsWithCorrespondences(), and mrpt::slam::CMetricMapBuilderICP::saveCurrentEstimationToImage().
void COccupancyGridMap2D::getAsImageFiltered | ( | img::CImage & | img, |
bool | verticalFlip = false , |
||
bool | forceRGB = false |
||
) | const |
Returns the grid as a 8-bit graylevel image, where each pixel is a cell (output image is RGB only if forceRGB is true) - This method filters the image for easy feature detection If "tricolor" is true, only three gray levels will appear in the image: gray for unobserved cells, and black/white for occupied/empty cells respectively.
Definition at line 132 of file COccupancyGridMap2D_getAs.cpp.
References mrpt::img::CImage::filterGaussian(), mrpt::img::CImage::filterMedian(), and mrpt::round().
Referenced by saveAsBitmapFileWithLandmarks(), and mrpt::slam::COccupancyGridMapFeatureExtractor::uncached_extractFeatures().
void COccupancyGridMap2D::getAsPointCloud | ( | mrpt::maps::CSimplePointsMap & | pm, |
const float | occup_threshold = 0.5f |
||
) | const |
Get a point cloud with all (border) occupied cells as points.
Definition at line 188 of file COccupancyGridMap2D_getAs.cpp.
References mrpt::maps::CMetricMap::clear(), mrpt::maps::CPointsMap::insertPoint(), and mrpt::maps::CSimplePointsMap::reserve().
|
inlinevirtualinherited |
If the map is a simple points map or it's a multi-metric map that contains EXACTLY one simple points map, return it.
Otherwise, return nullptr
Reimplemented in mrpt::maps::CMultiMetricMap, and mrpt::maps::CSimplePointsMap.
Definition at line 295 of file CMetricMap.h.
Referenced by mrpt::maps::CPointsMap::compute3DMatchingRatio(), and mrpt::maps::CMetricMap::getAsSimplePointsMap().
|
inlineinherited |
Definition at line 299 of file CMetricMap.h.
References mrpt::maps::CMetricMap::getAsSimplePointsMap().
|
inline |
Reads a cell in the "basis" maps.Used for Voronoi calculation.
Definition at line 424 of file COccupancyGridMap2D.h.
References ASSERT_ABOVEEQ_, ASSERT_BELOWEQ_, mrpt::containers::CDynamicGrid< T >::cellByIndex(), mrpt::containers::CDynamicGrid< T >::getSizeX(), mrpt::containers::CDynamicGrid< T >::getSizeY(), and m_basis_map.
|
inline |
Return the auxiliary "basis" map built while building the Voronoi diagram.
Definition at line 703 of file COccupancyGridMap2D.h.
References m_basis_map.
|
inline |
Read the real valued [0,1] contents of a cell, given its index.
Definition at line 357 of file COccupancyGridMap2D.h.
References l2p(), map, size_x, and size_y.
Referenced by mrpt::nav::PlannerSimple2D::computePath(), mrpt::maps::CRandomFieldGridMap2D::exist_relation_between2cells(), getPos(), mrpt::maps::CRandomFieldGridMap2D::internal_clear(), isStaticCell(), mrpt::slam::CMonteCarloLocalization2D::resetUniformFreeSpace(), subSample(), and TEST().
|
inlineprotected |
Read the real valued [0,1] contents of a cell, given its index.
Definition at line 136 of file COccupancyGridMap2D.h.
References l2p(), map, and size_x.
|
inlinestatic |
Definition at line 57 of file COccupancyGridMap2D.h.
|
inlineprivate |
Returns a byte with the occupancy of the 8 sorrounding cells.
cx | The cell index |
cy | The cell index |
Definition at line 493 of file COccupancyGridMap2D_voronoi.cpp.
|
inline |
Read the real valued [0,1] contents of a cell, given its coordinates.
Definition at line 394 of file COccupancyGridMap2D.h.
References getCell(), x2idx(), and y2idx().
Referenced by computePathCost(), and TEST().
|
inline |
Read-only access to the raw cell contents (cells are in log-odd units)
Definition at line 206 of file COccupancyGridMap2D.h.
References map.
|
inline |
Returns the resolution of the grid map.
Definition at line 286 of file COccupancyGridMap2D.h.
References resolution.
Referenced by mrpt::slam::CGridMapAligner::AlignPDF_robustMatch(), mrpt::nav::PlannerSimple2D::computePath(), internal_clear(), mrpt::maps::CRandomFieldGridMap2D::internal_clear(), mrpt::maps::CMultiMetricMapPDF::rebuildAverageMap(), mrpt::slam::CMonteCarloLocalization2D::resetUniformFreeSpace(), TEST(), mrpt::ros1bridge::toROS(), and mrpt::slam::COccupancyGridMapFeatureExtractor::uncached_extractFeatures().
|
inline |
Access to a "row": mainly used for drawing grid as a bitmap efficiently, do not use it normally.
Definition at line 369 of file COccupancyGridMap2D.h.
References map, size_x, and size_y.
Referenced by mrpt::slam::CGridMapAligner::AlignPDF_correlation(), mrpt::ros1bridge::fromROS(), and mrpt::ros1bridge::toROS().
|
inline |
Access to a "row": mainly used for drawing grid as a bitmap efficiently, do not use it normally.
Definition at line 379 of file COccupancyGridMap2D.h.
|
overridevirtual |
Returns information about the class of an object in runtime.
Reimplemented from mrpt::maps::CMetricMap.
|
static |
|
inline |
Returns the horizontal size of grid map in cells count.
Definition at line 274 of file COccupancyGridMap2D.h.
References size_x.
Referenced by mrpt::slam::CGridMapAligner::AlignPDF_correlation(), mrpt::slam::CGridMapAligner::AlignPDF_robustMatch(), mrpt::nav::PlannerSimple2D::computePath(), mrpt::maps::CRandomFieldGridMap2D::exist_relation_between2cells(), mrpt::maps::CRandomFieldGridMap2D::internal_clear(), mrpt::slam::CMonteCarloLocalization2D::resetUniformFreeSpace(), TEST(), and mrpt::ros1bridge::toROS().
|
inline |
Returns the vertical size of grid map in cells count.
Definition at line 276 of file COccupancyGridMap2D.h.
References size_y.
Referenced by mrpt::slam::CGridMapAligner::AlignPDF_correlation(), mrpt::slam::CGridMapAligner::AlignPDF_robustMatch(), mrpt::nav::PlannerSimple2D::computePath(), mrpt::maps::CRandomFieldGridMap2D::exist_relation_between2cells(), mrpt::maps::CRandomFieldGridMap2D::internal_clear(), mrpt::slam::CMonteCarloLocalization2D::resetUniformFreeSpace(), TEST(), and mrpt::ros1bridge::toROS().
|
inline |
Reads a the clearance of a cell (in centimeters), after building the Voronoi diagram with buildVoronoiDiagram.
Definition at line 673 of file COccupancyGridMap2D.h.
References ASSERT_ABOVEEQ_, ASSERT_BELOWEQ_, mrpt::containers::CDynamicGrid< T >::cellByIndex(), mrpt::containers::CDynamicGrid< T >::getSizeX(), mrpt::containers::CDynamicGrid< T >::getSizeY(), and m_voronoi_diagram.
|
inline |
Return the Voronoi diagram; each cell contains the distance to its closer obstacle, or 0 if not part of the Voronoi diagram.
Definition at line 711 of file COccupancyGridMap2D.h.
References m_voronoi_diagram.
|
inline |
Returns the "x" coordinate of right side of grid map.
Definition at line 280 of file COccupancyGridMap2D.h.
References x_max.
Referenced by mrpt::slam::CGridMapAligner::AlignPDF_correlation(), mrpt::nav::PlannerSimple2D::computePath(), and mrpt::maps::CRandomFieldGridMap2D::internal_clear().
|
inline |
Returns the "x" coordinate of left side of grid map.
Definition at line 278 of file COccupancyGridMap2D.h.
References x_min.
Referenced by mrpt::slam::CGridMapAligner::AlignPDF_correlation(), mrpt::nav::PlannerSimple2D::computePath(), mrpt::maps::CRandomFieldGridMap2D::internal_clear(), mrpt::maps::CMultiMetricMapPDF::rebuildAverageMap(), mrpt::slam::CMonteCarloLocalization2D::resetUniformFreeSpace(), mrpt::ros1bridge::toROS(), and mrpt::slam::COccupancyGridMapFeatureExtractor::uncached_extractFeatures().
|
inline |
Returns the "y" coordinate of bottom side of grid map.
Definition at line 284 of file COccupancyGridMap2D.h.
References y_max.
Referenced by mrpt::nav::PlannerSimple2D::computePath(), and mrpt::maps::CRandomFieldGridMap2D::internal_clear().
|
inline |
Returns the "y" coordinate of top side of grid map.
Definition at line 282 of file COccupancyGridMap2D.h.
References y_min.
Referenced by mrpt::nav::PlannerSimple2D::computePath(), mrpt::maps::CRandomFieldGridMap2D::internal_clear(), mrpt::slam::CMonteCarloLocalization2D::resetUniformFreeSpace(), mrpt::ros1bridge::toROS(), and mrpt::slam::COccupancyGridMapFeatureExtractor::uncached_extractFeatures().
|
inlinestaticprotected |
Entropy computation internal function:
Definition at line 119 of file COccupancyGridMap2D.h.
Referenced by computeEntropy(), and updateCell().
|
inlineprotectedinherited |
Can be called by a derived class before preparing an event for publishing with publishEvent to determine if there is no one subscribed, so it can save the wasted time preparing an event that will be not read.
Definition at line 53 of file CObservable.h.
References mrpt::system::CObservable::m_subscribers.
Referenced by mrpt::gui::CWindowDialog::OnMouseDown(), mrpt::gui::CWindowDialog::OnMouseMove(), mrpt::gui::CWindowDialog::OnResize(), mrpt::opengl::COpenGLViewport::render(), and mrpt::opengl::COpenGLViewport::renderNormalSceneMode().
|
inline |
Transform a cell index into a coordinate value.
Definition at line 307 of file COccupancyGridMap2D.h.
References resolution, and x_min.
Referenced by mrpt::slam::CGridMapAligner::AlignPDF_correlation(), mrpt::nav::PlannerSimple2D::computePath(), determineMatching2D(), and mrpt::slam::CMonteCarloLocalization2D::resetUniformFreeSpace().
|
inline |
Definition at line 311 of file COccupancyGridMap2D.h.
References resolution, and y_min.
Referenced by mrpt::slam::CGridMapAligner::AlignPDF_correlation(), mrpt::nav::PlannerSimple2D::computePath(), determineMatching2D(), and mrpt::slam::CMonteCarloLocalization2D::resetUniformFreeSpace().
|
inherited |
Insert the observation information into this map.
This method must be implemented in derived classes. See: Maps and observations compatibility matrix
obs | The observation |
robotPose | The 3D pose of the robot mobile base in the map reference system, or NULL (default) if you want to use the origin. |
Definition at line 93 of file CMetricMap.cpp.
Referenced by ICPTests::align2scans(), computeObservationLikelihood_likelihoodField_Thrun(), mrpt::graphslam::deciders::CRangeScanOps< GRAPH_T >::getICPEdge(), CAngularObservationMesh_fnctr::operator()(), TEST(), and mrpt::graphslam::CGraphSlamEngine< GRAPH_T >::updateMapVisualization().
|
inherited |
A wrapper for smart pointers, just calls the non-smart pointer version.
See: Maps and observations compatibility matrix
Definition at line 107 of file CMetricMap.cpp.
References MRPT_END, MRPT_START, and THROW_EXCEPTION.
Referenced by mrpt::slam::CMetricMapBuilderICP::processObservation().
|
overrideprivate |
Returns true if this map is able to compute a sensible likelihood function for this observation (i.e.
an occupancy grid map cannot with an image).
obs | The observation. |
Definition at line 958 of file COccupancyGridMap2D_likelihood.cpp.
References mrpt::maps::COccupancyGridMap2D::TInsertionOptions::horizontalTolerance, insertionOptions, IS_CLASS, mrpt::maps::COccupancyGridMap2D::TInsertionOptions::mapAltitude, and mrpt::maps::COccupancyGridMap2D::TInsertionOptions::useMapAltitude.
|
overrideprotectedvirtual |
Clear the map: It set all cells to their default occupancy value (0.5), without changing the resolution (the grid extension is reset to the default values).
Implements mrpt::maps::CMetricMap.
Definition at line 419 of file COccupancyGridMap2D_common.cpp.
References getResolution(), m_likelihoodCacheOutDated, and setSize().
|
overrideprivatevirtual |
Internal method called by computeObservationLikelihood()
Implements mrpt::maps::CMetricMap.
Definition at line 32 of file COccupancyGridMap2D_likelihood.cpp.
References IS_CLASS.
|
static |
Definition at line 71 of file COccupancyGridMap2D_common.cpp.
References mrpt::maps::COccupancyGridMap2D::TMapDefinition::insertionOpts, mrpt::maps::COccupancyGridMap2D::TMapDefinition::likelihoodOpts, mrpt::maps::COccupancyGridMap2D::TMapDefinition::max_x, mrpt::maps::COccupancyGridMap2D::TMapDefinition::max_y, mrpt::maps::COccupancyGridMap2D::TMapDefinition::min_x, mrpt::maps::COccupancyGridMap2D::TMapDefinition::min_y, and mrpt::maps::COccupancyGridMap2D::TMapDefinition::resolution.
|
overrideprotectedvirtual |
Insert the observation information into this map.
obs | The observation |
robotPose | The 3D pose of the robot mobile base in the map reference system, or nullptr (default) if you want to use CPose2D(0,0,deg) |
After successfull execution, "lastObservationInsertionInfo" is updated.
Implements mrpt::maps::CMetricMap.
Definition at line 43 of file COccupancyGridMap2D_insert.cpp.
References ASSERT_, TLocalPoint::cx, TLocalPoint::cy, mrpt::d2f(), for(), FRBITS, mrpt::poses::CPoseOrPoint< DERIVEDCLASS, DIM >::getHomogeneousMatrixVal(), IS_CLASS, mrpt::max3(), mrpt::min3(), mrpt_alloca, mrpt_alloca_free, MRPT_CHECK_NORMAL_NUMBER, mrpt::poses::CPose2D::phi(), R, mrpt::round(), TLocalPoint::x, mrpt::poses::CPoseOrPoint< DERIVEDCLASS, DIM >::x(), TLocalPoint::y, and mrpt::poses::CPoseOrPoint< DERIVEDCLASS, DIM >::y().
|
overridevirtual |
Returns true upon map construction or after calling clear(), the return changes to false upon successful insertObservation() or any other method to load data in the map.
Implements mrpt::maps::CMetricMap.
Definition at line 727 of file COccupancyGridMap2D_common.cpp.
References m_is_empty.
|
inline |
Definition at line 405 of file COccupancyGridMap2D.h.
References getCell().
Referenced by isStaticPos().
|
inline |
Returns "true" if cell is "static", i.e.if its occupancy is below a given threshold.
Definition at line 401 of file COccupancyGridMap2D.h.
References isStaticCell(), x2idx(), and y2idx().
|
inlinestatic |
Scales an integer representation of the log-odd into a real valued probability in [0,1], using p=exp(l)/(1+exp(l))
Definition at line 329 of file COccupancyGridMap2D.h.
References get_logodd_lut().
Referenced by computeEntropy(), getCell(), getCell_nocheck(), and updateCell().
|
inlinestatic |
Scales an integer representation of the log-odd into a linear scale [0,255], using p=exp(l)/(1+exp(l))
Definition at line 335 of file COccupancyGridMap2D.h.
References get_logodd_lut().
void COccupancyGridMap2D::laserScanSimulator | ( | mrpt::obs::CObservation2DRangeScan & | inout_Scan, |
const mrpt::poses::CPose2D & | robotPose, | ||
float | threshold = 0.6f , |
||
size_t | N = 361 , |
||
float | noiseStd = 0 , |
||
unsigned int | decimation = 1 , |
||
float | angleNoiseStd = mrpt::DEG2RAD(.0) |
||
) | const |
Simulates a laser range scan into the current grid map.
The simulated scan is stored in a CObservation2DRangeScan object, which is also used to pass some parameters: all previously stored characteristics (as aperture,...) are taken into account for simulation. Only a few more parameters are needed. Additive gaussian noise can be optionally added to the simulated scan.
inout_Scan | [IN/OUT] This must be filled with desired parameters before calling, and will contain the scan samples on return. |
robotPose | [IN] The robot pose in this map coordinates. Recall that sensor pose relative to this robot pose must be specified in the observation object. |
threshold | [IN] The minimum occupancy threshold to consider a cell to be occupied (Default: 0.5f) |
N | [IN] The count of range scan "rays", by default to 361. |
noiseStd | [IN] The standard deviation of measurement noise. If not desired, set to 0. |
decimation | [IN] The rays that will be simulated are at indexes: 0, D, 2D, 3D, ... Default is D=1 |
angleNoiseStd | [IN] The sigma of an optional Gaussian noise added to the angles at which ranges are measured (in radians). |
Definition at line 31 of file COccupancyGridMap2D_simulate.cpp.
References mrpt::obs::CObservation2DRangeScan::aperture, ASSERT_, mrpt::obs::CObservation2DRangeScan::maxRange, MRPT_END, MRPT_START, mrpt::poses::CPose2D::phi(), mrpt::obs::CObservation2DRangeScan::resizeScan(), mrpt::obs::CObservation2DRangeScan::rightToLeft, mrpt::obs::CObservation2DRangeScan::sensorPose, mrpt::obs::CObservation2DRangeScan::setScanRange(), mrpt::obs::CObservation2DRangeScan::setScanRangeValidity(), mrpt::poses::CPoseOrPoint< DERIVEDCLASS, DIM >::x(), and mrpt::poses::CPoseOrPoint< DERIVEDCLASS, DIM >::y().
Referenced by run_rnav_test().
void COccupancyGridMap2D::laserScanSimulatorWithUncertainty | ( | const TLaserSimulUncertaintyParams & | in_params, |
COccupancyGridMap2D::TLaserSimulUncertaintyResult & | out_results | ||
) | const |
Like laserScanSimulatorWithUncertainty() (see it for a discussion of most parameters) but taking into account the robot pose uncertainty and generating a vector of predicted variances for each ray.
Range uncertainty includes both, sensor noise and large non-linear effects caused by borders and discontinuities in the environment as seen from different robot poses.
in_params | [IN] Input settings. See TLaserSimulUncertaintyParams |
in_params | [OUT] Output range + uncertainty. |
Definition at line 241 of file COccupancyGridMap2D_simulate.cpp.
References mrpt::obs::CObservation2DRangeScan::aperture, mrpt::maps::COccupancyGridMap2D::TLaserSimulUncertaintyParams::aperture, mrpt::math::CMatrixDynamic< T >::asEigen(), mrpt::poses::CPoseOrPoint< DERIVEDCLASS, DIM >::asVectorVal(), mrpt::poses::CPosePDFGaussian::cov, func_laserSimul_callback(), TFunctorLaserSimulData::grid, mrpt::obs::CObservation2DRangeScan::maxRange, mrpt::maps::COccupancyGridMap2D::TLaserSimulUncertaintyParams::maxRange, mrpt::maps::COccupancyGridMap2D::TLaserSimulUncertaintyParams::MC_samples, mrpt::poses::CPosePDFGaussian::mean, mrpt::maps::COccupancyGridMap2D::TLaserSimulUncertaintyParams::method, mrpt::maps::COccupancyGridMap2D::TLaserSimulUncertaintyParams::nRays, TFunctorLaserSimulData::params, mrpt::obs::CObservation2DRangeScanWithUncertainty::rangeScan, mrpt::obs::CObservation2DRangeScanWithUncertainty::rangesCovar, mrpt::obs::CObservation2DRangeScanWithUncertainty::rangesMean, mrpt::obs::CObservation2DRangeScan::resizeScan(), mrpt::obs::CObservation2DRangeScan::rightToLeft, mrpt::maps::COccupancyGridMap2D::TLaserSimulUncertaintyParams::rightToLeft, mrpt::maps::COccupancyGridMap2D::TLaserSimulUncertaintyParams::robotPose, mrpt::maps::COccupancyGridMap2D::TLaserSimulUncertaintyResult::scanWithUncert, mrpt::obs::CObservation2DRangeScan::sensorPose, mrpt::maps::COccupancyGridMap2D::TLaserSimulUncertaintyParams::sensorPose, mrpt::obs::CObservation2DRangeScan::setScanRange(), mrpt::obs::CObservation2DRangeScan::setScanRangeValidity(), mrpt::math::transform_gaussian_montecarlo(), mrpt::math::transform_gaussian_unscented(), mrpt::maps::COccupancyGridMap2D::TLaserSimulUncertaintyParams::UT_alpha, mrpt::maps::COccupancyGridMap2D::TLaserSimulUncertaintyParams::UT_beta, and mrpt::maps::COccupancyGridMap2D::TLaserSimulUncertaintyParams::UT_kappa.
bool COccupancyGridMap2D::loadFromBitmap | ( | const mrpt::img::CImage & | img, |
float | resolution, | ||
const mrpt::math::TPoint2D & | origin = mrpt::math::TPoint2D( std::numeric_limits<double>::max(), std::numeric_limits<double>::max()) |
||
) |
Load the gridmap from a image in a file (the format can be any supported by CImage::loadFromFile).
See loadFromBitmapFile() for the meaning of parameters
Definition at line 282 of file COccupancyGridMap2D_io.cpp.
References mrpt::img::CImage::getAsFloat(), mrpt::img::CImage::getHeight(), mrpt::img::CImage::getWidth(), MRPT_END, and MRPT_START.
bool COccupancyGridMap2D::loadFromBitmapFile | ( | const std::string & | file, |
float | resolution, | ||
const mrpt::math::TPoint2D & | origin = mrpt::math::TPoint2D( std::numeric_limits<double>::max(), std::numeric_limits<double>::max()) |
||
) |
Load the gridmap from a image in a file (the format can be any supported by CImage::loadFromFile).
file | The file to be loaded. |
resolution | The size of a pixel (cell), in meters. Recall cells are always squared, so just a dimension is needed. |
origin | The x (0=first, increases left to right on the image) and y (0=first, increases BOTTOM upwards on the image) coordinates for the pixel which will be taken at the origin of map coordinates (0,0). (Default=center of the image) |
Definition at line 268 of file COccupancyGridMap2D_io.cpp.
References MRPT_END, and MRPT_START.
Referenced by mrpt::maps::CRandomFieldGridMap2D::internal_clear().
|
inherited |
Load the map contents from a CSimpleMap object, erasing all previous content of the map.
This is done invoking insertObservation()
for each observation at the mean 3D robot pose of each pose-observations pair in the CSimpleMap object.
std::exception | Some internal steps in invoked methods can raise exceptions on invalid parameters, etc... |
Definition at line 36 of file CMetricMap.cpp.
References ASSERTMSG_, mrpt::containers::clear(), mrpt::maps::CSimpleMap::get(), and mrpt::maps::CSimpleMap::size().
Referenced by mrpt::apps::MonteCarloLocalization_Base::do_pf_localization(), mrpt::maps::CMetricMap::loadFromSimpleMap(), mrpt::apps::CGridMapAlignerApp::run(), and run_test_pf_localization().
|
inlineinherited |
!
This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
Definition at line 106 of file CMetricMap.h.
References mrpt::maps::CMetricMap::loadFromProbabilisticPosesAndObservations().
Referenced by mrpt::maps::CRandomFieldGridMap2D::internal_clear(), and mrpt::ros1bridge::MapHdl::loadMap().
|
static |
Returns default map definition initializer.
See * mrpt::maps::TMetricMapInitializer
Definition at line 32 of file COccupancyGridMap2D_common.cpp.
|
overrideprotectedvirtual |
See base class.
Reimplemented from mrpt::maps::CMetricMap.
Definition at line 1369 of file COccupancyGridMap2D_insert.cpp.
|
inlinestatic |
Scales a real valued probability in [0,1] to an integer representation of: log(p)-log(1-p) in the valid range of cellType.
Definition at line 341 of file COccupancyGridMap2D.h.
References get_logodd_lut().
Referenced by determineMatching2D(), fill(), setCell(), setCell_nocheck(), setSize(), subSample(), and updateCell().
|
protectedinherited |
Called when you want this object to emit an event to all the observers currently subscribed to this object.
Definition at line 57 of file CObservable.cpp.
References MRPT_END, and MRPT_START.
Referenced by mrpt::gui::CDisplayWindow3D::internal_emitGrabImageEvent(), mrpt::gui::CWindowDialog::OnChar(), mrpt::gui::CWindowDialog::OnClose(), mrpt::gui::C3DWindowDialog::OnClose(), mrpt::gui::CWindowDialog::OnMouseDown(), mrpt::gui::CWindowDialog::OnMouseMove(), mrpt::gui::CWindowDialog::OnResize(), mrpt::gui::C3DWindowDialog::OnResize(), mrpt::opengl::COpenGLViewport::render(), and mrpt::opengl::COpenGLViewport::renderNormalSceneMode().
|
noexcept |
Change the size of gridmap, maintaining previous contents.
new_x_min | The "x" coordinates of new left most side of grid. |
new_x_max | The "x" coordinates of new right most side of grid. |
new_y_min | The "y" coordinates of new top most side of grid. |
new_y_max | The "y" coordinates of new bottom most side of grid. |
new_cells_default_value | The value of the new cells, tipically 0.5. |
additionalMargin | If set to true (default), an additional margin of a few meters will be added to the grid, ONLY if the new coordinates are larger than current ones. |
Definition at line 196 of file COccupancyGridMap2D_common.cpp.
References mrpt::system::os::memcpy(), and mrpt::round().
Referenced by mrpt::maps::CMultiMetricMapPDF::getCurrentJointEntropy(), and mrpt::maps::CMultiMetricMapPDF::rebuildAverageMap().
bool COccupancyGridMap2D::saveAsBitmapFile | ( | const std::string & | file | ) | const |
Saves the gridmap as a graphical file (BMP,PNG,...).
The format will be derived from the file extension (see CImage::saveToFile )
Definition at line 36 of file COccupancyGridMap2D_io.cpp.
References MRPT_END, and MRPT_START.
Referenced by mrpt::maps::CRandomFieldGridMap2D::internal_clear().
|
inline |
Saves the gridmap as a graphical bitmap file, 8 bit gray scale, 1 pixel is 1 cell, and with an overlay of landmarks.
Definition at line 976 of file COccupancyGridMap2D.h.
References mrpt::img::TColor::black(), mrpt::format(), getAsImageFiltered(), MRPT_END, MRPT_START, size_y, x2idx(), and y2idx().
|
static |
Saves a composite image with two gridmaps and lines representing a set of correspondences between them.
Definition at line 332 of file COccupancyGridMap2D_io.cpp.
References mrpt::img::CH_RGB, mrpt::img::CImage::drawImage(), mrpt::img::CCanvas::filledRectangle(), getAsImage(), mrpt::img::CImage::getHeight(), mrpt::random::getRandomGenerator(), mrpt::img::CImage::getWidth(), mrpt::img::CImage::line(), MRPT_END, MRPT_START, mrpt::img::CCanvas::rectangle(), mrpt::img::CImage::saveToFile(), x2idx(), and y2idx().
|
static |
Saves a composite image with two gridmaps and numbers for the correspondences between them.
/
for (i=0;i<n;i++) { lineColor = ((unsigned long)RandomUni(0,255.0f)) + (((unsigned long)RandomUni(0,255.0f)) << 8 ) + (((unsigned long)RandomUni(0,255.0f)) << 16 );
emf.line( m1->x2idx( corrs[i].this_x ), Ay1+ly1-1- m1->y2idx( corrs[i].this_y ), lx1+1+ m2->x2idx( corrs[i].other_x ), Ay2+ly2-1-m2->y2idx( corrs[i].other_y ), lineColor); } // i /
Definition at line 417 of file COccupancyGridMap2D_io.cpp.
References getAsImage(), mrpt::img::CImage::getHeight(), mrpt::img::CImage::getWidth(), MRPT_END, MRPT_START, mrpt::system::os::sprintf(), x2idx(), and y2idx().
|
overridevirtual |
This virtual method saves the map to a file "filNamePrefix"+< some_file_extension >, as an image or in any other applicable way (Notice that other methods to save the map may be implemented in classes implementing this virtual interface).
Implements mrpt::maps::CMetricMap.
Definition at line 523 of file COccupancyGridMap2D_io.cpp.
References mrpt::math::MATRIX_FORMAT_FIXED, and mrpt::math::MatrixVectorBase< Scalar, Derived >::saveToTextFile().
Referenced by mrpt::slam::COccupancyGridMapFeatureExtractor::uncached_extractFeatures().
|
overrideprotectedvirtual |
Pure virtual method for reading (deserializing) from an abstract archive.
Users don't call this method directly. Instead, use stream >> object;
.
in | The input binary stream where the object data must read from. |
version | The version of the object stored in the stream: use this version number in your code to know how to read the incoming data. |
std::exception | On any I/O error |
Implements mrpt::serialization::CSerializable.
Definition at line 110 of file COccupancyGridMap2D_io.cpp.
References ASSERT_, MRPT_THROW_UNKNOWN_SERIALIZATION_VERSION, mrpt::serialization::CArchive::ReadBuffer(), and mrpt::serialization::CArchive::ReadBufferFixEndianness().
|
inlineprotectedvirtualinherited |
Virtual method for reading (deserializing) from an abstract schema based archive.
Definition at line 74 of file CSerializable.h.
References mrpt::serialization::CSerializable::GetRuntimeClass(), and THROW_EXCEPTION.
|
overrideprotectedvirtual |
Must return the current versioning number of the object.
Start in zero for new classes, and increments each time there is a change in the stored format.
Implements mrpt::serialization::CSerializable.
Definition at line 53 of file COccupancyGridMap2D_io.cpp.
|
overrideprotectedvirtual |
Pure virtual method for writing (serializing) to an abstract archive.
Users don't call this method directly. Instead, use stream << object;
.
std::exception | On any I/O error |
Implements mrpt::serialization::CSerializable.
Definition at line 54 of file COccupancyGridMap2D_io.cpp.
|
inlineprotectedvirtualinherited |
Virtual method for writing (serializing) to an abstract schema based archive.
Definition at line 64 of file CSerializable.h.
References mrpt::serialization::CSerializable::GetRuntimeClass(), and THROW_EXCEPTION.
|
inline |
Change a cell in the "basis" maps.Used for Voronoi calculation.
Definition at line 411 of file COccupancyGridMap2D.h.
References ASSERT_ABOVEEQ_, ASSERT_BELOWEQ_, mrpt::containers::CDynamicGrid< T >::cellByIndex(), mrpt::containers::CDynamicGrid< T >::getSizeX(), mrpt::containers::CDynamicGrid< T >::getSizeY(), and m_basis_map.
|
inline |
Change the contents [0,1] of a cell, given its index.
Definition at line 346 of file COccupancyGridMap2D.h.
References map, p2l(), size_x, and size_y.
Referenced by run_rnav_test(), and setPos().
|
inlineprotected |
Change the contents [0,1] of a cell, given its index.
Definition at line 130 of file COccupancyGridMap2D.h.
References map, p2l(), and size_x.
|
inline |
Change the contents [0,1] of a cell, given its coordinates.
Definition at line 388 of file COccupancyGridMap2D.h.
References setCell(), x2idx(), and y2idx().
|
inlineprotected |
Changes a cell by its absolute index (Do not use it normally)
Definition at line 141 of file COccupancyGridMap2D.h.
void COccupancyGridMap2D::setSize | ( | float | x_min, |
float | x_max, | ||
float | y_min, | ||
float | y_max, | ||
float | resolution, | ||
float | default_value = 0.5f |
||
) |
Change the size of gridmap, erasing all its previous contents.
x_min | The "x" coordinates of left most side of grid. |
x_max | The "x" coordinates of right most side of grid. |
y_min | The "y" coordinates of top most side of grid. |
y_max | The "y" coordinates of bottom most side of grid. |
resolution | The new size of cells. |
default_value | The value of cells, tipically 0.5. |
Definition at line 140 of file COccupancyGridMap2D_common.cpp.
References ASSERT_, ASSERT_ABOVE_, ASSERT_ABOVEEQ_, ASSERT_BELOWEQ_, mrpt::containers::CDynamicGrid< T >::clear(), freeMap(), m_basis_map, m_is_empty, m_likelihoodCacheOutDated, m_voronoi_diagram, map, MRPT_END, MRPT_START, p2l(), resolution, mrpt::round(), size_x, size_y, x_max, x_min, y_max, and y_min.
Referenced by mrpt::slam::CGridMapAligner::AlignPDF_correlation(), COccupancyGridMap2D(), mrpt::ros1bridge::fromROS(), internal_clear(), mrpt::maps::CMultiMetricMapPDF::rebuildAverageMap(), run_rnav_test(), and subSample().
|
inlineprotected |
Used to set the clearance of a cell, while building the Voronoi diagram.
Definition at line 688 of file COccupancyGridMap2D.h.
References ASSERT_ABOVEEQ_, ASSERT_BELOWEQ_, mrpt::containers::CDynamicGrid< T >::cellByIndex(), mrpt::containers::CDynamicGrid< T >::getSizeX(), mrpt::containers::CDynamicGrid< T >::getSizeY(), and m_voronoi_diagram.
void COccupancyGridMap2D::simulateScanRay | ( | const double | x, |
const double | y, | ||
const double | angle_direction, | ||
float & | out_range, | ||
bool & | out_valid, | ||
const double | max_range_meters, | ||
const float | threshold_free = 0.4f , |
||
const double | noiseStd = .0 , |
||
const double | angleNoiseStd = .0 |
||
) | const |
Simulate just one "ray" in the grid map.
This method is used internally to sonarSimulator and laserScanSimulator.
Definition at line 110 of file COccupancyGridMap2D_simulate.cpp.
References mrpt::random::CRandomGenerator::drawGaussian1D_normalized(), mrpt::random::getRandomGenerator(), int_x2idx, int_y2idx, INTPRECNUMBIT, and mrpt::round().
Referenced by func_laserSimul_callback().
void COccupancyGridMap2D::sonarSimulator | ( | mrpt::obs::CObservationRange & | inout_observation, |
const mrpt::poses::CPose2D & | robotPose, | ||
float | threshold = 0.5f , |
||
float | rangeNoiseStd = 0.f , |
||
float | angleNoiseStd = mrpt::DEG2RAD(0.f) |
||
) | const |
Simulates the observations of a sonar rig into the current grid map.
The simulated ranges are stored in a CObservationRange object, which is also used to pass in some needed parameters, as the poses of the sonar sensors onto the mobile robot.
inout_observation | [IN/OUT] This must be filled with desired parameters before calling, and will contain the simulated ranges on return. |
robotPose | [IN] The robot pose in this map coordinates. Recall that sensor pose relative to this robot pose must be specified in the observation object. |
threshold | [IN] The minimum occupancy threshold to consider a cell to be occupied (Default: 0.5f) |
rangeNoiseStd | [IN] The standard deviation of measurement noise. If not desired, set to 0. |
angleNoiseStd | [IN] The sigma of an optional Gaussian noise added to the angles at which ranges are measured (in radians). |
Definition at line 70 of file COccupancyGridMap2D_simulate.cpp.
References ASSERT_, mrpt::obs::CObservationRange::begin(), mrpt::obs::CObservationRange::end(), mrpt::obs::CObservationRange::maxSensorDistance, mrpt::poses::CPose2D::phi(), mrpt::round(), mrpt::obs::CObservationRange::sensorConeApperture, mrpt::poses::CPoseOrPoint< DERIVEDCLASS, DIM >::x(), and mrpt::poses::CPoseOrPoint< DERIVEDCLASS, DIM >::y().
|
virtualinherited |
Returns the square distance from the 2D point (x0,y0) to the closest correspondence in the map.
Reimplemented in mrpt::maps::CPointsMap.
Definition at line 153 of file CMetricMap.cpp.
References MRPT_END, MRPT_START, and THROW_EXCEPTION.
void COccupancyGridMap2D::subSample | ( | int | downRatio | ) |
Performs a downsampling of the gridmap, by a given factor: resolution/=ratio.
Definition at line 482 of file COccupancyGridMap2D_common.cpp.
References ASSERT_, getCell(), map, p2l(), resolution, mrpt::round(), setSize(), x_max, x_min, y_max, and y_min.
void COccupancyGridMap2D::updateCell | ( | int | x, |
int | y, | ||
float | v | ||
) |
Performs the Bayesian fusion of a new observation of a cell.
Definition at line 440 of file COccupancyGridMap2D_common.cpp.
References mrpt::maps::COccupancyGridMap2D::TUpdateCellsInfoChangeOnly::cellsUpdated, mrpt::maps::COccupancyGridMap2D::TUpdateCellsInfoChangeOnly::enabled, H(), mrpt::maps::COccupancyGridMap2D::TUpdateCellsInfoChangeOnly::I_change, l2p(), map, MAX_H, OCCGRID_CELLTYPE_MAX, OCCGRID_CELLTYPE_MIN, p2l(), size_x, size_y, and updateInfoChangeOnly.
|
inlinestaticinherited |
Performs the Bayesian fusion of a new observation of a cell, without checking for grid limits nor updateInfoChangeOnly.
This method increases the "free-ness" of a cell, managing possible saturation.
x | Cell index in X axis. |
y | Cell index in Y axis. |
logodd_obs | Observation of the cell, in log-odd form as transformed by p2l. |
thres | This must be CELLTYPE_MAX-logodd_obs |
Definition at line 84 of file CLogOddsGridMap2D.h.
|
inlinestaticinherited |
Performs the Bayesian fusion of a new observation of a cell, without checking for grid limits nor updateInfoChangeOnly.
This method increases the "free-ness" of a cell, managing possible saturation.
x | Cell index in X axis. |
y | Cell index in Y axis. |
logodd_obs | Observation of the cell, in log-odd form as transformed by p2l. |
thres | This must be CELLTYPE_MAX-logodd_obs |
Definition at line 106 of file CLogOddsGridMap2D.h.
|
inlinestaticinherited |
Performs the Bayesian fusion of a new observation of a cell, without checking for grid limits nor updateInfoChangeOnly.
This method increases the "occupancy-ness" of a cell, managing possible saturation.
x | Cell index in X axis. |
y | Cell index in Y axis. |
logodd_obs | Observation of the cell, in log-odd form as transformed by p2l. |
thres | This must be CELLTYPE_MIN+logodd_obs |
Definition at line 43 of file CLogOddsGridMap2D.h.
|
inlinestaticinherited |
Performs the Bayesian fusion of a new observation of a cell, without checking for grid limits nor updateInfoChangeOnly.
This method increases the "occupancy-ness" of a cell, managing possible saturation.
theCell | The cell to modify |
logodd_obs | Observation of the cell, in log-odd form as transformed by p2l. |
thres | This must be CELLTYPE_MIN+logodd_obs |
Definition at line 64 of file CLogOddsGridMap2D.h.
|
inlinevirtualinherited |
Introduces a pure virtual method responsible for writing to a mxArray
Matlab object, typically a MATLAB struct
whose contents are documented in each derived class.
mxArray
(caller is responsible of memory freeing) or nullptr is class does not support conversion to MATLAB. Definition at line 90 of file CSerializable.h.
|
inline |
Transform a coordinate value into a cell index.
Definition at line 288 of file COccupancyGridMap2D.h.
References resolution, and x_min.
Referenced by mrpt::nav::PlannerSimple2D::computePath(), determineMatching2D(), getPos(), isStaticPos(), mrpt::slam::CMonteCarloLocalization2D::resetUniformFreeSpace(), run_rnav_test(), saveAsBitmapFileWithLandmarks(), saveAsBitmapTwoMapsWithCorrespondences(), saveAsEMFTwoMapsWithCorrespondences(), and setPos().
|
inline |
Definition at line 297 of file COccupancyGridMap2D.h.
References resolution, and x_min.
|
inline |
Transform a coordinate value into a cell index, using a diferent "x_min" value.
Definition at line 318 of file COccupancyGridMap2D.h.
References resolution.
|
inline |
Definition at line 292 of file COccupancyGridMap2D.h.
References resolution, and y_min.
Referenced by mrpt::nav::PlannerSimple2D::computePath(), determineMatching2D(), getPos(), isStaticPos(), mrpt::slam::CMonteCarloLocalization2D::resetUniformFreeSpace(), run_rnav_test(), saveAsBitmapFileWithLandmarks(), saveAsBitmapTwoMapsWithCorrespondences(), saveAsEMFTwoMapsWithCorrespondences(), and setPos().
|
inline |
Definition at line 301 of file COccupancyGridMap2D.h.
References resolution, and y_min.
|
inline |
Definition at line 322 of file COccupancyGridMap2D.h.
References resolution.
|
friend |
Definition at line 75 of file COccupancyGridMap2D.h.
|
friend |
Definition at line 76 of file COccupancyGridMap2D.h.
|
static |
Definition at line 57 of file COccupancyGridMap2D.h.
struct mrpt::maps::COccupancyGridMap2D::TCriticalPointsList mrpt::maps::COccupancyGridMap2D::CriticalPointsList |
|
private |
Used to store the 8 possible movements from a cell to the sorrounding ones.Filled in the constructor.
Definition at line 1141 of file COccupancyGridMap2D.h.
|
private |
Definition at line 1141 of file COccupancyGridMap2D.h.
|
staticprotected |
Internally used to speed-up entropy calculation.
Definition at line 127 of file COccupancyGridMap2D.h.
Referenced by computeEntropy().
|
inherited |
Common params to all maps.
Definition at line 274 of file CMetricMap.h.
Referenced by mrpt::maps::internal::TMetricMapTypesRegistry::factoryMapObjectFromDefinition(), mrpt::maps::CWirelessPowerGridMap2D::getAs3DObject(), mrpt::maps::CGasConcentrationGridMap2D::getAs3DObject(), mrpt::maps::CReflectivityGridMap2D::getAs3DObject(), mrpt::maps::CHeightGridMap2D::getAs3DObject(), mrpt::maps::CRandomFieldGridMap2D::getAs3DObject(), mrpt::maps::CLandmarksMap::getAs3DObject(), mrpt::maps::CPointsMap::getAs3DObject(), mrpt::maps::CGasConcentrationGridMap2D::serializeFrom(), mrpt::maps::CWirelessPowerGridMap2D::serializeFrom(), mrpt::maps::CReflectivityGridMap2D::serializeFrom(), mrpt::maps::CHeightGridMap2D::serializeFrom(), mrpt::maps::CWirelessPowerGridMap2D::serializeTo(), mrpt::maps::CGasConcentrationGridMap2D::serializeTo(), mrpt::maps::CReflectivityGridMap2D::serializeTo(), and mrpt::maps::CHeightGridMap2D::serializeTo().
TInsertionOptions mrpt::maps::COccupancyGridMap2D::insertionOptions |
With this struct options are provided to the observation insertion process.
Definition at line 530 of file COccupancyGridMap2D.h.
Referenced by internal_canComputeObservationLikelihood(), and mrpt::slam::CMetricMapBuilderICP::processObservation().
mrpt::maps::COccupancyGridMap2D::TLikelihoodOptions mrpt::maps::COccupancyGridMap2D::likelihoodOptions |
class mrpt::maps::COccupancyGridMap2D::TLikelihoodOutput mrpt::maps::COccupancyGridMap2D::likelihoodOutputs |
|
protected |
Used for Voronoi calculation.Same struct as "map", but contains a "0" if not a basis point.
Definition at line 100 of file COccupancyGridMap2D.h.
Referenced by copyMapContentFrom(), freeMap(), getBasisCell(), getBasisMap(), setBasisCell(), and setSize().
|
protected |
True upon construction; used by isEmpty()
Definition at line 109 of file COccupancyGridMap2D.h.
Referenced by copyMapContentFrom(), freeMap(), isEmpty(), and setSize().
|
protected |
Definition at line 96 of file COccupancyGridMap2D.h.
Referenced by copyMapContentFrom(), fill(), freeMap(), internal_clear(), and setSize().
|
static |
ID used to initialize class registration (just ignore it)
Definition at line 1157 of file COccupancyGridMap2D.h.
|
protected |
Used to store the Voronoi diagram.
Contains the distance of each cell to its closer obstacles in 1/100th distance units (i.e. in centimeters), or 0 if not into the Voronoi diagram
Definition at line 106 of file COccupancyGridMap2D.h.
Referenced by copyMapContentFrom(), freeMap(), getVoroniClearance(), getVoronoiDiagram(), setSize(), and setVoroniClearance().
|
protected |
Store of cell occupancy values.
Order: row by row, from left to right
Definition at line 84 of file COccupancyGridMap2D.h.
Referenced by computeEntropy(), copyMapContentFrom(), determineMatching2D(), fill(), freeMap(), getCell(), getCell_nocheck(), getRawMap(), getRow(), mrpt::maps::CMultiMetricMapPDF::rebuildAverageMap(), setCell(), setCell_nocheck(), setRawCell(), setSize(), subSample(), and updateCell().
|
static |
Definition at line 67 of file COccupancyGridMap2D.h.
Referenced by updateCell().
|
static |
Discrete to float conversion factors: The min/max values of the integer cell type, eg.
[0,255] or [0,65535]
Definition at line 65 of file COccupancyGridMap2D.h.
Referenced by updateCell().
|
protected |
Auxiliary variables to speed up the computation of observation likelihood values for LF method among others, at a high cost in memory (see TLikelihoodOptions::enableLikelihoodCache).
Definition at line 95 of file COccupancyGridMap2D.h.
|
static |
(Default:1.0) Can be set to <1 if a more fine raytracing is needed in sonarSimulator() and laserScanSimulator(), or >1 to speed it up.
Definition at line 72 of file COccupancyGridMap2D.h.
|
protected |
Cell size, i.e.
resolution of the grid map.
Definition at line 90 of file COccupancyGridMap2D.h.
Referenced by computeEntropy(), computePathCost(), copyMapContentFrom(), determineMatching2D(), getArea(), getResolution(), idx2x(), idx2y(), setSize(), subSample(), x2idx(), and y2idx().
|
staticprotected |
Definition at line 57 of file COccupancyGridMap2D.h.
|
protected |
The size of the grid in cells.
Definition at line 86 of file COccupancyGridMap2D.h.
Referenced by copyMapContentFrom(), determineMatching2D(), freeMap(), getArea(), getCell(), getCell_nocheck(), getRow(), getSizeX(), setCell(), setCell_nocheck(), setRawCell(), setSize(), and updateCell().
|
protected |
Definition at line 86 of file COccupancyGridMap2D.h.
Referenced by copyMapContentFrom(), determineMatching2D(), freeMap(), getArea(), getCell(), getRow(), getSizeY(), saveAsBitmapFileWithLandmarks(), setCell(), setRawCell(), setSize(), and updateCell().
struct mrpt::maps::COccupancyGridMap2D::TUpdateCellsInfoChangeOnly mrpt::maps::COccupancyGridMap2D::updateInfoChangeOnly |
Referenced by updateCell().
|
protected |
The free-cells threshold used to compute the Voronoi diagram.
Definition at line 115 of file COccupancyGridMap2D.h.
|
protected |
Definition at line 88 of file COccupancyGridMap2D.h.
Referenced by copyMapContentFrom(), determineMatching2D(), getXMax(), setSize(), and subSample().
|
protected |
The limits of the grid in "units" (meters)
Definition at line 88 of file COccupancyGridMap2D.h.
Referenced by copyMapContentFrom(), getXMin(), idx2x(), setSize(), subSample(), and x2idx().
|
protected |
Definition at line 88 of file COccupancyGridMap2D.h.
Referenced by copyMapContentFrom(), determineMatching2D(), getYMax(), setSize(), and subSample().
|
protected |
Definition at line 88 of file COccupancyGridMap2D.h.
Referenced by copyMapContentFrom(), determineMatching2D(), getYMin(), idx2y(), setSize(), subSample(), and y2idx().
Page generated by Doxygen 1.8.14 for MRPT 2.0.1 Git: 0fef1a6d7 Fri Apr 3 23:00:21 2020 +0200 at vie abr 3 23:20:28 CEST 2020 |