IMP logo
IMP Reference Guide  2.16.0
The Integrative Modeling Platform
pdb2density.py
1 ## \example em/pdb2density.py
2 # A simple example showing how to simulate density from a protein.
3 # IMP uses a Gaussian smoothing kernel. see SampledDensityMap::resample
4 # for documentation.
5 #
6 
7 import IMP.em
8 import IMP.core
9 import IMP.atom
10 import sys
11 
12 IMP.setup_from_argv(sys.argv, "pdb2density")
13 
14 m = IMP.Model()
15 # read protein
17 mh = IMP.atom.read_pdb(IMP.em.get_example_path("input.pdb"), m, sel)
18 # add radius info to each atom, otherwise the resampling would fail.
20 ps = IMP.core.get_leaves(mh)
21 # decide on resolution and spacing you would like to simulate to
22 resolution = 10.
23 apix = 1.5
24 dmap = IMP.em.particles2density(ps, resolution, apix)
25 # write out the map in the favorite format (xplor, mrc, em and spider are
26 # supported)
27 IMP.em.write_map(dmap, "example.mrc", IMP.em.MRCReaderWriter())
Strings setup_from_argv(const Strings &argv, std::string description, std::string positional_description, int num_positional)
void add_radii(Hierarchy d, const ForceFieldParameters *ffp=get_all_atom_CHARMM_parameters(), FloatKey radius_key=FloatKey("radius"))
SampledDensityMap * 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.
GenericHierarchies get_leaves(Hierarchy mhd)
Get all the leaves of the bit of hierarchy.
void read_pdb(TextInput input, int model, Hierarchy h)
Class for storing model, its restraints, constraints, and particles.
Definition: Model.h:73
Basic utilities for handling cryo-electron microscopy 3D density maps.
std::string get_example_path(std::string file_name)
Return the full path to one of this module's example files.
Basic functionality that is expected to be used by a wide variety of IMP users.
Select all non-water non-alternative ATOM and HETATM records.
Definition: pdb.h:287
Functionality for loading, creating, manipulating and scoring atomic structures.