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

an object which makes rigid-body transformation with translational and rotational perturbation More...

#include <Jump.hh>

Public Types

typedef numeric::xyzVector< RealVector
 
typedef numeric::xyzMatrix< RealMatrix
 

Public Member Functions

 Jump ()
 construction More...
 
 Jump (const RT &src_rt)
 constructor with only RT More...
 
 Jump (Stub const &stub1, Stub const &stub2)
 get RT from two stubs and ZERO rb_delta More...
 
 Jump (const Jump &)=default
 copy constructor More...
 
Jumpoperator= (Jump const &)=default
 copy operator More...
 
void from_stubs (Stub const &stub1, Stub const &stub2)
 get a jump from two stubs More...
 
void from_bond_cst (utility::vector1< Vector > &atoms, utility::vector1< Real > const &csts)
 get a jump from a bond cst definition More...
 
void make_jump (Stub const &stub1, Stub &stub2) const
 make a jump from stub1 to stub2 More...
 
bool ortho_check () const
 check RT's orthogonality More...
 
bool nonzero_deltas () const
 check whether there is rb_delta not being transformed. More...
 
void reset ()
 reset RT, rb_delta and rb_center More...
 
void random_trans (const float dist_in)
 translate along a randomly chosen vector by dist_in More...
 
utility::vector1< Realgaussian_move (int const dir, float const trans_mag, float const rot_mag)
 Given the desired magnitude of the translation and rotation components, applies Gaussian perturbation to this jump. Return the move that was applied. More...
 
void gaussian_move_single_rb (int const dir, float const mag, int rb)
 
void rotation_by_matrix (Stub const &stub, Vector const &center, Matrix const &matrix)
 make a rotation "matrix" about the center "center" More...
 
void rotation_by_axis (Stub const &stub, Vector const &axis, Vector const &center, float const alpha)
 make a rotation of alpha degrees around axix and center More...
 
void translation_along_axis (Stub const &stub, Vector const &axis, float const dist)
 make a translation along axis by dist More...
 
void identity_transform ()
 reset to identity matrix, 0 translation More...
 
void reverse ()
 change the direction of a jump, e.g, upstream stub becomes downstream stub now. More...
 
Jump reversed () const
 Return a new jump that is the inverse transformation. More...
 
Matrix const & get_rotation () const
 get rotation matrix More...
 
Vector const & get_translation () const
 get translation vector More...
 
void set_rotation (Matrix const &R_in)
 set rotation matrix More...
 
void set_translation (Vector const &t)
 set translation vector More...
 
void fold_in_rb_deltas ()
 transform small changes in translation and rotation into jump More...
 
RT rt () const
 return the RT modeled by this jump --> makes new Jump, calls fold_in_rb_delte and returns RT More...
 
utility::vector1< Real > const & get_rb_delta (int const dir) const
 get rb_delta by direction More...
 
Real get_rb_delta (int const rb_no, int const dir) const
 get rb_delta by direction and rb_number More...
 
void set_rb_delta (int const rb_no, int const dir, Real const value)
 set rb_delta by direction and rb_number More...
 
void set_rb_deltas (int const dir, utility::vector1< Real > const &)
 set rb_deltas by direction More...
 
Vector const & get_rb_center (int const dir) const
 get rb_center by direction More...
 
void set_rb_center (const int dir, Stub const &stub, Vector const &center)
 set rb_center by direction More...
 
bool get_invert_upstream () const
 
bool get_invert_downstream () const
 
void set_invert (bool upstream, bool downstream)
 
bool operator== (Jump const &) const
 
bool operator!= (Jump const &other) const
 

Static Public Attributes

static const int TRANS_X = 1
 
static const int TRANS_Y = 2
 
static const int TRANS_Z = 3
 
static const int ROT_X = 4
 
static const int ROT_Y = 5
 
static const int ROT_Z = 6
 
static Matrix mirror_z_transform = numeric::xyzMatrix< core::Real >::rows( 1.,0.,0., 0.,1.,0., 0.,0.,-1. )
 

Private Member Functions

int rb_index (int const dir) const
 get the lookup index into rb_delta and rb_center from the direction dir. More...
 

Private Attributes

RT rt_
 translation and rotation More...
 
utility::vector1< utility::vector1< Real > > rb_delta
 changes to translation and rotation More...
 
utility::vector1< Vectorrb_center
 the center around which the rotation is performed More...
 
bool invert_upstream_
 are the upstream or downstream residues an inverted coordinate system (Z <- -X cross Y) More...
 
bool invert_downstream_
 

Friends

std::ostream & operator<< (std::ostream &os, const Jump &jump)
 stream output operator More...
 
std::istream & operator>> (std::istream &is, Jump &jump)
 stream input operator More...
 
Real distance (Jump const &a_in, Jump const &b_in)
 RT root squared deviation. More...
 

Detailed Description

an object which makes rigid-body transformation with translational and rotational perturbation

See AtomTree overview and concepts for details.

Member Typedef Documentation

◆ Matrix

typedef numeric::xyzMatrix< Real > core::kinematics::Jump::Matrix

◆ Vector

typedef numeric::xyzVector< Real > core::kinematics::Jump::Vector

Constructor & Destructor Documentation

◆ Jump() [1/4]

core::kinematics::Jump::Jump ( )
inline

construction

◆ Jump() [2/4]

core::kinematics::Jump::Jump ( const RT src_rt)
inline

constructor with only RT

◆ Jump() [3/4]

core::kinematics::Jump::Jump ( Stub const &  stub1,
Stub const &  stub2 
)
inline

get RT from two stubs and ZERO rb_delta

◆ Jump() [4/4]

core::kinematics::Jump::Jump ( const Jump )
default

copy constructor

Member Function Documentation

◆ fold_in_rb_deltas()

void core::kinematics::Jump::fold_in_rb_deltas ( )

◆ from_bond_cst()

void core::kinematics::Jump::from_bond_cst ( utility::vector1< Vector > &  atoms,
utility::vector1< Real > const &  csts 
)

get a jump from a bond cst definition

this function translate a virtual bond type constraints into jump transformation. in atoms vector, A1, A2, A3, B1, B2, B3, and A1-B1 forms the virtual bond, like A3-A2-A1—B1-B2-B3 in csts vector, disAB, angleA and B, didedral A, AB and B. B1, B2 and B3 are moved given cst definition.

References distance, core::kinematics::RT::from_stubs(), core::kinematics::Stub::is_orthogonal(), rb_delta, rt_, core::kinematics::Stub::spherical(), and core::kinematics::ZERO().

Referenced by protocols::jumping::ResiduePairJump::build_cst_conformer_jumps().

◆ from_stubs()

void core::kinematics::Jump::from_stubs ( Stub const &  stub1,
Stub const &  stub2 
)

get a jump from two stubs

Note
: we dont reset rb_center!!!!!!!!!

References core::kinematics::RT::from_stubs(), rb_delta, rt_, and core::kinematics::ZERO().

Referenced by core::kinematics::tree::JumpAtom::update_internal_coords().

◆ gaussian_move()

utility::vector1< Real > core::kinematics::Jump::gaussian_move ( int const  dir,
float const  trans_mag,
float const  rot_mag 
)

Given the desired magnitude of the translation and rotation components, applies Gaussian perturbation to this jump. Return the move that was applied.

clear existing rb_delta first if any and then apply the gaussian move Return the move that was applied.

References fold_in_rb_deltas(), get_rb_delta(), core::scoring::rg, and set_rb_delta().

Referenced by protocols::enzdes::PredesignPerturbMover::apply(), protocols::ligand_docking::RigidSearchMover::apply(), protocols::rigid::RigidBodyPerturbMover::apply(), protocols::rigid::RigidBodyPerturbNoCenterMover::apply(), and protocols::rigid::gaussian_jump_move().

◆ gaussian_move_single_rb()

void core::kinematics::Jump::gaussian_move_single_rb ( int const  dir,
float const  mag,
int  rb 
)

◆ get_invert_downstream()

bool core::kinematics::Jump::get_invert_downstream ( ) const
inline

◆ get_invert_upstream()

bool core::kinematics::Jump::get_invert_upstream ( ) const
inline

◆ get_rb_center()

Jump::Vector const & core::kinematics::Jump::get_rb_center ( int const  dir) const
inline

get rb_center by direction

References rb_center, and rb_index().

Referenced by core::kinematics::tree::JumpAtom::get_dof_axis_and_end_pos().

◆ get_rb_delta() [1/2]

const utility::vector1< Real > & core::kinematics::Jump::get_rb_delta ( int const  dir) const
inline

◆ get_rb_delta() [2/2]

Real core::kinematics::Jump::get_rb_delta ( int const  rb_no,
int const  dir 
) const
inline

get rb_delta by direction and rb_number

References rb_delta, and rb_index().

◆ get_rotation()

Jump::Matrix const & core::kinematics::Jump::get_rotation ( ) const
inline

◆ get_translation()

Jump::Vector const & core::kinematics::Jump::get_translation ( ) const
inline

get translation vector

References core::kinematics::RT::get_translation(), and rt_.

Referenced by protocols::multistage_rosetta_scripts::cluster::metrics::JumpMetric::analyze(), protocols::rigid::UniformRigidBodyMover::apply(), protocols::cryst::UpdateCrystInfo::apply(), protocols::rigid::RigidBodyDofRandomizeMover::apply(), protocols::rigid::RigidBodyDofTransMover::apply(), protocols::matdes::SymDofMover::apply(), protocols::simple_moves::PeriodicBoxMover::apply(), core::kinematics::inverse::calculate_new_jump(), protocols::simple_moves::PeriodicBoxMover::change_volume_move(), protocols::matdes::GetRBDOFValues::compute(), protocols::simple_moves::PeriodicBoxMover::dump_ASU(), core::import_pose::atom_tree_diffs::dump_atom_tree_diff(), core::scoring::loop_graph::evaluator::SixDTransRotPotential::evaluate(), core::scoring::loop_graph::evaluator::SixDTransRotPotential::get_derivative(), protocols::ligand_docking::ga_ligand_dock::GALigandDock::make_starting_pose_for_virtual_screening(), protocols::recces::sampler::rna::MC_RNA_OneJump::operator++(), operator==(), protocols::rna::denovo::output::RNA_FragmentMonteCarloOutputter::output_jump_information(), protocols::cryst::DockLatticeMover::perturb_lattice(), protocols::simple_moves::PeriodicBoxMover::perturb_molecule_move(), protocols::rna::denovo::movers::RNA_DeNovoMasterMover::randomize_rigid_body_orientations(), protocols::forge::remodel::RemodelLoopMover::repeat_generation_with_additional_residue(), protocols::features::PoseConformationFeatures::report_features_implementation(), protocols::simple_moves::PeriodicBoxMover::report_thermodynamics(), protocols::environment::ClientMover::sandboxed_copy(), protocols::loop_grower::SheetPositions::SheetPositions(), protocols::rna::denovo::coarse::MultipleDomainMover::slide_back_to_origin(), protocols::rna::denovo::coarse::MultipleDomainMover::try_to_slide_into_contact(), and protocols::ligand_docking::ga_ligand_dock::LigandConformer::update_conf().

◆ identity_transform()

void core::kinematics::Jump::identity_transform ( )

reset to identity matrix, 0 translation

References core::kinematics::RT::identity_transform(), rb_center, rb_delta, rt_, and core::kinematics::ZERO().

◆ make_jump()

void core::kinematics::Jump::make_jump ( Stub const &  stub1,
Stub stub2 
) const

◆ nonzero_deltas()

bool core::kinematics::Jump::nonzero_deltas ( ) const

check whether there is rb_delta not being transformed.

References rb_delta, and core::simple_metrics::metrics::sum.

Referenced by rt().

◆ operator!=()

bool core::kinematics::Jump::operator!= ( Jump const &  other) const
inline

References operator==().

◆ operator=()

Jump& core::kinematics::Jump::operator= ( Jump const &  )
default

copy operator

◆ operator==()

bool core::kinematics::Jump::operator== ( Jump const &  other) const

References get_rotation(), and get_translation().

Referenced by operator!=().

◆ ortho_check()

bool core::kinematics::Jump::ortho_check ( ) const

check RT's orthogonality

References core::kinematics::RT::ortho_check(), and rt_.

Referenced by core::conformation::Conformation::set_jump().

◆ random_trans()

void core::kinematics::Jump::random_trans ( const float  dist_in)

translate along a randomly chosen vector by dist_in

choose a random unit vector as the direction of translation, by selecting its spherical coordinates phi(0-180) and theta(0-360) then move dist_in along that direction

References fold_in_rb_deltas(), core::scoring::rg, rt_, and core::kinematics::RT::set_translation().

Referenced by protocols::topology_broker::BasicJumpClaimer::initialize_dofs().

◆ rb_index()

int core::kinematics::Jump::rb_index ( int const  dir) const
inlineprivate

get the lookup index into rb_delta and rb_center from the direction dir.

dir should be 1 or -1

Referenced by get_rb_center(), get_rb_delta(), set_rb_delta(), and set_rb_deltas().

◆ reset()

void core::kinematics::Jump::reset ( )

◆ reverse()

void core::kinematics::Jump::reverse ( )

◆ reversed()

Jump core::kinematics::Jump::reversed ( ) const

Return a new jump that is the inverse transformation.

References fold_in_rb_deltas(), and reverse().

Referenced by protocols::motifs::Motif::Motif().

◆ rotation_by_axis()

void core::kinematics::Jump::rotation_by_axis ( Stub const &  stub,
Vector const &  axis,
Vector const &  center,
float const  alpha 
)

make a rotation of alpha degrees around axix and center

stub should be the Stub of the upstream jump. this does a rotation by alpha degrees around the specified the axis and center, both of which are written in xyz frame

Parameters
alphadegrees

References core::conformation::membrane::center, fold_in_rb_deltas(), and rotation_by_matrix().

Referenced by protocols::nonlocal::HelixRotate::apply(), protocols::protein_interface_design::movers::SpinMover::apply(), protocols::rigid::WholeBodyRotationMover::apply(), protocols::rigid::RigidBodySpinMover::apply(), protocols::rigid::RigidBodyDeterministicSpinMover::apply(), protocols::denovo_design::movers::RotateSegmentMover::apply(), and protocols::rigid::RigidBodyTiltMover::tilt().

◆ rotation_by_matrix()

void core::kinematics::Jump::rotation_by_matrix ( Stub const &  stub,
Vector const &  center,
Matrix const &  matrix 
)

◆ rt()

RT core::kinematics::Jump::rt ( ) const

◆ set_invert()

void core::kinematics::Jump::set_invert ( bool  upstream,
bool  downstream 
)

◆ set_rb_center()

void core::kinematics::Jump::set_rb_center ( const int  dir,
Stub const &  stub,
Vector const &  center 
)

set rb_center by direction

stub defines the local reference frame for the corresponding direction of the jump (forward == folding direction, backward = reverse) by which to transform center
center is an xyz position in the absolute (protein) reference frame
the column vectors of stub are unit vectors defined in the absolute reference frame( the frame in which center is defined)
so left-multiplying by stub.M.transposed() puts an absolute ref frame vector (like center) into the local jump reference frame
dir is the direction through the jump, 1 == forward == folding order, -1 == backward == reverse folding order
if dir == 1, the stub should be the stub for the downstream jump atom, (eg returned by "folder.downstream_jump_stub( jump_number )" or "downstream_jump_atom.get_stub()" ) and the center is the center of rotation for the downstream segment (ie the segment which is built later during folding). The center determines how the forward rb_deltas are interpreted, and is the center for rigid-body rotation during minimization.
if dir == -1, the stub should be the stub for the upstream jump atom, (eg returned by folder.upstream_jump_stub( jump_number ) ) and the center is the center of rotation for the upstream segment. The center determines how the reverse rb_deltas are interpreted, eg in perturbation moves when we modify the reverse rb-deltas directly, eg in the call "Jump::gaussian_move(-1,tmag,rmag)"

References core::conformation::membrane::center, fold_in_rb_deltas(), core::kinematics::Stub::M, rb_center, and core::kinematics::Stub::v.

Referenced by protocols::enzdes::PredesignPerturbMover::apply(), protocols::ligand_docking::RigidSearchMover::apply(), protocols::rigid::RigidBodyPerturbMover::apply(), protocols::rigid::RigidBodyRandomizeMover::apply(), protocols::rigid::RigidBodySpinMover::apply(), protocols::rigid::RigidBodyDeterministicSpinMover::apply(), protocols::ligand_docking::LigandDockProtocol::optimize_orientation3(), and protocols::rigid::RigidBodyTiltMover::tilt().

◆ set_rb_delta()

void core::kinematics::Jump::set_rb_delta ( int const  rb_no,
int const  dir,
Real const  value 
)

◆ set_rb_deltas()

void core::kinematics::Jump::set_rb_deltas ( int const  dir,
utility::vector1< Real > const &  rb_deltas 
)

set rb_deltas by direction

References rb_delta, and rb_index().

Referenced by protocols::rigid::RigidBodyPerturbMover::apply().

◆ set_rotation()

void core::kinematics::Jump::set_rotation ( Matrix const &  R_in)

◆ set_translation()

void core::kinematics::Jump::set_translation ( Vector const &  t)

set translation vector

References fold_in_rb_deltas(), rt_, core::kinematics::RT::set_translation(), and protocols::hybridization::t.

Referenced by protocols::rigid::UniformRigidBodyMover::apply(), protocols::dna::SeparateDnaFromNonDna::apply(), protocols::rigid::RigidBodyDofRandomizeMover::apply(), protocols::rigid::RigidBodyDofTransMover::apply(), protocols::rigid::RigidBodyDofPerturbMover::apply(), protocols::matdes::SymDofMover::apply(), protocols::loop_grower::LoopPartialSolution::apply_sheets(), core::kinematics::inverse::calculate_new_jump(), protocols::simple_moves::PeriodicBoxMover::change_volume_move(), protocols::toolbox::sample_around::do_xy_scan(), protocols::ligand_docking::ga_ligand_dock::GALigandDock::make_starting_pose_for_virtual_screening(), protocols::recces::sampler::rna::MC_RNA_OneJump::operator++(), protocols::cryst::DockLatticeMover::perturb_lattice(), core::import_pose::atom_tree_diffs::pose_from_atom_tree_diff(), protocols::stepwise::modeler::protein::StepWiseProteinBackboneSampler::prepare_ghost_pose(), protocols::rna::denovo::movers::RNA_DeNovoMasterMover::randomize_rigid_body_orientations(), protocols::forge::remodel::RemodelLoopMover::repeat_generation_with_additional_residue(), protocols::stepwise::monte_carlo::mover::AddMover::setup_initial_jump(), protocols::rna::denovo::coarse::MultipleDomainMover::slide_back_to_origin(), protocols::stepwise::modeler::split_pose(), protocols::ligand_docking::ga_ligand_dock::LigandConformer::to_pose(), and protocols::rna::denovo::coarse::MultipleDomainMover::try_to_slide_into_contact().

◆ translation_along_axis()

void core::kinematics::Jump::translation_along_axis ( Stub const &  stub,
Vector const &  axis,
float const  dist 
)

Friends And Related Function Documentation

◆ distance

Real distance ( Jump const &  a_in,
Jump const &  b_in 
)
friend

RT root squared deviation.

Referenced by from_bond_cst().

◆ operator<<

std::ostream& operator<< ( std::ostream &  os,
const Jump jump 
)
friend

stream output operator

◆ operator>>

std::istream& operator>> ( std::istream &  is,
Jump jump 
)
friend

stream input operator

Member Data Documentation

◆ invert_downstream_

bool core::kinematics::Jump::invert_downstream_
private

◆ invert_upstream_

bool core::kinematics::Jump::invert_upstream_
private

are the upstream or downstream residues an inverted coordinate system (Z <- -X cross Y)

Referenced by get_invert_upstream(), make_jump(), reset(), and set_invert().

◆ mirror_z_transform

numeric::xyzMatrix< core::Real > core::kinematics::Jump::mirror_z_transform = numeric::xyzMatrix< core::Real >::rows( 1.,0.,0., 0.,1.,0., 0.,0.,-1. )
static

◆ rb_center

utility::vector1< Vector > core::kinematics::Jump::rb_center
private

the center around which the rotation is performed

3x2 table, for each of the two folding directions, the rotation center is written in the local frame of the downstream stub.

Referenced by fold_in_rb_deltas(), get_rb_center(), identity_transform(), make_jump(), reset(), and set_rb_center().

◆ rb_delta

utility::vector1< utility::vector1<Real> > core::kinematics::Jump::rb_delta
private

changes to translation and rotation

6x2 table, for each of the two folding directions, the first three are for translations along xyz, and the next three are for rotatation around xyz axes.

Referenced by fold_in_rb_deltas(), from_bond_cst(), from_stubs(), get_rb_delta(), identity_transform(), make_jump(), nonzero_deltas(), reset(), set_rb_delta(), and set_rb_deltas().

◆ ROT_X

const int core::kinematics::Jump::ROT_X = 4
static

◆ ROT_Y

const int core::kinematics::Jump::ROT_Y = 5
static

◆ ROT_Z

const int core::kinematics::Jump::ROT_Z = 6
static

◆ rt_

RT core::kinematics::Jump::rt_
private

◆ TRANS_X

const int core::kinematics::Jump::TRANS_X = 1
static

◆ TRANS_Y

const int core::kinematics::Jump::TRANS_Y = 2
static

◆ TRANS_Z

const int core::kinematics::Jump::TRANS_Z = 3
static

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