IMP logo
IMP Reference Guide  develop.d97d4ead1f,2024/11/21
The Integrative Modeling Platform
IMP::bff::PathMap Class Reference

Class to search path on grids. More...

#include <IMP/bff/PathMap.h>

+ Inheritance diagram for IMP::bff::PathMap:

Detailed Description

Class to search path on grids.

Definition at line 45 of file PathMap.h.

Public Member Functions

 PathMap (PathMapHeader &header, std::string name="PathMap%1%", IMP::em::KernelType kt=IMP::em::BINARIZED_SPHERE, float resolution=-1.0)
 Constructs a PathMap object. More...
 
void fill_sphere (IMP::algebra::Vector3D r0, double radius, double value, bool inverse=true)
 Change the value of a density inside or outside of a sphere. More...
 
void find_path (long path_begin_idx, long path_end_idx=-1, int heuristic_mode=0)
 Finds a path between two indices in the path map. More...
 
void find_path_astar (long path_begin_idx, long path_end_idx=-1)
 Finds the shortest path between two indices using the A* algorithm. More...
 
void find_path_dijkstra (long path_begin_idx, long path_end_idx=-1)
 Finds the shortest path between two nodes using Dijkstra's algorithm. More...
 
int get_dim_index_by_voxel (long index, int dim)
 Get index of voxel in an axis. More...
 
std::vector< int > get_neighbor_idx_offsets (double neighbor_radius=-1)
 
const PathMapHeaderget_path_map_header () const
 Returns a read-only pointer to the header of the map. More...
 
PathMapHeaderget_path_map_header_writable ()
 Returns a pointer to the header of the map in a writable version. More...
 
void get_tile_values (float **output, int *nx, int *ny, int *nz, int value_type=PM_TILE_COST, std::pair< float, float > bounds=std::pair< float, float >({std::numeric_limits< float >::min(), std::numeric_limits< float >::max()}), const std::string &feature_name="")
 Retrieves the values of the tiles in the path map. More...
 
std::vector< PathMapTile > & get_tiles ()
 
std::vector
< IMP::algebra::Vector4D
get_xyz_density ()
 Get the XYZ density of the path map. This function returns a vector of IMP::algebra::Vector4D objects representing the XYZ density of the path map. More...
 
void get_xyz_density (double **output, int *n_output1, int *n_output2)
 Get the XYZ density of the path map. This function returns the XYZ density of the path map as a 2D array. More...
 
void resize (unsigned int nvox)
 Resizes the PathMap object. More...
 
void sample_obstacles (double extra_radius=0.0)
 Resamples the obstacles in the path map. More...
 
void set_data (double *input, int n_input, float obstacle_threshold=-1, bool binarize=true, float obstacle_penalty=TILE_PENALTY_DEFAULT)
 Sets the data for the path map. More...
 
void set_path_map_header (PathMapHeader &path_map_header, float resolution=-1.0)
 Set the path map header. More...
 
void update_tiles (float obstacle_threshold=-1.0, bool binarize=true, float obstacle_penalty=TILE_PENALTY_DEFAULT, bool reset_tile_edges=true)
 Updates the tiles in the path map. More...
 
- Public Member Functions inherited from IMP::em::SampledDensityMap
 SampledDensityMap (KernelType kt=GAUSSIAN)
 Creates a new density map for sampled map. More...
 
 SampledDensityMap (const DensityHeader &header, KernelType kt=GAUSSIAN)
 The size of the map is determined by the header and the data is allocated. More...
 
 SampledDensityMap (const ParticlesTemp &ps, double resolution, double voxel_size, IMP::FloatKey mass_key=IMP::atom::Mass::get_mass_key(), int sig_cutoff=3, KernelType kt=GAUSSIAN)
 Generate a sampled density map from the particles. More...
 
KernelParametersget_kernel_params ()
 
float get_minimum_resampled_value ()
 Get minimum density value between voxels that correspond to particles. More...
 
const Particlesget_sampled_particles () const
 
FloatKey get_weight_key () const
 
const core::XYZRsget_xyzr_particles () const
 
void project (const ParticlesTemp &ps, int x_margin, int y_margin, int z_margin, algebra::Vector3D shift=algebra::Vector3D(0., 0., 0.), FloatKey mass_key=atom::Mass::get_mass_key())
 Project particles on the grid by their mass value. More...
 
virtual void resample ()
 Resample beads on an EM grid. More...
 
void set_particles (const ParticlesTemp &ps, IMP::FloatKey mass_key=IMP::atom::Mass::get_mass_key())
 setting particles in case they were not set by the constructor More...
 
void update_resolution (Float res)
 Update the simulation resolution. More...
 
- Public Member Functions inherited from IMP::em::DensityMap
 DensityMap (std::string name="DensityMap%1%")
 
 DensityMap (const DensityHeader &header, std::string name="DensityMap%1%")
 Construct a density map as instructed in the input header. More...
 
void add (const DensityMap *other)
 Sums two grids; the result is kept in the map. More...
 
void add (Float d)
 Add a number to every voxel in the map. More...
 
void calc_all_voxel2loc ()
 Calculates the coordinates that correspond to all voxels. More...
 
double calcRMS ()
 Calculates RMSD and mean of a map values as stored in the header. More...
 
void copy_map (const DensityMap *other)
 Copy map into this map. More...
 
algebra::Vector3D get_centroid (double threshold=0.0) const
 Get the centroid of all voxels with density above a given threshold. More...
 
DensityMapget_cropped (float threshold)
 Create and return a new cropped map. More...
 
DensityMapget_cropped (Particles ps, double distance, bool inverse=false, bool keep_map_dimensions=false)
 Create and return a new cropped map based on particles. More...
 
DensityMapget_cropped (const algebra::BoundingBox3D &bb)
 Create and return a new cropped map with the given bounding box extent. More...
 
int get_dim_index_by_location (const algebra::Vector3D &v, int ind) const
 Calculate dimension index of a given location. More...
 
int get_dim_index_by_location (float loc_val, int ind) const
 
const DensityHeaderget_header () const
 Returns a read-only pointer to the header of the map. More...
 
DensityHeaderget_header_writable ()
 Returns a pointer to the header of the map in a writable version. More...
 
algebra::Vector3D get_location_by_voxel (long index) const
 Calculate the location of a given voxel. More...
 
float get_location_in_dim_by_voxel (long index, int dim) const
 Calculate the location of a given voxel in a given dimension. More...
 
std::string get_locations_string (float t)
 Print the locations of all voxels with value above a given threshold. More...
 
double get_max_value () const
 Returns the value of the voxel with the highest density. More...
 
float get_maximum_value_in_xy_plane (int z_ind)
 Get the maximum value in a XY plane indicated by a Z index. More...
 
float get_maximum_value_in_xz_plane (int y_ind)
 Get the maximum value in a XZ plane indicated by a Y index. More...
 
float get_maximum_value_in_yz_plane (int x_ind)
 Get the maximum value in a YZ plane indicated by a X index. More...
 
double get_min_value () const
 Returns the value of the voxel with the lowest density. More...
 
long get_number_of_voxels () const
 Get the number of map voxels. More...
 
algebra::Vector3D get_origin () const
 
bool get_rms_calculated () const
 
Float get_spacing () const
 Get the voxel size of the map. More...
 
algebra::Vector3D get_top () const
 
virtual std::string get_type_name () const override
 
double get_value (float x, float y, float z) const
 Gets the value of the voxel located at (x,y,z) More...
 
double get_value (const algebra::Vector3D &point) const
 
double get_value (long index) const
 Gets the value of the voxel at a given index. More...
 
virtual ::IMP::VersionInfo get_version_info () const override
 Get information about the module and version of the object. More...
 
long get_voxel_by_location (float x, float y, float z) const
 Calculate the voxel of a given location. More...
 
long get_voxel_by_location (const algebra::Vector3D &v) const
 Calculate the voxel of a given location. More...
 
bool is_normalized () const
 
bool is_part_of_volume (float x, float y, float z) const
 Checks whether a given point is in the grid. More...
 
bool is_part_of_volume (const algebra::Vector3D &v) const
 Checks whether a given point is in the grid. More...
 
bool is_xyz_ind_part_of_volume (int ix, int iy, int iz) const
 
int lower_voxel_shift (double loc, double kdist, double orig, int ndim) const
 
void multiply (float factor)
 Multiply in place each voxel in the map by the input factor. More...
 
void pad (int nx, int ny, int nz, float val=0.0)
 Increase the dimension of the map. More...
 
DensityMappad_margin (int mrg_x, int mrg_y, int mrg_z, float val=0.0)
 Create a new padded map. More...
 
void pick_max (const DensityMap *other)
 Pick the max value between two corresponding voxels between two maps. More...
 
void reset_data (float value=0.0)
 Set the density voxels to some value and reset the management flags. More...
 
bool same_dimensions (const DensityMap *other) const
 Checks if two maps have the same dimensions. More...
 
bool same_origin (const DensityMap *other) const
 Checks if two maps have the same origin. More...
 
bool same_voxel_size (const DensityMap *other) const
 Checks if two maps have the same voxel size. More...
 
void set_origin (float x, float y, float z)
 Sets the origin of the header. More...
 
void set_origin (const IMP::algebra::Vector3D &v)
 
void set_value (long index, double value)
 Set the value of the voxel at a given index. More...
 
void set_value (float x, float y, float z, double value)
 Set the value of the voxel at given coordinates. More...
 
void set_void_map (int nx, int ny, int nz)
 Set the map dimension and reset all voxels to 0. More...
 
void std_normalize ()
 Normalize the density voxels according to standard deviation (stdv). More...
 
void update_voxel_size (float new_apix)
 Update the voxel size of the map. More...
 
int upper_voxel_shift (double loc, double kdist, double orig, int ndim) const
 
long xyz_ind2voxel (int x, int y, int z) const
 Calculate the voxel of a given xyz index. More...
 
- Public Member Functions inherited from IMP::Object
virtual void clear_caches ()
 
CheckLevel get_check_level () const
 
LogLevel get_log_level () const
 
void set_check_level (CheckLevel l)
 
void set_log_level (LogLevel l)
 Set the logging level used in this object. More...
 
void set_was_used (bool tf) const
 
void show (std::ostream &out=std::cout) const
 
const std::string & get_name () const
 
void set_name (std::string name)
 

Protected Member Functions

std::vector< PathMapTileEdge > & get_edges (int tile_idx)
 
- Protected Member Functions inherited from IMP::em::SampledDensityMap
IMP::algebra::BoundingBoxD< 3 > calculate_particles_bounding_box (const Particles &ps)
 Calculate the parameters of the particles bounding box. More...
 
void determine_grid_size (double resolution, double voxel_size, int sig_cutoff)
 
void set_header (const algebra::Vector3D &lower_bound, const algebra::Vector3D &upper_bound, double maxradius, double resolution, double voxel_size, int sig_offset)
 
void set_neighbor_mask (float radius)
 
- Protected Member Functions inherited from IMP::em::DensityMap
void allocated_data ()
 
void float2real (float *f_data, boost::scoped_array< double > &r_data)
 
void real2float (double *r_data, boost::scoped_array< float > &f_data)
 
void reset_all_voxel2loc ()
 
void update_header ()
 update the header values – still in work More...
 
- Protected Member Functions inherited from IMP::Object
 Object (std::string name)
 Construct an object with the given name. More...
 
virtual void do_destroy ()
 

Protected Attributes

std::vector< int > offsets_
 
PathMapHeader pathMapHeader_
 
std::vector< PathMapTiletiles
 
- Protected Attributes inherited from IMP::em::SampledDensityMap
KernelParameters kernel_params_
 kernel handling More...
 
KernelType kt_
 
Particles ps_
 
FloatKey weight_key_
 
FloatKey x_key_
 
core::XYZRs xyzr_
 
FloatKey y_key_
 
FloatKey z_key_
 
- Protected Attributes inherited from IMP::em::DensityMap
boost::scoped_array< double > data_
 
bool data_allocated_
 
DensityHeader header_
 
bool loc_calculated_
 true if the locations have already been computed More...
 
bool normalized_
 
bool rms_calculated_
 
boost::scoped_array< float > x_loc_
 
boost::scoped_array< float > y_loc_
 
boost::scoped_array< float > z_loc_
 

Related Functions

(Note that these are not member functions.)

std::vector< float > get_tile_values (int value_type=PM_TILE_COST, std::pair< float, float > bounds=std::pair< float, float >({std::numeric_limits< float >::min(), std::numeric_limits< float >::max()}), const std::string &feature_name="")
 Get the values of all tiles. More...
 

Constructor & Destructor Documentation

IMP::bff::PathMap::PathMap ( PathMapHeader header,
std::string  name = "PathMap%1%",
IMP::em::KernelType  kt = IMP::em::BINARIZED_SPHERE,
float  resolution = -1.0 
)
explicit

Constructs a PathMap object.

Parameters
headerThe PathMapHeader object.
nameThe name of the PathMap.
ktThe kernel type.
resolutionThe resolution of the PathMap.

Member Function Documentation

void IMP::bff::PathMap::fill_sphere ( IMP::algebra::Vector3D  r0,
double  radius,
double  value,
bool  inverse = true 
)

Change the value of a density inside or outside of a sphere.

Changes the value of the density inside or outside of a sphere. The density is used in the path-search (i.e., the accessible volume calculation

Parameters
r0location of the sphere
radiusradius of the sphere
valuevalue inside or outside of the sphere
inverseif set to true (default) the values outside of the sphere are modified. If false the values inside of the sphere are modified.
void IMP::bff::PathMap::find_path ( long  path_begin_idx,
long  path_end_idx = -1,
int  heuristic_mode = 0 
)

Finds a path between two indices in the path map.

This function finds a path between the specified path begin index and path end index in the path map. If the path end index is not specified, the function will find a path from the path begin index to the last index in the path map. The heuristic mode parameter determines the heuristic function to be used for path finding.

Parameters
path_begin_idxThe index of the path begin point in the path map.
path_end_idxThe index of the path end point in the path map. Default value is -1, which means the last index in the path map.
heuristic_modeThe mode of the heuristic function to be used for path finding. Default value is 0.
Returns
void
void IMP::bff::PathMap::find_path_astar ( long  path_begin_idx,
long  path_end_idx = -1 
)

Finds the shortest path between two indices using the A* algorithm.

This function uses the A* algorithm to find the shortest path between two indices in the path map. The path is stored in the path member variable.

Parameters
path_begin_idxThe index of the starting point of the path.
path_end_idxThe index of the ending point of the path. If not provided, the function will use the last index in the path map.
void IMP::bff::PathMap::find_path_dijkstra ( long  path_begin_idx,
long  path_end_idx = -1 
)

Finds the shortest path between two nodes using Dijkstra's algorithm.

This function finds the shortest path between two nodes in the path map using Dijkstra's algorithm. The path is calculated from the node at index path_begin_idx to the node at index path_end_idx. If path_end_idx is not provided, the function will calculate the path to the last node in the map.

Parameters
path_begin_idxThe index of the starting node.
path_end_idxThe index of the ending node (optional).
int IMP::bff::PathMap::get_dim_index_by_voxel ( long  index,
int  dim 
)

Get index of voxel in an axis.

Returns the index of a voxel on a grid in a certain dimension.

Parameters
indexvoxel index
dimdimension
Returns
index of voxel in axis of dimension
std::vector<int> IMP::bff::PathMap::get_neighbor_idx_offsets ( double  neighbor_radius = -1)

Returns a vector of neighbor index offsets within a given radius.

Parameters
neighbor_radiusThe radius within which to find neighbors. If negative, the radius is obtained from the path map header.
Returns
A vector of neighbor index offsets, where each offset consists of three integers (z, y, x) representing the relative position of the neighbor, and one integer representing the tile offset. The vector also includes the edge cost between the current tile and the neighbor.

Definition at line 121 of file PathMap.h.

+ Here is the call graph for this function:

const PathMapHeader* IMP::bff::PathMap::get_path_map_header ( ) const

Returns a read-only pointer to the header of the map.

Returns
const PathMapHeader* A read-only pointer to the header of the map.

Definition at line 177 of file PathMap.h.

PathMapHeader* IMP::bff::PathMap::get_path_map_header_writable ( )

Returns a pointer to the header of the map in a writable version.

Returns
A pointer to the header of the map in a writable version.

Definition at line 186 of file PathMap.h.

void IMP::bff::PathMap::get_tile_values ( float **  output,
int *  nx,
int *  ny,
int *  nz,
int  value_type = PM_TILE_COST,
std::pair< float, float >  bounds = std::pair< float, float >({std::numeric_limits< float >::min(), std::numeric_limits< float >::max()}),
const std::string &  feature_name = "" 
)

Retrieves the values of the tiles in the path map.

This function retrieves the values of the tiles in the path map and stores them in the output array.

Parameters
outputA pointer to a 2D array of floats where the tile values will be stored.
nxA pointer to an integer that will store the number of tiles in the x-direction.
nyA pointer to an integer that will store the number of tiles in the y-direction.
nzA pointer to an integer that will store the number of tiles in the z-direction.
value_typeThe type of value to retrieve for each tile. Defaults to PM_TILE_COST.
boundsA pair of floats representing the lower and upper bounds for the tile values. Defaults to the minimum and maximum float values.
feature_nameThe name of the feature for which to retrieve the tile values. Defaults to an empty string.
std::vector<PathMapTile>& IMP::bff::PathMap::get_tiles ( )

Values of tiles

Returns
vector of all tiles in the accessible volume
std::vector<IMP::algebra::Vector4D> IMP::bff::PathMap::get_xyz_density ( )

Get the XYZ density of the path map. This function returns a vector of IMP::algebra::Vector4D objects representing the XYZ density of the path map.

Returns
std::vector The XYZ density of the path map.
void IMP::bff::PathMap::get_xyz_density ( double **  output,
int *  n_output1,
int *  n_output2 
)

Get the XYZ density of the path map. This function returns the XYZ density of the path map as a 2D array.

Parameters
outputA pointer to a 2D array to store the XYZ density.
n_output1A pointer to an integer to store the number of rows in the output array.
n_output2A pointer to an integer to store the number of columns in the output array.
void IMP::bff::PathMap::resize ( unsigned int  nvox)

Resizes the PathMap object.

This function resizes the PathMap object to accommodate the specified number of voxels.

Parameters
nvoxThe number of voxels to resize the PathMap to.
void IMP::bff::PathMap::sample_obstacles ( double  extra_radius = 0.0)

Resamples the obstacles in the path map.

This function resamples the obstacles in the path map, updating their positions and sizes.

Parameters
extra_radiusThe extra radius to add to the obstacles (optional, default is 0.0).
void IMP::bff::PathMap::set_data ( double *  input,
int  n_input,
float  obstacle_threshold = -1,
bool  binarize = true,
float  obstacle_penalty = TILE_PENALTY_DEFAULT 
)

Sets the data for the path map.

This function sets the input data for the path map. The input data is a 1D array of doubles representing the map.

Parameters
inputPointer to the input data array.
n_inputNumber of elements in the input data array.
obstacle_thresholdThe threshold value for considering a cell as an obstacle. Default value is -1.
binarizeFlag indicating whether to binarize the input data. Default value is true.
obstacle_penaltyThe penalty value for obstacle cells. Default value is TILE_PENALTY_DEFAULT.
Note
The input data array should be of size n_input.
If binarize is set to true, the input data will be converted to binary values based on the obstacle_threshold.
The obstacle_penalty is used to assign penalty values to obstacle cells in the path map.
void IMP::bff::PathMap::set_path_map_header ( PathMapHeader path_map_header,
float  resolution = -1.0 
)

Set the path map header.

This function sets the path map header for the path map.

Parameters
path_map_headerThe path map header to set.
resolutionThe resolution of the path map. Default value is -1.0.
void IMP::bff::PathMap::update_tiles ( float  obstacle_threshold = -1.0,
bool  binarize = true,
float  obstacle_penalty = TILE_PENALTY_DEFAULT,
bool  reset_tile_edges = true 
)

Updates the tiles in the path map.

This function updates the tiles in the path map based on the given parameters.

Parameters
obstacle_thresholdThe threshold value for considering a cell as an obstacle. Default value is -1.0.
binarizeA flag indicating whether to binarize the path map. Default value is true.
obstacle_penaltyThe penalty value for obstacle cells. Default value is TILE_PENALTY_DEFAULT.
reset_tile_edgesA flag indicating whether to reset the edges of the tiles. Default value is true.

Friends And Related Function Documentation

std::vector<float> get_tile_values ( int  value_type = PM_TILE_COST,
std::pair< float, float >  bounds = std::pair< float, float >({std::numeric_limits< float >::min(), std::numeric_limits< float >::max()}),
const std::string &  feature_name = "" 
)
related

Get the values of all tiles.

A tile in a path map contains information on the penalty for visiting a tile, the cost of a path from the origin of a path search to the tile, the density of the tile, and other user-defined information.

When getting information from a tile, the returned values can be cropped to a specified range.

Parameters
value_typeSpecifies the type of the returned information (see: PathMapTileOutputs). Depending on the value type, the output can be the penalty for visiting the tile, the total cost of a path to the tile, or the density of the tile. Additional user-defined content can also be accessed.
boundsBound for cropping the output values.
feature_nameName of a feature when accessing additional information.
Returns
A vector of values for the specified parameters.

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