IMP logo

IMP::em Namespace Reference


Detailed Description

This module allows density maps to be used to generate restraints.

Author:
Keren Lasker, Javier Velazquez-Muriel, Friedrich Foerster
Version:
1.0
Overview:
This module provides basic utilities for handling 2D and 3D density maps. The main functionalities are: (1) Reading and writing various density map formts such as XPLOR, MRC, EM and SPIDER. (2) Simulating density maps of particles, supports particles of any radiuii and mass. (3) Calculating cross-correlation scores bewteen a density map and a set of particles.
Examples
Examples can be found on the IMP.em examples page.
License:
LGPL. This library is free software; you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation; either version 2 of the License, or (at your option) any later version.
Publications:
  • Daniel Russel, Keren Lasker, Ben Webb, Dina Schneidman, Javier Velazquez-Muriel, Andrej Sali, “Integrative assembly modeling using IMP”, submitted, 2010. This paper provides an overview of the key concepts in IMP and how to apply them to biological problems.
  • Frank Alber, Friedrich Foerster, Dmitry Korkin, Maya Topf, Andrej Sali, “Integrating diverse data for structure determination of macromolecular assemblies”, Annual Review of Biochemistry, 2008. This paper provides a review of the integrative structure determination methodology and various data sources that can be used.


Data Structures

class  CoarseCC
 Responsible for performing coarse fitting between two density objects. More...
class  CoarseCCatIntervals
 Cross correlation coefficient calculator. More...
class  CoarseConvolution
 Convolutes two grids. More...
class  DensityHeader
class  DensityMap
 Class for handling density maps. More...
class  DensityMaps
class  DensityMapsTemp
 rotate a grid More...
class  FilterByThreshold
class  FitRestraint
 Calculate score based on fit to EM map. More...
class  Image
 Template class for managing 2D Electron Microscopy images in IMP. More...
class  ImageHeader
class  KernelParameters
class  MapFilterByThreshold
 Class to filter by threshold (only DensityMap). More...
class  MRCHeader
 Class to deal with the header of MRC files. More...
class  MRCReaderWriter
class  RadiusDependentKernelParameters
 Calculates kernel parameters as a function of a specific radius. More...
class  SampledDensityMap
 Class for sampling a density map from particles. More...
class  SpiderHeader
 Header for Spider images. IMP-EM is designed to be compatible with it. More...
class  SpiderImageReaderWriter
class  SpiderMapReaderWriter
 Class to read EM maps (3D) in Spider and Xmipp formats. More...
class  SurfaceShellDensityMap
 The class repersents a molecule as shells of distance from the surface. More...
class  Volume
 Template class for managing 3D Electron Microscopy volumes in IMP. More...
class  Voxel

Typedefs

typedef double emreal
typedef Decorators< Voxel,
Particles
Voxels

Functions

template<typename T >
void add_noise (T &v, double op1, double op2, const String &mode="uniform", double df=3)
 Add noise to actual values of the image.
Float approximate_molecular_mass (DensityMap *m, Float threshold)
void byte_swap (unsigned char *ch, int n_array)
 Swaps the byte order in an array of 32-bit ints.
Float compute_fitting_score (Particles &ps, DensityMap *em_map, FloatKey rad_key=core::XYZR::get_default_radius_key(), FloatKey wei_key=atom::Mass::get_mass_key())
 Compute fitting scores for a given set of rigid transformations.
FittingSolutions compute_fitting_scores (Particles &ps, DensityMap *em_map, const FloatKey &rad_key, const FloatKey &wei_key, const std::vector< IMP::algebra::Transformation3D > &transformations, bool fast_version=false)
 Compute fitting scores for a given set of rigid transformations.
Particles density2particles (DensityMap &dmap, Float threshold, Model *m)
 Converts a density grid to a set of paritlces.
void DensityHeader_to_ImageHeader (const DensityHeader &header, ImageHeader &h)
double EXP (float y)
algebra::BoundingBoxD< 3 > get_bounding_box (const DensityMap *m)
algebra::BoundingBoxD< 3 > get_bounding_box (DensityMap *m, Float threshold)
std::string get_data_path (std::string file_name)
 Return the path to installed data for this module.
double get_density (DensityMap *m, const algebra::VectorD< 3 > &v)
std::string get_example_path (std::string file_name)
 Return the path to installed example data for this module.
int get_machine_stamp ()
std::string get_module_name ()
const VersionInfoget_module_version_info ()
long get_number_of_particles_outside_of_the_density (DensityMap *dmap, const Particles &ps)
 Get the number of particles that are outside of the density.
DensityMapget_resampled (DensityMap *in, double scaling)
DensityMapget_transformed (DensityMap *in, const algebra::Transformation3D &tr)
DensityMapget_transformed (DensityMap *in, const algebra::Transformation3D &tr, double threshold)
void ImageHeader_to_DensityHeader (const ImageHeader &h, DensityHeader &header)
int is_bigendian ()
 Returns true if this machine is big endian.
FittingSolutions local_rigid_fitting (core::RigidBody &rb, const FloatKey &radius_key, const FloatKey &weight_key, DensityMap *dmap, OptimizerState *display_log=NULL, Int number_of_optimization_runs=5, Int number_of_mc_steps=10, Int number_of_cg_steps=100, Float max_translation=2., Float max_rotation=5.)
 Local rigid fitting of a rigid body.
FittingSolutions local_rigid_fitting_around_point (core::RigidBody &rb, const FloatKey &radius_key, const FloatKey &weight_key, DensityMap *dmap, const algebra::VectorD< 3 > &anchor_centroid, OptimizerState *display_log, Int number_of_optimization_runs=5, Int number_of_mc_steps=10, Int number_of_cg_steps=100, Float max_translation=2., Float max_rotation=5.)
 Local rigid fitting of a rigid body around a center point.
FittingSolutions local_rigid_fitting_around_points (core::RigidBody &rb, const FloatKey &rad_key, const FloatKey &wei_key, DensityMap *dmap, const std::vector< algebra::VectorD< 3 > > &anchor_centroids, OptimizerState *display_log, Int number_of_optimization_runs=5, Int number_of_mc_steps=10, Int number_of_cg_steps=100, Float max_translation=2., Float max_rotation=5.)
 Local rigid fitting of a rigid body around a set of center points.
FittingSolutions local_rigid_fitting_grid_search (Particles &ps, const FloatKey &rad_key, const FloatKey &wei_key, DensityMap *dmap, Int max_voxels_translation=2, Int translation_step=1, Float max_angle_in_radians=0.174, Int number_of_rotations=100)
 Local grid search rigid fitting.
SampledDensityMapparticles2density (Particles &ps, Float resolution, Float apix, int sig_cutoff=3, const FloatKey &rad_key=IMP::core::XYZR::get_default_radius_key(), const FloatKey &weight_key=IMP::atom::Mass::get_mass_key())
 Resample a set of particles into a density grid.
void project_given_direction (DensityMap &map, IMP::algebra::Matrix2D< double > &m2, const int Ydim, const int Xdim, IMP::algebra::VectorD< 3 > &direction, const IMP::algebra::VectorD< 3 > &shift, const double equality_tolerance)
template<typename T >
void project_given_direction1 (IMP::algebra::Matrix3D< T > &m3, IMP::algebra::Matrix2D< T > &m2, const int Ydim, const int Xdim, IMP::algebra::VectorD< 3 > &direction, const IMP::algebra::VectorD< 3 > &shift, const double equality_tolerance)
void project_given_rotation (DensityMap &map, IMP::algebra::Matrix2D< double > &m2, const int Ydim, const int Xdim, const IMP::algebra::Rotation3D &Rot, const IMP::algebra::VectorD< 3 > &shift, const double equality_tolerance)
template<typename T >
void project_given_rotation1 (IMP::algebra::Matrix3D< T > &m3, IMP::algebra::Matrix2D< double > &m2, const int Ydim, const int Xdim, const IMP::algebra::Rotation3D &Rot, const IMP::algebra::VectorD< 3 > &shift, const double equality_tolerance)
DensityMapread_map (const char *filename)
DensityMapread_map (const char *filename, MapReaderWriter &reader)
void write_map (DensityMap *m, const char *filename, MapReaderWriter &writer)

Variables

const double EPS = 1e-16
 Image = _Image
 ImageReaderWriter = _ImageReaderWriter
 SpiderImageReaderWriter = _SpiderImageReaderWriter
 Volume = _Volume

Function Documentation

template<typename T >
void IMP::em::add_noise ( T &  v,
double  op1,
double  op2,
const String &  mode = "uniform",
double  df = 3 
)

Add noise to actual values of the image.

Adds random noise to the image array. Supported distributions:

  • uniform distribution, giving the range (lower, upper). DEFAULT
  • gaussian distribution, giving the mean and the standard deviation
     add_noise(v1,0, 1);
     // uniform distribution between 0 and 1
    
     v1.add_noise(0, 1, "uniform");
     // the same
    
     v1.add_noise(0, 1, "gaussian");
     // gaussian distribution with 0 mean and stddev=1
    
    Note:
    Tested with MultiArray, Matrix2D and Matrix3D

Float IMP::em::approximate_molecular_mass ( DensityMap *  m,
Float  threshold 
)

Parameters:
[in] m a density map
[in] threshold consider volume of only voxels above this threshold
Note:
The method assumes 1.21 cubic A per dalton (Harpaz 1994).

Float IMP::em::compute_fitting_score ( Particles &  ps,
DensityMap *  em_map,
FloatKey  rad_key = core::XYZR::get_default_radius_key(),
FloatKey  wei_key = atom::Mass::get_mass_key() 
)

Compute fitting scores for a given set of rigid transformations.

Score how well a set of particles fit a map

Parameters:
[in] ps The particles to be fitted
[in] em_map The density map to fit to
[in] rad_key The raidus key of the particles in the rigid body
[in] wei_key The weight key of the particles in the rigid body
Note:
the function assumes the density map holds its density

FittingSolutions IMP::em::compute_fitting_scores ( Particles &  ps,
DensityMap *  em_map,
const FloatKey &  rad_key,
const FloatKey &  wei_key,
const std::vector< IMP::algebra::Transformation3D > &  transformations,
bool  fast_version = false 
)

Compute fitting scores for a given set of rigid transformations.

Score how well a set of particles fit to the map in various rigid transformations.

Parameters:
[in] ps The particles to be fitted (treated rigid)
[in] em_map The density map to fit to
[in] rad_key The raidus key of the particles in the rigid body
[in] wei_key The weight key of the particles in the rigid body
[in] fast_version If fast_version is used the sampled density map of the input particles (ps) is not resampled for each transformation but instead the corresponding grid is rotated. This option significantly improves the running times but the returned scores are less accurate
[in] transformations A set of rigid transformations
Returns:
The scored fitting solutions
Note:
the function assumes the density map holds its density

Particles IMP::em::density2particles ( DensityMap &  dmap,
Float  threshold,
Model *  m 
)

Converts a density grid to a set of paritlces.

Each such particle will be have xyz attributes and a density_val attribute of type Float.

Parameters:
[in] dmap the density map
[in] threshold only voxels with density above the given threshold will be converted to particles
[in] m model to store the new particles
Returns:
particles corresponding to all voxels above the threshold

void IMP::em::DensityHeader_to_ImageHeader ( const DensityHeader &  header,
ImageHeader &  h 
)

Function to transfer the (compatible) information content from DensityHeader to ImageHeader

algebra::BoundingBoxD<3> IMP::em::get_bounding_box ( DensityMap *  m,
Float  threshold 
)

Parameters:
[in] m a density map
[in] threshold find the boudning box for voxels with value above the threshold

std::string IMP::em::get_data_path ( std::string  file_name  ) 

Return the path to installed data for this module.

Each module has its own data directory, so be sure to use the version of this function in the correct module. To read the data file "data_library" that was placed in the data directory of module "mymodule", do something like

    std::ifstream in(IMP::mymodule::get_data_path("data_library"));
This will ensure that the code works when IMP is installed or used via the tools/imppy.sh script.

double get_density ( DensityMap *  m,
const algebra::VectorD< 3 > &  v 
)

Return the value for the density map, m, at point v, interpolating linearly from the sample values. The resulting function is C0 over R3.

std::string IMP::em::get_example_path ( std::string  file_name  ) 

Return the path to installed example data for this module.

Each module has its own example directory, so be sure to use the version of this function in the correct module. For example to read the file example_protein.pdb located in the examples directory of the IMP::atom module, do

    IMP::atom::read_pdb(IMP::atom::get_example_path("example_protein.pdb", model));
This will ensure that the code works when IMP is installed or used via the tools/imppy.sh script.

int IMP::em::get_machine_stamp (  ) 

Returns a CCP4 convention machine stamp: 0x11110000 for big endian, or 0x44440000 for little endian

long IMP::em::get_number_of_particles_outside_of_the_density ( DensityMap *  dmap,
const Particles &  ps 
)

Get the number of particles that are outside of the density.

/note the function assumes that all of the particles have XYZ coordinates

DensityMap* IMP::em::get_resampled ( DensityMap *  in,
double  scaling 
)

Get a resampled version of the map. The spacing is multiplied by scaling. That means, scaling values greater than 1 increase the voxel size.

DensityMap * get_transformed ( DensityMap *  in,
const algebra::Transformation3D &  tr 
)

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.

DensityMap * get_transformed ( DensityMap *  in,
const algebra::Transformation3D &  tr,
double  threshold 
)

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()).

void IMP::em::ImageHeader_to_DensityHeader ( const ImageHeader &  h,
DensityHeader &  header 
)

Function to transfer the (compatible) information content from ImageHeader to DensityHeader

FittingSolutions IMP::em::local_rigid_fitting ( core::RigidBody &  rb,
const FloatKey &  radius_key,
const FloatKey &  weight_key,
DensityMap *  dmap,
OptimizerState *  display_log = NULL,
Int  number_of_optimization_runs = 5,
Int  number_of_mc_steps = 10,
Int  number_of_cg_steps = 100,
Float  max_translation = 2.,
Float  max_rotation = 5. 
)

Local rigid fitting of a rigid body.

Fit a set of particles to a density map around their centroid. The fitting is assessed using the cross-correaltion score. The optimization is a standard MC/CG procedure. The function returns a list of solutions sortedo the cross-correlation score.

Note:
The transformations are the rigid-body tranformation (with respect to an internal coordinate system). To get the actual delta transformation from the original placement of the rigid body, use the operator/ with a reference trasnforamtion outside of this function.

The returned cross-correlation score is 1-cc, as we usually want to minimize a scroing function. Thus a score of 1 means no-correlation and a score of 0. is perfect correlation.

The input rigid body should be also IMP::atom::Hierarchy

Parameters:
[in] rb The rigid body to fit
[in] radius_key The raidus key of the particles in the rigid body
[in] weight_key The weight key of the particles in the rigid body
[in] dmap The density map to fit to
[in] display_log If provided, then intermediate states in during the optimization procedure are printed
[in] number_of_optimization_runs number of Monte Carlo optimizations
[in] number_of_mc_steps number of steps in a Monte Carlo optimization
[in] number_of_cg_steps number of Conjugate Gradients steps in a Monte Carlo step
[in] max_translation maximum translation step in a MC optimization step
[in] max_rotation maximum rotation step in a single MC optimization step
Returns:
the refined fitting solutions

FittingSolutions IMP::em::local_rigid_fitting_around_point ( core::RigidBody &  rb,
const FloatKey &  radius_key,
const FloatKey &  weight_key,
DensityMap *  dmap,
const algebra::VectorD< 3 > &  anchor_centroid,
OptimizerState *  display_log,
Int  number_of_optimization_runs = 5,
Int  number_of_mc_steps = 10,
Int  number_of_cg_steps = 100,
Float  max_translation = 2.,
Float  max_rotation = 5. 
)

Local rigid fitting of a rigid body around a center point.

Fit a set of particles to a density map around an anchor point. The fitting is assessed using the cross-correaltion score. The optimization is a standard MC/CG procedure. The function returns a list of solutions sortedo the cross-correlation score.

Note:
The transformations are the rigid-body tranformation (with respect to an internal coordinate system). To get the actual delta transformation from the original placement of the rigid body, use the operator/ with a reference trasnforamtion outside of this function.

The returned cross-correlation score is 1-cc, as we usually want to minimize a scroing function. Thus a score of 1 means no-correlation and a score of 0. is perfect correlation.

The input rigid body should be also IMP::atom::Hierarchy

Parameters:
[in] rb The rigid body to fit
[in] radius_key The raidus key of the particles in the rigid body
[in] weight_key The weight key of the particles in the rigid body
[in] dmap The density map to fit to
[in] anchor_centroid The point to fit the particles around
[in] display_log If provided, then intermediate states in during the optimization procedure are printed
[in] number_of_optimization_runs number of Monte Carlo optimizations
[in] number_of_mc_steps number of steps in a Monte Carlo optimization
[in] number_of_cg_steps number of Conjugate Gradients steps in a Monte Carlo step
[in] max_translation maximum translation step in a MC optimization step
[in] max_rotation maximum rotation step in a single MC optimization step
Returns:
the refined fitting solutions

FittingSolutions IMP::em::local_rigid_fitting_around_points ( core::RigidBody &  rb,
const FloatKey &  rad_key,
const FloatKey &  wei_key,
DensityMap *  dmap,
const std::vector< algebra::VectorD< 3 > > &  anchor_centroids,
OptimizerState *  display_log,
Int  number_of_optimization_runs = 5,
Int  number_of_mc_steps = 10,
Int  number_of_cg_steps = 100,
Float  max_translation = 2.,
Float  max_rotation = 5. 
)

Local rigid fitting of a rigid body around a set of center points.

Fit a set of particles to a density map around each of the input points. The function apply local_rigid_fitting_around_point around each center.

Note:
The input rigid body should be also IMP::atom::Hierarchy
Parameters:
[in] rb The rigid body to fit
[in] rad_key The raidus key of the particles in the rigid body
[in] wei_key The weight key of the particles in the rigid body
[in] dmap The density map to fit to
[in] anchor_centroids The points to fit the particles around
[in] display_log If provided, then intermediate states in during the optimization procedure are printed
[in] number_of_optimization_runs number of Monte Carlo optimizations
[in] number_of_mc_steps number of steps in a Monte Carlo optimization
[in] number_of_cg_steps number of Conjugate Gradients steps in a Monte Carlo step
[in] max_translation maximum translation step in a MC optimization step
[in] max_rotation maximum rotation step in a single MC optimization step
Returns:
the refined fitting solutions

FittingSolutions IMP::em::local_rigid_fitting_grid_search ( Particles &  ps,
const FloatKey &  rad_key,
const FloatKey &  wei_key,
DensityMap *  dmap,
Int  max_voxels_translation = 2,
Int  translation_step = 1,
Float  max_angle_in_radians = 0.174,
Int  number_of_rotations = 100 
)

Local grid search rigid fitting.

Fit a set of particles to a density map around their centroid. The fitting is assessed using the cross-correaltion score. The optimization is a grid search

Note:
The transformations are not clustered.

The returned cross-correlation score is 1-cc, as we usually want to minimize a scroing function. Thus a score of 1 means no-correlation and a score of 0. is perfect correlation.

Parameters:
[in] ps The particles to be fitted (treated rigid)
[in] rad_key The raidus key of the particles in the rigid body
[in] wei_key The weight key of the particles in the rigid body
[in] dmap The density map to fit to
[in] max_voxels_translation Sample translations within -max_voxel_translation to max_voxel_translation
[in] translation_step The translation sampling step
[in] max_angle_in_radians Sample rotations with +- max_angle_in_radians around the current rotation
[in] number_of_rotations The number of rotations to sample
Returns:
the refiend fitting solutions

SampledDensityMap * particles2density ( Particles &  ps,
Float  resolution,
Float  apix,
int  sig_cutoff = 3,
const FloatKey &  rad_key = IMP::core::XYZR::get_default_radius_key(),
const FloatKey &  weight_key = IMP::atom::Mass::get_mass_key() 
)

Resample a set of particles into a density grid.

Each such particle should be have xyz radius and weight attributes

Parameters:
[in] ps the particles to sample
[in] resolution the resolution of the new sampled map
[in] apix the voxel size of the sampled map
[in] sig_cutoff sigma cutoff used in sampling
[in] rad_key the radius attribute key of the particles
[in] weight_key the weight attribute key of the particles
Returns:
the sampled density grid

void IMP::em::project_given_direction ( DensityMap &  map,
IMP::algebra::Matrix2D< double > &  m2,
const int  Ydim,
const int  Xdim,
IMP::algebra::VectorD< 3 > &  direction,
const IMP::algebra::VectorD< 3 > &  shift,
const double  equality_tolerance 
)

Projects a given DensityMap into a 2D matrix given a direction and a shift vector. The direction is used to build a rotation that ultimately describes the projection plane. This rotation builds the coordinate system for the projection. The projection plane is the XY plane (Z=0) of this new coordinate system. The additional shift vector is applied in the coordinate system for the projection. The projection model is: p = Rot * r + v where: p - coordinates of a point in the projection coordinate system Rot - Rotation3D applied to a point r of the universal coordinate system employed for the DensityMap r - coordinates of a point r in the of the universal coordinate system of DensityMap v - shift vector applied to p in the coordinate system of the projection.

Parameters:
[in] map A DensityMap of values to project.
[out] m2 A Matrix2D of floats to store the projection.
[in] Ydim size in rows for the Matrix2D that is to contain the projection
[in] Xdim size in columns for the Matrix2D that is to contain the projection
[in] direction vector containing the direction of projection desired
[in] shift Shift vector applied to p in the coordinate system of the projection.
[in] equality_tolerance tolerance allowed to consider a value in the direction as zero.
Note:
The function assumes that the matrices are given and stored with the (z,y,x) convention for 3D and (y,x) for 2D. But it expects and operates all the vector3D using the (x,y,z) convention.

void IMP::em::project_given_rotation ( DensityMap &  map,
IMP::algebra::Matrix2D< double > &  m2,
const int  Ydim,
const int  Xdim,
const IMP::algebra::Rotation3D Rot,
const IMP::algebra::VectorD< 3 > &  shift,
const double  equality_tolerance 
)

Projects a given DensityMap into a 2D matrix given the euler angles (ZYZ) and a shift vector. The euler angles are used to build a rotation that ultimately describes the projection. This rotation builds the coordinate system for the projection. The projection plane is the XY plane (Z=0) of this new coordinate system. The additional shift vector is applied in the coordinate system for the projection. The projection model is: p = Rot * r + v where: p - coordinates of a point in the projection coordinate system Rot - Rotation3D applied to a point r of the universal coordinate system employed for the DensityMap r - coordinates of a point r in the of the universal coordinate system of DensityMap v - shift vector applied to p in the coordinate system of the projection.

Parameters:
[in] map A DensityMap of values to project.
[out] m2 A Matrix2D of floats to store the projection.
[in] Ydim size in rows for the Matrix2D that is to contain the projection
[in] Xdim size in columns for the Matrix2D that is to contain the projection
[in] Rot The rotation that is applied before projection along z
[in] shift Shift vector applied to p in the coordinate system of the projection.
[in] equality_tolerance tolerance allowed to consider a value in the direction as zero.
Note:
The function assumes that the matrices are given and stored with the (z,y,x) convention for 3D and (y,x) for 2D. But it expects and operates all the vector3D using the (x,y,z) convention.

template<typename T >
void IMP::em::project_given_rotation1 ( IMP::algebra::Matrix3D< T > &  m3,
IMP::algebra::Matrix2D< double > &  m2,
const int  Ydim,
const int  Xdim,
const IMP::algebra::Rotation3D Rot,
const IMP::algebra::VectorD< 3 > &  shift,
const double  equality_tolerance 
)

Projects a given 3D matrix into a 2D matrix given a direction shift vector. The direction is used to build a rotation that ultimately describes the projection plane. This rotation builds the coordinate system for the projection. The projection plane is the XY plane (Z=0) of this new coordinate system. The additional shift vector is applied in the coordinate system for the projection. The projection model is: p = Rot * r + v where: p - coordinates of a point in the projection coordinate system Rot - Rotation3D applied to a point r of the universal coordinate system employed for the Matrix3D r - coordinates of a point r in the of the universal coordinate system of Matrix3D v - shift vector applied to p in the coordinate system of the projection.

Parameters:
[in] m3 A matrix3D of values to project.
[out] m2 A Matrix2D of floats to store the projection.
[in] Ydim size in rows for the Matrix2D that is to contain the projection
[in] Xdim size in columns for the Matrix2D that is to contain the projection
[in] Rot the rotation to apply before projection along z
[in] shift Shift vector applied to p in the coordinate system of the projection.
[in] equality_tolerance tolerance allowed to consider a value in the direction as zero.
Note:
The function assumes that the matrices are given and stored with the (z,y,x) convention for 3D and (y,x) for 2D. But it expects and operates all the vector3D using the (x,y,z) convention.

DensityMap * read_map ( const char *  filename  ) 

Read a density map from a file and return it. Guess the file type from the file name. The file formats supported are:

  • .mrc for MRC files

DensityMap * read_map ( const char *  filename,
MapReaderWriter &  reader 
)

Read a density map from a file and return it.

void write_map ( DensityMap *  m,
const char *  filename,
MapReaderWriter &  writer 
)

Write a density map to a file.


Generated on Mon Mar 8 23:09:01 2010 for IMP by doxygen 1.5.8