IMP Reference Guide
2.10.0
The Integrative Modeling Platform
|
The npctransport module is a module for simulating transport through the NPC.
Authors: Barak Raveh, Daniel Russel
License:
Publications: Timney*, Raveh* et al., JCB 2016
Building from source in a nutshell (see https://integrativemodeling.org/latest/doc/manual/installation.html for general IMP installation instruction):
Versions: 4.5
The fg_simulation
command line tool can be used to run simulations of NPC FG repeat domains using this module. For an example of its use, see our 2018 study in Nature.
Functions | |
void | add_hierarchies_with_sites (RMF::FileHandle fh, const atom::Hierarchies &hs) |
void | add_hierarchies_with_sites (RMF::NodeHandle fh, const atom::Hierarchies &hs) |
void | add_hierarchy_with_sites (RMF::FileHandle fh, atom::Hierarchy hs) |
void | add_test_sites (RMF::FileHandle fh, core::ParticleType t, double display_radius, algebra::Vector3Ds sites) |
void | add_test_sites (RMF::FileHandle fh, core::ParticleType t, algebra::Sphere3Ds sites) |
int | assign_ranges (std::string input_config_fname, std::string output_assignment_fname, unsigned int work_unit, bool show_steps, boost::uint64_t random_seed) |
void | copy_FGs_coordinates (SimulationData const *src_sd, SimulationData *trg_sd) |
copy coordinates of src_sd to trg_sd for FG repeats only More... | |
void | copy_hierarchy_reference_frame_recursive (Particle *src_p, Particle *trg_p) |
void | copy_particle_reference_frame_if_applicable (Particle *src_p, Particle *trg_p) |
boost::timer | create_boost_timer () |
FGChain * | create_fg_chain (IMP::npctransport::SimulationData *sd, atom::Hierarchy parent, const ::npctransport_proto::Assignment_FGAssignment &fg_data, display::Color c) |
atom::Hierarchies | create_hierarchies_with_sites (RMF::FileConstHandle fh, Model *m) |
double | do_evaluate_index (algebra::Sphere3D &d_xyzr0, algebra::Sphere3D &d_xyzr1, DerivativeAccumulator *da, const algebra::Vector3D &delta, double delta_length, double x0, double k) |
double | do_evaluate_index_harmonic (Model *m, const ParticleIndexPair &pp, DerivativeAccumulator *da, const algebra::Vector3D &delta, double delta_length, double x0, double k) |
void | do_main_loop (SimulationData *sd, const RestraintsTemp &init_restraints) |
unsigned int | find_or_add_fg_bead_of_type (::npctransport_proto::Statistics *s, IMP::core::ParticleType pt) |
unsigned int | find_or_add_fg_chain_of_type (::npctransport_proto::Statistics *s, IMP::core::ParticleType pt) |
unsigned int | find_or_add_floater_of_type (::npctransport_proto::Statistics *s, IMP::core::ParticleType pt) |
unsigned int | find_or_add_interaction_of_type (::npctransport_proto::Statistics *s, IMP::npctransport::InteractionType it) |
internal_avro::ValidSchema | get_avro_data_file_schema () |
double | get_close_pairs_range (double max_range, double max_range_factor) |
double | get_close_pairs_range (const ::npctransport_proto::Assignment &config) |
int | get_dump_interval_in_frames (const ::npctransport_proto::Assignment &config, double time_step) |
FGChain * | get_fg_chain (atom::Hierarchy root) |
FGChain * | get_fg_chain (Particle *p_root) |
int | get_frames_from_ns (double ns, double time_step) |
algebra::Vector3D | get_global_from_local_v3 (Particle *p, const algebra::Vector3D &local) |
unsigned int | get_maximal_number_of_unordered_pairs (ParticlesTemp const &ps0, ParticlesTemp const &ps1) |
template<class t_ordered_set > | |
boost::tuple< unsigned int, unsigned int > | get_n_lost_and_gained (t_ordered_set old, t_ordered_set cur) |
ParticlesTemp | get_non_optimizable_particles (ParticlesTemp const &particles) |
int | get_number_of_frames (const ::npctransport_proto::Assignment &config, double time_step) |
int | get_number_of_work_units (std::string configuration_file) |
ParticlesTemp | get_optimizable_particles (ParticlesTemp const &particles) |
int | get_output_statistics_interval_in_frames (const ::npctransport_proto::Assignment &assign, double time_step, double default_value_ns=1.0) |
ParticleIndexes | get_particle_indexes (ParticlesTemp const &particles) |
void | get_protobuf_configuration_from_text (std::string config_txt, std::string config_pb) |
WeakObjectKey | get_simulation_data_key () |
algebra::Vector3Ds | get_spheres_centers (algebra::Sphere3Ds const &spheres) |
template<typename V3iter > | |
algebra::Sphere3Ds | get_spheres_from_vectors (V3iter first, V3iter last, double radius) |
convert vectors to spheres of passed radius More... | |
algebra::Sphere3Ds | get_spheres_from_vectors (algebra::Vector3Ds const &vs, double radius) |
convert vectors to spheres of passed radius More... | |
int | get_statistics_interval_in_frames (const ::npctransport_proto::Assignment &assign, double time_step, double default_value_ns=0.1) |
double | get_time_step (double max_d_factor, double max_k, double min_radius, double min_range, double max_trans_relative_to_radius=0.1, double time_step_factor=1.0) |
double | get_time_step (const ::npctransport_proto::Assignment &config, double max_trans_relative_to_radius=0.1) |
void | inflate_floater (SimulationData *sd, const std::string floater_name, const float new_radius) |
inflate floater of specified type to new_radius More... | |
void | initialize_positions (SimulationData *sd, const RestraintsTemp &extra_restraints=RestraintsTemp(), bool debug=false, double short_init_factor=1.0, bool is_disable_randomize=false, bool are_fgs_pre_initialized=false) |
void | link_hierarchies_with_sites (RMF::FileConstHandle fh, const atom::Hierarchies &hs) |
bool | load_output_protobuf (std::string output_fname,::npctransport_proto::Output &output) |
void | load_pb_conformation (const ::npctransport_proto::Conformation &conformation, IMP::SingletonContainerAdaptor beads, boost::unordered_map< core::ParticleType, algebra::Sphere3Ds > &sites) |
InteractionType | make_ordered_interaction_type (IMP::core::ParticleType t0, IMP::core::ParticleType t1) |
convenience method for creating and interaction type in swig More... | |
InteractionType | make_unordered_interaction_type (IMP::core::ParticleType t0, IMP::core::ParticleType t1) |
template<class t_value > | |
std::pair< t_value, t_value > | make_unordered_pair (t_value v0, t_value v1) |
Canonize such that v0>=v1 so order doesn't matter. More... | |
IMP::ParticleIndexPair | make_unordered_particle_index_pair (IMP::ParticleIndex pi0, IMP::ParticleIndex pi1) |
ParticleIndexPair | make_unordered_particle_index_pair (ParticleIndexPair pip) |
template<class ParticlesList , class BoundingVolume > | |
void | randomize_particles (const ParticlesList &ps, const BoundingVolume &bv) |
template<class RigidBody , class BoundingVolume > | |
void | randomize_rigid_body (RigidBody rbi, const BoundingVolume &bv) |
void | remove_Nup42 (SimulationData *sd) |
remove nup42 and its anchors (also from obstacles) More... | |
void | reset_box_size (SimulationData *sd, double box_size) |
change box size sd to specified box size and update output file More... | |
void | save_pb_conformation (IMP::SingletonContainerAdaptor beads, const boost::unordered_map< core::ParticleType, algebra::Sphere3Ds > &sites,::npctransport_proto::Conformation *conformation) |
void | show_ranges (std::string fname) |
IMP::npctransport::SimulationData * | startup (int argc, char *argv[]) |
void | write_geometry (const ParticlesTemp &kaps, const algebra::Sphere3Ds &kap_sites, const ParticlesTemp &chains, const algebra::Sphere3Ds &chain_sites, const RestraintsTemp &rs, display::Writer *out) |
Variables | |
const double | FS_IN_NS = 1.0E+6 |
const double | HALF_SQRT_MAX_DOUBLE = 0.5 * std::sqrt( MAX_DOUBLE ) |
const double | MAX_DOUBLE = std::numeric_limits< double >::max() |
Standard module functions | |
All | |
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 IMP::Vector< Avro2PBReader > IMP::npctransport::Avro2PBReaders |
Pass or store a set of Avro2PBReader .
Definition at line 90 of file Avro2PBReader.h.
typedef IMP::Vector<IMP::Pointer< BipartitePairsStatisticsOptimizerState > > IMP::npctransport::BipartitePairsStatisticsOptimizerStates |
A vector of reference-counting object pointers.
Definition at line 260 of file BipartitePairsStatisticsOptimizerState.h.
typedef IMP::Vector<IMP::WeakPointer< BipartitePairsStatisticsOptimizerState > > IMP::npctransport::BipartitePairsStatisticsOptimizerStatesTemp |
A vector of weak (non reference-counting) pointers to specified objects.
Definition at line 260 of file BipartitePairsStatisticsOptimizerState.h.
typedef IMP::Vector<IMP::Pointer< BodyStatisticsOptimizerState > > IMP::npctransport::BodyStatisticsOptimizerStates |
A vector of reference-counting object pointers.
Definition at line 70 of file BodyStatisticsOptimizerState.h.
typedef IMP::Vector<IMP::WeakPointer< BodyStatisticsOptimizerState > > IMP::npctransport::BodyStatisticsOptimizerStatesTemp |
A vector of weak (non reference-counting) pointers to specified objects.
Definition at line 70 of file BodyStatisticsOptimizerState.h.
typedef IMP::Vector<IMP::Pointer< ChainStatisticsOptimizerState > > IMP::npctransport::ChainStatisticsOptimizerStates |
A vector of reference-counting object pointers.
Definition at line 115 of file ChainStatisticsOptimizerState.h.
typedef IMP::Vector<IMP::WeakPointer< ChainStatisticsOptimizerState > > IMP::npctransport::ChainStatisticsOptimizerStatesTemp |
A vector of weak (non reference-counting) pointers to specified objects.
Definition at line 115 of file ChainStatisticsOptimizerState.h.
typedef IMP::Vector<IMP::Pointer< FGChain > > IMP::npctransport::FGChains |
typedef IMP::Vector<IMP::Pointer< GlobalStatisticsOptimizerState > > IMP::npctransport::GlobalStatisticsOptimizerStates |
A vector of reference-counting object pointers.
Definition at line 52 of file GlobalStatisticsOptimizerState.h.
typedef IMP::Vector<IMP::WeakPointer< GlobalStatisticsOptimizerState > > IMP::npctransport::GlobalStatisticsOptimizerStatesTemp |
A vector of weak (non reference-counting) pointers to specified objects.
Definition at line 52 of file GlobalStatisticsOptimizerState.h.
typedef IMP::Vector<IMP::Pointer< HarmonicSpringSingletonScore > > IMP::npctransport::HarmonicSpringSingletonScores |
A vector of reference-counting object pointers.
Definition at line 152 of file HarmonicSpringSingletonScore.h.
typedef IMP::Vector<IMP::WeakPointer< HarmonicSpringSingletonScore > > IMP::npctransport::HarmonicSpringSingletonScoresTemp |
A vector of weak (non reference-counting) pointers to specified objects.
Definition at line 152 of file HarmonicSpringSingletonScore.h.
typedef IMP::Vector<IMP::Pointer< HarmonicWellPairScore > > IMP::npctransport::HarmonicWellPairScores |
A vector of reference-counting object pointers.
Definition at line 114 of file harmonic_distance_pair_scores.h.
typedef IMP::Vector<IMP::WeakPointer< HarmonicWellPairScore > > IMP::npctransport::HarmonicWellPairScoresTemp |
A vector of weak (non reference-counting) pointers to specified objects.
Definition at line 114 of file harmonic_distance_pair_scores.h.
typedef std::pair<IMP::core::ParticleType, IMP::core::ParticleType> IMP::npctransport::InteractionType |
an interaction that involves particles of two types
Definition at line 21 of file typedefs.h.
Pass or store a set of InteractionType .
Definition at line 22 of file typedefs.h.
A vector of reference-counting object pointers.
Definition at line 85 of file functor_linear_distance_pair_scores_typedefs.h.
typedef IMP::Vector<IMP::WeakPointer< LinearInteraction > > IMP::npctransport::LinearInteractionsTemp |
A vector of weak (non reference-counting) pointers to specified objects.
Definition at line 85 of file functor_linear_distance_pair_scores_typedefs.h.
A vector of reference-counting object pointers.
Definition at line 468 of file linear_distance_pair_scores.h.
typedef IMP::Vector<IMP::WeakPointer< LinearWellPairScore > > IMP::npctransport::LinearWellPairScoresTemp |
A vector of weak (non reference-counting) pointers to specified objects.
Definition at line 468 of file linear_distance_pair_scores.h.
typedef IMP::Vector<IMP::Pointer< ParticleTransportStatisticsOptimizerState > > IMP::npctransport::ParticleTransportStatisticsOptimizerStates |
A vector of reference-counting object pointers.
Definition at line 113 of file ParticleTransportStatisticsOptimizerState.h.
typedef IMP::Vector<IMP::WeakPointer< ParticleTransportStatisticsOptimizerState > > IMP::npctransport::ParticleTransportStatisticsOptimizerStatesTemp |
A vector of weak (non reference-counting) pointers to specified objects.
Definition at line 113 of file ParticleTransportStatisticsOptimizerState.h.
Pass or store a set of SitesPairScoreParameters .
Definition at line 97 of file SitesPairScoreParameters.h.
void IMP::npctransport::add_hierarchies_with_sites | ( | RMF::FileHandle | fh, |
const atom::Hierarchies & | hs | ||
) |
Functions for adding a hierarchy to an RMF file or linking an RMF file to an existing hierarchy, including support for particles with sites.
These functions are practically identical to the add / link methods in modules/rmf/include/hierarchy_io.h, such as IMP::rmf::link_hierarchies(), except here NPC particle sites support is included.Add objects to the file.
void IMP::npctransport::add_hierarchies_with_sites | ( | RMF::NodeHandle | fh, |
const atom::Hierarchies & | hs | ||
) |
Add objects to the file under the specified node.
void IMP::npctransport::add_hierarchy_with_sites | ( | RMF::FileHandle | fh, |
atom::Hierarchy | hs | ||
) |
Add a single HierarchyWithSites object to the RMF file.
void IMP::npctransport::add_test_sites | ( | RMF::FileHandle | fh, |
core::ParticleType | t, | ||
double | display_radius, | ||
algebra::Vector3Ds | sites | ||
) |
for testing - adds the list of sites with specified radius, to be associated with particle type t. The file handle fh relies on this list only if it doesn't have particles with simulation data keys
void IMP::npctransport::add_test_sites | ( | RMF::FileHandle | fh, |
core::ParticleType | t, | ||
algebra::Sphere3Ds | sites | ||
) |
for testing - adds the list of sites with specified display radius, to be associated with particle type t. The file handle fh relies on this list only if it doesn't have particles with simulation data keys
int IMP::npctransport::assign_ranges | ( | std::string | input_config_fname, |
std::string | output_assignment_fname, | ||
unsigned int | work_unit, | ||
bool | show_steps, | ||
boost::uint64_t | random_seed | ||
) |
reads the protobuf message in [input_config_fname], which may contain submessages with ranges of values (indicated by presence of .upper and .lower fields, with .steps possible steps for each such field). The output is a message with the [work_unit]'th possible combination of these ranges, to the file [output_assignment_fname].
Note: the range values are enumerated as if they lie on a grid with log-evenly distributed axis-aligned grid points, using the .base field as the log base for each ranged field, such that e.g. iterating over the range [1..8] with 3 steps and base 2 will be enumerated as (1,4,8)
fname | input_config_fname configuration file name |
output | output_assignment_fname assignment file name |
work_unit | the index of combination of range values to be used. If the total of possible combinations of all fields with ranges is k, it is guaranteed that iterating over work_unit between 0..k-1 will enumerate over all possible combinations, and that work_unit and (work_unit % k) will return the same output for the same input. |
show_steps | show the steps that occur |
random_seed | the random seed used to initialize the IMP random number generator for this simulation |
IMP::ValueException | if any of the values in the configuration file are in conflict (e.g., simulation time and maximal number of frames) |
void IMP::npctransport::copy_FGs_coordinates | ( | SimulationData const * | src_sd, |
SimulationData * | trg_sd | ||
) |
copy coordinates of src_sd to trg_sd for FG repeats only
void IMP::npctransport::copy_hierarchy_reference_frame_recursive | ( | Particle * | src_p, |
Particle * | trg_p | ||
) |
Copy XYZ coordinates or RigidBody reference frame from src_pi to trg_pi if applicable, and if src_pi and trg_pi are an atom hierarchy, proceed recursively to their children. If so, assumes identical topology of hierarchies for src_pi and trg_pi
void IMP::npctransport::copy_particle_reference_frame_if_applicable | ( | Particle * | src_p, |
Particle * | trg_p | ||
) |
Copy XYZ coordinates or RigidBody reference frame from src_p to trg_p, if it is decorated with XYZ or RigidBody. Do nothing otherwise.
FGChain* IMP::npctransport::create_fg_chain | ( | IMP::npctransport::SimulationData * | sd, |
atom::Hierarchy | parent, | ||
const ::npctransport_proto::Assignment_FGAssignment & | fg_data, | ||
display::Color | c | ||
) |
Create a chain particle hierarchy, to be owned by the model of sd, with restraint bonding consecutive particles added to sd, according to the parameters specified in fg_data.
Notes:
The type of the chain root is specified in fg_data.type(). Individual particles may have different types if fg_data.type_suffix_list is non-empty. In this case, the type of particle i is suffixed with fg_data.type_suffix_list[i], or remains fg_data.type() if the suffix is "".
The rest length between two consecutive chain beads is fg_data.radius() * 2.0 * fg_data.rest_length_factor() and the spring constant is the simulation backbone_k parameter.
If fg_data.is_tamd() is true, created a TAMD hierarchy, otherwise a simple parent + beads structure. In the TAMD case, the custom restraint are added to sd->get_scoring() and the tamd images are added to sd->root()
[in,out] | sd | the simulation data whose model is associated with the new chain. The chain is also added to the simulation data scoring object. |
parent | parent hierarchy to which chain is added | |
[in] | fg_data | data about the FG chain |
[in] | c | color of chain particles |
atom::Hierarchies IMP::npctransport::create_hierarchies_with_sites | ( | RMF::FileConstHandle | fh, |
Model * | m | ||
) |
Create HierarchyWithSites objects from the RMF file.
double IMP::npctransport::do_evaluate_index | ( | algebra::Sphere3D & | d_xyzr0, |
algebra::Sphere3D & | d_xyzr1, | ||
DerivativeAccumulator * | da, | ||
const algebra::Vector3D & | delta, | ||
double | delta_length, | ||
double | x0, | ||
double | k | ||
) |
Evaluates a linear pair potential, such that score = k * (delta_length - x0) is returned. Also, updates the derivative vectors of the particles if da is not null.
[out] | d_xyzr0 | pointer to derivative of particle 0 |
[out] | d_xyzr1 | pointer to derivative of particle 1 |
[in,out] | da | accumulator for score derivatives to be updated |
[in] | delta | a vector that represents the displacement between the two particles |
delta_length | the cached length of delta, assumed correct, and required for faster calculation | |
x0 | resting distance (where score = 0) | |
k | score linear coefficient. Note that the score is attractive for a positive k, and repulsive for a negative k (assuming the lower the score the better). |
Definition at line 41 of file linear_distance_pair_scores.h.
double IMP::npctransport::do_evaluate_index_harmonic | ( | Model * | m, |
const ParticleIndexPair & | pp, | ||
DerivativeAccumulator * | da, | ||
const algebra::Vector3D & | delta, | ||
double | delta_length, | ||
double | x0, | ||
double | k | ||
) |
Evaluates a linear pair potential, such that force = k * (delta_length - x0) is returned. Also, updates the derivative vectors of the particles in the model m.
[out] | m | a model for the particles, in which particle derivatives are updated |
[in] | pp | a pair of particle indices for fast access through internal model methods |
[in,out] | da | accumulator for score derivatives to be updated |
[in] | delta | a vector from pp[1] to pp[0] |
delta_length | the cached length of delta, assumed correct, and required for faster calculation | |
x0 | resting distance (where score = 0) | |
k | score linear coefficient. Note that the score is attractive for a positive k, and repulsive for a negative k (assuming the lower the score the better). |
Definition at line 41 of file harmonic_distance_pair_scores.h.
void IMP::npctransport::do_main_loop | ( | SimulationData * | sd, |
const RestraintsTemp & | init_restraints | ||
) |
Run simulation using preconstructed SimulationData object sd.
sd | SimulationData object to optimize |
init_restraints | ad-hoc restraints during initialization only |
unsigned int IMP::npctransport::find_or_add_fg_bead_of_type | ( | ::npctransport_proto::Statistics * | s, |
IMP::core::ParticleType | pt | ||
) |
finds the index of s.fg_beads() whose type equals pt.get_string() if it does not exist, add it to s
s | the statistics message for searching the fg |
pt | the type of fg to look for |
unsigned int IMP::npctransport::find_or_add_fg_chain_of_type | ( | ::npctransport_proto::Statistics * | s, |
IMP::core::ParticleType | pt | ||
) |
finds the index of s.fgs() whose type equals pt.get_string() if it does not exist, add it to s
s | the statistics message for searching the fg |
pt | the type of fg to look for |
unsigned int IMP::npctransport::find_or_add_floater_of_type | ( | ::npctransport_proto::Statistics * | s, |
IMP::core::ParticleType | pt | ||
) |
finds the index of s->floaters() whose type equals pt.get_string() if it does not exist, add it to s
s | the statistics message for searching t |
pt | the type to look for |
unsigned int IMP::npctransport::find_or_add_interaction_of_type | ( | ::npctransport_proto::Statistics * | s, |
IMP::npctransport::InteractionType | it | ||
) |
finds the index of s.interactions() whose type0 and type1 particle type equivalents are equale to it. If it does not exist, add it to s
s | the statistics message for searching t |
it | the type to look for |
double IMP::npctransport::get_close_pairs_range | ( | const ::npctransport_proto::Assignment & | config | ) |
returns an upper bound on the contact range (maximal sphere-sphere distance that enables a contact) between any two particles in the system
std::string IMP::npctransport::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
This will ensure that the code works both when IMP is installed or if used via the setup_environment.sh
script.
int IMP::npctransport::get_dump_interval_in_frames | ( | const ::npctransport_proto::Assignment & | config, |
double | time_step | ||
) |
computes the number of frames for the specified dump_interval, which is normalized by the time step if the interval is specified in ns
config | the simulation parameters |
time_step | the time step in femtoseconds |
std::string IMP::npctransport::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
This will ensure that the code works both when IMP is installed or if used via the setup_environment.sh
script.
FGChain* IMP::npctransport::get_fg_chain | ( | atom::Hierarchy | root | ) |
gets a newly allocated chain structure from a root of an FG nup (by adding its ordered leaves)
FGChain* IMP::npctransport::get_fg_chain | ( | Particle * | p_root | ) |
gets a newly allocated chain structure from a root of an FG nup (by adding its ordered leaves)
p_root | a particle that is assumed to be Hierarchy decorated |
int IMP::npctransport::get_frames_from_ns | ( | double | ns, |
double | time_step | ||
) |
computes the number of frames that approximate a certain number of nanoseconds using the specified time_step per step in fs
ns | time in ns |
time_step | step in fs per frame |
algebra::Vector3D IMP::npctransport::get_global_from_local_v3 | ( | Particle * | p, |
const algebra::Vector3D & | local | ||
) |
unsigned int IMP::npctransport::get_maximal_number_of_unordered_pairs | ( | ParticlesTemp const & | ps0, |
ParticlesTemp const & | ps1 | ||
) |
gets the maximal theoretical number of unordered pairs of different particles between two sets of particles (note that the calculation is not entirely trivial since ps0 and ps1 may be overlapping sets while the pairs are unordered)
Definition at line 172 of file util.h.
boost::tuple<unsigned int, unsigned int> IMP::npctransport::get_n_lost_and_gained | ( | t_ordered_set | old, |
t_ordered_set | cur | ||
) |
returns |old/core| and |cur/old|, i.e., the number of items lost from old, and the ones gained in cur.
ParticlesTemp IMP::npctransport::get_non_optimizable_particles | ( | ParticlesTemp const & | particles | ) |
returns particles with non-optimizable coordinates from particles
int IMP::npctransport::get_number_of_frames | ( | const ::npctransport_proto::Assignment & | config, |
double | time_step | ||
) |
computes the number of frames needed to achieve simulation time required in the configuration time, with time step computed from config.
config | the simulation parameters |
time_step | the time step in femtoseconds |
ValueException | if maximum number of frames specified in config is exceeded |
ParticlesTemp IMP::npctransport::get_optimizable_particles | ( | ParticlesTemp const & | particles | ) |
returns particles with optimizable coordinates from particles
int IMP::npctransport::get_output_statistics_interval_in_frames | ( | const ::npctransport_proto::Assignment & | assign, |
double | time_step, | ||
double | default_value_ns = 1.0 |
||
) |
computes the maximal interval for outputting statistics to output file (and to dump order params) in frames, for the output_statistics_interval_ns param in assignment, based on specified time_step
assign | the simulation assignment parameters |
time_step | the time step in femtoseconds |
default_value_ns | interval in ns to use if none specified in assignment parameters |
ParticleIndexes IMP::npctransport::get_particle_indexes | ( | ParticlesTemp const & | particles | ) |
returns particle indexes from a list of particles
void IMP::npctransport::get_protobuf_configuration_from_text | ( | std::string | config_txt, |
std::string | config_pb | ||
) |
Converts protobuf configuration file config_txt (which is in pretty protobuf textual output format) to binary protobuf format (in file config_pb)
config_txt | the input textual protobuf config file |
config_pb | the output binary protobuf file |
algebra::Sphere3Ds IMP::npctransport::get_spheres_from_vectors | ( | V3iter | first, |
V3iter | last, | ||
double | radius | ||
) |
algebra::Sphere3Ds IMP::npctransport::get_spheres_from_vectors | ( | algebra::Vector3Ds const & | vs, |
double | radius | ||
) |
convert vectors to spheres of passed radius
int IMP::npctransport::get_statistics_interval_in_frames | ( | const ::npctransport_proto::Assignment & | assign, |
double | time_step, | ||
double | default_value_ns = 0.1 |
||
) |
computes the number of frames for the specified statistics_interval, normalized by the time step if the interval is specified in ns
assign | the simulation assignment parameters |
time_step | the time step in femtoseconds |
default_value_ns | interval in ns to use if none specified in assignment parameters |
double IMP::npctransport::get_time_step | ( | double | max_d_factor, |
double | max_k, | ||
double | min_radius, | ||
double | min_range, | ||
double | max_trans_relative_to_radius = 0.1 , |
||
double | time_step_factor = 1.0 |
||
) |
Computes the time step size that is required for a stable simulation. Formally, first computes the time step size that restricts the estimated translation size of any particle, at any given simulation step to [max_trans_relative_to_radius * min_radius], given the simulation parameters max_d_factor and max_k. This base time step is then multiplied by time_step_factor.
max_d_factor | the maximal diffusion factor of any particle in the system, which factors the Einstein diffusion coefficient) |
max_k | the maximal force applied on any particle in the simulation |
min_radius | the minimal radius of any interaction in the system |
min_range | the minimal range of any interaction |
max_trans_relative_to_radius | the maximal estimated translation allowed for any particle as fraction of min_radius (before factoring by time_step_factor) |
time_step_factor | multiply final time step in this factor |
double IMP::npctransport::get_time_step | ( | const ::npctransport_proto::Assignment & | config, |
double | max_trans_relative_to_radius = 0.1 |
||
) |
computes the time step size that is required for a stable simulation, where the translation of any particle at any simulation time step is restricted, based on the simulation parameters in config.
config | the simulation parameters used to compute the time step size |
max_trans_relative_to_radius | the maximal estimated translation allowed for any particle as fraction of its radius (before factoring by time step factor, specified in config) |
void IMP::npctransport::inflate_floater | ( | SimulationData * | sd, |
const std::string | floater_name, | ||
const float | new_radius | ||
) |
inflate floater of specified type to new_radius
void IMP::npctransport::initialize_positions | ( | SimulationData * | sd, |
const RestraintsTemp & | extra_restraints = RestraintsTemp() , |
||
bool | debug = false , |
||
double | short_init_factor = 1.0 , |
||
bool | is_disable_randomize = false , |
||
bool | are_fgs_pre_initialized = false |
||
) |
pre-optimize the positions of all diffusing particles in 'sd' whose coordinates are optimizable, using only chain restraints, excluded volumes and bounding box restraints (but not interaction restraints).
The initialization data is dumped to the RMF file sd->get_rmf_file_name()
using dump interval sd->get_rmf_dump_interval_frames() * 100
, or every frame in case that debug
is true.
sd | the SimulationData object containing diffusing particles |
extra_restraints | a list of additional ad-hoc restraints that will be used only throughout initialization |
debug | if true, the initialization will dump much more output (e.g. every frame to RMF file) |
short_init_factor | a factor between >0 and 1 for decreasing the number of optimization cycles at each round |
is_disable_randomize | if true, do not initially randomize particle positions, essentially performing an extended relaxation from the starting coordinates |
are_fgs_pre_initialized | if true, do not try to pre-optimize FGs before adding diffusers |
void IMP::npctransport::link_hierarchies_with_sites | ( | RMF::FileConstHandle | fh, |
const atom::Hierarchies & | hs | ||
) |
Link HierarchyWithSites objects with the RMF file, possibly overwriting an existing link for loading from the file. This does not alter the object, but will affect the behavior of functions like load_frame() and save_frame(). See IMP::rmf::link_hierarchies() for more details. The only difference is the addition of particle sites support.
bool IMP::npctransport::load_output_protobuf | ( | std::string | output_fname, |
::npctransport_proto::Output & | output | ||
) |
load file output_fname into protobuf output object output return true if succesful
void IMP::npctransport::load_pb_conformation | ( | const ::npctransport_proto::Conformation & | conformation, |
IMP::SingletonContainerAdaptor | beads, | ||
boost::unordered_map< core::ParticleType, algebra::Sphere3Ds > & | sites | ||
) |
Loads a protobuf conformation into the diffusers and sites
conformation | the saved conformation protobuf message |
beads | corresponding diffusers to be updated |
sites | a map of sites for each diffuser particle type to be updated |
InteractionType IMP::npctransport::make_ordered_interaction_type | ( | IMP::core::ParticleType | t0, |
IMP::core::ParticleType | t1 | ||
) |
convenience method for creating and interaction type in swig
Definition at line 26 of file typedefs.h.
InteractionType IMP::npctransport::make_unordered_interaction_type | ( | IMP::core::ParticleType | t0, |
IMP::core::ParticleType | t1 | ||
) |
an interaction type canonized to (t0,t1) s.t. t0 >= t1, so order doesn't matter
Definition at line 35 of file typedefs.h.
std::pair<t_value, t_value> IMP::npctransport::make_unordered_pair | ( | t_value | v0, |
t_value | v1 | ||
) |
IMP::ParticleIndexPair IMP::npctransport::make_unordered_particle_index_pair | ( | IMP::ParticleIndex | pi0, |
IMP::ParticleIndex | pi1 | ||
) |
ParticleIndexPair canonized to (pi0',pi1') s.t. pi0' >= pi1' so order doesn't matter
Definition at line 45 of file typedefs.h.
ParticleIndexPair IMP::npctransport::make_unordered_particle_index_pair | ( | ParticleIndexPair | pip | ) |
ParticleIndexPair canonized to (pi0',pi1') s.t. pi0' >= pi1' so order doesn't matter
Definition at line 55 of file typedefs.h.
void IMP::npctransport::randomize_particles | ( | const ParticlesList & | ps, |
const BoundingVolume & | bv | ||
) |
Randomize the positions of a set of particles within a bounding volume. Rigid bodies have their orientation randomized too.
Definition at line 35 of file randomize_particles.h.
void IMP::npctransport::remove_Nup42 | ( | SimulationData * | sd | ) |
remove nup42 and its anchors (also from obstacles)
void IMP::npctransport::reset_box_size | ( | SimulationData * | sd, |
double | box_size | ||
) |
change box size sd to specified box size and update output file
void IMP::npctransport::save_pb_conformation | ( | IMP::SingletonContainerAdaptor | beads, |
const boost::unordered_map< core::ParticleType, algebra::Sphere3Ds > & | sites, | ||
::npctransport_proto::Conformation * | conformation | ||
) |
Saves a protobuf conformation from the diffusers and sites
beads | beads to save |
sites | a map of sites for each diffuser particle type to be saved |
conformation | the conformation protobuf message to be save |
IMP::npctransport::SimulationData* IMP::npctransport::startup | ( | int | argc, |
char * | argv[] | ||
) |
initialize and return a simulation data object based on program command line parameters
IMP::IOException | if there was any IO problem |