Package rosetta :: Package core :: Package pose :: Module _core_pose_
[hide private]
[frames] | no frames]

Module _core_pose_

Classes [hide private]
  CrystInfo
core/pose/CrystInfo.hh:27
  CrystInfoAP
  CrystInfoCAP
  MiniPose
lightweight version of the pose with stuff I need.
  MiniPoseAP
  MiniPoseCAP
  PDBInfo
maintains pdb residue & atom information inside a Pose Upon creation of new residue records, e.g.
  PDBInfoAP
  PDBInfoCAP
  PDBPoseMap
PDBPoseMap can be queried with PDB information (chain, sequence position) and returns a pose's resid position.
  PDBPoseMapAP
  PDBPoseMapCAP
  Pose
A molecular system including residues, kinematics, and energies The Pose class represents a molecular system (protein-dna-ligand...) as a container of Rosetta Residue objects together with a Conformation object that defines how internal coordinate changes propagate through the system and an Energies object that stores information from the last energy evaluation.
  PoseAP
  PoseCAP
  PoseCoordPickMode
core/pose/xyzStripeHashPose.fwd.hh:19
  PosePyObserver
  Py_basic_datacache_DataCache_CacheableObserver
  RemarkInfo
core/pose/Remarks.hh:32
  RemarkInfoAP
  RemarkInfoCAP
  Remarks
core/pose/Remarks.hh:50
  RemarksAP
  RemarksCAP
  SingletonBase_T_core_chemical_ChemicalManager_T
SingletonBase is meant to serve as a base class for singleton classes in Rosetta handling the initialization of the singleton in a thread-safe way.
  SingletonBase_T_core_chemical_ChemicalManager_TAP
  SingletonBase_T_core_chemical_ChemicalManager_TCAP
  UnrecognizedAtomRecord
info about an atom in a unrecognized res (not in pose, but we want to remember it)
  UnrecognizedAtomRecordAP
  UnrecognizedAtomRecordCAP
  xyzStripeHashPose
core/pose/xyzStripeHashPose.hh:31
  xyzStripeHashPoseAP
  xyzStripeHashPoseCAP
Functions [hide private]
None :
QQQ_PosePyObserverTesterFunction()
C++ signature :...
 
addVirtualResAsRoot(...)
addVirtualResAsRoot( (xyzVector_Real)xyz, (Pose)pose) -> None : Adds a virtual residue to the end of the pose at the specified location.
 
add_comment(...)
add_comment( (Pose)pose, (str)key, (str)val) -> None : Adds a key-value pair to the STRING_MAP in the Pose DataCache.
 
add_lower_terminus_type_to_pose_residue(...)
add_lower_terminus_type_to_pose_residue( (Pose)pose, (int)seqpos) -> None : core/pose/util.hh:436
 
add_score_line_string(...)
add_score_line_string( (Pose)pose, (str)key, (str)val) -> None : Dumps a pdb with comments at end of file Sets a PDB-style REMARK entry in the Pose.
 
add_upper_terminus_type_to_pose_residue(...)
add_upper_terminus_type_to_pose_residue( (Pose)pose, (int)seqpos) -> None : core/pose/util.hh:443
 
add_variant_type_to_pose_residue(...)
add_variant_type_to_pose_residue( (Pose)pose, (VariantType)variant_type, (int)seqpos) -> None : Construct a variant of an existing pose residue.
 
add_variant_type_to_residue(...)
add_variant_type_to_residue( (Residue)old_rsd, (VariantType)variant_type, (Pose)pose) -> Residue : Construct a variant of an existing residue.
 
annotated_to_oneletter_sequence(...)
annotated_to_oneletter_sequence( (str)annotated_seq) -> str : Returns the oneletter_sequence that corresponds to the given annotated sequence.
 
append_pose_to_pose(...)
append_pose_to_pose( (Pose)pose1, (Pose)pose2 [, (bool)new_chain=True]) -> None : Append residues of pose2 to pose1.
 
append_subpose_to_pose(...)
append_subpose_to_pose( (Pose)pose1, (Pose)pose2, (int)start_res, (int)end_res [, (bool)new_chain=True]) -> None : Append specified residues of pose2 to pose1.
 
atom_id_to_named_atom_id(...)
atom_id_to_named_atom_id( (AtomID)atom_id, (Pose)pose) -> NamedAtomID : core/pose/util.hh:353
 
canonical_atom_count(...)
canonical_atom_count( (Pose)pose) -> int : count the number of canonical amino acid atoms in the pose
 
canonical_residue_count(...)
canonical_residue_count( (Pose)pose) -> int : count the number of canonical residues in the pose
 
center_of_mass(...)
center_of_mass( (Pose)pose, (int)start, (int)stop) -> xyzVector_Real : core/pose/util.hh:669
 
clearPoseExtraScore(...)
clearPoseExtraScore( (Pose)pose, (str)name) -> None : core/pose/util.hh:170
 
clearPoseExtraScores(...)
clearPoseExtraScores( (Pose)pose) -> None : core/pose/util.hh:174
 
compare_atom_coordinates(...)
compare_atom_coordinates( (Pose)lhs, (Pose)rhs [, (int)n_dec_places=3]) -> bool : this function compares pose atom coordinates for equality; it is not the == operator because it does not compare all pose data.
 
compare_binary_protein_silent_struct(...)
compare_binary_protein_silent_struct( (Pose)lhs, (Pose)rhs) -> bool : this function compares poses for equality up to the information stored in the binary protein silent struct format.
 
conf2pdb_chain(...)
conf2pdb_chain( (Pose)pose) -> object : get Conformation chain -> PDBInfo chain mapping @remarks Any chains whose PDBInfo chain records are marked entirely as PDBInfo::empty_record() will be mapped to that character.
 
convert_from_std_map(...)
convert_from_std_map( (object)atom_map, (Pose)pose) -> AtomID_Map_T_core_id_AtomID_T : core/pose/util.hh:687
 
correctly_add_cutpoint_variants(...)
correctly_add_cutpoint_variants( (Pose)pose, (int)cutpoint_res [, (bool)check_fold_tree=True]) -> None : core/pose/util.hh:696
 
create_subpose(...)
create_subpose( (Pose)src, (vector1_Size)positions, (FoldTree)f, (Pose)pose) -> None : Create a subpose of the src pose.
 
delete_comment(...)
delete_comment( (Pose)pose, (str)key) -> None : Deletes the entry in the STRING_MAP associated with the given key.
 
dump_comment_pdb(...)
dump_comment_pdb( (str)file_name, (Pose)pose) -> None : dumps pose+ comments to pdb file
 
energy_from_pose(...)
energy_from_pose( (Pose)pose, (str)sc_type) -> float : core/pose/util.hh:387
 
extract_tag_from_pose(...)
extract_tag_from_pose( (Pose)pose) -> str : Returns a string giving the pose's tag if there is such a thing or "UnknownTag" otherwise.
 
getPoseExtraScore(...)
getPoseExtraScore( (Pose)pose, (str)name, (str)value) -> bool : core/pose/util.hh:181
 
get_all_comments(...)
get_all_comments( (Pose)pose) -> map_string_string : Gets a map< string, string > representing comments about the Pose in the form of key-value pairs.
 
get_all_score_line_strings(...)
get_all_score_line_strings( (Pose)pose) -> map_string_string : Gets a map< string, string > representing score_line_strings about the Pose in the form of key-value pairs.
 
get_center_of_mass(...)
get_center_of_mass( (Pose)pose) -> xyzVector_Real : Get center of mass of a pose.
 
get_chain_from_chain_id(...)
get_chain_from_chain_id( (int)chain_id, (Pose)pose) -> str : core/pose/util.hh:508
 
get_chain_from_jump_id(...)
get_chain_from_jump_id( (int)jump_id, (Pose)pose) -> str : core/pose/util.hh:526
 
get_chain_id_from_chain(...)
get_chain_id_from_chain( (str)chain, (Pose)pose) -> int : core/pose/util.hh:496
 
get_chain_id_from_jump_id(...)
get_chain_id_from_jump_id( (int)jump_id, (Pose)pose) -> int : core/pose/util.hh:523
 
get_chain_ids_from_chain(...)
get_chain_ids_from_chain( (str)chain, (Pose)pose) -> vector1_Size : core/pose/util.hh:502
 
get_chain_ids_from_chains(...)
get_chain_ids_from_chains( (vector1_string)chains, (Pose)pose) -> vector1_Size : core/pose/util.hh:505
 
get_chain_residues(...)
get_chain_residues( (Pose)pose, (int)chain_id) -> object : core/pose/util.hh:529
 
get_comment(...)
get_comment( (Pose)pose, (str)key, (str)val) -> bool : Attempts to access the entry in the STRING_MAP associated with the given key.
 
get_hash_excluding_chain(...)
get_hash_excluding_chain( (str)chain, (Pose)pose) -> int : core/pose/util.hh:572
 
get_hash_from_chain(...)
get_hash_from_chain( (str)chain, (Pose)pose) -> int : core/pose/util.hh:569
 
get_jump_id_from_chain(...)
get_jump_id_from_chain( (str)chain, (Pose)pose) -> int : core/pose/util.hh:514
 
get_jump_id_from_chain_id(...)
get_jump_id_from_chain_id( (int)chain_id, (Pose)pose) -> int : core/pose/util.hh:488
 
get_jump_ids_from_chain(...)
get_jump_ids_from_chain( (str)chain, (Pose)pose) -> vector1_Size : core/pose/util.hh:520
 
get_jump_ids_from_chain_ids(...)
get_jump_ids_from_chain_ids( (set_Size)chain_ids, (Pose)pose) -> set_Size : core/pose/util.hh:485
 
get_resnum(...)
get_resnum( (Tag)tag_ptr, (Pose)pose [, (str)prefix='']) -> int : core/pose/selection.hh:37
 
get_resnum_list(...)
get_resnum_list( (str)str, (Pose)pose) -> set_Size : returns a resnum list directly from a string
 
get_resnum_list_ordered(...)
get_resnum_list_ordered( (str)str, (Pose)pose) -> vector1_Size : returns a resnum list directly from a string, preserving order
 
get_score_line_string(...)
get_score_line_string( (Pose)pose, (str)key, (str)val) -> bool : core/pose/util.hh:250
 
get_sequence_len(...)
get_sequence_len( (str)sequence_in) -> int : Get the real length of a annotated sequence
 
get_sha1_hash_excluding_chain(...)
get_sha1_hash_excluding_chain( (str)chain, (Pose)pose) -> str : core/pose/util.hh:574
 
hasPoseExtraScore(...)
hasPoseExtraScore( (Pose)pose, (str)name) -> bool : core/pose/util.hh:159
 
has_chain(...)
has_chain( (int)chain_id, (Pose)pose) -> bool : core/pose/util.hh:482
 
initialize_atomid_map(...)
initialize_atomid_map( (AtomID_Map_T_core_id_AtomID_T)atom_map, (Conformation)conformation, (AtomID)value) -> None : Initialize an AtomID_Map for a given Conformation using a specified fill value
 
initialize_atomid_map_AtomID(...)
initialize_atomid_map_AtomID( (AtomID_Map_T_core_id_AtomID_T)atom_map, (Conformation)conformation, (AtomID)value) -> None : core/pose/util.tmpl.hh:379
 
initialize_disulfide_bonds(...)
initialize_disulfide_bonds( (Pose)pose) -> None : detect and fix disulfide bonds
 
is_ideal_pose(...)
is_ideal_pose( (Pose)pose) -> bool : Returns true if the <pose> geometry is ideal [in] pose The Pose to check.
 
is_ideal_position(...)
is_ideal_position( (int)seqpos, (Pose)pose) -> bool : Returns true if the <pose> geometry is ideal in position <seqpos> [in] pose The Pose to check.
 
is_lower_terminus(...)
is_lower_terminus( (Pose)pose, (int)resid) -> bool : checks to see if this is a lower chain ending more intelligently than just checking residue variants
 
is_position_conserved_residue(...)
is_position_conserved_residue( (Pose)pose, (int)residue) -> bool : Returns true if <residue> is positionally conserved, false otherwise
 
is_upper_terminus(...)
is_upper_terminus( (Pose)pose, (int)resid) -> bool : checks to see if this is a lower chain ending more intelligently than just checking residue variants
 
jumps_from_pose(...)
jumps_from_pose( (Pose)pose, (set_int)jumps) -> None : Retrieves jump information from <pose>, storing the result in <jumps>.
 
make_atom_map(...)
make_atom_map( (Pose)p, (PoseCoordPickMode)m) -> AtomID_Map_T_double_T : core/pose/xyzStripeHashPose.hh:29
 
make_pose_from_saccharide_sequence(...)
make_pose_from_saccharide_sequence( (Pose)pose, (str)sequence [, (str)type_set_name='fa_standard' [, (bool)auto_termini=True]]) -> None : Create a Pose from an annotated, linear, IUPAC polysaccharide sequence <sequence> with residue type set name <type_set_name> and store it in <pose>.
 
make_pose_from_sequence(...)
make_pose_from_sequence( (Pose)pose, (str)sequence, (str)type_set_name [, (bool)auto_termini=True]) -> None : Creates a Pose from the annotated protein sequence <sequence> with the desired <type_set_name> and stores it in <pose> : any existing data in <pose> is cleared, auto_termini mark position 1, last_residue with lower, upper termini; default true
 
mass(...)
mass( (int)begin, (int)end, (Pose)pose) -> float : core/pose/util.hh:566
 
named_atom_id_to_atom_id(...)
named_atom_id_to_atom_id( (NamedAtomID)named_atom_id, (Pose)pose [, (bool)raise_exception=True]) -> AtomID : core/pose/util.hh:360
 
named_stub_id_to_stub_id(...)
named_stub_id_to_stub_id( (NamedStubID)named_stub_id, (Pose)pose) -> StubID : core/pose/util.hh:372
 
noncanonical_atom_count(...)
noncanonical_atom_count( (Pose)pose) -> int : count the number of non-canonical amino acids in thepose
 
noncanonical_chi_count(...)
noncanonical_chi_count( (Pose)pose) -> int : count the number of non-canonical chi angles in the pose
 
noncanonical_residue_count(...)
noncanonical_residue_count( (Pose)pose) -> int : count the number of non-canonical residues in the pose
 
nres_protein(...)
nres_protein( (Pose)pose) -> int : Number of protein residues in the pose No virtuals, membrane residues or embedding residues counted
 
num_atoms(...)
num_atoms( (int)begin, (int)end, (Pose)pose) -> int : core/pose/util.hh:541
 
num_chi_angles(...)
num_chi_angles( (int)begin, (int)end, (Pose)pose) -> int : core/pose/util.hh:559
 
num_hbond_acceptors(...)
num_hbond_acceptors( (int)begin, (int)end, (Pose)pose) -> int : core/pose/util.hh:547
 
num_hbond_donors(...)
num_hbond_donors( (int)begin, (int)end, (Pose)pose) -> int : core/pose/util.hh:553
 
num_heavy_atoms(...)
num_heavy_atoms( (int)begin, (int)end, (Pose)pose) -> int : core/pose/util.hh:535
 
parse_resnum(...)
parse_resnum( (str)resnum, (Pose)pose) -> int : Extracts a residue number from a string.
 
parse_sequence(...)
parse_sequence( (str)sequence_in, (vector1_string)fullname_list, (std_vector_Size)oneletter_to_fullname_index, (str)one_letter_sequence) -> None : Parse the input annotated sequence
 
partition_pose_by_jump(...)
partition_pose_by_jump( (Pose)src, (int)jump_number, (Pose)partner1, (Pose)partner2) -> None : core/pose/util.hh:125
 
pdbslice(...)
pdbslice( (Pose)pose, (vector1_Size)slice_res) -> None : Create a subpose of the src pose -- figures out a reasonable fold tree.
 
pose_from_saccharide_sequence(...)
pose_from_saccharide_sequence( (str)sequence [, (str)type_set_name='fa_standard' [, (bool)auto_termini=True]]) -> Pose : Return a Pose from an annotated, linear, IUPAC polysaccharide sequence <sequence> with residue type set name <type_set_name>.
 
pose_max_nbr_radius(...)
pose_max_nbr_radius( (Pose)pose) -> float : returns a Distance
 
pose_residue_is_terminal(...)
pose_residue_is_terminal( (Pose)pose, (int)resid) -> bool : returns true if the given residue in the pose is a chain ending or has upper/lower terminal variants
 
read_comment_pdb(...)
read_comment_pdb( (str)file_name, (Pose)pose) -> None : Reads the comments from the pdb file and adds it into comments
 
read_psipred_ss2_file(...)
read_psipred_ss2_file( (Pose)pose) -> vector1_char : core/pose/util.hh:141
 
remove_ligand_canonical_residues(...)
remove_ligand_canonical_residues( (Pose)pose) -> None : this function removes all residues with both UPPER and LOWER terminus types.
 
remove_lower_terminus_type_from_pose_residue(...)
remove_lower_terminus_type_from_pose_residue( (Pose)pose, (int)seqpos) -> None : core/pose/util.hh:450
 
remove_nonprotein_residues(...)
remove_nonprotein_residues( (Pose)pose) -> None : this function removes all residues from the pose which are not protein residues.
 
remove_upper_terminus_type_from_pose_residue(...)
remove_upper_terminus_type_from_pose_residue( (Pose)pose, (int)seqpos) -> None : core/pose/util.hh:457
 
remove_variant_type_from_pose_residue(...)
remove_variant_type_from_pose_residue( (Pose)pose, (VariantType)variant_type, (int)seqpos) -> None : Construct a non-variant of an existing pose residue.
 
remove_variant_type_from_residue(...)
remove_variant_type_from_residue( (Residue)old_rsd, (VariantType)variant_type, (Pose)pose) -> Residue : Remove variant from an existing residue.
 
remove_virtual_residues(...)
remove_virtual_residues( (Pose)pose) -> None : Removes all virtual residues from <pose>
 
renumber_pdbinfo_based_on_conf_chains(...)
renumber_pdbinfo_based_on_conf_chains( (Pose)pose [, (bool)fix_chains=True [, (bool)start_from_existing_numbering=True [, (bool)keep_insertion_codes=False [, (bool)rotate_chain_ids=False]]]]) -> bool : renumber PDBInfo based on Conformation chains; each chain starts from 1 [in,out] pose The Pose to modify.
 
replace_pose_residue_copying_existing_coordinates(...)
replace_pose_residue_copying_existing_coordinates( (Pose)pose, (int)seqpos, (ResidueType)new_rsd_type) -> None : core/pose/util.hh:405
 
residue_center_of_mass(...)
residue_center_of_mass( (Pose)pose, (int)start, (int)stop) -> int : core/pose/util.hh:676
 
residue_types_from_saccharide_sequence(...)
residue_types_from_saccharide_sequence( (str)sequence, (ResidueTypeSet)residue_set) -> object : Return a list of carbohydrate ResidueTypes corresponding to an annotated, linear, IUPAC polysaccharide sequence.
 
residue_types_from_sequence(...)
residue_types_from_sequence( (str)sequence_in, (ResidueTypeSet)residue_set [, (bool)auto_termini=True]) -> object : return a list of ResidueTypes corresponding to an annotated protein sequence [in] sequence_in an annotated sequence [in] residue_set the desired residue set [in] auto_termini mark position 1, last_residue with lower, upper termini; default true
 
return_nearest_residue(...)
return_nearest_residue( (Pose)pose, (int)begin, (int)end, (xyzVector_Real)center) -> int : core/pose/util.hh:684
 
sequence_map_from_pdbinfo(...)
sequence_map_from_pdbinfo( (Pose)first, (Pose)second) -> SequenceMapping : Create a sequence map of first pose onto the second, matching the PDBInfo If the PDBInfo of either Pose is missing or invalid, do a simple sequence alignment matching.
 
setPoseExtraScore(...)
setPoseExtraScore( (Pose)pose, (str)name, (str)value) -> None : core/pose/util.hh:187
 
set_ss_from_phipsi(...)
set_ss_from_phipsi( (Pose)pose) -> None : Analyzes <pose> residue phi/psi sets and guesses the secondary structure, ideally dssp should be used for that
 
setup_dof_mask_from_move_map(...)
setup_dof_mask_from_move_map( (MoveMap)mm, (Pose)pose, (object)dof_mask) -> None : convert from allow-bb/allow-chi MoveMap to simple DOF_ID boolean mask needed by the minimizer
 
setup_dof_to_torsion_map(...)
setup_dof_to_torsion_map( (Pose)pose, (object)dof_map) -> None : set up a map to look up TORSION_ID by DOF_ID (Map[DOF_ID] = TORISION_ID)
 
sort_pose_by_score(...)
sort_pose_by_score( (Pose)pose1, (Pose)pose2) -> bool : core/pose/util.hh:379
 
stub_id_to_named_stub_id(...)
stub_id_to_named_stub_id( (StubID)stub_id, (Pose)pose) -> NamedStubID : core/pose/util.hh:366
 
swap_transform(...)
swap_transform( (int)jump_num, (RT)xform, (Pose)pose) -> None : Updates the rigid-body transform of the specified jump in <pose>
 
tag_from_pose(...)
tag_from_pose( (Pose)pose) -> str : ////////////////////////////////////////////////////////////////
 
tag_into_pose(...)
tag_into_pose( (Pose)pose, (str)tag) -> None : core/pose/util.hh:376
 
total_energy_from_pose(...)
total_energy_from_pose( (Pose)pose) -> float : core/pose/util.hh:389
 
transfer_jumps(...)
transfer_jumps( (Pose)srcpose, (Pose)tgtpose) -> None : core/pose/util.hh:398
 
transfer_phi_psi(...)
transfer_phi_psi( (Pose)srcpose, (Pose)tgtpose) -> None : core/pose/util.hh:395
Variables [hide private]
  PoseCoordPickMode_ALL = rosetta.core.pose._core_pose_.PoseCoor...
  PoseCoordPickMode_BB = rosetta.core.pose._core_pose_.PoseCoord...
  PoseCoordPickMode_BNP = rosetta.core.pose._core_pose_.PoseCoor...
  PoseCoordPickMode_CA = rosetta.core.pose._core_pose_.PoseCoord...
  PoseCoordPickMode_CB = rosetta.core.pose._core_pose_.PoseCoord...
  PoseCoordPickMode_CB_else_CA = rosetta.core.pose._core_pose_.P...
  PoseCoordPickMode_CBorCA = rosetta.core.pose._core_pose_.PoseC...
  PoseCoordPickMode_HVY = rosetta.core.pose._core_pose_.PoseCoor...
  PoseCoordPickMode_HVY_IF_NP = rosetta.core.pose._core_pose_.Po...
  PoseCoordPickMode_NBR = rosetta.core.pose._core_pose_.PoseCoor...
  PoseCoordPickMode_NUL = rosetta.core.pose._core_pose_.PoseCoor...
  PoseCoordPickMode_N_CA_C = rosetta.core.pose._core_pose_.PoseC...
  PoseCoordPickMode_N_CA_C_CB = rosetta.core.pose._core_pose_.Po...
  PoseCoordPickMode_N_C_O = rosetta.core.pose._core_pose_.PoseCo...
  __package__ = None
Function Details [hide private]

QQQ_PosePyObserverTesterFunction()

 
    C++ signature :
        void QQQ_PosePyObserverTesterFunction()

Returns: None :

addVirtualResAsRoot(...)

 

addVirtualResAsRoot( (xyzVector_Real)xyz, (Pose)pose) -> None :
    Adds a virtual residue to the end of the pose at the specified location.
    Roots the structure on this residue.
    

    C++ signature :
        void addVirtualResAsRoot(numeric::xyzVector<double>,core::pose::Pose {lvalue})

addVirtualResAsRoot( (Pose)pose) -> None :
    Adds a VRT res to the end of the pose at the center of mass.
    Reroots the structure on this res.
    

    C++ signature :
        void addVirtualResAsRoot(core::pose::Pose {lvalue})

add_comment(...)

 

add_comment( (Pose)pose, (str)key, (str)val) -> None :
    Adds a key-value pair to the STRING_MAP in the Pose DataCache. If
    there is no STRING_MAP in the DataCache, one is created.
    

    C++ signature :
        void add_comment(core::pose::Pose {lvalue},std::string,std::string)

add_lower_terminus_type_to_pose_residue(...)

 

add_lower_terminus_type_to_pose_residue( (Pose)pose, (int)seqpos) -> None :
    core/pose/util.hh:436

    C++ signature :
        void add_lower_terminus_type_to_pose_residue(core::pose::Pose {lvalue},unsigned long)

add_score_line_string(...)

 

add_score_line_string( (Pose)pose, (str)key, (str)val) -> None :
    Dumps a pdb with comments at end of file
    Sets a PDB-style REMARK entry in the Pose.
    This is different from a comment in its interpretation by the
    silent-file output machinery. A REMARK is written on its own separate line
    in the output silent-file, while a comment is written as part of the Pose
    SCORE: lines.
    

    C++ signature :
        void add_score_line_string(core::pose::Pose {lvalue},std::string,std::string)

add_upper_terminus_type_to_pose_residue(...)

 

add_upper_terminus_type_to_pose_residue( (Pose)pose, (int)seqpos) -> None :
    core/pose/util.hh:443

    C++ signature :
        void add_upper_terminus_type_to_pose_residue(core::pose::Pose {lvalue},unsigned long)

add_variant_type_to_pose_residue(...)

 

add_variant_type_to_pose_residue( (Pose)pose, (VariantType)variant_type, (int)seqpos) -> None :
    Construct a variant of an existing pose residue.
    

    C++ signature :
        void add_variant_type_to_pose_residue(core::pose::Pose {lvalue},core::chemical::VariantType,unsigned long)

add_variant_type_to_residue(...)

 

add_variant_type_to_residue( (Residue)old_rsd, (VariantType)variant_type, (Pose)pose) -> Residue :
    Construct a variant of an existing residue.
    

    C++ signature :
        boost::shared_ptr<core::conformation::Residue> add_variant_type_to_residue(core::conformation::Residue,core::chemical::VariantType,core::pose::Pose)

annotated_to_oneletter_sequence(...)

 

annotated_to_oneletter_sequence( (str)annotated_seq) -> str :
    Returns the oneletter_sequence that corresponds to the given
    annotated sequence.
    

    C++ signature :
        std::string annotated_to_oneletter_sequence(std::string)

append_pose_to_pose(...)

 

append_pose_to_pose( (Pose)pose1, (Pose)pose2 [, (bool)new_chain=True]) -> None :
    Append residues of pose2 to pose1.
    

    C++ signature :
        void append_pose_to_pose(core::pose::Pose {lvalue},core::pose::Pose [,bool=True])

append_subpose_to_pose(...)

 

append_subpose_to_pose( (Pose)pose1, (Pose)pose2, (int)start_res, (int)end_res [, (bool)new_chain=True]) -> None :
    Append specified residues of pose2 to pose1.
    

    C++ signature :
        void append_subpose_to_pose(core::pose::Pose {lvalue},core::pose::Pose,unsigned long,unsigned long [,bool=True])

atom_id_to_named_atom_id(...)

 

atom_id_to_named_atom_id( (AtomID)atom_id, (Pose)pose) -> NamedAtomID :
    core/pose/util.hh:353

    C++ signature :
        core::id::NamedAtomID atom_id_to_named_atom_id(core::id::AtomID,core::pose::Pose)

canonical_atom_count(...)

 

canonical_atom_count( (Pose)pose) -> int :
    count the number of canonical amino acid atoms in the pose
    

    C++ signature :
        unsigned long canonical_atom_count(core::pose::Pose)

canonical_residue_count(...)

 

canonical_residue_count( (Pose)pose) -> int :
    count the number of canonical residues in the pose
    

    C++ signature :
        unsigned long canonical_residue_count(core::pose::Pose)

center_of_mass(...)

 

center_of_mass( (Pose)pose, (int)start, (int)stop) -> xyzVector_Real :
    core/pose/util.hh:669

    C++ signature :
        numeric::xyzVector<double> center_of_mass(core::pose::Pose,int,int)

clearPoseExtraScore(...)

 

clearPoseExtraScore( (Pose)pose, (str)name) -> None :
    core/pose/util.hh:170

    C++ signature :
        void clearPoseExtraScore(core::pose::Pose {lvalue},std::string)

clearPoseExtraScores(...)

 

clearPoseExtraScores( (Pose)pose) -> None :
    core/pose/util.hh:174

    C++ signature :
        void clearPoseExtraScores(core::pose::Pose {lvalue})

compare_atom_coordinates(...)

 

compare_atom_coordinates( (Pose)lhs, (Pose)rhs [, (int)n_dec_places=3]) -> bool :
    this function compares pose atom coordinates for equality; it is not the == operator because it does not compare all pose data.
    Steven Lewis smlewi@gmail.com
    [in] lhs one pose to compare
    [in] rhs one pose to compare
    [in] n_dec_places number of decimal places to compare for the coordinates (remember == doesn't work for float); defaults to 3 which is PDB accuracy
    

    C++ signature :
        bool compare_atom_coordinates(core::pose::Pose,core::pose::Pose [,unsigned long=3])

compare_binary_protein_silent_struct(...)

 

compare_binary_protein_silent_struct( (Pose)lhs, (Pose)rhs) -> bool :
    this function compares poses for equality up to the
    information stored in the binary protein silent struct format.
    

    C++ signature :
        bool compare_binary_protein_silent_struct(core::pose::Pose,core::pose::Pose)

conf2pdb_chain(...)

 

conf2pdb_chain( (Pose)pose) -> object :
    get Conformation chain -> PDBInfo chain mapping
    @remarks Any chains whose PDBInfo chain records are marked entirely as
     PDBInfo::empty_record() will be mapped to that character.  Note that
     Conformation -> PDBInfo is always unique, but the reverse may not be true.
    the mapping if PDBInfo available and chains appear consistent,
     otherwise returns an empty mapping
    

    C++ signature :
        std::map<int, char, std::less<int>, std::allocator<std::pair<int const, char> > > conf2pdb_chain(core::pose::Pose)

convert_from_std_map(...)

 

convert_from_std_map( (object)atom_map, (Pose)pose) -> AtomID_Map_T_core_id_AtomID_T :
    core/pose/util.hh:687

    C++ signature :
        core::id::AtomID_Map<core::id::AtomID> convert_from_std_map(std::map<core::id::AtomID, core::id::AtomID, std::less<core::id::AtomID>, std::allocator<std::pair<core::id::AtomID const, core::id::AtomID> > >,core::pose::Pose)

correctly_add_cutpoint_variants(...)

 

correctly_add_cutpoint_variants( (Pose)pose, (int)cutpoint_res [, (bool)check_fold_tree=True]) -> None :
    core/pose/util.hh:696

    C++ signature :
        void correctly_add_cutpoint_variants(core::pose::Pose {lvalue},unsigned long [,bool=True])

correctly_add_cutpoint_variants( (Pose)pose) -> None :
    Add cutpoint variants to all residues annotated as cutpoints in the pose.
    

    C++ signature :
        void correctly_add_cutpoint_variants(core::pose::Pose {lvalue})

create_subpose(...)

 

create_subpose( (Pose)src, (vector1_Size)positions, (FoldTree)f, (Pose)pose) -> None :
    Create a subpose of the src pose.  PDBInfo is set as NULL.
    

    C++ signature :
        void create_subpose(core::pose::Pose,utility::vector1<unsigned long, std::allocator<unsigned long> >,core::kinematics::FoldTree,core::pose::Pose {lvalue})

delete_comment(...)

 

delete_comment( (Pose)pose, (str)key) -> None :
    Deletes the entry in the STRING_MAP associated with the
    given key.
    

    C++ signature :
        void delete_comment(core::pose::Pose {lvalue},std::string)

dump_comment_pdb(...)

 

dump_comment_pdb( (str)file_name, (Pose)pose) -> None :
    dumps pose+ comments to pdb file
    

    C++ signature :
        void dump_comment_pdb(std::string,core::pose::Pose)

energy_from_pose(...)

 

energy_from_pose( (Pose)pose, (str)sc_type) -> float :
    core/pose/util.hh:387

    C++ signature :
        double energy_from_pose(core::pose::Pose,std::string)

energy_from_pose( (Pose)pose, (ScoreType)sc_type) -> float :
    core/pose/util.hh:383

    C++ signature :
        double energy_from_pose(core::pose::Pose,core::scoring::ScoreType)

extract_tag_from_pose(...)

 

extract_tag_from_pose( (Pose)pose) -> str :
    Returns a string giving the pose's tag if there is such a thing or "UnknownTag" otherwise.
    

    C++ signature :
        std::string extract_tag_from_pose(core::pose::Pose {lvalue})

getPoseExtraScore(...)

 

getPoseExtraScore( (Pose)pose, (str)name, (str)value) -> bool :
    core/pose/util.hh:181

    C++ signature :
        bool getPoseExtraScore(core::pose::Pose,std::string,std::string {lvalue})

getPoseExtraScore( (Pose)pose, (str)name) -> float :
    core/pose/util.hh:154

    C++ signature :
        double getPoseExtraScore(core::pose::Pose,std::string)

getPoseExtraScore( (Pose)pose, (str)name, (float)value) -> bool :
    getters/setters for things in the Pose DataCache
    

    C++ signature :
        bool getPoseExtraScore(core::pose::Pose,std::string,double {lvalue})

get_all_comments(...)

 

get_all_comments( (Pose)pose) -> map_string_string :
    Gets a map< string, string > representing comments about the Pose in
    the form of key-value pairs.
    

    C++ signature :
        std::map<std::string, std::string, std::less<std::string>, std::allocator<std::pair<std::string const, std::string> > > get_all_comments(core::pose::Pose)

get_all_score_line_strings(...)

 

get_all_score_line_strings( (Pose)pose) -> map_string_string :
    Gets a map< string, string > representing score_line_strings about the Pose in
    the form of key-value pairs.
    

    C++ signature :
        std::map<std::string, std::string, std::less<std::string>, std::allocator<std::pair<std::string const, std::string> > > get_all_score_line_strings(core::pose::Pose)

get_center_of_mass(...)

 

get_center_of_mass( (Pose)pose) -> xyzVector_Real :
    Get center of mass of a pose.
    

    C++ signature :
        numeric::xyzVector<double> get_center_of_mass(core::pose::Pose)

get_chain_from_chain_id(...)

 

get_chain_from_chain_id( (int)chain_id, (Pose)pose) -> str :
    core/pose/util.hh:508

    C++ signature :
        char get_chain_from_chain_id(unsigned long,core::pose::Pose)

get_chain_from_jump_id(...)

 

get_chain_from_jump_id( (int)jump_id, (Pose)pose) -> str :
    core/pose/util.hh:526

    C++ signature :
        char get_chain_from_jump_id(unsigned long,core::pose::Pose)

get_chain_id_from_chain(...)

 

get_chain_id_from_chain( (str)chain, (Pose)pose) -> int :
    core/pose/util.hh:496

    C++ signature :
        unsigned long get_chain_id_from_chain(char,core::pose::Pose)

get_chain_id_from_chain( (str)chain, (Pose)pose) -> int :
    core/pose/util.hh:493

    C++ signature :
        unsigned long get_chain_id_from_chain(std::string,core::pose::Pose)

get_chain_id_from_jump_id(...)

 

get_chain_id_from_jump_id( (int)jump_id, (Pose)pose) -> int :
    core/pose/util.hh:523

    C++ signature :
        unsigned long get_chain_id_from_jump_id(unsigned long,core::pose::Pose)

get_chain_ids_from_chain(...)

 

get_chain_ids_from_chain( (str)chain, (Pose)pose) -> vector1_Size :
    core/pose/util.hh:502

    C++ signature :
        utility::vector1<unsigned long, std::allocator<unsigned long> > get_chain_ids_from_chain(char,core::pose::Pose)

get_chain_ids_from_chain( (str)chain, (Pose)pose) -> vector1_Size :
    core/pose/util.hh:499

    C++ signature :
        utility::vector1<unsigned long, std::allocator<unsigned long> > get_chain_ids_from_chain(std::string,core::pose::Pose)

get_chain_ids_from_chains(...)

 

get_chain_ids_from_chains( (vector1_string)chains, (Pose)pose) -> vector1_Size :
    core/pose/util.hh:505

    C++ signature :
        utility::vector1<unsigned long, std::allocator<unsigned long> > get_chain_ids_from_chains(utility::vector1<std::string, std::allocator<std::string> >,core::pose::Pose)

get_chain_residues(...)

 

get_chain_residues( (Pose)pose, (int)chain_id) -> object :
    core/pose/util.hh:529

    C++ signature :
        utility::vector1<boost::shared_ptr<core::conformation::Residue const>, std::allocator<boost::shared_ptr<core::conformation::Residue const> > > get_chain_residues(core::pose::Pose,unsigned long)

get_comment(...)

 

get_comment( (Pose)pose, (str)key, (str)val) -> bool :
    Attempts to access the entry in the STRING_MAP associated with the
    given key. If an entry for the key exists, the value associated with the key
    is put into val, and this function returns true. Otherwise, this function
    returns false and val left unmodified.
    

    C++ signature :
        bool get_comment(core::pose::Pose,std::string,std::string {lvalue})

get_hash_excluding_chain(...)

 

get_hash_excluding_chain( (str)chain, (Pose)pose) -> int :
    core/pose/util.hh:572

    C++ signature :
        unsigned long get_hash_excluding_chain(char,core::pose::Pose)

get_hash_from_chain(...)

 

get_hash_from_chain( (str)chain, (Pose)pose) -> int :
    core/pose/util.hh:569

    C++ signature :
        unsigned long get_hash_from_chain(char,core::pose::Pose)

get_jump_id_from_chain(...)

 

get_jump_id_from_chain( (str)chain, (Pose)pose) -> int :
    core/pose/util.hh:514

    C++ signature :
        unsigned long get_jump_id_from_chain(char,core::pose::Pose)

get_jump_id_from_chain( (str)chain, (Pose)pose) -> int :
    core/pose/util.hh:511

    C++ signature :
        unsigned long get_jump_id_from_chain(std::string,core::pose::Pose)

get_jump_id_from_chain_id(...)

 

get_jump_id_from_chain_id( (int)chain_id, (Pose)pose) -> int :
    core/pose/util.hh:488

    C++ signature :
        unsigned long get_jump_id_from_chain_id(unsigned long,core::pose::Pose)

get_jump_ids_from_chain(...)

 

get_jump_ids_from_chain( (str)chain, (Pose)pose) -> vector1_Size :
    core/pose/util.hh:520

    C++ signature :
        utility::vector1<unsigned long, std::allocator<unsigned long> > get_jump_ids_from_chain(std::string,core::pose::Pose)

get_jump_ids_from_chain( (str)chain, (Pose)pose) -> vector1_Size :
    core/pose/util.hh:517

    C++ signature :
        utility::vector1<unsigned long, std::allocator<unsigned long> > get_jump_ids_from_chain(char,core::pose::Pose)

get_jump_ids_from_chain_ids(...)

 

get_jump_ids_from_chain_ids( (set_Size)chain_ids, (Pose)pose) -> set_Size :
    core/pose/util.hh:485

    C++ signature :
        std::set<unsigned long, std::less<unsigned long>, std::allocator<unsigned long> > get_jump_ids_from_chain_ids(std::set<unsigned long, std::less<unsigned long>, std::allocator<unsigned long> >,core::pose::Pose)

get_resnum(...)

 

get_resnum( (Tag)tag_ptr, (Pose)pose [, (str)prefix='']) -> int :
    core/pose/selection.hh:37

    C++ signature :
        unsigned long get_resnum(boost::shared_ptr<utility::tag::Tag const>,core::pose::Pose [,std::string=''])

get_resnum_list(...)

 

get_resnum_list( (str)str, (Pose)pose) -> set_Size :
    returns a resnum list directly from a string
    

    C++ signature :
        std::set<unsigned long, std::less<unsigned long>, std::allocator<unsigned long> > get_resnum_list(std::string,core::pose::Pose)

get_resnum_list( (Tag)tag_ptr, (str)tag, (Pose)pose) -> vector1_Size :
    Extracts a list of residue numbers from a tag.
    The tag should contain a comma-separated list of numbers, in either
      pdb or rosetta format (@see parse_resnum for details)
    

    C++ signature :
        utility::vector1<unsigned long, std::allocator<unsigned long> > get_resnum_list(boost::shared_ptr<utility::tag::Tag const>,std::string,core::pose::Pose)

get_resnum_list_ordered(...)

 

get_resnum_list_ordered( (str)str, (Pose)pose) -> vector1_Size :
    returns a resnum list directly from a string, preserving order
    

    C++ signature :
        utility::vector1<unsigned long, std::allocator<unsigned long> > get_resnum_list_ordered(std::string,core::pose::Pose)

get_score_line_string(...)

 

get_score_line_string( (Pose)pose, (str)key, (str)val) -> bool :
    core/pose/util.hh:250

    C++ signature :
        bool get_score_line_string(core::pose::Pose,std::string,std::string {lvalue})

get_sequence_len(...)

 

get_sequence_len( (str)sequence_in) -> int :
    Get the real length of a annotated sequence
    

    C++ signature :
        unsigned long get_sequence_len(std::string)

get_sha1_hash_excluding_chain(...)

 

get_sha1_hash_excluding_chain( (str)chain, (Pose)pose) -> str :
    core/pose/util.hh:574

    C++ signature :
        std::string get_sha1_hash_excluding_chain(char,core::pose::Pose)

hasPoseExtraScore(...)

 

hasPoseExtraScore( (Pose)pose, (str)name) -> bool :
    core/pose/util.hh:159

    C++ signature :
        bool hasPoseExtraScore(core::pose::Pose,std::string)

has_chain(...)

 

has_chain( (int)chain_id, (Pose)pose) -> bool :
    core/pose/util.hh:482

    C++ signature :
        bool has_chain(unsigned long,core::pose::Pose)

has_chain( (str)chain, (Pose)pose) -> bool :
    core/pose/util.hh:479

    C++ signature :
        bool has_chain(char,core::pose::Pose)

has_chain( (str)chain, (Pose)pose) -> bool :
    core/pose/util.hh:476

    C++ signature :
        bool has_chain(std::string,core::pose::Pose)

initialize_atomid_map(...)

 

initialize_atomid_map( (AtomID_Map_T_core_id_AtomID_T)atom_map, (Conformation)conformation, (AtomID)value) -> None :
    Initialize an AtomID_Map for a given Conformation using a specified fill value
    

    C++ signature :
        void initialize_atomid_map(core::id::AtomID_Map<core::id::AtomID> {lvalue},core::conformation::Conformation,core::id::AtomID)

initialize_atomid_map( (AtomID_Map_T_core_id_AtomID_T)atom_map, (Conformation)conformation) -> None :
    Initialize an AtomID_Map for a given Conformation using the AtomID_Map's current default fill values
    

    C++ signature :
        void initialize_atomid_map(core::id::AtomID_Map<core::id::AtomID> {lvalue},core::conformation::Conformation)

initialize_atomid_map( (AtomID_Map_T_core_id_AtomID_T)atom_map, (Pose)pose, (AtomID)value) -> None :
    Initialize an AtomID_Map for a given Pose using a specified fill value
    

    C++ signature :
        void initialize_atomid_map(core::id::AtomID_Map<core::id::AtomID> {lvalue},core::pose::Pose,core::id::AtomID)

initialize_atomid_map( (AtomID_Map_T_core_id_AtomID_T)atom_map, (Pose)pose) -> None :
    Initialize an AtomID_Map for a given Pose using the AtomID_Map's current default fill values
    

    C++ signature :
        void initialize_atomid_map(core::id::AtomID_Map<core::id::AtomID> {lvalue},core::pose::Pose)

initialize_atomid_map_AtomID(...)

 

initialize_atomid_map_AtomID( (AtomID_Map_T_core_id_AtomID_T)atom_map, (Conformation)conformation, (AtomID)value) -> None :
    core/pose/util.tmpl.hh:379

    C++ signature :
        void initialize_atomid_map_AtomID(core::id::AtomID_Map<core::id::AtomID> {lvalue},core::conformation::Conformation,core::id::AtomID)

initialize_atomid_map_AtomID( (AtomID_Map_T_core_id_AtomID_T)atom_map, (Conformation)conformation) -> None :
    core/pose/util.tmpl.hh:376

    C++ signature :
        void initialize_atomid_map_AtomID(core::id::AtomID_Map<core::id::AtomID> {lvalue},core::conformation::Conformation)

initialize_atomid_map_AtomID( (AtomID_Map_T_core_id_AtomID_T)atom_map, (Pose)pose, (AtomID)value) -> None :
    core/pose/util.tmpl.hh:372

    C++ signature :
        void initialize_atomid_map_AtomID(core::id::AtomID_Map<core::id::AtomID> {lvalue},core::pose::Pose,core::id::AtomID)

initialize_atomid_map_AtomID( (AtomID_Map_T_core_id_AtomID_T)atom_map, (Pose)pose) -> None :
    core/pose/util.tmpl.hh:369

    C++ signature :
        void initialize_atomid_map_AtomID(core::id::AtomID_Map<core::id::AtomID> {lvalue},core::pose::Pose)

initialize_disulfide_bonds(...)

 

initialize_disulfide_bonds( (Pose)pose) -> None :
    detect and fix disulfide bonds
    

    C++ signature :
        void initialize_disulfide_bonds(core::pose::Pose {lvalue})

is_ideal_pose(...)

 

is_ideal_pose( (Pose)pose) -> bool :
    Returns true if the  <pose>  geometry is ideal
    [in] pose The Pose to check.
    true if all pose positions have ideal bond lengths and angles
     up to some very small epsilon
    

    C++ signature :
        bool is_ideal_pose(core::pose::Pose)

is_ideal_position(...)

 

is_ideal_position( (int)seqpos, (Pose)pose) -> bool :
    Returns true if the  <pose> geometry is ideal in position  <seqpos>
    [in] pose The Pose to check.
    true if position seqpos has ideal bond lengths and angles
     up to some very small epsilon
    

    C++ signature :
        bool is_ideal_position(unsigned long,core::pose::Pose)

is_lower_terminus(...)

 

is_lower_terminus( (Pose)pose, (int)resid) -> bool :
    checks to see if this is a lower chain ending more intelligently than just checking residue variants
    

    C++ signature :
        bool is_lower_terminus(core::pose::Pose,unsigned long)

is_position_conserved_residue(...)

 

is_position_conserved_residue( (Pose)pose, (int)residue) -> bool :
    Returns true if <residue> is positionally conserved, false otherwise
    

    C++ signature :
        bool is_position_conserved_residue(core::pose::Pose,unsigned long)

is_upper_terminus(...)

 

is_upper_terminus( (Pose)pose, (int)resid) -> bool :
    checks to see if this is a lower chain ending more intelligently than just checking residue variants
    

    C++ signature :
        bool is_upper_terminus(core::pose::Pose,unsigned long)

jumps_from_pose(...)

 

jumps_from_pose( (Pose)pose, (set_int)jumps) -> None :
    Retrieves jump information from <pose>, storing the result in <jumps>.
    Jumps are keyed by their jump id.
    

    C++ signature :
        void jumps_from_pose(core::pose::Pose,std::set<int, std::less<int>, std::allocator<int> >*)

make_atom_map(...)

 

make_atom_map( (Pose)p, (PoseCoordPickMode)m) -> AtomID_Map_T_double_T :
    core/pose/xyzStripeHashPose.hh:29

    C++ signature :
        core::id::AtomID_Map<double> make_atom_map(core::pose::Pose,core::pose::PoseCoordPickMode)

make_pose_from_saccharide_sequence(...)

 

make_pose_from_saccharide_sequence( (Pose)pose, (str)sequence [, (str)type_set_name='fa_standard' [, (bool)auto_termini=True]]) -> None :
    Create a Pose from an annotated, linear, IUPAC polysaccharide sequence <sequence> with residue type set name
    <type_set_name> and store it in <pose>.
    

    C++ signature :
        void make_pose_from_saccharide_sequence(core::pose::Pose {lvalue},std::string [,std::string='fa_standard' [,bool=True]])

make_pose_from_saccharide_sequence( (Pose)pose, (str)sequence, (ResidueTypeSet)residue_set [, (bool)auto_termini=True]) -> None :
    Create a Pose from an annotated, linear, IUPAC polysaccharide sequence <sequence> with ResidueTypeSet
    <residue_set> and store it in <pose>.
    

    C++ signature :
        void make_pose_from_saccharide_sequence(core::pose::Pose {lvalue},std::string,core::chemical::ResidueTypeSet [,bool=True])

make_pose_from_sequence(...)

 

make_pose_from_sequence( (Pose)pose, (str)sequence, (str)type_set_name [, (bool)auto_termini=True]) -> None :
    Creates a Pose from the annotated protein sequence  <sequence>
    with the desired  <type_set_name>  and stores it in  <pose>
    : any existing data in  <pose>  is cleared, auto_termini
    mark position 1, last_residue with lower, upper termini; default true
    

    C++ signature :
        void make_pose_from_sequence(core::pose::Pose {lvalue},std::string,std::string [,bool=True])

make_pose_from_sequence( (Pose)pose, (str)sequence, (ResidueTypeSet)residue_set [, (bool)auto_termini=True]) -> None :
    core/pose/annotated_sequence.hh:85

    C++ signature :
        void make_pose_from_sequence(core::pose::Pose {lvalue},std::string,core::chemical::ResidueTypeSet [,bool=True])

make_pose_from_sequence( (Pose)pose, (object)requested_types [, (bool)auto_termini=True]) -> None :
    Creates a Pose from the annotated protein sequence  <sequence>
    with ResidueTypeSet  <residue_set>  and stores it in  <pose>
    : any existing data in  <pose>  is cleared, auto_termini
    mark position 1, last_residue with lower, upper termini; default true
    example(s):
        make_pose_from_sequence(pose,"THANKSEVAN","fa_standard")
    See also:
        Pose
        PDBInfo
        pose_from_pdb
        pose_from_rcsb
        pose_from_sequence
    

    C++ signature :
        void make_pose_from_sequence(core::pose::Pose {lvalue},utility::vector1<boost::shared_ptr<core::chemical::ResidueType const>, std::allocator<boost::shared_ptr<core::chemical::ResidueType const> > > [,bool=True])

mass(...)

 

mass( (int)begin, (int)end, (Pose)pose) -> float :
    core/pose/util.hh:566

    C++ signature :
        double mass(unsigned long,unsigned long,core::pose::Pose)

named_atom_id_to_atom_id(...)

 

named_atom_id_to_atom_id( (NamedAtomID)named_atom_id, (Pose)pose [, (bool)raise_exception=True]) -> AtomID :
    core/pose/util.hh:360

    C++ signature :
        core::id::AtomID named_atom_id_to_atom_id(core::id::NamedAtomID,core::pose::Pose [,bool=True])

named_stub_id_to_stub_id(...)

 

named_stub_id_to_stub_id( (NamedStubID)named_stub_id, (Pose)pose) -> StubID :
    core/pose/util.hh:372

    C++ signature :
        core::id::StubID named_stub_id_to_stub_id(core::id::NamedStubID,core::pose::Pose)

noncanonical_atom_count(...)

 

noncanonical_atom_count( (Pose)pose) -> int :
    count the number of non-canonical amino acids in thepose
    

    C++ signature :
        unsigned long noncanonical_atom_count(core::pose::Pose)

noncanonical_chi_count(...)

 

noncanonical_chi_count( (Pose)pose) -> int :
    count the number of non-canonical chi angles in the pose
    

    C++ signature :
        unsigned long noncanonical_chi_count(core::pose::Pose)

noncanonical_residue_count(...)

 

noncanonical_residue_count( (Pose)pose) -> int :
    count the number of non-canonical residues in the pose
    

    C++ signature :
        unsigned long noncanonical_residue_count(core::pose::Pose)

nres_protein(...)

 

nres_protein( (Pose)pose) -> int :
    Number of protein residues in the pose
    No virtuals, membrane residues or embedding residues counted
    

    C++ signature :
        unsigned long nres_protein(core::pose::Pose)

num_atoms(...)

 

num_atoms( (int)begin, (int)end, (Pose)pose) -> int :
    core/pose/util.hh:541

    C++ signature :
        unsigned long num_atoms(unsigned long,unsigned long,core::pose::Pose)

num_chi_angles(...)

 

num_chi_angles( (int)begin, (int)end, (Pose)pose) -> int :
    core/pose/util.hh:559

    C++ signature :
        unsigned long num_chi_angles(unsigned long,unsigned long,core::pose::Pose)

num_hbond_acceptors(...)

 

num_hbond_acceptors( (int)begin, (int)end, (Pose)pose) -> int :
    core/pose/util.hh:547

    C++ signature :
        unsigned long num_hbond_acceptors(unsigned long,unsigned long,core::pose::Pose)

num_hbond_donors(...)

 

num_hbond_donors( (int)begin, (int)end, (Pose)pose) -> int :
    core/pose/util.hh:553

    C++ signature :
        unsigned long num_hbond_donors(unsigned long,unsigned long,core::pose::Pose)

num_heavy_atoms(...)

 

num_heavy_atoms( (int)begin, (int)end, (Pose)pose) -> int :
    core/pose/util.hh:535

    C++ signature :
        unsigned long num_heavy_atoms(unsigned long,unsigned long,core::pose::Pose)

parse_resnum(...)

 

parse_resnum( (str)resnum, (Pose)pose) -> int :
    Extracts a residue number from a string.
    @detail Recognizes two forms of numbering:
      - Rosetta residue numbers (numbered sequentially from 1 to the last residue
        in the pose). These have the form [0-9]+
      - PDB numbers. These have the form [0-9]+[A-Z], where the trailing letter
        is the chain ID.
    the rosetta residue number for the string, or 0 upon an error
    

    C++ signature :
        unsigned long parse_resnum(std::string,core::pose::Pose)

parse_sequence(...)

 

parse_sequence( (str)sequence_in, (vector1_string)fullname_list, (std_vector_Size)oneletter_to_fullname_index, (str)one_letter_sequence) -> None :
    Parse the input annotated sequence
    

    C++ signature :
        void parse_sequence(std::string,utility::vector1<std::string, std::allocator<std::string> > {lvalue},std::vector<unsigned long, std::allocator<unsigned long> > {lvalue},std::string {lvalue})

partition_pose_by_jump(...)

 

partition_pose_by_jump( (Pose)src, (int)jump_number, (Pose)partner1, (Pose)partner2) -> None :
    core/pose/util.hh:125

    C++ signature :
        void partition_pose_by_jump(core::pose::Pose,int,core::pose::Pose {lvalue},core::pose::Pose {lvalue})

pdbslice(...)

 

pdbslice( (Pose)pose, (vector1_Size)slice_res) -> None :
    Create a subpose of the src pose -- figures out a reasonable fold tree.
    

    C++ signature :
        void pdbslice(core::pose::Pose {lvalue},utility::vector1<unsigned long, std::allocator<unsigned long> >)

pdbslice( (Pose)new_pose, (Pose)pose, (vector1_Size)slice_res) -> None :
    Create a subpose of the src pose -- figures out a reasonable fold tree.
    

    C++ signature :
        void pdbslice(core::pose::Pose {lvalue},core::pose::Pose,utility::vector1<unsigned long, std::allocator<unsigned long> >)

pose_from_saccharide_sequence(...)

 

pose_from_saccharide_sequence( (str)sequence [, (str)type_set_name='fa_standard' [, (bool)auto_termini=True]]) -> Pose :
    Return a Pose from an annotated, linear, IUPAC polysaccharide sequence <sequence> with residue type set name
    <type_set_name>.
    

    C++ signature :
        boost::shared_ptr<core::pose::Pose> pose_from_saccharide_sequence(std::string [,std::string='fa_standard' [,bool=True]])

pose_max_nbr_radius(...)

 

pose_max_nbr_radius( (Pose)pose) -> float :
    returns a Distance
    

    C++ signature :
        double pose_max_nbr_radius(core::pose::Pose)

pose_residue_is_terminal(...)

 

pose_residue_is_terminal( (Pose)pose, (int)resid) -> bool :
    returns true if the given residue in the pose is a chain ending or has upper/lower terminal variants
    

    C++ signature :
        bool pose_residue_is_terminal(core::pose::Pose,unsigned long)

read_comment_pdb(...)

 

read_comment_pdb( (str)file_name, (Pose)pose) -> None :
    Reads the comments from the pdb file and adds it into comments
    

    C++ signature :
        void read_comment_pdb(std::string,core::pose::Pose {lvalue})

read_psipred_ss2_file(...)

 

read_psipred_ss2_file( (Pose)pose) -> vector1_char :
    core/pose/util.hh:141

    C++ signature :
        utility::vector1<char, std::allocator<char> > read_psipred_ss2_file(core::pose::Pose)

remove_ligand_canonical_residues(...)

 

remove_ligand_canonical_residues( (Pose)pose) -> None :
    this function removes all residues with both UPPER and LOWER terminus types.  This is intended for removing ligands that are canonical residues.
    

    C++ signature :
        void remove_ligand_canonical_residues(core::pose::Pose {lvalue})

remove_lower_terminus_type_from_pose_residue(...)

 

remove_lower_terminus_type_from_pose_residue( (Pose)pose, (int)seqpos) -> None :
    core/pose/util.hh:450

    C++ signature :
        void remove_lower_terminus_type_from_pose_residue(core::pose::Pose {lvalue},unsigned long)

remove_nonprotein_residues(...)

 

remove_nonprotein_residues( (Pose)pose) -> None :
    this function removes all residues from the pose which are not protein residues.  This removal includes, but is not limited to, metals, DNA, RNA, and ligands.  It will NOT remove ligands which are canonical residues (for example, if a protein binds an alanine monomer, the monomer will be untouched).
    

    C++ signature :
        void remove_nonprotein_residues(core::pose::Pose {lvalue})

remove_upper_terminus_type_from_pose_residue(...)

 

remove_upper_terminus_type_from_pose_residue( (Pose)pose, (int)seqpos) -> None :
    core/pose/util.hh:457

    C++ signature :
        void remove_upper_terminus_type_from_pose_residue(core::pose::Pose {lvalue},unsigned long)

remove_variant_type_from_pose_residue(...)

 

remove_variant_type_from_pose_residue( (Pose)pose, (VariantType)variant_type, (int)seqpos) -> None :
    Construct a non-variant of an existing pose residue.
    

    C++ signature :
        void remove_variant_type_from_pose_residue(core::pose::Pose {lvalue},core::chemical::VariantType,unsigned long)

remove_variant_type_from_residue(...)

 

remove_variant_type_from_residue( (Residue)old_rsd, (VariantType)variant_type, (Pose)pose) -> Residue :
    Remove variant from an existing residue.
    

    C++ signature :
        boost::shared_ptr<core::conformation::Residue> remove_variant_type_from_residue(core::conformation::Residue,core::chemical::VariantType,core::pose::Pose)

remove_virtual_residues(...)

 

remove_virtual_residues( (Pose)pose) -> None :
    Removes all virtual residues from <pose>
    

    C++ signature :
        void remove_virtual_residues(core::pose::Pose*)

renumber_pdbinfo_based_on_conf_chains(...)

 

renumber_pdbinfo_based_on_conf_chains( (Pose)pose [, (bool)fix_chains=True [, (bool)start_from_existing_numbering=True [, (bool)keep_insertion_codes=False [, (bool)rotate_chain_ids=False]]]]) -> bool :
    renumber PDBInfo based on Conformation chains; each chain starts from 1
    [in,out] pose The Pose to modify.
    [in] fix_chains If true, the procedure will attempt to fix any empty record
     characters it finds in the PDBInfo. (default true)
    [in] start_from_existing_numbering If true, will attempt to start each
     chain from the existing numbering in the PDBInfo.  E.g. if the first residue
     of chain 2 in the Conformation is 27, then the renumbering of the chain in
     PDBInfo will start from 27. (default true)
    [in] keep_insertion_codes If true, will maintain insertion codes and
     will not increment the pdb residue numbering for those residues.  This means
     new numbering with insertion codes will only reflect properly if the
     old numbering included the base numbering of the insertion code residues,
     i.e. 100 100A 100B and not just 100A 100B (with 100 never appearing).
     (default false)
    [in] rotate_chain_ids If true, allows support for more than 26 pdb chains
     by rotating [A,Z] continuously.  WARNING: This will break the assumption
     made by the PDBPoseMap that each pdb chain id is unique, so make sure you
     are not using the PDBPoseMap feature downstream in your code path without
     corrections! (default false)
    @remarks If fixing chains and there is only one chain and the PDBInfo exists
     but all records are marked as empty, will renumber and set the PDBInfo chain
     to 'A'.
    true if renumbering successful, false otherwise
    

    C++ signature :
        bool renumber_pdbinfo_based_on_conf_chains(core::pose::Pose {lvalue} [,bool=True [,bool=True [,bool=False [,bool=False]]]])

replace_pose_residue_copying_existing_coordinates(...)

 

replace_pose_residue_copying_existing_coordinates( (Pose)pose, (int)seqpos, (ResidueType)new_rsd_type) -> None :
    core/pose/util.hh:405

    C++ signature :
        void replace_pose_residue_copying_existing_coordinates(core::pose::Pose {lvalue},unsigned long,core::chemical::ResidueType)

residue_center_of_mass(...)

 

residue_center_of_mass( (Pose)pose, (int)start, (int)stop) -> int :
    core/pose/util.hh:676

    C++ signature :
        int residue_center_of_mass(core::pose::Pose,int,int)

residue_types_from_saccharide_sequence(...)

 

residue_types_from_saccharide_sequence( (str)sequence, (ResidueTypeSet)residue_set) -> object :
    Return a list of carbohydrate ResidueTypes corresponding to an annotated, linear, IUPAC polysaccharide
    sequence.
    

    C++ signature :
        utility::vector1<boost::shared_ptr<core::chemical::ResidueType const>, std::allocator<boost::shared_ptr<core::chemical::ResidueType const> > > residue_types_from_saccharide_sequence(std::string,core::chemical::ResidueTypeSet)

residue_types_from_sequence(...)

 

residue_types_from_sequence( (str)sequence_in, (ResidueTypeSet)residue_set [, (bool)auto_termini=True]) -> object :
    return a list of ResidueTypes corresponding to an annotated protein sequence
    [in] sequence_in an annotated sequence
    [in] residue_set the desired residue set
    [in] auto_termini mark position 1, last_residue with lower, upper termini; default true
    

    C++ signature :
        utility::vector1<boost::shared_ptr<core::chemical::ResidueType const>, std::allocator<boost::shared_ptr<core::chemical::ResidueType const> > > residue_types_from_sequence(std::string,core::chemical::ResidueTypeSet [,bool=True])

return_nearest_residue(...)

 

return_nearest_residue( (Pose)pose, (int)begin, (int)end, (xyzVector_Real)center) -> int :
    core/pose/util.hh:684

    C++ signature :
        int return_nearest_residue(core::pose::Pose,int,int,numeric::xyzVector<double>)

sequence_map_from_pdbinfo(...)

 

sequence_map_from_pdbinfo( (Pose)first, (Pose)second) -> SequenceMapping :
    Create a sequence map of first pose onto the second, matching the PDBInfo
       If the PDBInfo of either Pose is missing or invalid, do a simple sequence alignment matching.
    

    C++ signature :
        core::id::SequenceMapping sequence_map_from_pdbinfo(core::pose::Pose,core::pose::Pose)

setPoseExtraScore(...)

 

setPoseExtraScore( (Pose)pose, (str)name, (str)value) -> None :
    core/pose/util.hh:187

    C++ signature :
        void setPoseExtraScore(core::pose::Pose {lvalue},std::string,std::string)

setPoseExtraScore( (Pose)pose, (str)name, (float)value) -> None :
    core/pose/util.hh:165

    C++ signature :
        void setPoseExtraScore(core::pose::Pose {lvalue},std::string,double)

set_ss_from_phipsi(...)

 

set_ss_from_phipsi( (Pose)pose) -> None :
    Analyzes  <pose>  residue phi/psi sets and guesses the secondary
    structure, ideally dssp should be used for that
    

    C++ signature :
        void set_ss_from_phipsi(core::pose::Pose {lvalue})

setup_dof_mask_from_move_map(...)

 

setup_dof_mask_from_move_map( (MoveMap)mm, (Pose)pose, (object)dof_mask) -> None :
    convert from allow-bb/allow-chi MoveMap to simple DOF_ID boolean mask needed by the minimizer
    

    C++ signature :
        void setup_dof_mask_from_move_map(core::kinematics::MoveMap,core::pose::Pose,core::id::DOF_ID_Map<bool> {lvalue})

setup_dof_to_torsion_map(...)

 

setup_dof_to_torsion_map( (Pose)pose, (object)dof_map) -> None :
    set up a map to look up TORSION_ID by DOF_ID (Map[DOF_ID] = TORISION_ID)
    

    C++ signature :
        void setup_dof_to_torsion_map(core::pose::Pose,core::id::DOF_ID_Map<core::id::TorsionID> {lvalue})

sort_pose_by_score(...)

 

sort_pose_by_score( (Pose)pose1, (Pose)pose2) -> bool :
    core/pose/util.hh:379

    C++ signature :
        bool sort_pose_by_score(boost::shared_ptr<core::pose::Pose>,boost::shared_ptr<core::pose::Pose>)

stub_id_to_named_stub_id(...)

 

stub_id_to_named_stub_id( (StubID)stub_id, (Pose)pose) -> NamedStubID :
    core/pose/util.hh:366

    C++ signature :
        core::id::NamedStubID stub_id_to_named_stub_id(core::id::StubID,core::pose::Pose)

swap_transform(...)

 

swap_transform( (int)jump_num, (RT)xform, (Pose)pose) -> None :
    Updates the rigid-body transform of the specified jump in <pose>
    

    C++ signature :
        void swap_transform(unsigned long,core::kinematics::RT,core::pose::Pose*)

tag_from_pose(...)

 

tag_from_pose( (Pose)pose) -> str :
    ////////////////////////////////////////////////////////////////
    

    C++ signature :
        std::string tag_from_pose(core::pose::Pose)

tag_into_pose(...)

 

tag_into_pose( (Pose)pose, (str)tag) -> None :
    core/pose/util.hh:376

    C++ signature :
        void tag_into_pose(core::pose::Pose {lvalue},std::string)

total_energy_from_pose(...)

 

total_energy_from_pose( (Pose)pose) -> float :
    core/pose/util.hh:389

    C++ signature :
        double total_energy_from_pose(core::pose::Pose)

transfer_jumps(...)

 

transfer_jumps( (Pose)srcpose, (Pose)tgtpose) -> None :
    core/pose/util.hh:398

    C++ signature :
        void transfer_jumps(core::pose::Pose,core::pose::Pose {lvalue})

transfer_phi_psi(...)

 

transfer_phi_psi( (Pose)srcpose, (Pose)tgtpose) -> None :
    core/pose/util.hh:395

    C++ signature :
        void transfer_phi_psi(core::pose::Pose,core::pose::Pose {lvalue})

transfer_phi_psi( (Pose)srcpose, (Pose)tgtpose, (int)ir, (int)jr) -> None :
    core/pose/util.hh:392

    C++ signature :
        void transfer_phi_psi(core::pose::Pose,core::pose::Pose {lvalue},unsigned long,unsigned long)


Variables Details [hide private]

PoseCoordPickMode_ALL

Value:
rosetta.core.pose._core_pose_.PoseCoordPickMode.PoseCoordPickMode_ALL

PoseCoordPickMode_BB

Value:
rosetta.core.pose._core_pose_.PoseCoordPickMode.PoseCoordPickMode_BB

PoseCoordPickMode_BNP

Value:
rosetta.core.pose._core_pose_.PoseCoordPickMode.PoseCoordPickMode_BNP

PoseCoordPickMode_CA

Value:
rosetta.core.pose._core_pose_.PoseCoordPickMode.PoseCoordPickMode_CA

PoseCoordPickMode_CB

Value:
rosetta.core.pose._core_pose_.PoseCoordPickMode.PoseCoordPickMode_CB

PoseCoordPickMode_CB_else_CA

Value:
rosetta.core.pose._core_pose_.PoseCoordPickMode.PoseCoordPickMode_CB_e\
lse_CA

PoseCoordPickMode_CBorCA

Value:
rosetta.core.pose._core_pose_.PoseCoordPickMode.PoseCoordPickMode_CBor\
CA

PoseCoordPickMode_HVY

Value:
rosetta.core.pose._core_pose_.PoseCoordPickMode.PoseCoordPickMode_HVY

PoseCoordPickMode_HVY_IF_NP

Value:
rosetta.core.pose._core_pose_.PoseCoordPickMode.PoseCoordPickMode_HVY_\
IF_NP

PoseCoordPickMode_NBR

Value:
rosetta.core.pose._core_pose_.PoseCoordPickMode.PoseCoordPickMode_NBR

PoseCoordPickMode_NUL

Value:
rosetta.core.pose._core_pose_.PoseCoordPickMode.PoseCoordPickMode_NUL

PoseCoordPickMode_N_CA_C

Value:
rosetta.core.pose._core_pose_.PoseCoordPickMode.PoseCoordPickMode_N_CA\
_C

PoseCoordPickMode_N_CA_C_CB

Value:
rosetta.core.pose._core_pose_.PoseCoordPickMode.PoseCoordPickMode_N_CA\
_C_CB

PoseCoordPickMode_N_C_O

Value:
rosetta.core.pose._core_pose_.PoseCoordPickMode.PoseCoordPickMode_N_C_\
O