IMP logo
IMP Reference Guide  2.10.0
The Integrative Modeling Platform
IMP::em::DensityMap Class Reference

Class for handling density maps. More...

#include <IMP/em/DensityMap.h>

+ Inheritance diagram for IMP::em::DensityMap:

Detailed Description

Class for handling density maps.

Note
The location of a voxel is its center. That is important for sampling function as well as for functions like get_location_in_dim_by_voxel.

Definition at line 93 of file DensityMap.h.

Public Member Functions

 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...
 
emreal 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 (emreal 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 (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...
 
emreal 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...
 
emreal 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
 
emreal get_value (float x, float y, float z) const
 Gets the value of the voxel located at (x,y,z) More...
 
emreal get_value (const algebra::Vector3D &point) const
 
emreal get_value (long index) const
 Gets the value of the voxel at a given index. More...
 
virtual ::IMP::VersionInfo get_version_info () const
 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 (emreal loc, emreal kdist, emreal 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, emreal value)
 Set the value of the voxel at a given index. More...
 
void set_value (float x, float y, float z, emreal 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 (emreal loc, emreal kdist, emreal 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

void allocated_data ()
 
void float2real (float *f_data, boost::scoped_array< emreal > &r_data)
 
void real2float (emreal *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

boost::scoped_array< emreal > 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_
 

Friends

DensityMapread_map (std::string filename, MapReaderWriter *reader)
 
void write_map (DensityMap *m, std::string filename, MapReaderWriter *writer)
 

Related Functions

(Note that these are not member functions.)

void add_to_map (DensityMap *dm, const Particles &pis)
 Rasterize the particles into an existing density map. More...
 
DensityMapcreate_density_map (algebra::DenseGrid3D< float > &grid)
 Return a density map with the values taken from the grid. More...
 
DensityMapcreate_density_map (algebra::DenseGrid3D< double > &grid)
 Return a density map with the values taken from the grid. More...
 
DensityMapget_binarized_interior (DensityMap *dmap)
 Return a binarized map with 1 for voxels that are internal in the input map. More...
 
double get_density (const DensityMap *m, const algebra::Vector3D &v)
 Get density value at point v, interpolating linearly from the sample values. More...
 
DensityMapget_resampled (DensityMap *input, double scaling)
 Get a resampled version of the map. More...
 
DensityMapget_transformed (const DensityMap *input, const algebra::Transformation3D &tr, double threshold)
 Return a new density map containing a rotated version of the old one. More...
 
DensityMapget_transformed (DensityMap *input, const algebra::Transformation3D &tr)
 Return a new density map containing a rotated version of the old one. More...
 
void get_transformed_into (const DensityMap *source, const algebra::Transformation3D &tr, DensityMap *into, bool calc_rms=true)
 Rotate a density map into another map. More...
 
DensityMapread_map (std::string filename, MapReaderWriter *reader)
 Read a density map from a file (using the given reader) and return it. More...
 
DensityMapread_map (std::string filename)
 Read a density map from a file and return it. More...
 
void write_map (DensityMap *m, std::string filename, MapReaderWriter *writer)
 Write a density map to a file using the given writer. More...
 
void write_map (DensityMap *m, std::string filename)
 Write a density map to a file. More...
 

Constructor & Destructor Documentation

IMP::em::DensityMap::DensityMap ( const DensityHeader header,
std::string  name = "DensityMap%1%" 
)

Construct a density map as instructed in the input header.

Member Function Documentation

void IMP::em::DensityMap::add ( const DensityMap other)

Sums two grids; the result is kept in the map.

Parameters
[in]otherthe other map
Note
The shared extend is summed
The two maps should have the same voxelsize and other should be contained within this map
void IMP::em::DensityMap::add ( Float  d)

Add a number to every voxel in the map.

Parameters
[in]dValue to add
void IMP::em::DensityMap::calc_all_voxel2loc ( )

Calculates the coordinates that correspond to all voxels.

emreal IMP::em::DensityMap::calcRMS ( )

Calculates RMSD and mean of a map values as stored in the header.

The header stores whether the map is normalized.

void IMP::em::DensityMap::copy_map ( const DensityMap other)

Copy map into this map.

algebra::Vector3D IMP::em::DensityMap::get_centroid ( emreal  threshold = 0.0) const

Get the centroid of all voxels with density above a given threshold.

Parameters
[in]thresholdthe input threshold
DensityMap* IMP::em::DensityMap::get_cropped ( float  threshold)

Create and return a new cropped map.

The margins are determined to be the bounding box with density values below the input

Parameters
[in]thresholdused for cropping
Returns
the new cropped map.
DensityMap* IMP::em::DensityMap::get_cropped ( const algebra::BoundingBox3D bb)

Create and return a new cropped map with the given bounding box extent.

Parameters
[in]bbthe bounding box
Returns
the new cropped map.
Note
If the input bounding box is larger than the density box, it is snapped to the right size.
int IMP::em::DensityMap::get_dim_index_by_location ( const algebra::Vector3D v,
int  ind 
) const

Calculate dimension index of a given location.

Parameters
[in]vThe position (in angstroms)
[in]inddimension index (X:0,Y:1 or Z:2)
const DensityHeader* IMP::em::DensityMap::get_header ( ) const

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

Definition at line 243 of file DensityMap.h.

DensityHeader* IMP::em::DensityMap::get_header_writable ( )

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

Definition at line 246 of file DensityMap.h.

algebra::Vector3D IMP::em::DensityMap::get_location_by_voxel ( long  index) const

Calculate the location of a given voxel.

Parameters
[in]indexThe voxel index
Returns
the location (x,y,z) (in angstroms) of a given voxel.

Definition at line 168 of file DensityMap.h.

float IMP::em::DensityMap::get_location_in_dim_by_voxel ( long  index,
int  dim 
) const

Calculate the location of a given voxel in a given dimension.

Parameters
[in]indexThe voxel index
[in]dimThe dimension of interest (x:=0,y:=1,z:=2)
Returns
the location (x,y,z) (in angstroms) of a given voxel. If the index is not part of the map, the function returns -1.
std::string IMP::em::DensityMap::get_locations_string ( float  t)

Print the locations of all voxels with value above a given threshold.

emreal IMP::em::DensityMap::get_max_value ( ) const

Returns the value of the voxel with the highest density.

float IMP::em::DensityMap::get_maximum_value_in_xy_plane ( int  z_ind)

Get the maximum value in a XY plane indicated by a Z index.

float IMP::em::DensityMap::get_maximum_value_in_xz_plane ( int  y_ind)

Get the maximum value in a XZ plane indicated by a Y index.

float IMP::em::DensityMap::get_maximum_value_in_yz_plane ( int  x_ind)

Get the maximum value in a YZ plane indicated by a X index.

emreal IMP::em::DensityMap::get_min_value ( ) const

Returns the value of the voxel with the lowest density.

long IMP::em::DensityMap::get_number_of_voxels ( ) const

Get the number of map voxels.

Float IMP::em::DensityMap::get_spacing ( ) const

Get the voxel size of the map.

Note
Use update_voxel_size() to set the spacing value.

Definition at line 407 of file DensityMap.h.

emreal IMP::em::DensityMap::get_value ( float  x,
float  y,
float  z 
) const

Gets the value of the voxel located at (x,y,z)

Parameters
[in]xThe position (in angstroms) of the x coordinate
[in]yThe position (in angstroms) of the y coordinate
[in]zThe position (in angstroms) of the z coordinate
Returns
the value of the voxel located at (x,y,z)
Exceptions
IndexExceptionThe point is not covered by the grid.
Note
the value is not interpolated between this and neighboring voxels. For that, see get_density().
emreal IMP::em::DensityMap::get_value ( long  index) const

Gets the value of the voxel at a given index.

Parameters
[in]indexvoxel number in physical sense, NOT logical
virtual ::IMP::VersionInfo IMP::em::DensityMap::get_version_info ( ) const
virtual

Get information about the module and version of the object.

Reimplemented from IMP::Object.

Definition at line 412 of file DensityMap.h.

long IMP::em::DensityMap::get_voxel_by_location ( float  x,
float  y,
float  z 
) const

Calculate the voxel of a given location.

Parameters
[in]xThe position (in angstroms) of the x coordinate
[in]yThe position (in angstroms) of the y coordinate
[in]zThe position (in angstroms) of the z coordinate
Returns
the voxel index of a given position. If the position is out of the boundaries of the map, the function returns -1.
long IMP::em::DensityMap::get_voxel_by_location ( const algebra::Vector3D v) const

Calculate the voxel of a given location.

Parameters
[in]vThe position (in angstroms)
Returns
the voxel index of a given position. If the position is out of the boundaries of the map, the function returns -1.

Definition at line 155 of file DensityMap.h.

bool IMP::em::DensityMap::is_part_of_volume ( float  x,
float  y,
float  z 
) const

Checks whether a given point is in the grid.

Parameters
[in]xThe position (in angstroms) of the x coordinate
[in]yThe position (in angstroms) of the y coordinate
[in]zThe position (in angstroms) of the z coordinate
Returns
true if the point is part of the grid, false otherwise.
bool IMP::em::DensityMap::is_part_of_volume ( const algebra::Vector3D v) const

Checks whether a given point is in the grid.

Parameters
[in]vThe position (in angstroms)
Returns
true if the point is part of the grid, false otherwise.

Definition at line 190 of file DensityMap.h.

void IMP::em::DensityMap::multiply ( float  factor)

Multiply in place each voxel in the map by the input factor.

The result is kept in the map.

Parameters
[in]factorthe multiplication factor
void IMP::em::DensityMap::pad ( int  nx,
int  ny,
int  nz,
float  val = 0.0 
)

Increase the dimension of the map.

The function pads zeroes to the right-upper section on the map. The original content of the map will be in the lower XYZ part of the map.

Parameters
[in]nxthe number of voxels on the X axis
[in]nythe number of voxels on the Y axis
[in]nzthe number of voxels on the Z axis
[in]valall additional voxels will have this value
Exceptions
UsageExceptionif the input numbers of x/y/z voxels are smaller than those currently in the map
DensityMap* IMP::em::DensityMap::pad_margin ( int  mrg_x,
int  mrg_y,
int  mrg_z,
float  val = 0.0 
)

Create a new padded map.

Given this map of size [nx,ny,nz], the new map is of size [2*mrg_x+nx,2*mrg_y+ny,2*mrg_z+nz]. The new map will consist of the values of the old map, with a padded margin on all sides.

Parameters
[in]mrg_xnumber of margin voxels to add on both right and left on the X axis
[in]mrg_ynumber of margin voxels to add on both right and left on the Y axis
[in]mrg_znumber of margin voxels to add on both right and left on the Z axis
[in]valignored
Returns
the new padded map.
Exceptions
UsageExceptionif the input numbers of x/y/z voxels are smaller than those currently in the map
void IMP::em::DensityMap::pick_max ( const DensityMap other)

Pick the max value between two corresponding voxels between two maps.

The result is kept in the map.

Parameters
[in]otherthe other map
Note
The two maps should have the same voxelsize and the same dimensions
void IMP::em::DensityMap::reset_data ( float  value = 0.0)

Set the density voxels to some value and reset the management flags.

Parameters
[in]valueall of the density voxels will have this value
bool IMP::em::DensityMap::same_dimensions ( const DensityMap other) const

Checks if two maps have the same dimensions.

Parameters
[in]otherthe map to compare with
Returns
true if the two maps have the same dimensions
bool IMP::em::DensityMap::same_origin ( const DensityMap other) const

Checks if two maps have the same origin.

Parameters
[in]otherthe map to compare with
Returns
true if the two maps have the same origin
bool IMP::em::DensityMap::same_voxel_size ( const DensityMap other) const

Checks if two maps have the same voxel size.

Parameters
[in]otherthe map to compare with
Returns
true if the two maps have the same voxel size
void IMP::em::DensityMap::set_origin ( float  x,
float  y,
float  z 
)

Sets the origin of the header.

Parameters
xthe new x (angstroms)
ythe new y (angstroms)
zthe new z (angstroms)
void IMP::em::DensityMap::set_value ( long  index,
emreal  value 
)

Set the value of the voxel at a given index.

Parameters
[in]indexvoxel number in physical sense, NOT logical
[in]valuevalue
void IMP::em::DensityMap::set_value ( float  x,
float  y,
float  z,
emreal  value 
)

Set the value of the voxel at given coordinates.

void IMP::em::DensityMap::set_void_map ( int  nx,
int  ny,
int  nz 
)

Set the map dimension and reset all voxels to 0.

Parameters
[in]nxx-dimension (voxels)
[in]nyy-dimension (voxels)
[in]nzz-dimension (voxels)
Note
the origin and spacing remain unchanged
void IMP::em::DensityMap::std_normalize ( )

Normalize the density voxels according to standard deviation (stdv).

The mean is subtracted from the map, which is then divided by the stdv. The normalization flag is set to avoid repeated computation.

void IMP::em::DensityMap::update_header ( )
protected

update the header values – still in work

void IMP::em::DensityMap::update_voxel_size ( float  new_apix)

Update the voxel size of the map.

long IMP::em::DensityMap::xyz_ind2voxel ( int  x,
int  y,
int  z 
) const

Calculate the voxel of a given xyz index.

Parameters
[in]xThe voxel index on the x axis of the grid
[in]yThe voxel index on the y axis of the grid
[in]zThe voxel index on the z axis of the grid
Returns
the voxel index.

Definition at line 137 of file DensityMap.h.

Friends And Related Function Documentation

void add_to_map ( DensityMap dm,
const Particles pis 
)
related

Rasterize the particles into an existing density map.

DensityMap * create_density_map ( algebra::DenseGrid3D< float > &  grid)
related

Return a density map with the values taken from the grid.

DensityMap * create_density_map ( algebra::DenseGrid3D< double > &  grid)
related

Return a density map with the values taken from the grid.

DensityMap * get_binarized_interior ( DensityMap dmap)
related

Return a binarized map with 1 for voxels that are internal in the input map.

double get_density ( const DensityMap m,
const algebra::Vector3D v 
)
related

Get density value at point v, interpolating linearly from the sample values.

The resulting function is C0 over R3.

DensityMap * get_resampled ( DensityMap input,
double  scaling 
)
related

Get a resampled version of the map.

The spacing is multiplied by scaling, i.e. scaling values greater than 1 increase the voxel size.

DensityMap * get_transformed ( const DensityMap input,
const algebra::Transformation3D tr,
double  threshold 
)
related

Return a new density map containing a rotated version of the old one.

Only voxels whose value is above threshold are considered when computing the bounding box of the new map (set IMP::em::get_bounding_box()).

DensityMap * get_transformed ( DensityMap input,
const algebra::Transformation3D tr 
)
related

Return a new density map containing a rotated version of the old one.

The dimension of the new map is the same as the old one.

void get_transformed_into ( const DensityMap source,
const algebra::Transformation3D tr,
DensityMap into,
bool  calc_rms = true 
)
related

Rotate a density map into another map.

Parameters
[in]sourcethe map to transform
[in]trtransform the from density map by this transformation
[out]intothe map to transform into
[in]calc_rmsif true RMS is calculated on the transformed map
DensityMap * read_map ( std::string  filename,
MapReaderWriter reader 
)
related

Read a density map from a file (using the given reader) and return it.

DensityMap * read_map ( std::string  filename)
related

Read a density map from a file and return it.

Guess the file type from the file name. The file formats supported are:

  • .mrc
  • .em
  • .vol
  • .xplor
void write_map ( DensityMap m,
std::string  filename,
MapReaderWriter writer 
)
related

Write a density map to a file using the given writer.

void write_map ( DensityMap m,
std::string  filename 
)
related

Write a density map to a file.

Guess the file type from the file name. The file formats supported are:

  • .mrc
  • .em
  • .vol
  • .xplor

Member Data Documentation

bool IMP::em::DensityMap::loc_calculated_
protected

true if the locations have already been computed

Definition at line 458 of file DensityMap.h.

boost::scoped_array<float> IMP::em::DensityMap::x_loc_
protected

Locations (centers) for each of the voxels of the map (they are precomputed and each one is of size nvox, where nvox is the size of the map)

Definition at line 456 of file DensityMap.h.


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