IMP
2.2.0
The Integrative Modeling Platform
|
See IMP.em for more information.
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 Microsocopy images in IMP. More... | |
class | KernelParameters |
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 |
class | PCAFitRestraint |
Calculate score based on fit to EM map. More... | |
class | RadiusDependentDistanceMask |
class | RadiusDependentKernelParameters |
Calculates kernel parameters as a function of a specific radius. More... | |
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 repersents a molecule as shells of distance from the surface. More... | |
class | Voxel |
class | XplorReaderWriter |
Typedefs | |
typedef IMP::base::Vector < IMP::base::Pointer < DensityMap > > | DensityMaps |
rotate a grid More... | |
typedef IMP::base::Vector < IMP::base::WeakPointer < DensityMap > > | DensityMapsTemp |
typedef double | emreal |
typedef IMP::base::Vector < FittingSolutions > | FittingSolutionsList |
typedef IMP::base::Vector < KernelParameters > | KernelParametersList |
typedef IMP::base::Vector < RadiusDependentKernelParameters > | RadiusDependentKernelParametersList |
typedef IMP::base::Vector < IMP::base::Pointer < SampledDensityMap > > | SampledDensityMaps |
typedef IMP::base::Vector < IMP::base::WeakPointer < SampledDensityMap > > | SampledDensityMapsTemp |
typedef IMP::base::Vector < IMP::base::Pointer < SurfaceShellDensityMap > > | SurfaceShellDensityMaps |
typedef IMP::base::Vector < IMP::base::WeakPointer < SurfaceShellDensityMap > > | SurfaceShellDensityMapsTemp |
typedef IMP::base::Vector< Voxel > | Voxels |
Enumerations | |
enum | KernelType { GAUSSIAN, BINARIZED_SPHERE, SPHERE } |
Functions | |
void | add_to_map (DensityMap *dm, const kernel::Particles &pis) |
Float | approximate_molecular_mass (DensityMap *m, Float threshold) |
DensityMap * | binarize (DensityMap *orig_map, float threshold, bool reverse=false) |
Float | calculate_intersection_score (const SurfaceShellDensityMap *d1, const SurfaceShellDensityMap *d2) |
Float | compute_fitting_score (const kernel::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 kernel::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) |
DensityHeader | create_density_header (const algebra::BoundingBoxD< 3 > &bb, float spacing) |
Create a header from a bounding box 3D. More... | |
DensityMap * | create_density_map (const DensityMap *other) |
create a copy of another map | |
DensityMap * | create_density_map (const algebra::BoundingBox3D &bb, double spacing) |
Create an empty density map from a boudning box. | |
DensityMap * | create_density_map (int nx, int ny, int nz, double spacing) |
Create an empty density map. | |
template<class S , class V , class E > | |
DensityMap * | create_density_map (const IMP::algebra::GridD< 3, S, V, E > &arg) |
DensityMap * | create_density_map (const algebra::GridD< 3, algebra::DenseGridStorageD< 3, float >, float > &grid) |
kernel::Particles | density2particles (DensityMap *dmap, Float threshold, kernel::Model *m, int step=1) |
Converts a density grid to a set of paritlces. More... | |
algebra::Vector3Ds | density2vectors (DensityMap *dmap, Float threshold) |
Converts a density grid to a set of paritlces. More... | |
void | DensityHeader_to_ImageHeader (const DensityHeader &header, ImageHeader &h) |
double | EXP (float y) |
DensityMap * | get_binarized_interior (DensityMap *dmap) |
Return a binaries density map with 1 for voxels that are internal. | |
algebra::BoundingBoxD< 3 > | get_bounding_box (const DensityMap *m, Float threshold) |
algebra::BoundingBoxD< 3 > | get_bounding_box (const DensityMap *m) |
std::string | get_data_path (std::string file_name) |
Return the full path to installed data. More... | |
double | get_density (const DensityMap *m, const algebra::Vector3D &v) |
std::string | get_example_path (std::string file_name) |
Return the path to installed example data for this module. More... | |
algebra::GridD < 3, algebra::DenseGridStorageD < 3, float >, float > | get_grid (DensityMap *in_map) |
bool | get_interiors_intersect (const DensityMap *d1, const DensityMap *d2) |
DensityMap * | get_max_map (DensityMaps maps) |
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 kernel::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 kernel::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 kernel::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... | |
DensityMap * | get_resampled (DensityMap *input, double scaling) |
Get a resampled version of the map. More... | |
DensityMap * | 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. More... | |
DensityMap * | get_segment (DensityMap *map_to_segment, algebra::Vector3Ds vecs, float dist) |
Get a segment of the map covered by the input points. | |
DensityMap * | get_segment_by_masking (DensityMap *map_to_segment, DensityMap *mask, float mas_threshold) |
Get a segment of the map covered by another map. | |
double | get_sum (const DensityMap *m1) |
Return the sum of all voxels. | |
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... | |
DensityMap * | get_threshold_map (const DensityMap *orig_map, float threshold) |
DensityMap * | get_transformed (const DensityMap *input, const algebra::Transformation3D &tr, double threshold) |
DensityMap * | get_transformed (DensityMap *input, const algebra::Transformation3D &tr) |
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... | |
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) |
DensityMap * | interpolate_map (DensityMap *in_map, double new_spacing) |
FittingSolutions | local_rigid_fitting (kernel::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 (kernel::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 (kernel::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 kernel::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... | |
DensityMap * | mask_and_norm (em::DensityMap *dmap, em::DensityMap *mask) |
DensityMap * | multiply (const DensityMap *m1, const DensityMap *m2) |
SampledDensityMap * | particles2binarized_density (const kernel::ParticlesTemp &ps, Float resolution, Float apix, int sig_cutoff=3, const FloatKey &weight_key=IMP::atom::Mass::get_mass_key()) |
SampledDensityMap * | particles2density (const kernel::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... | |
SurfaceShellDensityMap * | particles2surface (const kernel::ParticlesTemp &ps, Float apix, const FloatKey &weight_key=IMP::atom::Mass::get_mass_key()) |
Resample a set of particles into a density grid. More... | |
DensityMap * | read_map (std::string filename, MapReaderWriter *reader) |
DensityMap * | read_map (std::string filename) |
void | write_map (DensityMap *m, std::string filename, MapReaderWriter *writer) |
void | write_map (DensityMap *m, std::string filename) |
def | write_pca_cmm |
Write out principal components to a file in Chimera Marker format. | |
Variables | |
const double | EPS = 1e-16 |
/param[in] orig_dens the density map to rotate /param[in] trans the transformation
Definition at line 508 of file DensityMap.h.
Pass a set of objects. See DensityMap
Definition at line 508 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 139 of file KernelParameters.h.
typedef IMP::base::Vector< RadiusDependentKernelParameters > IMP::em::RadiusDependentKernelParametersList |
Pass or store a set of RadiusDependentKernelParameters .
Definition at line 71 of file KernelParameters.h.
Store a set of objects.
Definition at line 131 of file SampledDensityMap.h.
typedef IMP::base::Vector<IMP::base::WeakPointer< SampledDensityMap > > IMP::em::SampledDensityMapsTemp |
Pass a set of objects. See SampledDensityMap
Definition at line 131 of file SampledDensityMap.h.
typedef IMP::base::Vector<IMP::base::Pointer< SurfaceShellDensityMap > > IMP::em::SurfaceShellDensityMaps |
Store a set of objects.
Definition at line 93 of file SurfaceShellDensityMap.h.
typedef IMP::base::Vector<IMP::base::WeakPointer< SurfaceShellDensityMap > > IMP::em::SurfaceShellDensityMapsTemp |
Pass a set of objects. See SurfaceShellDensityMap
Definition at line 93 of file SurfaceShellDensityMap.h.
void IMP::em::add_to_map | ( | DensityMap * | dm, |
const kernel::Particles & | pis | ||
) |
Rasterize the particles into an existing density map.
Float IMP::em::approximate_molecular_mass | ( | DensityMap * | m, |
Float | threshold | ||
) |
[in] | m | a density map |
[in] | threshold | consider volume of only voxels above this threshold |
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 thoes above
[in] | orig_map | the map to binarize |
[in] | threshold | values below the threshold are set to 0 and 1 otherwise |
[in] | reverse | if true values above the threshold are set to 0 and 1 otherwise |
Float IMP::em::compute_fitting_score | ( | const kernel::ParticlesTemp & | ps, |
DensityMap * | em_map, | ||
FloatKey | wei_key = atom::Mass::get_mass_key() |
||
) |
Scores how well a set of particles fit a map
[in] | ps | The particles to be fitted |
[in] | em_map | The density map to fit to |
[in] | wei_key | The weight key of the particles in the rigid body |
FittingSolutions IMP::em::compute_fitting_scores | ( | const kernel::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() |
||
) |
Score how well a set of particles fit to the map in various rigid transformations.
[in] | ps | The particles to be fitted (treated rigid) |
[in] | em_map | The density map to fit to |
[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 |
[in] | fast_version | if true the density map of each transformation is interpolated |
[in] | local_score | if true a local cross correlation score is used |
FittingSolutions IMP::em::compute_fitting_scores | ( | DensityMap * | em_map, |
core::RigidBody | rb, | ||
Refiner * | refiner, | ||
const algebra::Transformation3Ds & | transformations | ||
) |
Score how well a rigid body fits to the map
[in] | em_map | The density map to fit to |
[in] | rb | The rigid body |
[in] | refiner | The rigid body refiner |
[in] | transformations | Transformations of the rigid body |
Definition at line 272 of file rigid_fitting.h.
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 | ||
) |
See DensityHeader
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 637 of file DensityMap.h.
DensityMap* IMP::em::create_density_map | ( | const algebra::GridD< 3, algebra::DenseGridStorageD< 3, float >, float > & | grid | ) |
Return a density map with the values taken from the grid.
kernel::Particles IMP::em::density2particles | ( | DensityMap * | dmap, |
Float | threshold, | ||
kernel::Model * | m, | ||
int | step = 1 |
||
) |
Each such particle will be have xyz attributes and a density_val attribute of type Float.
[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 |
[in] | step | sample every X steps in each direction |
algebra::Vector3Ds IMP::em::density2vectors | ( | DensityMap * | dmap, |
Float | threshold | ||
) |
Each such particle will be have xyz attributes and a density_val attribute of type Float.
[in] | dmap | the density map |
[in] | threshold | only voxels with density above the given threshold will be converted to particles |
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 | ||
) |
[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 | ) |
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
This will ensure that the code works when IMP
is installed or used via the setup_environment.sh
script.
double IMP::em::get_density | ( | const DensityMap * | m, |
const algebra::Vector3D & | 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. See DensityMap
std::string IMP::em::get_example_path | ( | std::string | file_name | ) |
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
This will ensure that the code works when IMP
is installed or used via the setup_environment.sh
script.
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 as well as the same bounding box.
DensityMap* IMP::em::get_max_map | ( | DensityMaps | maps | ) |
Return a density map for which each voxel is the maximum value from the input densities.
Float IMP::em::get_molecular_mass_at_threshold | ( | DensityMap * | m, |
Float | threshold, | ||
atom::ProteinDensityReference | ref = atom::HARPAZ |
||
) |
Compute an approximate molecular mass for the set of voxels with intensity under a given threshold
[in] | m | a density map |
[in] | threshold | only voxels above this threshold will be considered |
[in] | ref | the protein density reference to use in the computation. The default protein density for this computation is HARPAZ |
long IMP::em::get_number_of_particles_outside_of_the_density | ( | DensityMap * | dmap, |
const kernel::Particles & | ps, | ||
const IMP::algebra::Transformation3D & | t = IMP::algebra::get_identity_transformation_3d() , |
||
float | thr = 0.0 |
||
) |
/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 kernel::Particles & | ps, | ||
const IMP::algebra::Transformation3Ds & | transformations, | ||
float | thr = 0.0 |
||
) |
/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 kernel::Particles & | ps, | ||
float | smoothing_radius = 3. , |
||
const IMP::algebra::Transformation3D & | t = IMP::algebra::get_identity_transformation_3d() , |
||
float | thr = 0.0 |
||
) |
/note the function assumes that all of the particles have XYZ coordinates
DensityMap* IMP::em::get_resampled | ( | DensityMap * | input, |
double | scaling | ||
) |
The spacing is multiplied by scaling. That means, scaling values greater than 1 increase the voxel size. See DensityMap
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 | ||
) |
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 (only voxels with intensity greater than the threshold are considered)
[in] | m | a density map |
[in] | desired_mass | (in Da) |
[in] | ref | the 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 (i.e, the set of voxels with intensity greater than the threshold occupies that volume)
[in] | m | a 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
[in] | orig_map | the map to binarize |
[in] | threshold | values below the threshold are set to 0 and 1 otherwise |
DensityMap* IMP::em::get_transformed | ( | const DensityMap * | input, |
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()). See DensityMap
DensityMap* IMP::em::get_transformed | ( | DensityMap * | input, |
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. See DensityMap
void IMP::em::get_transformed_into | ( | const DensityMap * | source, |
const algebra::Transformation3D & | tr, | ||
DensityMap * | into, | ||
bool | calc_rms = true |
||
) |
[in] | source | the map to transform |
[in] | tr | transform the from density map by this transformation |
[out] | into | the map to tranform into |
[in] | calc_rms | if true RMS is calculated on the transformed map See DensityMap |
Float IMP::em::get_volume_at_threshold | ( | DensityMap * | m, |
Float | threshold | ||
) |
Compute an approximate volume for the set of voxels with intensity under a given threshold
[in] | m | a density map |
[in] | threshold | consider volume of only voxels above this 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 and interpolate input data accordinly
FittingSolutions IMP::em::local_rigid_fitting | ( | kernel::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 |
||
) |
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.
[in] | p | The root of the hierarchy to fit |
[in] | refiner | The refiner to get the leaves of the particle |
[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 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 radians in a single MC optimization step |
[in] | fast | if true the density map of the rigid body is not resampled but transformed at each iteration of the optimization |
Definition at line 168 of file rigid_fitting.h.
FittingSolutions IMP::em::local_rigid_fitting_around_point | ( | kernel::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 |
||
) |
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.
[in] | p | The rigid body to fit |
[in] | refiner | Refiner to yield rigid body members |
[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 |
[in] | fast | if true the density map of the rigid body is not resampled but transformed at each iteration of the optimization |
FittingSolutions IMP::em::local_rigid_fitting_around_points | ( | kernel::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 |
||
) |
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.
[in] | p | The rigid body to fit |
[in] | refiner | Refiner to yield rigid body members |
[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 |
FittingSolutions IMP::em::local_rigid_fitting_grid_search | ( | const kernel::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 |
||
) |
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
[in] | ps | The particles to be fitted (treated rigid) |
[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 |
DensityMap* IMP::em::mask_and_norm | ( | em::DensityMap * | dmap, |
em::DensityMap * | mask | ||
) |
Return a masked density , and normalize the output map within the masked region
[in] | dmap | the density map to mask |
[in] | mask | the mask |
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 kernel::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
[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] | weight_key | the weight attribute key of the particles |
Definition at line 80 of file converters.h.
SampledDensityMap* IMP::em::particles2density | ( | const kernel::ParticlesTemp & | ps, |
Float | resolution, | ||
Float | apix, | ||
int | sig_cutoff = 3 , |
||
const FloatKey & | weight_key = IMP::atom::Mass::get_mass_key() |
||
) |
Each such particle should have xyz, radius and weight attributes
[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] | weight_key | the weight attribute key of the particles |
SurfaceShellDensityMap* IMP::em::particles2surface | ( | const kernel::ParticlesTemp & | ps, |
Float | apix, | ||
const FloatKey & | weight_key = IMP::atom::Mass::get_mass_key() |
||
) |
Each such particle should have xyz radius and weight attributes
[in] | ps | the particles to sample |
[in] | apix | the voxel size of the surface map |
[in] | weight_key | the weight attribute key of the particles |
DensityMap* IMP::em::read_map | ( | std::string | filename, |
MapReaderWriter * | reader | ||
) |
Read a density map from a file and return it. See DensityMap
DensityMap* IMP::em::read_map | ( | std::string | filename | ) |
Read a density map from a file and return it. Guess the file type from the file name. The file formats supported are:
void IMP::em::write_map | ( | DensityMap * | m, |
std::string | filename, | ||
MapReaderWriter * | writer | ||
) |
Write a density map to a file. See DensityMap
void IMP::em::write_map | ( | DensityMap * | m, |
std::string | filename | ||
) |
Write a density map to a file. Guess the file type from the file name. The file formats supported are: