IMP logo
IMP Reference Guide  2.6.0
The Integrative Modeling Platform
IMP::em Namespace Reference

Basic utilities for handling cryo-electron microscopy 3D density maps. More...

Detailed Description

Basic utilities for handling cryo-electron microscopy 3D density maps.

The main functionalities are:

  1. Reading and writing various density map formats such as XPLOR, MRC, EM and SPIDER.
  2. Simulating density maps of particles, supports particles of any radii and mass.
  3. Calculating cross-correlation scores between a density map and a set of particles.
  4. Implements several Restraints.

The restraints vary based on how accurately the fit to the density is scored, which is usually related to the evaluation speed. Pick more accurate Restraints for higher resolution maps and models. Below is the Restraints list sorted from the fastest and most simple to the slowest and most accurate.

We also provide a number of command line tools to handle electron microscopy density maps and associated files:

estimate_threshold_from_molecular_mass

Estimate EM density threshold.

map2pca

Write out density map principal components in cmm format.

mol2pca

Write out protein principal components in cmm format.

resample_density

Resample a density map.

simulate_density_from_pdb

Samples a protein into a simulated 3D density map.

view_density_header

Display the header information for a density map.

Info

Author(s): Keren Lasker, Javier Velázquez-Muriel, Friedrich Förster, Daniel Russel, Dina Schneidman

Maintainer: benmwebb

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:

Classes

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  DensityFillingRestraint
 Calculate score based on fit to EM map. More...
 
class  DensityHeader
 
class  DensityMap
 Class for handling density maps. More...
 
class  DistanceMask
 Calculates and stores a distance mask. More...
 
class  EMReaderWriter
 
class  EnvelopeFitRestraint
 A restraint for envelope-based scoring of particles in the density map. More...
 
class  EnvelopePenetrationRestraint
 Calculate score based on fit to EM map. More...
 
class  EnvelopeScore
 class for envelope based scoring using MapDistanceTransform More...
 
class  FitRestraint
 Calculate score based on fit to EM map. More...
 
class  FittingSolutions
 A simple list of fitting solutions. More...
 
class  HighDensityEmbedding
 
class  ImageHeader
 Class to deal with the header of Electron Microscopy images in IMP. More...
 
class  KernelParameters
 Calculates and stores Gaussian kernel parameters. More...
 
class  MapDistanceTransform
 Class for getting the distance from the map envelope. More...
 
class  MapReaderWriter
 The base class to handle reading and writing of density maps. More...
 
class  MRCReaderWriter
 
class  PCAAligner
 Fast alignment of points to a density map using principal components. More...
 
class  PCAFitRestraint
 Calculate score based on fit to EM map. More...
 
class  RadiusDependentDistanceMask
 
class  SampledDensityMap
 Class for sampling a density map from particles. More...
 
struct  SpiderHeader
 Header for Spider images. IMP-EM is designed to be compatible with it. More...
 
class  SpiderMapReaderWriter
 Class to read EM maps (3D) in Spider and Xmipp formats. More...
 
class  SurfaceShellDensityMap
 The class represents a molecule as shells of distance from the surface. More...
 
class  Voxel
 
class  XplorReaderWriter
 

Typedefs

typedef IMP::Vector
< IMP::Pointer< DensityMap > > 
DensityMaps
 
typedef IMP::Vector
< IMP::WeakPointer< DensityMap > > 
DensityMapsTemp
 
typedef double emreal
 
typedef IMP::Vector
< FittingSolutions
FittingSolutionsList
 
typedef IMP::Vector
< KernelParameters
KernelParametersList
 
typedef IMP::Vector
< IMP::Pointer
< SampledDensityMap > > 
SampledDensityMaps
 
typedef IMP::Vector
< IMP::WeakPointer
< SampledDensityMap > > 
SampledDensityMapsTemp
 
typedef IMP::Vector
< IMP::Pointer
< SurfaceShellDensityMap > > 
SurfaceShellDensityMaps
 
typedef IMP::Vector
< IMP::WeakPointer
< SurfaceShellDensityMap > > 
SurfaceShellDensityMapsTemp
 
typedef IMP::Vector< VoxelVoxels
 

Enumerations

enum  KernelType { GAUSSIAN, BINARIZED_SPHERE, SPHERE }
 

Functions

Float approximate_molecular_mass (DensityMap *m, Float threshold)
 Estimate the molecular mass from a map. More...
 
DensityMapbinarize (DensityMap *orig_map, float threshold, bool reverse=false)
 Return a map with 0 for all voxels below the threshold and 1 for those above. More...
 
Float calculate_intersection_score (const SurfaceShellDensityMap *d1, const SurfaceShellDensityMap *d2)
 
Float compute_fitting_score (const ParticlesTemp &ps, DensityMap *em_map, FloatKey wei_key=atom::Mass::get_mass_key())
 Compute fitting scores for a given set of rigid transformations. More...
 
FittingSolutions compute_fitting_scores (const ParticlesTemp &ps, DensityMap *em_map, const algebra::Transformation3Ds &transformations, bool fast_version=false, bool local_score=false, const FloatKey &wei_key=atom::Mass::get_mass_key())
 Compute fitting scores for a given set of rigid transformations. More...
 
FittingSolutions compute_fitting_scores (DensityMap *em_map, core::RigidBody rb, Refiner *refiner, const algebra::Transformation3Ds &transformations)
 Compute fitting scores for a given set of rigid transformations. More...
 
double convolute (const DensityMap *m1, const DensityMap *m2)
 Return a convolution between density maps m1 and m2. More...
 
DensityHeader create_density_header (const algebra::BoundingBoxD< 3 > &bb, float spacing)
 Create a header from a bounding box in 3D. More...
 
DensityMapcreate_density_map (const DensityMap *other)
 Create a copy of another map. More...
 
DensityMapcreate_density_map (const algebra::BoundingBox3D &bb, double spacing)
 Create an empty density map from a bounding box. More...
 
DensityMapcreate_density_map (int nx, int ny, int nz, double spacing)
 Create an empty density map. More...
 
template<class S , class V , class E >
DensityMapcreate_density_map (const IMP::algebra::GridD< 3, S, V, E > &arg)
 Create a density map from an arbitrary IMP::algebra::GridD. More...
 
Particles density2particles (DensityMap *dmap, Float threshold, Model *m, int step=1)
 Convert a density grid to a set of particles. More...
 
algebra::Vector3Ds density2vectors (DensityMap *dmap, Float threshold)
 Convert a density grid to a set of vectors. More...
 
void DensityHeader_to_ImageHeader (const DensityHeader &header, ImageHeader &h)
 
double EXP (float y)
 
algebra::BoundingBoxD< 3 > get_bounding_box (const DensityMap *m, Float threshold)
 Get the bounding box for a map. More...
 
algebra::BoundingBoxD< 3 > get_bounding_box (const DensityMap *m)
 
algebra::GridD
< 3, algebra::DenseGridStorageD
< 3, float >, float > 
get_grid (DensityMap *in_map)
 Return a dense grid containing the voxels of the passed density map. More...
 
bool get_interiors_intersect (const DensityMap *d1, const DensityMap *d2)
 
DensityMapget_max_map (DensityMaps maps)
 Return a map where each voxel is the maximum value from the input maps. More...
 
Float get_molecular_mass_at_threshold (DensityMap *m, Float threshold, atom::ProteinDensityReference ref=atom::HARPAZ)
 Compute an approximate molecular mass. More...
 
long get_number_of_particles_outside_of_the_density (DensityMap *dmap, const Particles &ps, const IMP::algebra::Transformation3D &t=IMP::algebra::get_identity_transformation_3d(), float thr=0.0)
 Get the number of particles that are outside of the density. More...
 
Ints get_numbers_of_particles_outside_of_the_density (DensityMap *dmap, const Particles &ps, const IMP::algebra::Transformation3Ds &transformations, float thr=0.0)
 Get numbers of particles (mult transforms) that are outside the density. More...
 
double get_percentage_of_voxels_covered_by_particles (DensityMap *dmap, const Particles &ps, float smoothing_radius=3., const IMP::algebra::Transformation3D &t=IMP::algebra::get_identity_transformation_3d(), float thr=0.0)
 Get the number of density voxels that are not covered by particles. More...
 
DensityMapget_segment (DensityMap *map_to_segment, int nx_start, int nx_end, int ny_start, int ny_end, int nz_start, int nz_end)
 Get a segment of the map according to xyz indexes. More...
 
DensityMapget_segment (DensityMap *map_to_segment, algebra::Vector3Ds vecs, float dist)
 Get a segment of the map covered by the input points. More...
 
DensityMapget_segment_by_masking (DensityMap *map_to_segment, DensityMap *mask, float mas_threshold)
 Get a segment of the map covered by another map. More...
 
double get_sum (const DensityMap *m1)
 Return the sum of all voxels. More...
 
Float get_threshold_for_approximate_mass (DensityMap *m, Float desired_mass, atom::ProteinDensityReference ref=atom::HARPAZ)
 Computes the threshold to consider in an EM map to get a desired mass. More...
 
Float get_threshold_for_approximate_volume (DensityMap *m, Float desired_volume)
 Computes the threshold consider in an EM map to get a desired volume. More...
 
DensityMapget_threshold_map (const DensityMap *orig_map, float threshold)
 Return a map with 0 for all voxels below the threshold. More...
 
void get_transformed_into2 (const DensityMap *source, const algebra::Transformation3D &tr, DensityMap *into)
 
Float get_volume_at_threshold (DensityMap *m, Float threshold)
 Compute an approximate volume. More...
 
void ImageHeader_to_DensityHeader (const ImageHeader &h, DensityHeader &header)
 
DensityMapinterpolate_map (DensityMap *in_map, double new_spacing)
 Return a new map with an updated spacing. More...
 
FittingSolutions local_rigid_fitting (Particle *p, Refiner *refiner, const FloatKey &weight_key, DensityMap *dmap, OptimizerStates 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=.3, bool fast=true)
 Local rigid fitting of a rigid body. More...
 
FittingSolutions local_rigid_fitting_around_point (Particle *p, Refiner *refiner, const FloatKey &weight_key, DensityMap *dmap, const algebra::Vector3D &anchor_centroid, OptimizerStates 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=.3, bool fast=false)
 Local rigid fitting of a rigid body around a center point. More...
 
FittingSolutions local_rigid_fitting_around_points (Particle *p, Refiner *refiner, const FloatKey &wei_key, DensityMap *dmap, const algebra::Vector3Ds &anchor_centroids, OptimizerStates 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=.3)
 Local rigid fitting of a rigid body around a set of center points. More...
 
FittingSolutions local_rigid_fitting_grid_search (const ParticlesTemp &ps, 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. More...
 
DensityMapmask_and_norm (em::DensityMap *dmap, em::DensityMap *mask)
 Return a masked density, and normalize the output map within the mask region. More...
 
DensityMapmultiply (const DensityMap *m1, const DensityMap *m2)
 Return a density map for which voxel i contains the result of m1[i]*m2[i]. More...
 
SampledDensityMapparticles2binarized_density (const ParticlesTemp &ps, Float resolution, Float apix, int sig_cutoff=3, const FloatKey &weight_key=IMP::atom::Mass::get_mass_key())
 
SampledDensityMapparticles2density (const ParticlesTemp &ps, Float resolution, Float apix, int sig_cutoff=3, const FloatKey &weight_key=IMP::atom::Mass::get_mass_key())
 Resample a set of particles into a density grid. More...
 
SurfaceShellDensityMapparticles2surface (const ParticlesTemp &ps, Float apix, const FloatKey &weight_key=IMP::atom::Mass::get_mass_key())
 Resample a set of particles into a density grid. More...
 
def write_pca_cmm
 Write out principal components to a file in Chimera Marker format. More...
 

Variables

const double EPS = 1e-16
 

Standard module functions

All IMP modules have a set of standard functions to help get information about the module and about files associated with the module.

std::string get_module_version ()
 
std::string get_module_name ()
 
std::string get_data_path (std::string file_name)
 Return the full path to one of this module's data files. More...
 
std::string get_example_path (std::string file_name)
 Return the full path to one of this module's example files. More...
 

Typedef Documentation

Store a set of objects.

Definition at line 498 of file DensityMap.h.

Pass a set of objects.

See Also
DensityMap

Definition at line 498 of file DensityMap.h.

Pass or store a set of FittingSolutions .

Definition at line 103 of file rigid_fitting.h.

Pass or store a set of KernelParameters .

Definition at line 97 of file KernelParameters.h.

Store a set of objects.

Definition at line 136 of file SampledDensityMap.h.

Pass a set of objects.

See Also
SampledDensityMap

Definition at line 136 of file SampledDensityMap.h.

Store a set of objects.

Definition at line 102 of file SurfaceShellDensityMap.h.

Function Documentation

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

Estimate the molecular mass from a map.

Parameters
[in]ma density map
[in]thresholdconsider volume of only voxels above this threshold
Note
The method assumes 1.21 cubic Å per Dalton (Harpaz 1994).
DensityMap* IMP::em::binarize ( DensityMap *  orig_map,
float  threshold,
bool  reverse = false 
)

Return a map with 0 for all voxels below the threshold and 1 for those above.

Parameters
[in]orig_mapthe map to binarize
[in]thresholdvalues below the threshold are set to 0 and 1 otherwise
[in]reverseif true values above the threshold are set to 0 and 1 otherwise
Float IMP::em::compute_fitting_score ( const ParticlesTemp &  ps,
DensityMap *  em_map,
FloatKey  wei_key = atom::Mass::get_mass_key() 
)

Compute fitting scores for a given set of rigid transformations.

Scores how well a set of particles fit a map

Parameters
[in]psThe particles to be fitted
[in]em_mapThe density map to fit to
[in]wei_keyThe 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 ( const ParticlesTemp &  ps,
DensityMap *  em_map,
const algebra::Transformation3Ds &  transformations,
bool  fast_version = false,
bool  local_score = false,
const 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 to the map in various rigid transformations.

Parameters
[in]psThe particles to be fitted (treated rigid)
[in]em_mapThe density map to fit to
[in]wei_keyThe weight key of the particles in the rigid body
[in]fast_versionIf 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]transformationsA set of rigid transformations
[in]fast_versionif true the density map of each transformation is interpolated
[in]local_scoreif true a local cross correlation score is used
Returns
The scored fitting solutions
Note
the function assumes the density map holds its density
FittingSolutions IMP::em::compute_fitting_scores ( DensityMap *  em_map,
core::RigidBody  rb,
Refiner *  refiner,
const algebra::Transformation3Ds &  transformations 
)

Compute fitting scores for a given set of rigid transformations.

Score how well a rigid body fits to the map

Parameters
[in]em_mapThe density map to fit to
[in]rbThe rigid body
[in]refinerThe rigid body refiner
[in]transformationsTransformations of the rigid body
Returns
The scored fitting solutions
Note
the function assumes the density map holds its density

Definition at line 272 of file rigid_fitting.h.

+ Here is the call graph for this function:

double IMP::em::convolute ( const DensityMap *  m1,
const DensityMap *  m2 
)

Return a convolution between density maps m1 and m2.

The function assumes m1 and m2 are of the same dimensions.

DensityHeader IMP::em::create_density_header ( const algebra::BoundingBoxD< 3 > &  bb,
float  spacing 
)

Create a header from a bounding box in 3D.

See Also
DensityHeader
DensityMap* IMP::em::create_density_map ( const DensityMap *  other)

Create a copy of another map.

DensityMap* IMP::em::create_density_map ( const algebra::BoundingBox3D &  bb,
double  spacing 
)

Create an empty density map from a bounding box.

DensityMap* IMP::em::create_density_map ( int  nx,
int  ny,
int  nz,
double  spacing 
)

Create an empty density map.

template<class S , class V , class E >
DensityMap* IMP::em::create_density_map ( const IMP::algebra::GridD< 3, S, V, E > &  arg)

Create a density map from an arbitrary IMP::algebra::GridD.

Definition at line 626 of file DensityMap.h.

+ Here is the call graph for this function:

Particles IMP::em::density2particles ( DensityMap *  dmap,
Float  threshold,
Model *  m,
int  step = 1 
)

Convert a density grid to a set of particles.

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

Parameters
[in]dmapthe density map
[in]thresholdonly voxels with density above the given threshold will be converted to particles
[in]mmodel to store the new particles
[in]stepsample every X steps in each direction
Returns
particles corresponding to all voxels above the threshold
algebra::Vector3Ds IMP::em::density2vectors ( DensityMap *  dmap,
Float  threshold 
)

Convert a density grid to a set of vectors.

Parameters
[in]dmapthe density map
[in]thresholdonly voxels with density above the given threshold will be converted to vectors
Returns
Vector3Ds corresponding to the positions of 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 ( const DensityMap *  m,
Float  threshold 
)

Get the bounding box for a map.

Parameters
[in]ma density map
[in]thresholdfind the bounding box for voxels with value above the threshold
std::string IMP::em::get_data_path ( std::string  file_name)

Return the full path to one of this module's data files.

To read the data file "data_library" that was placed in the data directory of this module, do something like

std::ifstream in(IMP::em::get_data_path("data_library"));

This will ensure that the code works both when IMP is installed or if used via the setup_environment.sh script.

Note
Each module has its own data directory, so be sure to use this function from the correct module.
std::string IMP::em::get_example_path ( std::string  file_name)

Return the full path to one of this module's example files.

To read the example file "example_protein.pdb" that was placed in the examples directory of this module, do something like

std::ifstream in(IMP::em::get_example_path("example_protein.pdb"));

This will ensure that the code works both when IMP is installed or if used via the setup_environment.sh script.

Note
Each module has its own example directory, so be sure to use this function from the correct module.
algebra::GridD<3, algebra::DenseGridStorageD<3, float>, float> IMP::em::get_grid ( DensityMap *  in_map)

Return a dense grid containing the voxels of the passed density map.

The returned grid has the same bounding box as the map.

DensityMap* IMP::em::get_max_map ( DensityMaps  maps)

Return a map where each voxel is the maximum value from the input maps.

Note
all input maps should have the same extent
Float IMP::em::get_molecular_mass_at_threshold ( DensityMap *  m,
Float  threshold,
atom::ProteinDensityReference  ref = atom::HARPAZ 
)

Compute an approximate molecular mass.

Compute an approximate molecular mass for the set of voxels with intensity under a given threshold

Parameters
[in]ma density map
[in]thresholdonly voxels above this threshold will be considered
[in]refthe protein density reference to use in the computation. The default protein density for this computation is HARPAZ
Returns
an approximate molecular mass for the set of voxels with intensity under the provided threshold (mass in Da)
long IMP::em::get_number_of_particles_outside_of_the_density ( DensityMap *  dmap,
const Particles &  ps,
const IMP::algebra::Transformation3D t = IMP::algebra::get_identity_transformation_3d(),
float  thr = 0.0 
)

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

Note
the function assumes that all of the particles have XYZ coordinates
Ints IMP::em::get_numbers_of_particles_outside_of_the_density ( DensityMap *  dmap,
const Particles &  ps,
const IMP::algebra::Transformation3Ds transformations,
float  thr = 0.0 
)

Get numbers of particles (mult transforms) that are outside the density.

Note
the function assumes that all of the particles have XYZ coordinates
double IMP::em::get_percentage_of_voxels_covered_by_particles ( DensityMap *  dmap,
const Particles &  ps,
float  smoothing_radius = 3.,
const IMP::algebra::Transformation3D t = IMP::algebra::get_identity_transformation_3d(),
float  thr = 0.0 
)

Get the number of density voxels that are not covered by particles.

Note
the function assumes that all of the particles have XYZ coordinates
DensityMap* IMP::em::get_segment ( DensityMap *  map_to_segment,
int  nx_start,
int  nx_end,
int  ny_start,
int  ny_end,
int  nz_start,
int  nz_end 
)

Get a segment of the map according to xyz indexes.

Note
the output map will cover the region [[nx_start,nx_end],[]ny_start,ny_end,[nz_start,nz_end]]
DensityMap* IMP::em::get_segment ( DensityMap *  map_to_segment,
algebra::Vector3Ds  vecs,
float  dist 
)

Get a segment of the map covered by the input points.

DensityMap* IMP::em::get_segment_by_masking ( DensityMap *  map_to_segment,
DensityMap *  mask,
float  mas_threshold 
)

Get a segment of the map covered by another map.

double IMP::em::get_sum ( const DensityMap *  m1)

Return the sum of all voxels.

Float IMP::em::get_threshold_for_approximate_mass ( DensityMap *  m,
Float  desired_mass,
atom::ProteinDensityReference  ref = atom::HARPAZ 
)

Computes the threshold to consider in an EM map to get a desired mass.

Computes the threshold to consider in an EM map to get a desired mass (only voxels with intensity greater than the threshold are considered)

Parameters
[in]ma density map
[in]desired_mass(in Da)
[in]refthe protein density reference to use in the computation. The default protein density for this computation is HARPAZ
Float IMP::em::get_threshold_for_approximate_volume ( DensityMap *  m,
Float  desired_volume 
)

Computes the threshold consider in an EM map to get a desired volume.

Computes the threshold consider in an EM map to get a desired volume (i.e, the set of voxels with intensity greater than the threshold occupies that volume)

Parameters
[in]ma density map
[in]desired_volume(in A^3)
DensityMap* IMP::em::get_threshold_map ( const DensityMap *  orig_map,
float  threshold 
)

Return a map with 0 for all voxels below the threshold.

Parameters
[in]orig_mapthe map to binarize
[in]thresholdvalues below the threshold are set to 0 and 1 otherwise
Float IMP::em::get_volume_at_threshold ( DensityMap *  m,
Float  threshold 
)

Compute an approximate volume.

Compute an approximate volume for the set of voxels with intensity under a given threshold

Parameters
[in]ma density map
[in]thresholdconsider volume of only voxels above this threshold
Returns
a volume for the set of voxels with intensity under the provided threshold
void IMP::em::ImageHeader_to_DensityHeader ( const ImageHeader &  h,
DensityHeader &  header 
)

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

DensityMap* IMP::em::interpolate_map ( DensityMap *  in_map,
double  new_spacing 
)

Return a new map with an updated spacing.

Input data is interpolated accordingly.

FittingSolutions IMP::em::local_rigid_fitting ( Particle *  p,
Refiner *  refiner,
const FloatKey &  weight_key,
DensityMap *  dmap,
OptimizerStates  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 = .3,
bool  fast = true 
)

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-correlation score. The optimization is a standard MC/CG procedure. The function returns a list of solutions sorted by the cross-correlation score.

Note
The returned cross-correlation score is 1-cc, as we usually want to minimize a scoring 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]pThe root of the hierarchy to fit
[in]refinerThe refiner to get the leaves of the particle
[in]weight_keyThe weight key of the particles in the rigid body
[in]dmapThe density map to fit to
[in]display_logIf provided, then intermediate states during the optimization procedure are printed
[in]number_of_optimization_runsnumber of Monte Carlo optimizations
[in]number_of_mc_stepsnumber of steps in a Monte Carlo optimization
[in]number_of_cg_stepsnumber of Conjugate Gradients steps in a Monte Carlo step
[in]max_translationmaximum translation step in a MC optimization step
[in]max_rotationmaximum rotation step in radians in a single MC optimization step
[in]fastif true the density map of the rigid body is not resampled but transformed at each iteration of the optimization
Returns
the refined fitting solutions

Definition at line 168 of file rigid_fitting.h.

+ Here is the call graph for this function:

FittingSolutions IMP::em::local_rigid_fitting_around_point ( Particle *  p,
Refiner *  refiner,
const FloatKey &  weight_key,
DensityMap *  dmap,
const algebra::Vector3D &  anchor_centroid,
OptimizerStates  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 = .3,
bool  fast = false 
)

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-correlation score. The optimization is a standard MC/CG procedure. The function returns a list of solutions sorted by the cross-correlation score.

Note
The returned cross-correlation score is 1-cc, as we usually want to minimize a scoring 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]pThe rigid body to fit
[in]refinerRefiner to yield rigid body members
[in]weight_keyThe weight key of the particles in the rigid body
[in]dmapThe density map to fit to
[in]anchor_centroidThe point to fit the particles around
[in]display_logIf provided, then intermediate states in during the optimization procedure are printed
[in]number_of_optimization_runsnumber of Monte Carlo optimizations
[in]number_of_mc_stepsnumber of steps in a Monte Carlo optimization
[in]number_of_cg_stepsnumber of Conjugate Gradients steps in a Monte Carlo step
[in]max_translationmaximum translation step in a MC optimization step
[in]max_rotationmaximum rotation step in a single MC optimization step
[in]fastif true the density map of the rigid body is not resampled but transformed at each iteration of the optimization
Returns
the refined fitting solutions
FittingSolutions IMP::em::local_rigid_fitting_around_points ( Particle *  p,
Refiner *  refiner,
const FloatKey &  wei_key,
DensityMap *  dmap,
const algebra::Vector3Ds &  anchor_centroids,
OptimizerStates  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 = .3 
)

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]pThe rigid body to fit
[in]refinerRefiner to yield rigid body members
[in]wei_keyThe weight key of the particles in the rigid body
[in]dmapThe density map to fit to
[in]anchor_centroidsThe points to fit the particles around
[in]display_logIf provided, then intermediate states in during the optimization procedure are printed
[in]number_of_optimization_runsnumber of Monte Carlo optimizations
[in]number_of_mc_stepsnumber of steps in a Monte Carlo optimization
[in]number_of_cg_stepsnumber of Conjugate Gradients steps in a Monte Carlo step
[in]max_translationmaximum translation step in a MC optimization step
[in]max_rotationmaximum rotation step in a single MC optimization step
Returns
the refined fitting solutions
FittingSolutions IMP::em::local_rigid_fitting_grid_search ( const ParticlesTemp &  ps,
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-correlation 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 scoring function. Thus a score of 1 means no-correlation and a score of 0. is perfect correlation.
Parameters
[in]psThe particles to be fitted (treated rigid)
[in]wei_keyThe weight key of the particles in the rigid body
[in]dmapThe density map to fit to
[in]max_voxels_translationSample translations within -max_voxel_translation to max_voxel_translation
[in]translation_stepThe translation sampling step
[in]max_angle_in_radiansSample rotations with +- max_angle_in_radians around the current rotation
[in]number_of_rotationsThe number of rotations to sample
Returns
the refined fitting solutions
DensityMap* IMP::em::mask_and_norm ( em::DensityMap *  dmap,
em::DensityMap *  mask 
)

Return a masked density, and normalize the output map within the mask region.

Parameters
[in]dmapthe density map to mask
[in]maskthe mask
Returns
the masked and normalized map
DensityMap* IMP::em::multiply ( const DensityMap *  m1,
const DensityMap *  m2 
)

Return a density map for which voxel i contains the result of m1[i]*m2[i].

The function assumes m1 and m2 are of the same dimensions.

SampledDensityMap* IMP::em::particles2binarized_density ( const ParticlesTemp &  ps,
Float  resolution,
Float  apix,
int  sig_cutoff = 3,
const FloatKey &  weight_key = IMP::atom::Mass::get_mass_key() 
)

Resample a set of particles into a binarized density map 1 for voxels containing particles and 0 otherwise Each such particle should have xyz radius and weight attributes

Parameters
[in]psthe particles to sample
[in]resolutionthe resolution of the new sampled map
[in]apixthe voxel size of the sampled map
[in]sig_cutoffsigma cutoff used in sampling
[in]weight_keythe weight attribute key of the particles
Returns
the sampled density grid
See Also
SampledDensityMap

Definition at line 76 of file converters.h.

+ Here is the call graph for this function:

SampledDensityMap* IMP::em::particles2density ( const ParticlesTemp &  ps,
Float  resolution,
Float  apix,
int  sig_cutoff = 3,
const FloatKey &  weight_key = IMP::atom::Mass::get_mass_key() 
)

Resample a set of particles into a density grid.

Each such particle should have xyz, radius and weight attributes

Parameters
[in]psthe particles to sample
[in]resolutionthe resolution of the new sampled map
[in]apixthe voxel size of the sampled map
[in]sig_cutoffsigma cutoff used in sampling
[in]weight_keythe weight attribute key of the particles
Returns
the sampled density grid
See Also
SampledDensityMap
SurfaceShellDensityMap* IMP::em::particles2surface ( const ParticlesTemp &  ps,
Float  apix,
const FloatKey &  weight_key = IMP::atom::Mass::get_mass_key() 
)

Resample a set of particles into a density grid.

Each such particle should have xyz radius and weight attributes

Parameters
[in]psthe particles to sample
[in]apixthe voxel size of the surface map
[in]weight_keythe weight attribute key of the particles
Returns
the surface grid
See Also
SampledDensityMap
def IMP.em.write_pca_cmm (   pca,
  fh 
)

Write out principal components to a file in Chimera Marker format.

Note
This function is only available in Python.

Definition at line 2866 of file em/__init__.py.