PSS - Protein Structural Statistics - object-oriented classes and methods |
PSS::Atom - Protein Structural Statistics - Atom objects
use PSS::Atom;
$atom = PSS::Atom->new(ID => '1', NAME => 'C', ALT => ' ', RESNAME => 'GLY', CHAIN => 'A', RESID => '1', INS => ' ', COORD => [$x,$y,$z], OCCUP => '1.00', BFACTOR => '1.00');
print $atom->print_pdb_atom;
This module contains methods to describe atom objects. An atom is defined by its id, name, alternative position id, residue name, chain id, residue id, insertion code, coordinates, occupancy factor, and B-factor. These attributes are commonly used in molecular file formats like PDB.
These attributes are accessed or modified through methods defined below.
In general, to get the value of a property, use $object->method
without
any parameters. To set the value, use $object->method($new_value)
.
Methods are also provided to read and print atom description in PDB format, and to compute geometrical values.
Create a new PSS::Atom object with the specified attributes.
$atom = PSS::Atom->new(ID=>$id, NAME=>$name, ALT=>$alt, RESNAME=>$resname, CHAIN=>$chain, RESID=>$resid, INS=>$ins, COORD=>\@coord, OCCUP=>$occup, BFACTOR=>$bfactor);
Create a deep copy of a PSS::Atom object.
$newatom = PSS::Atom->clone($atom);
Get or set the identifier of the atom.
$id = $atom->id(); $atom->id($newid);
Get or set the name of the atom.
$name = $atom->name(); $atom->name($newname);
Get or set the alternative position identifier of the atom.
$alt = $atom->alt(); $atom->alt($newalt);
Get or set the name (3 letters) of the residue that contains the atom.
$resname = $atom->resname(); $atom->resname($newresname);
Get or set the name (1 letter) of the residue that contains the atom.
$resname1 = $atom->resname1(); $atom->resname1($newresname1);
Get or set the chain identifier of the molecule that contains the atom.
$chain = $atom->chain(); $atom->chain($newchain);
Get or set the identifier of the residue that contains the atom.
$resid = $atom->resid(); $atom->resid($newresid);
Get or set the insertion code of the residue that contains the atom.
$ins = $atom->ins(); $atom->ins($newins);
Get or set the Cartesian coordinates of the atom, as a 3 members array.
@coord = $atom->coord(); $atom->coord(@newcoord);
Get or set the x Cartesian coordinate of the atom.
$x = $atom->x(); $atom->x($newx);
Get or set the y Cartesian coordinate of the atom.
$y = $atom->y(); $atom->y($newy);
Get or set the z Cartesian coordinate of the atom.
$z = $atom->z(); $atom->z($newz);
Get or set the occupancy factor of the atom.
$occup = $atom->occup(); $atom->occup($newoccup);
Get or set the B factor of the atom.
$bfactor = $atom->bfactor(); $atom->bfactor($newbfactor);
Read a PDB file line and create a PSS::Atom object.
$atom = PSS::Atom->read_pdb_atom($pdbline);
Return the distance between 2 atoms.
$dist = PSS::Atom->dist_atom($atom1,$atom2);
Return the value in degrees of the angle formed by 3 atoms.
$angle = PSS::Atom->angle_atom($atom1,$atom2,$atom3);
Return the value in degrees of the dihedral angle formed by 4 atoms.
$dihedral = PSS::Atom->dihedral_atom($atom1,$atom2,$atom3,$atom4);
Return a string containing a PDB file format line describing the atom.
$pdbline = $atom->print_pdb_atom();
PSS::Residue - Protein Structural Statistics - Residue objects
use PSS::Residue;
@atoms = PSS::Residue->resatoms($resname3);
$resname1 = PSS::Residue->resname3to1($resname3);
$resname3 = PSS::Residue->resname1to3($resname1);
This module contains methods to describe residue objects. Residues are
regular groups of atoms in polymer molecules. A residue is defined by
its name, its identifier, and insertion code. These attributes are
accessed or modified through methods defined below. In general, to
get the value of a property, use $object->method
without any
parameters. To set the value, use $object->method($new_value)
.
Class methods to handle residue name conversions and to show atoms contained by a given residue are also provided.
Create a new PSS::Residue object with the specified attributes.
$residue = PSS::Residue->new(RESNAME=>$resname, RESID=>$resid, INS=>$ins);
RESNAME is a 3 letter name.
Create a deep copy of a PSS::Residue object.
$newresidue = PSS::Residue->clone($residue);
Get or set the 3 letter name of the residue.
$resname = $residue->resname(); $residue->resname($newresname);
Get or set the identifier of the residue.
$resid = $residue->resid(); $residue->resid($newresid);
Get or set the insertion code of the residue.
$ins = $residue->ins(); $residue->ins($newins);
Get the 1 letter name corresponding to the 3 letter name of the residue.
$resname1 = $residue->resname1();
Return a list of atoms contained in a given 3 letter residue type.
@resatoms = PSS::Residue->resatoms($resname);
Return the number of sidechain dihedral angles of a given 3 letter residue type.
$nbchis = PSS::Residue->nbchis($resname);
Return a list of atoms defining the chi_$number sidechain dihedral angle of a given 3 letter residue type.
@chiatoms = PSS::Residue->chiatoms($resname,$number);
Return the 3 letter name of the residue corresponding to the given 1 letter name.
$resname3 = PSS::Residue->resname1to3($resname1);
Return the 1 letter name of the residue corresponding to the given 3 letter name.
$resname1 = PSS::Residue->resname3to1($resname3);
PSS::Sequence - Protein Structural Statistics - Sequence objects
use PSS::Sequence;
$sequence = PSS::Sequence->read_pdb_chain($pdbname,$pdbchain);
@pir_file_lines = $sequence->print_pir_sequence;
This module contains methods to describe sequence objects. A sequence is the list of residues of a polymer molecule. A sequence is defined by its name, and a list of PSS::Residue objects. Indexations of Residue objects is added to facilitate direct access to these elements within the Sequence object.
These attributes are accessed or modified through methods defined below.
In general, to get the value of a property, use $object->method
without
any parameters. To set the value, use $object->method($new_value)
.
Methods are also provided to read a Sequence object from PDB, PIR, or FASTA format files, and to print a Sequence object in PIR or FASTA formats.
Create a new PSS::Sequence object with the specified attributes.
$molecule = PSS::Sequence->new(NAME=>$name, RESIDUES=>\@residues, RESBYID=>\%resbyid, HASINS=>$hasins);
@residues is an array containing PSS::Residue objects.
%resbyid is a hash of PSS::Residue objects, indexed by "RESID,INS" attributes.
$hasins is a flag (0/1) indicating the presence of residues with insertion codes.
Create a deep copy of a PSS::Sequence object.
$newsequence = PSS::Sequence->clone($sequence);
Read a PDB file chain and create a PSS::Sequence object.
$sequence = PSS::Sequence->read_pdb_chain($pdbname,$pdbchain);
Get or set the name of the sequence.
$name = $sequence->name(); $sequence->name($newname);
Get an array of PSS::Residue objects, corresponding to the residues contained in the sequence.
@residues = $sequence->residues();
Get a hash of PSS::Residue objects, corresponding to the residues contained in the sequence, indexed by PSS::Residue "RESID,INS" attributes. If a "$resid,$ins" value is passed, the matching PSS::Residue object is returned.
%resbyid = $sequence->resbyid(); $residue = $sequence->resbyid("$resid,$ins");
Get or set the flag (0/1) indicating the presence of residues with insertion codes.
$hasins = $sequence->hasins(); $sequence->hasins($newhasins);
Return the number of residues in the sequence.
$nbresidues = $sequence->nbresidues();
Return the PSS::Residue object corresponding to the first residue of the sequence.
$firstres = $sequence->firstres();
Return the PSS::Residue object corresponding to the last residue of the sequence.
$lastres = $sequence->lastres();
Read lines describing a sequence in PIR format and create a PSS::Sequence object.
$sequence = PSS::Sequence->read_pir_sequence(\@pir_file_lines);
Read lines describing a sequence in FASTA format and create a PSS::Sequence object.
$sequence = PSS::Sequence->read_fasta_sequence(\@fasta_file_lines);
Return an array, containing the lines of a file in PIR format, describing the sequence.
@pir_file_lines = $sequence->print_pir_sequence();
A PIR file can be readily written with:
open(PIR,">$pirname"); for (@pir_file_lines) {print PIR $_};
Return an array, containing the lines of a file in FASTA format, describing the sequence.
@fasta_file_lines = $sequence->print_fasta_sequence();
A FASTA file can be readily written with:
open(FASTA,">$fastaname"); for (@fasta_file_lines) {print FASTA $_};
PSS::Alignment - Protein Structural Statistics - Alignment objects
use PSS::Alignment;
$alignment = PSS::Alignment->read_pir_alignment($pirname);
@pir_file_lines = $alignment->print_pir_alignment;
This module contains methods to describe alignment objects. An alignment is a collection of sequences of the same length. An alignment is defined by its name and a list of PSS::Sequence objects. Verification is made that sequences contained in the alignment have the same length.
These attributes are accessed or modified through methods defined below.
In general, to get the value of a property, use $object->method
without
any parameters. To set the value, use $object->method($new_value)
.
Methods are also provided to read an Alignment object from a PIR or FASTA format file, and to print an Alignment object in PIR, PAP, or FASTA formats.
Create a new PSS::Alignment object with the specified attributes.
$alignment = PSS::Alignment->new(NAME=>$name, SEQUENCES=>\@sequences);
@sequences is an array containing PSS::Sequence objects.
Create a deep copy of a PSS::Alignment object.
$newalignment = PSS::Alignment->clone($alignment);
Read a PIR alignment file and create a PSS::Alignment object.
$alignment = PSS::Alignment->read_pir_alignment($pirname);
Read a FASTA alignment file and create a PSS::Alignment object.
$alignment = PSS::Alignment->read_fasta_alignment($fastaname);
Get or set the name of the alignment.
$name = $alignment->name(); $alignment->name($newname);
Get an array of PSS::Sequence objects, corresponding to the sequences contained in the alignment.
@sequences = $alignment->sequences();
Return the number of sequences in the alignment.
$nbseq = $alignment->nbseq();
Return the length of the alignment.
$length = $alignment->length();
Return an array, containing the lines of an alignment file in PIR format.
@pir_file_lines = $alignment->print_pir_alignment();
A PIR file can be readily written with:
open(PIR,">$pirname"); for (@pir_file_lines) {print PIR $_};
Return an array, containing the lines of an alignment file in FASTA format.
@fasta_file_lines = $alignment->print_fasta_alignment();
A FASTA file can be readily written with:
open(FASTA,">$fastaname"); for (@fasta_file_lines) {print FASTA $_};
Return an array, containing the lines of an alignment file in PAP format.
@pap_file_lines = $alignment->print_pap_alignment();
A PAP file can be readily written with:
open(PAP,">$papname"); for (@pap_file_lines) {print PAP $_};
PSS::Molecule - Protein Structural Statistics - Molecule objects
use PSS::Molecule;
$molecule = PSS::Molecule->read_pdb_chain($pdbname,$pdbchain);
@pdb_file_lines = $molecule->print_pdb_molecule;
This module contains methods to describe molecule objects. A molecule is formed by atoms and can contain residues. A molecule is defined by its name, a list of PSS::Atom objects, and a list of PSS::Residue objects. Indexations of Atom and Residue objects is added to facilitate direct access to these elements within the Molecule object.
These attributes are accessed or modified through methods defined below.
In general, to get the value of a property, use $object->method
without
any parameters. To set the value, use $object->method($new_value)
.
Methods are also provided to read and print a Molecule object in PDB format, and to print the sequence of the molecule residues in PIR or FASTA formats.
Create a new PSS::Molecule object with the specified attributes.
$molecule = PSS::Molecule->new(NAME=>$name, ATOMS=>\@atoms, BYID=>\%byid, BYNAMERESID=>\%bynameresid, RESIDUES=>\@residues, RESBYID=>\%resbyid, HASINS=>$hasins);
@atoms is an array containing PSS::Atom objects.
%byid and %bynameresid are hashes of PSS::Atom objects, respectively indexed by ID and "NAME,RESID,INS" attributes.
@residues is an array containing PSS::Residue objects.
%resbyid is a hash of PSS::Residue objects, indexed by "RESID,INS" attributes.
$hasins is a flag (0/1) indicating the presence of residues with insertion codes.
Create a deep copy of a PSS::Molecule object.
$newmolecule = PSS::Molecule->clone($molecule);
Read a PDB file chain and create a PSS::Molecule object.
$molecule = PSS::Molecule->read_pdb_chain($pdbname,$pdbchain);
Get or set the name of the molecule.
$name = $molecule->name(); $molecule->name($newname);
Get an array of PSS::Atom objects, corresponding to the atoms contained in the molecule.
@atoms = $molecule->atoms();
Get a hash of PSS::Atom objects, corresponding to the atoms contained in the molecule, indexed by PSS::Atom ID attribute. If a $id value is passed, the matching PSS::Atom object is returned.
%byid = $molecule->byid(); $atom = $molecule->byid($id);
Get a hash of PSS::Atom objects, corresponding to the atoms contained in the molecule, indexed by PSS::Atom "NAME,RESID,INS" attributes. If a "$name,$resid,$ins" value is passed, the matching PSS::Atom object is returned.
%bynameresid = $molecule->bynameresid(); $atom = $molecule->bynameresid("$name,$resid,$ins");
Get an array of PSS::Residue objects, corresponding to the residues contained in the molecule.
@residues = $molecule->residues();
Get a hash of PSS::Residue objects, corresponding to the residues contained in the molecule, indexed by PSS::Residue "RESID,INS" attributes. If a "$resid,$ins" value is passed, the matching PSS::Residue object is returned.
%resbyid = $molecule->resbyid(); $residue = $molecule->resbyid("$resid,$ins");
Get or set the flag (0/1) indicating the presence of residues with insertion codes.
$name = $molecule->hasins(); $molecule->hasins($newhasins);
Return the number of atoms in the molecule.
$nbatoms = $molecule->nbatoms();
Return the number of residues in the molecule.
$nbresidues = $molecule->nbresidues();
Return the PSS::Residue object corresponding to the first residue of the molecule.
$firstres = $molecule->firstres();
Return the PSS::Residue object corresponding to the last residue of the molecule.
$lastres = $molecule->lastres();
Return an array, containing the lines of a file in PDB format, describing the molecule.
@pdb_file_lines = $molecule->print_pdb_molecule();
A PDB file can be readily written with:
open(PDB,">$pdbname"); for (@pdb_file_lines) {print PDB $_};
Return an array, containing the lines of a file in PIR format, describing the sequence of the molecule residues.
@pir_file_lines = $molecule->print_pir_sequence();
A PIR file can be readily written with:
open(PIR,">$pirname"); for (@pir_file_lines) {print PIR $_};
Return an array, containing the lines of a file in FASTA format, describing the sequence of the molecule residues.
@fasta_file_lines = $molecule->print_fasta_sequence();
A FASTA file can be readily written with:
open(FASTA,">$fastaname"); for (@fasta_file_lines) {print FASTA $_};
Return the value in degrees of the Phi dihedral angle of a residue determined by its "RESID,INS" attributes. The previous residue is also needed to calculate Phi. Verification is made that the two residues are contiguous.
$phi = $molecule->phi("$resid,$ins","$prevresid,$previns");
Return the value in degrees of the Psi dihedral angle of a residue determined by its "RESID,INS" attributes. The next residue is also needed to calculate Psi. Verification is made that the two residues are contiguous.
$psi = $molecule->psi("$resid,$ins","$nextresid,$nextins");
Return the value in degrees of the Omega dihedral angle of a residue determined by its "RESID,INS" attributes. The next residue is also needed to calculate Omega. Verification is made that the two residues are contiguous.
$omega = $molecule->omega("$resid,$ins","$nextresid,$nextins");
Return the value in degrees of a Chi dihedral angle of a residue determined by its "RESID,INS" attributes. The Chi angle number (from 1 to 4) is provided as first argument and the "RESID,INS" as second argument.
$chi1 = $molecule->chi(1,"$resid,$ins");
Return the value in degrees of the Chi1 dihedral angle of a residue determined by its "RESID,INS" attributes.
$chi1 = $molecule->chi1("$resid,$ins");
Return the value in degrees of the Chi2 dihedral angle of a residue determined by its "RESID,INS" attributes.
$chi2 = $molecule->chi2("$resid,$ins");
Return the value in degrees of the Chi3 dihedral angle of a residue determined by its "RESID,INS" attributes.
$chi3 = $molecule->chi2("$resid,$ins");
Return the value in degrees of the Chi4 dihedral angle of a residue determined by its "RESID,INS" attributes.
$chi4 = $molecule->chi2("$resid,$ins");
We appreciate if users report bugs, feature requests, improvements, or patches to <thomas.gaillard@polytechnique.edu>.
The project homepage is http://thomasgaillard.fr/pss.
Please cite the following article in a publication using PSS:
Protein structural statistics with PSS
Thomas Gaillard, Benjamin B.L. Schwarz, Yassmine Chebaro, Roland H. Stote, and Annick Dejaegere
J. Chem. Inf. Model. 2013, 53, 2471-2482.
DOI:10.1021/ci400233j
1.2
Thomas Gaillard <thomas.gaillard@polytechnique.edu>, with contributions from Benjamin Schwarz, Roland H. Stote, and Annick Dejaegere
Copyright (c) 2006-2016 Thomas Gaillard, Benjamin Schwarz, Roland H. Stote, Annick Dejaegere
This program is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; version 2 or any later version.
PSS - Protein Structural Statistics - object-oriented classes and methods |