Rosetta
Public Types | Public Member Functions | Private Member Functions | Private Attributes | List of all members
core::kinematics::AtomTree Class Reference

The AtomTree class is a wrapper for a tree of kinematic Atoms. More...

#include <AtomTree.hh>

Inheritance diagram for core::kinematics::AtomTree:
Inheritance graph
[legend]

Public Types

typedef id::AtomID AtomID
 
typedef id::AtomID_Mask AtomID_Mask
 
typedef id::DOF_Type DOF_Type
 
typedef id::DOF_ID DOF_ID
 
typedef id::DOF_ID_Mask DOF_ID_Mask
 
typedef id::StubID StubID
 
typedef std::vector< AtomIDAtomIDs
 
typedef std::map< AtomID, VectorFragXYZ
 
typedef std::map< StubID, RTFragRT
 
typedef tree::Atom Atom
 
typedef tree::AtomOP AtomOP
 
typedef tree::AtomCOP AtomCOP
 

Public Member Functions

 AtomTree (AtomPointer2D const &new_atom_pointer, bool const from_xyz=true)
 construction: from a tree of atoms; More...
 
 AtomTree ()
 default constructor More...
 
 ~AtomTree () override
 Destructor. More...
 
 AtomTree (AtomTree const &src)
 copy constructor More...
 
AtomTreeoperator= (AtomTree const &src)
 Copy assignment. If these two AtomTrees do not know anything about each other, then this function makes this AtomTree a complete copy of the src AtomTree, and then keeps track of the fact that the topology of the trees are identical by keeping a (raw) pointer to the src AtomTree and letting the src AtomTree keep a (raw) pointer to it. The next time this AtomTree's topology changes, or the next time the src AtomTree's topology changes, the "observer" connection between the two of them is severed – each AtomTree getting modified. This is not a thread-safe operation! If the AtomTree's have matching topologies, then this function simply iterates across all the atoms in both trees and copies the DOFs and coordinates from src into this. More...
 
void detached_copy (AtomTree const &src)
 Copy assignment that first invokes operator =, and then immediately severs the connection between the src and this AtomTrees. This function modifies the src AtomTree briefly while it executes, but then makes sure that there is no further (behind the sceens) communication between the two trees (e.g. that might happen when the topology of one of the trees changes, so that the "topological match to" status has to be severed). More...
 
void set_weak_pointer_to_self (AtomTreeCAP self_pointer)
 Weak-pointer setter. The object that instantiates an owning pointer to an AtomTree object must hand that AtomTree a weak pointer to itself so that the AtomTree may share that weak pointer with other AtomTrees. Such sharing allows for crucial speedups when copying between AtomTrees. If the object that instantiates this AtomTree does not provide it with a pointer-to-self, the AtomTree will still function, but it will not share its pointers properly. More...
 
Size size () const
 number of residues More...
 
bool has (AtomID const &id) const
 true only if AtomID is in the range of the map and its Atom pointer has value More...
 
void replace_residue_subtree (id::BondID const &incoming, utility::vector1< id::BondID > const &outgoing, AtomPointer1D const &new_atoms)
 
void promote_sameresidue_nonjump_child (AtomID const &parent_atom_id)
 Useful for guaranteeing that a stub remains within a single residue. More...
 
void delete_seqpos (Size const seqpos)
 Deletes atoms for seqpos. Does not renumber other atoms – need to call update_sequence_numbering for that designed for the simple case of 1 incoming connecting and 0 or 1 outgoing connection, where the desired behavior is to rewire the outgoing connection in place of seqpos' tree. More...
 
void update_sequence_numbering (Size const new_size, utility::vector1< int > const &old2new)
 updates the Atom's AtomID's and the atom_pointer array More...
 
void replace_tree (AtomPointer2D const &new_atom_pointer, bool const from_xyz=true)
 replaces the entire tree More...
 
void copy_coords (AtomTree const &src)
 copy the internal and xyz coords from src tree, which should have the same topology as me More...
 
void set_dof (DOF_ID const &id, Real const setting)
 set a specific DOF in the tree More...
 
void set_xyz (AtomID const &id, PointPosition const &xyz)
 set a specific atom xyz position More...
 
void batch_set_xyz (utility::vector1< AtomID > const &id, utility::vector1< PointPosition > const &xyz)
 simultaniously set several atom xyz positions More...
 
void set_jump (AtomID const &id, Jump const &jump)
 set a specific jump transform More...
 
void set_jump_now (AtomID const &id, Jump const &jump)
 set a specific jump transform and immediately refold downstream atoms More...
 
AtomID get_jump_atom_id (StubID const &stub_id1, StubID const &stub_id2, int &direction) const
 find the atom giving rise to the jump connecting two stubs. More...
 
AtomID set_stub_transform (StubID const &stub_id1, StubID const &stub_id2, RT const &target_rt)
 Set the transform between two stubs, returns the atomid of the jump atom which moved (for book-keeping) More...
 
RT get_stub_transform (StubID const &stub_id1, StubID const &stub_id2) const
 get the transform between two stubs More...
 
DOF_ID set_torsion_angle (AtomID const &atom1, AtomID const &atom2, AtomID const &atom3, AtomID const &atom4, Real const setting, bool const quiet=false)
 set a torsion angle "setting" to a specifc DOF derived by the four atoms More...
 
DOF_ID set_bond_angle (AtomID const &atom1, AtomID const &atom2, AtomID const &atom3, Real const setting)
 
DOF_ID set_bond_length (AtomID const &atom1, AtomID const &atom2, Real const setting)
 
void update_domain_map (DomainMap &domain_map, AtomID_Mask const &dof_moved, AtomID_Mask const &xyz_moved) const
 generates a "domain_map" defining the rigid body regions whose internal coords have not changed, according to the informaiton in the two bool Mask's update domain map from dof_moved and xyz_moved More...
 
void clear ()
 clear the content of an AtomTree object, delete root atom and all pointers More...
 
void insert_fragment (StubID const &instub_id, FragRT const &outstub_transforms, FragXYZ const &frag_xyz, utility::vector1< AtomID > &moving_atoms)
 
ResidueListIterator residue_xyz_change_list_begin () const
 The AtomTree provides to the Conformation object a list of residues whose xyz coordinates have changed. When the Conformation has finished reading off residues that have changed from the AtomTree, and has copied the coordinates of those residues into its conformation::Residue objects, it informs the AtomTree to reset this list by a call to mark_changed_residues_registered. More...
 
ResidueListIterator residue_xyz_change_list_end () const
 
void note_coordinate_change_registered () const
 The AtomTree provides a list of residues who's xyz coordinates have changed to the Conformation object. When the Conformation has finished reading off residues that have changed from the AtomTree, and has copied the coordinates of those residues into its conformation::Residue objects, it informs the AtomTree to reset this list by a call to mark_changed_residues_registered. More...
 
bool empty () const
 is there any atom in the tree yet? More...
 
AtomCOP root () const
 const-access to the root of the tree More...
 
Real dof (DOF_ID const &id) const
 get value of DOF( PHI, THETA, D, RB1, ....) given its DOF_ID More...
 
PointPosition const & xyz (AtomID const &id) const
 get xyz position of an atom given its AtomID More...
 
Atom const & atom (AtomID const &id) const
 retreive a kinematic Atom in the tree by its AtomID More...
 
Atom const & atom_dont_do_update (AtomID const &id) const
 retreive a kinematic Atom in the tree by its AtomID – no update! More...
 
Jump const & jump (AtomID const &id) const
 retrieve a Jump in the tree by that JumpAtom's AtomID. More...
 
DOF_ID torsion_angle_dof_id (AtomID const &atom1, AtomID const &atom2, AtomID const &atom3, AtomID const &atom4, bool const quiet=false) const
 a wrapper function to get the DOF_ID of a torsion angle given those four atoms which define this torsion More...
 
DOF_ID torsion_angle_dof_id (AtomID const &atom1_in_id, AtomID const &atom2_in_id, AtomID const &atom3_in_id, AtomID const &atom4_in_id, Real &offset, bool const quiet=false) const
 get the DOF_ID of a torsion angle given those four atoms which define this torsion More...
 
DOF_ID bond_angle_dof_id (AtomID const &atom1_in_id, AtomID const &atom2_in_id, AtomID const &atom3_in_id, Real &offset) const
 get the DOF_ID of a bond angle given those 3 atoms which define this torsion More...
 
DOF_ID bond_length_dof_id (AtomID const &atom1_in_id, AtomID const &atom2_in_id) const
 get the DOF_ID of a bond length given those 2 atoms which define this torsion More...
 
Real torsion_angle (AtomID const &atom1, AtomID const &atom2, AtomID const &atom3, AtomID const &atom4) const
 calculate torsion angle defined by four atoms in the atom tree More...
 
Real bond_angle (AtomID const &atom1, AtomID const &atom2, AtomID const &atom3) const
 
Real bond_length (AtomID const &atom1, AtomID const &atom2) const
 
void set_jump_atom_stub_id (StubID const &id)
 
Stub stub_from_id (StubID const &id) const
 
AtomTreeCAP topological_match_to () const
 Functions only necessary for unit tests. More...
 
utility::vector1< AtomTreeCAP > const & topological_observers () const
 For testing purposes only: report the list of observer AtomTrees that are topological copies of this tree. The fact that AtomTrees keep track of other atom trees is "private" in the sense that no other class needs to worry about it. However, to test that the topological match algorithm is working properly, this private data needs to be readable. Do not use this function outside of the unit tests. More...
 

Private Member Functions

void attach_topological_observer (AtomTreeCAP observer) const
 When an atom tree copies the topology of another atom tree, it must register itself as a topological observer of that other tree. When the other tree changes its topology, then the tree being observed must notify its observers that they are no longer a topological copy of this tree. An atom tree may only be the topological copy of a single other atom tree, though several atom trees may be copies of a single atom tree. More...
 
void notify_topological_change (AtomTreeCAP observee) const
 When an atom tree changes its topology, it must inform all of its observers that they are no longer the same topology as this tree. More...
 
void detach_topological_observer (AtomTreeCAP observer) const
 When an atom tree observing this tree decides it wants to become an observer of another tree, it must notify the tree that it formerly observed of this change. More...
 
void find_root_from_atom_pointer ()
 Deduce root_ from atom_pointer_ – look for atom with atom->parent() == 0. More...
 
void get_frag_atoms (StubID const &id, FragXYZ const &frag_xyz, AtomCOP &frag_atom, AtomCOP &nonfrag_atom) const
 
StubID get_frag_pseudo_stub_id (AtomID const &id, FragXYZ const &frag_xyz, bool &fail) const
 
Stub get_frag_local_stub (StubID const &stubid, FragXYZ const &frag_xyz, bool &fail) const
 
Vector get_frag_local_xyz (AtomID const &id, FragXYZ const &frag_xyz, bool &fail) const
 
Vector get_frag_descendant_local_xyz (AtomCOP atom, FragXYZ const &frag_xyz, bool &fail) const
 
Vector get_frag_parent_local_xyz (AtomCOP child, FragXYZ const &frag_xyz, bool &fail) const
 
void insert_single_fragment (StubID const &instub_id, FragRT const &outstub_transforms, FragXYZ const &frag_xyz, utility::vector1< AtomID > &moving_atoms)
 
void update_atom_ids_from_atom_pointer ()
 bookkeeping – set the Atoms' atomIDs from the atom_pointer_ array More...
 
AtomOP atom_pointer (AtomID const &id)
 Convenience when we want an Atom*. More...
 
AtomCOP atom_pointer (AtomID const &id) const
 Convenience when we want an Atom*. More...
 
void update_internal_coords () const
 Update the internal coordinates using the xyz (cartesian) coords. More...
 
void update_xyz_coords () const
 Update the xyz coordinates using the internal coords. More...
 
void set_new_topology ()
 Notify self of new tree topology Useful if we move to caching some things that depend on the tree. More...
 

Private Attributes

AtomTreeCAP this_weak_ptr_
 A weak pointer to self (this). More...
 
AtomOP root_
 Root Atom. More...
 
AtomPointer2D atom_pointer_
 Atom pointers map (map[AtomID] = AtomPointer) More...
 
bool internal_coords_need_updating_
 List of the jump atom ID's, excluding the root. Order matters (for movemap indexing) More...
 
bool xyz_coords_need_updating_
 XYZ coords out of date? More...
 
AtomTreeCAP topological_match_to_
 pointer to the atom tree this tree has an exact topological match to since that tree was the last tree copied from without subsequence topological modifications – or at most one modification when that tree copied this tree's topology More...
 
utility::vector1< AtomTreeCAPtopological_observers_
 pointers to all atom trees that are observing this tree's topology. On topological changes (including the destruction of this tree!), each of these trees have their topological_match_to_ pointers set to null and this list is cleared. More...
 
AtomDOFChangeSet dof_changeset_
 A list of the atoms that have had changed DOFs since the last refold. More...
 
ResidueCoordinateChangeListOP external_coordinate_residues_changed_
 A list of residues that have had xyz coordinate changes since the last time the owning Conformation object has asked for an update. More...
 

Detailed Description

The AtomTree class is a wrapper for a tree of kinematic Atoms.

The responsibilities of the class are:

See AtomTree overview and concepts for details.

Member Typedef Documentation

◆ Atom

◆ AtomCOP

◆ AtomID

◆ AtomID_Mask

◆ AtomIDs

◆ AtomOP

◆ DOF_ID

◆ DOF_ID_Mask

◆ DOF_Type

◆ FragRT

◆ FragXYZ

◆ StubID

Constructor & Destructor Documentation

◆ AtomTree() [1/3]

core::kinematics::AtomTree::AtomTree ( AtomPointer2D const &  new_atom_pointer,
bool const  from_xyz = true 
)

construction: from a tree of atoms;

this will claim the tree as our own. new_root has information about its children, and they have information about their children. From those atom positions, internal coordinates can be updated and atom pointers will be added into the AtomTree map.

Note
that we steal the atoms, ie it's incorporated into the AtomTree (by recording their pointers),not cloned

References external_coordinate_residues_changed_, replace_tree(), and core::id::AtomID_Map< T >::size().

◆ AtomTree() [2/3]

core::kinematics::AtomTree::AtomTree ( )

default constructor

◆ ~AtomTree()

core::kinematics::AtomTree::~AtomTree ( )
override

Destructor.

References clear().

◆ AtomTree() [3/3]

core::kinematics::AtomTree::AtomTree ( AtomTree const &  src)

copy constructor

copy ctor, uses operator= without this initialization, the destruction of this

Member Function Documentation

◆ atom()

tree::Atom const & core::kinematics::AtomTree::atom ( AtomID const &  id) const

retreive a kinematic Atom in the tree by its AtomID

References atom_pointer_, update_internal_coords(), and update_xyz_coords().

Referenced by core::scoring::MultipoleElecPotential::align_residue_multipole_axes(), protocols::simple_filters::JumpEvaluator::apply(), core::optimization::atom_tree_dfunc(), core::optimization::symmetry::atom_tree_dfunc(), protocols::backrub::BackrubSegment::bond_angle_atoms(), core::import_pose::RNA_JumpMover::check_forward_backward(), protocols::simple_moves::ConstrainToIdealMover::check_if_really_connected(), protocols::environment::collect_dofs(), protocols::simple_filters::SimpleHbondsToAtomFilter::compute(), protocols::loops::loop_closure::jacobi::JacobiLoopClosureMover::default_target_icoors(), delete_seqpos(), protocols::ligand_docking::ga_ligand_dock::GriddedAtomTreeMultifunc::dfunc(), core::energy_methods::LinearChainbreakEnergy::do_score_ovp(), protocols::backrub::BackrubMover::dof_id_ranges(), core::import_pose::atom_tree_diffs::dump_atom_tree_diff(), protocols::backrub::BackrubSegment::end_atoms1(), core::conformation::Conformation::fill_missing_atoms(), protocols::abinitio::abscript::fix_internal_coords_of_siblings(), protocols::topology_broker::fix_internal_coords_of_siblings(), get_frag_descendant_local_xyz(), get_frag_local_xyz(), get_frag_pseudo_stub_id(), core::pose::rna::initialize_atoms_for_which_we_need_new_dofs(), protocols::simple_moves::DOFHistogramRecorder::insert_dofs_by_residue(), insert_fragment(), insert_single_fragment(), core::conformation::is_ideal_position(), protocols::simple_filters::JumpEvaluator::JumpEvaluator(), protocols::branch_angle::BranchAngleOptimizer::optimize_angles(), protocols::branch_angle::BranchAngleOptimizer::overall_params(), protocols::rna::denovo::print_internal_coords(), protocols::backrub::BackrubMover::random_angle(), core::import_pose::atom_tree_diffs::rms_error_with_noise(), protocols::rna::movers::RNA_LoopCloser::rna_ccd_close(), protocols::backrub::BackrubMover::rotate_segment(), protocols::backrub::BackrubSegment::start_atoms1(), protocols::backrub::BackrubSegment::start_atoms2(), protocols::loops::loop_closure::jacobi::JacobiLoopClosureMover::store_target(), core::pack::scmin::ResidueAtomTreeCollection::update_atom_tree(), protocols::loops::loop_closure::jacobi::JacobiLoopClosureMover::update_current(), core::kinematics::jacobian::ModuleType1::update_screw_vectors(), protocols::rna::movers::ErraserMinimizerMover::vary_bond_geometry(), protocols::simple_moves::ConstrainToIdealMover::vary_bond_geometry(), and protocols::environment::claims::XYZClaim::yield_elements().

◆ atom_dont_do_update()

tree::Atom const & core::kinematics::AtomTree::atom_dont_do_update ( AtomID const &  id) const

◆ atom_pointer() [1/2]

AtomOP core::kinematics::AtomTree::atom_pointer ( AtomID const &  id)
inlineprivate

◆ atom_pointer() [2/2]

AtomCOP core::kinematics::AtomTree::atom_pointer ( AtomID const &  id) const
inlineprivate

Convenience when we want an Atom*.

References atom_pointer_.

◆ attach_topological_observer()

void core::kinematics::AtomTree::attach_topological_observer ( AtomTreeCAP  observer) const
private

When an atom tree copies the topology of another atom tree, it must register itself as a topological observer of that other tree. When the other tree changes its topology, then the tree being observed must notify its observers that they are no longer a topological copy of this tree. An atom tree may only be the topological copy of a single other atom tree, though several atom trees may be copies of a single atom tree.

References this_weak_ptr_, and topological_observers_.

Referenced by operator=().

◆ batch_set_xyz()

void core::kinematics::AtomTree::batch_set_xyz ( utility::vector1< AtomID > const &  id,
utility::vector1< PointPosition > const &  xyz 
)

simultaniously set several atom xyz positions

References atom_pointer_, internal_coords_need_updating_, and update_xyz_coords().

◆ bond_angle()

Real core::kinematics::AtomTree::bond_angle ( AtomID const &  atom1,
AtomID const &  atom2,
AtomID const &  atom3 
) const

◆ bond_angle_dof_id()

id::DOF_ID core::kinematics::AtomTree::bond_angle_dof_id ( AtomID const &  atom1_in_id,
AtomID const &  atom2_in_id,
AtomID const &  atom3_in_id,
Real offset 
) const

◆ bond_length()

Real core::kinematics::AtomTree::bond_length ( AtomID const &  atom1,
AtomID const &  atom2 
) const

◆ bond_length_dof_id()

id::DOF_ID core::kinematics::AtomTree::bond_length_dof_id ( AtomID const &  atom1_in_id,
AtomID const &  atom2_in_id 
) const

◆ clear()

void core::kinematics::AtomTree::clear ( )

clear the content of an AtomTree object, delete root atom and all pointers

This is done by releasing memories of all atoms including the root atom , and clearing atom_pointer map

References atom_pointer_, core::id::AtomID_Map< T >::clear(), external_coordinate_residues_changed_, root_, and set_new_topology().

Referenced by operator=(), replace_tree(), and ~AtomTree().

◆ copy_coords()

void core::kinematics::AtomTree::copy_coords ( AtomTree const &  src)

copy the internal and xyz coords from src tree, which should have the same topology as me

Note
this may trigger a coordinate update of src (access to root atom)

References dof_changeset_, external_coordinate_residues_changed_, internal_coords_need_updating_, root(), root_, and xyz_coords_need_updating_.

Referenced by operator=().

◆ delete_seqpos()

void core::kinematics::AtomTree::delete_seqpos ( Size const  seqpos)

Deletes atoms for seqpos. Does not renumber other atoms – need to call update_sequence_numbering for that designed for the simple case of 1 incoming connecting and 0 or 1 outgoing connection, where the desired behavior is to rewire the outgoing connection in place of seqpos' tree.

assumes one incoming and at most one outgoing Need fancier fxn for general case

References atom(), atom_pointer_, core::kinematics::tree::Atom::atoms_begin(), core::kinematics::tree::Atom::atoms_end(), core::kinematics::tree::Atom::parent(), core::id::AtomID_Map< T >::resize(), root(), set_new_topology(), size(), and update_sequence_numbering().

◆ detach_topological_observer()

void core::kinematics::AtomTree::detach_topological_observer ( AtomTreeCAP  observer_ap) const
private

When an atom tree observing this tree decides it wants to become an observer of another tree, it must notify the tree that it formerly observed of this change.

When an AtomTree which was observing this AtomTree changes its topology or is deleted, it must invoke this method on this AtomTree. When this happens, this AtomTree marks the observer's position in its list of observers as null – and it also resets the "topological_match_to_" pointer in the oberserver_ap AtomTree.

References topological_observers_.

◆ detached_copy()

void core::kinematics::AtomTree::detached_copy ( AtomTree const &  src)

Copy assignment that first invokes operator =, and then immediately severs the connection between the src and this AtomTrees. This function modifies the src AtomTree briefly while it executes, but then makes sure that there is no further (behind the sceens) communication between the two trees (e.g. that might happen when the topology of one of the trees changes, so that the "topological match to" status has to be severed).

This works by invoking this AtomTree's opeator = and then resetting the topological_match_to_ pointer by calling set_new_topology();

References set_new_topology().

◆ dof()

Real core::kinematics::AtomTree::dof ( DOF_ID const &  id) const

◆ empty()

bool core::kinematics::AtomTree::empty ( ) const
inline

is there any atom in the tree yet?

References root_.

◆ find_root_from_atom_pointer()

void core::kinematics::AtomTree::find_root_from_atom_pointer ( )
private

Deduce root_ from atom_pointer_ – look for atom with atom->parent() == 0.

References atom_pointer_, root_, and core::id::AtomID_Map< T >::size().

Referenced by replace_tree().

◆ get_frag_atoms()

void core::kinematics::AtomTree::get_frag_atoms ( StubID const &  id,
FragXYZ const &  frag_xyz,
AtomCOP frag_atom,
AtomCOP nonfrag_atom 
) const
private

◆ get_frag_descendant_local_xyz()

Vector core::kinematics::AtomTree::get_frag_descendant_local_xyz ( AtomCOP  atom,
FragXYZ const &  frag_xyz,
bool &  fail 
) const
private

◆ get_frag_local_stub()

Stub core::kinematics::AtomTree::get_frag_local_stub ( StubID const &  stubid,
FragXYZ const &  frag_xyz,
bool &  fail 
) const
private

◆ get_frag_local_xyz()

Vector core::kinematics::AtomTree::get_frag_local_xyz ( AtomID const &  id,
FragXYZ const &  frag_xyz,
bool &  fail 
) const
private

private helper for fragment insertion routines

id is either in frag or a child or a gchild or a parent

References atom(), atom_pointer(), core::kinematics::tree::Atom::child(), get_frag_descendant_local_xyz(), get_frag_parent_local_xyz(), core::kinematics::tree::Atom::n_children(), and core::kinematics::tree::Atom::parent().

Referenced by get_frag_local_stub().

◆ get_frag_parent_local_xyz()

Vector core::kinematics::AtomTree::get_frag_parent_local_xyz ( AtomCOP  child,
FragXYZ const &  frag_xyz,
bool &  fail 
) const
private

◆ get_frag_pseudo_stub_id()

id::StubID core::kinematics::AtomTree::get_frag_pseudo_stub_id ( AtomID const &  id,
FragXYZ const &  frag_xyz,
bool &  fail 
) const
private

id is a frag atom

look for two more nearby frag atoms to define a pseudo-stub for getting coords for parents or children of this atom

References atom(), atom_pointer(), core::kinematics::tree::Atom::get_nonjump_atom(), core::kinematics::tree::Atom::is_jump(), and core::kinematics::tree::Atom::parent().

Referenced by get_frag_descendant_local_xyz(), get_frag_parent_local_xyz(), and insert_fragment().

◆ get_jump_atom_id()

id::AtomID core::kinematics::AtomTree::get_jump_atom_id ( StubID const &  stub_id1,
StubID const &  stub_id2,
int &  direction 
) const

find the atom giving rise to the jump connecting two stubs.

Parameters
directionis set to 1 or -1 depending on the direction of the jump

References core::id::StubID::atom(), and atom_pointer().

Referenced by set_stub_transform().

◆ get_stub_transform()

RT core::kinematics::AtomTree::get_stub_transform ( StubID const &  stub_id1,
StubID const &  stub_id2 
) const

get the transform between two stubs

References stub_from_id().

◆ has()

bool core::kinematics::AtomTree::has ( AtomID const &  id) const
inline

◆ insert_fragment()

void core::kinematics::AtomTree::insert_fragment ( StubID const &  instub_id,
FragRT const &  outstub_transforms,
FragXYZ const &  frag_xyz,
utility::vector1< AtomID > &  moving_atoms 
)

◆ insert_single_fragment()

void core::kinematics::AtomTree::insert_single_fragment ( StubID const &  instub_id,
FragRT const &  outstub_transforms,
FragXYZ const &  frag_xyz,
utility::vector1< AtomID > &  moving_atoms 
)
private

◆ jump()

Jump const & core::kinematics::AtomTree::jump ( AtomID const &  id) const

retrieve a Jump in the tree by that JumpAtom's AtomID.

Note
will abort if a BondedAtom's AtomID is passed in.

References atom_pointer_, and update_internal_coords().

Referenced by set_jump(), and set_jump_now().

◆ note_coordinate_change_registered()

void core::kinematics::AtomTree::note_coordinate_change_registered ( ) const

The AtomTree provides a list of residues who's xyz coordinates have changed to the Conformation object. When the Conformation has finished reading off residues that have changed from the AtomTree, and has copied the coordinates of those residues into its conformation::Residue objects, it informs the AtomTree to reset this list by a call to mark_changed_residues_registered.

References external_coordinate_residues_changed_.

◆ notify_topological_change()

void core::kinematics::AtomTree::notify_topological_change ( AtomTreeCAP  observee) const
private

When an atom tree changes its topology, it must inform all of its observers that they are no longer the same topology as this tree.

The AtomTree being observed calls this notify_topological_change on all of the AtomeTrees that are observing it. The this object "detaches" itself from the observee in this function. In debug mode, this function ensures that this object was actually observing the observee – if it wasn't, then an internal error has occurred. The observee and the observer have gotten out of sync.

References topological_match_to_.

◆ operator=()

AtomTree & core::kinematics::AtomTree::operator= ( AtomTree const &  src)

Copy assignment. If these two AtomTrees do not know anything about each other, then this function makes this AtomTree a complete copy of the src AtomTree, and then keeps track of the fact that the topology of the trees are identical by keeping a (raw) pointer to the src AtomTree and letting the src AtomTree keep a (raw) pointer to it. The next time this AtomTree's topology changes, or the next time the src AtomTree's topology changes, the "observer" connection between the two of them is severed – each AtomTree getting modified. This is not a thread-safe operation! If the AtomTree's have matching topologies, then this function simply iterates across all the atoms in both trees and copies the DOFs and coordinates from src into this.

makes a complete copy of another AtomTree src by cloning the root_atom and all its offspring atoms. Also update atom_pointer map.

Note
clear the content of tree first to release memory before doing the assignment.

References atom_pointer_, attach_topological_observer(), clear(), copy_coords(), dof_changeset_, external_coordinate_residues_changed_, internal_coords_need_updating_, root_, this_weak_ptr_, topological_match_to(), topological_match_to_, and xyz_coords_need_updating_.

◆ promote_sameresidue_nonjump_child()

void core::kinematics::AtomTree::promote_sameresidue_nonjump_child ( AtomID const &  parent_atom_id)

Useful for guaranteeing that a stub remains within a single residue.

Useful for guaranteeing that a stub remains within a single residue

References atom_pointer(), internal_coords_need_updating_, set_new_topology(), core::kinematics::TR(), and update_xyz_coords().

Referenced by core::conformation::promote_sameresidue_child_of_jump_atom().

◆ replace_residue_subtree()

void core::kinematics::AtomTree::replace_residue_subtree ( id::BondID const &  incoming,
utility::vector1< id::BondID > const &  outgoing,
AtomPointer1D const &  new_atoms 
)
Note
This can also handle the case of inserting a residue, proved that we've already renumbered atomtree leaving an empty slot at seqpos.
If the new residue contains the root of the tree, incoming.atom1 should be AtomID::BOGUS_ATOM_ID()
Note that for all BondID's atom1 should be the parent and atom2 should be the child Couple of special cases: – appending a residue – inserting a residue – replacing the root residue – inserting a new root residue

References core::id::BondID::atom1, core::id::BondID::atom2, atom_pointer(), atom_pointer_, core::id::AtomID::atomno(), core::id::AtomID_Map< T >::clear(), internal_coords_need_updating_, root_, core::id::AtomID::rsd(), set_new_topology(), update_xyz_coords(), and core::id::AtomID::valid().

Referenced by core::conformation::replace_residue_in_atom_tree().

◆ replace_tree()

void core::kinematics::AtomTree::replace_tree ( AtomPointer2D const &  new_atom_pointer,
bool const  from_xyz = true 
)

replaces the entire tree

fill the AtomTree with a new tree of atoms by recording their pointers in the map. Sync internal and xyz coords.

References atom_pointer_, clear(), external_coordinate_residues_changed_, find_root_from_atom_pointer(), internal_coords_need_updating_, set_new_topology(), core::id::AtomID_Map< T >::size(), update_internal_coords(), update_xyz_coords(), and xyz_coords_need_updating_.

Referenced by AtomTree().

◆ residue_xyz_change_list_begin()

ResidueListIterator core::kinematics::AtomTree::residue_xyz_change_list_begin ( ) const

The AtomTree provides to the Conformation object a list of residues whose xyz coordinates have changed. When the Conformation has finished reading off residues that have changed from the AtomTree, and has copied the coordinates of those residues into its conformation::Residue objects, it informs the AtomTree to reset this list by a call to mark_changed_residues_registered.

The list of which residues have had coordinate changes is unknown until the DFS has completed. The DFS must be triggered before iterators are given to the Conformation object.

References external_coordinate_residues_changed_, update_xyz_coords(), and xyz_coords_need_updating_.

◆ residue_xyz_change_list_end()

ResidueListIterator core::kinematics::AtomTree::residue_xyz_change_list_end ( ) const

The list of which residues have had coordinate changes is unknown until the DFS has completed. The DFS must be triggered before iterators are given to the Conformation object.

References external_coordinate_residues_changed_, update_xyz_coords(), and xyz_coords_need_updating_.

◆ root()

AtomCOP core::kinematics::AtomTree::root ( ) const
inline

◆ set_bond_angle()

id::DOF_ID core::kinematics::AtomTree::set_bond_angle ( AtomID const &  atom1,
AtomID const &  atom2,
AtomID const &  atom3,
Real const  setting 
)

◆ set_bond_length()

id::DOF_ID core::kinematics::AtomTree::set_bond_length ( AtomID const &  atom1,
AtomID const &  atom2,
Real const  setting 
)

◆ set_dof()

void core::kinematics::AtomTree::set_dof ( DOF_ID const &  id,
Real const  setting 
)

◆ set_jump()

void core::kinematics::AtomTree::set_jump ( AtomID const &  id,
Jump const &  jump 
)

set a specific jump transform

Note
AtomID must point to a JumpAtom, otherwise it will die.

References atom_pointer_, dof_changeset_, jump(), update_internal_coords(), and xyz_coords_need_updating_.

Referenced by set_jump_now().

◆ set_jump_atom_stub_id()

void core::kinematics::AtomTree::set_jump_atom_stub_id ( StubID const &  id)

◆ set_jump_now()

void core::kinematics::AtomTree::set_jump_now ( AtomID const &  id,
Jump const &  jump 
)

set a specific jump transform and immediately refold downstream atoms

Note
DEPRECATED. AtomID must point to a JumpAtom, otherwise it will die. Hopefully, the new output-sensitive refold will make this method unnecessary

References jump(), and set_jump().

◆ set_new_topology()

void core::kinematics::AtomTree::set_new_topology ( )
private

Notify self of new tree topology Useful if we move to caching some things that depend on the tree.

Note
Probably need to go through and put more calls of this guy

If this tree changes its topology, then its no longer a match in topology to the tree it was copied from, nor are the trees that are a copy of it. Detach this as an observer of the AtomTree it is observing (if any) and inform all the trees observing this tree that the topology has changed before clearing the topological_observers_ array.

References this_weak_ptr_, topological_match_to(), topological_match_to_, and topological_observers_.

Referenced by clear(), delete_seqpos(), detached_copy(), promote_sameresidue_nonjump_child(), replace_residue_subtree(), and replace_tree().

◆ set_stub_transform()

id::AtomID core::kinematics::AtomTree::set_stub_transform ( StubID const &  stub_id1,
StubID const &  stub_id2,
RT const &  target_rt 
)

Set the transform between two stubs, returns the atomid of the jump atom which moved (for book-keeping)

Set the transform between two stubs Returns the atomid of the jump atom which moved.

Note
Requires that there be a jump in the atomtree between a pair of atoms in the stubs

References protocols::comparative_modeling::features::A, atom_pointer(), protocols::match::upstream::b, core::kinematics::tree::distance_squared(), external_coordinate_residues_changed_, core::kinematics::find_stub_transform(), get_jump_atom_id(), internal_coords_need_updating_, core::kinematics::RT::reverse(), and stub_from_id().

◆ set_torsion_angle()

id::DOF_ID core::kinematics::AtomTree::set_torsion_angle ( AtomID const &  atom1,
AtomID const &  atom2,
AtomID const &  atom3,
AtomID const &  atom4,
Real const  setting,
bool const  quiet = false 
)

set a torsion angle "setting" to a specifc DOF derived by the four atoms

It is possible that no DOF can be obtained given these four atoms or the torsion does not match exactly the DOF(PHI) angle. In the former case, a DOF_ID::BOGUS_DOF_ID() is returned as indicator and in the latter case, an offset value is deducted from input setting to set the real DOF value properly.

References set_dof(), torsion_angle_dof_id(), and core::kinematics::TR().

◆ set_weak_pointer_to_self()

void core::kinematics::AtomTree::set_weak_pointer_to_self ( AtomTreeCAP  self_pointer)

Weak-pointer setter. The object that instantiates an owning pointer to an AtomTree object must hand that AtomTree a weak pointer to itself so that the AtomTree may share that weak pointer with other AtomTrees. Such sharing allows for crucial speedups when copying between AtomTrees. If the object that instantiates this AtomTree does not provide it with a pointer-to-self, the AtomTree will still function, but it will not share its pointers properly.

References this_weak_ptr_.

◆ set_xyz()

void core::kinematics::AtomTree::set_xyz ( AtomID const &  id,
PointPosition const &  xyz 
)

◆ size()

Size core::kinematics::AtomTree::size ( ) const
inline

◆ stub_from_id()

Stub core::kinematics::AtomTree::stub_from_id ( StubID const &  id) const
inline

◆ topological_match_to()

AtomTreeCAP core::kinematics::AtomTree::topological_match_to ( ) const
inline

Functions only necessary for unit tests.

For testing purposes only: report the address of the AtomTree this tree is a topological copy of. The fact that AtomTrees keep track of other atom trees is "private" in the sense that no other class needs to worry about it. However, to test that the topological match algorithm is working properly, this private data needs to be readable. Do not use this function outside of the unit tests.

References topological_match_to_.

Referenced by operator=(), and set_new_topology().

◆ topological_observers()

utility::vector1< AtomTreeCAP > const& core::kinematics::AtomTree::topological_observers ( ) const
inline

For testing purposes only: report the list of observer AtomTrees that are topological copies of this tree. The fact that AtomTrees keep track of other atom trees is "private" in the sense that no other class needs to worry about it. However, to test that the topological match algorithm is working properly, this private data needs to be readable. Do not use this function outside of the unit tests.

References topological_observers_.

◆ torsion_angle()

Real core::kinematics::AtomTree::torsion_angle ( AtomID const &  atom1,
AtomID const &  atom2,
AtomID const &  atom3,
AtomID const &  atom4 
) const

calculate torsion angle defined by four atoms in the atom tree

Note
this is NOT calculated straight from atom positions.Instead, it is calculated from internal DOFs. If no DOF can be found from these four atoms, 0.0 will be returned and a warning message is printed.

References atom_pointer_, torsion_angle_dof_id(), core::kinematics::TR(), and update_internal_coords().

Referenced by protocols::chemically_conjugated_docking::UBQ_GTPaseMover::analyze_and_filter(), and protocols::simple_moves::TorsionDOFMover::apply().

◆ torsion_angle_dof_id() [1/2]

DOF_ID core::kinematics::AtomTree::torsion_angle_dof_id ( AtomID const &  atom1,
AtomID const &  atom2,
AtomID const &  atom3,
AtomID const &  atom4,
bool const  quiet = false 
) const
inline

◆ torsion_angle_dof_id() [2/2]

id::DOF_ID core::kinematics::AtomTree::torsion_angle_dof_id ( AtomID const &  atom1_in_id,
AtomID const &  atom2_in_id,
AtomID const &  atom3_in_id,
AtomID const &  atom4_in_id,
Real offset,
bool const  quiet = false 
) const

◆ update_atom_ids_from_atom_pointer()

void core::kinematics::AtomTree::update_atom_ids_from_atom_pointer ( )
private

bookkeeping – set the Atoms' atomIDs from the atom_pointer_ array

References atom_pointer_, core::id::AtomID_Map< T >::size(), and size().

Referenced by update_sequence_numbering().

◆ update_domain_map()

void core::kinematics::AtomTree::update_domain_map ( DomainMap domain_map,
AtomID_Mask const &  dof_moved,
AtomID_Mask const &  xyz_moved 
) const

generates a "domain_map" defining the rigid body regions whose internal coords have not changed, according to the informaiton in the two bool Mask's update domain map from dof_moved and xyz_moved

domain map is residue-based and indicate which residues stay relatively rigid to each other in one domain since last move and therefore their interaction energy does not have to be recalculated. This function calls atom->update_domain_map from the root atom of the tree.

References root_.

Referenced by core::optimization::CartesianMinimizerMap::setup(), core::optimization::MinimizerMap::setup(), and core::optimization::symmetry::SymMinimizerMap::SymMinimizerMap().

◆ update_internal_coords()

void core::kinematics::AtomTree::update_internal_coords ( ) const
private

Update the internal coordinates using the xyz (cartesian) coords.

Note
Updated by Vikram K. Mulligan (vmull.nosp@m.igan.nosp@m.@flat.nosp@m.iron.nosp@m.insti.nosp@m.tute.nosp@m..org) on 24 Sept 2020 to use ordinary iteration instead of recursion.

update_internal_coords would be called in situations where we want to ensure that the internal degrees are valid, e.g., when we are about to set a new internal DOF.

Note
private method, const to allow lazy updating of
see usage in AtomTree.cc
these guys are const because of the lazy updating scheme we are using: they need to be called from within const accessing functions.
Updated by Vikram K. Mulligan (vmull.nosp@m.igan.nosp@m.@flat.nosp@m.iron.nosp@m.insti.nosp@m.tute.nosp@m..org) on 24 Sept 2020 to use ordinary iteration instead of recursion.

References core::kinematics::default_stub, internal_coords_need_updating_, root_, and xyz_coords_need_updating_.

Referenced by atom(), dof(), jump(), replace_tree(), set_dof(), set_jump(), torsion_angle(), and torsion_angle_dof_id().

◆ update_sequence_numbering()

void core::kinematics::AtomTree::update_sequence_numbering ( Size const  new_size,
utility::vector1< int > const &  old2new 
)

◆ update_xyz_coords()

void core::kinematics::AtomTree::update_xyz_coords ( ) const
private

◆ xyz()

PointPosition const & core::kinematics::AtomTree::xyz ( AtomID const &  id) const

get xyz position of an atom given its AtomID

retrieve the xyz position of an atom in the tree by its AtomID

References atom_pointer_, and update_xyz_coords().

Referenced by get_frag_descendant_local_xyz(), get_frag_parent_local_xyz(), set_xyz(), and stub_from_id().

Member Data Documentation

◆ atom_pointer_

AtomPointer2D core::kinematics::AtomTree::atom_pointer_
private

◆ dof_changeset_

AtomDOFChangeSet core::kinematics::AtomTree::dof_changeset_
mutableprivate

A list of the atoms that have had changed DOFs since the last refold.

Referenced by copy_coords(), operator=(), set_dof(), set_jump(), and update_xyz_coords().

◆ external_coordinate_residues_changed_

ResidueCoordinateChangeListOP core::kinematics::AtomTree::external_coordinate_residues_changed_
private

A list of residues that have had xyz coordinate changes since the last time the owning Conformation object has asked for an update.

Referenced by AtomTree(), clear(), copy_coords(), insert_single_fragment(), note_coordinate_change_registered(), operator=(), replace_tree(), residue_xyz_change_list_begin(), residue_xyz_change_list_end(), set_stub_transform(), update_sequence_numbering(), and update_xyz_coords().

◆ internal_coords_need_updating_

bool core::kinematics::AtomTree::internal_coords_need_updating_
mutableprivate

List of the jump atom ID's, excluding the root. Order matters (for movemap indexing)

Internal coords out of date?

Referenced by batch_set_xyz(), copy_coords(), insert_single_fragment(), operator=(), promote_sameresidue_nonjump_child(), replace_residue_subtree(), replace_tree(), set_jump_atom_stub_id(), set_stub_transform(), set_xyz(), update_internal_coords(), and update_xyz_coords().

◆ root_

AtomOP core::kinematics::AtomTree::root_
private

◆ this_weak_ptr_

AtomTreeCAP core::kinematics::AtomTree::this_weak_ptr_
private

A weak pointer to self (this).

The weak pointer must be provided to the AtomTree immediately after creation.

Referenced by attach_topological_observer(), operator=(), set_new_topology(), and set_weak_pointer_to_self().

◆ topological_match_to_

AtomTreeCAP core::kinematics::AtomTree::topological_match_to_
mutableprivate

pointer to the atom tree this tree has an exact topological match to since that tree was the last tree copied from without subsequence topological modifications – or at most one modification when that tree copied this tree's topology

Referenced by notify_topological_change(), operator=(), set_new_topology(), and topological_match_to().

◆ topological_observers_

utility::vector1< AtomTreeCAP > core::kinematics::AtomTree::topological_observers_
mutableprivate

pointers to all atom trees that are observing this tree's topology. On topological changes (including the destruction of this tree!), each of these trees have their topological_match_to_ pointers set to null and this list is cleared.

Referenced by attach_topological_observer(), detach_topological_observer(), set_new_topology(), and topological_observers().

◆ xyz_coords_need_updating_

bool core::kinematics::AtomTree::xyz_coords_need_updating_
mutableprivate

The documentation for this class was generated from the following files: