jgromacs.data
Class Residue

java.lang.Object
  extended by jgromacs.data.Residue
All Implemented Interfaces:
java.lang.Cloneable

public class Residue
extends java.lang.Object
implements java.lang.Cloneable

Objects of this class represent a single residue


Constructor Summary
Residue()
          Constructs a new Residue object
Residue(int index, java.lang.String name, ResidueType residueType)
          Constructs a new Residue object of given index, name and residue type
Residue(int index, java.lang.String name, java.lang.String chainID, ResidueType residueType)
          Constructs a new Residue object of given index, name, chainID and residue type
Residue(java.lang.String name)
          Constructs a new Residue object of given name
 
Method Summary
 void addAtom(Atom atom)
          Adds a new atom to the residue
 java.lang.Object clone()
          Returns an identical Residue object
 double distanceAlphaCarbons(Residue other)
          Returns the Euclidean distance of alpha-carbon atoms of this amino acid an another
 double distanceClosest(Residue other)
          Returns the Euclidean distance of the closest atoms of this amino acid and another
 double distanceClosestHeavy(Residue other)
          Returns the Euclidean distance of the closest heavy atoms of this amino acid and another
 boolean equals(java.lang.Object other)
          Returns true if the two residues are identical
 java.lang.String get1LetterCode()
          Returns the 1 letter code of residue type
 java.lang.String get3LetterCode()
          Returns the 3 letter code of residue type
 PointList getAllAtomCoordinates()
          Returns the coordinates of all atoms
 Atom getAlphaCarbon()
          Returns the alpha carbon atom if any (otherwise returns null)
 Point3D getAlphaCarbonCoordinates()
          Returns the position of alpha-carbon atom
 Atom getAtom(int i)
          Returns atom #i of the residue
 Atom getAtomByIndex(int index)
          Returns the atom of given index
 Atom getAtomByName(java.lang.String name)
          Returns the atom of name
 IndexSet getAtomIndices()
          Returns the index set of all atoms
 java.util.ArrayList<Atom> getAtomsAsArrayList()
          Returns the list of atoms as an ArrayList object
 IndexSet getBackBoneAtomIndices()
          Returns the index set of backbone atoms
 Atom getBetaCarbon()
          Returns the beta carbon atom if any (otherwise returns null)
 Atom getCarbonylOxygen()
          Returns the carbonyl oxygen atom of amino acid
 java.lang.String getChainID()
          Returns the chain ID of residue
 java.lang.String getCombinedCode()
          Returns residue index and the 3 letter code of residue type
 Atom getCTerminalCarbon()
          Returns the C-terminal carbon atom of amino acid
 Atom getDeltaCarbon()
          Returns the delta carbon atom if any (otherwise returns null)
 Atom getEpsilonCarbon()
          Returns the epsilon carbon atom if any (otherwise returns null)
 Atom getGammaCarbon()
          Returns the gamma carbon atom if any (otherwise returns null)
 IndexSet getHeavyAtomIndices()
          Returns the index set of heavy atoms
 IndexSet getHydrogenAtomIndices()
          Returns the index set of hydrogen atoms
 int getIndex()
          Returns the index of residue
 IndexSet getMainChainAtomIndices()
          Returns the index set of main chain atoms
 IndexSet getMainChainPlusCbAtomIndices()
          Returns the index set of main chain atoms and beta carbon atom
 IndexSet getMainChainPlusHAtomIndices()
          Returns the index set of main chain atoms and hydrogen atoms
 java.lang.String getName()
          Returns the name of residue
 Atom getNTerminalNitrogen()
          Returns the N-terminal nitrogen atom of amino acid
 int getNumberOfAtoms()
          Returns the number of atoms in the residue
 ResidueType getResidueType()
          Returns the type of residue
 IndexSet getSideChainAtomIndices()
          Returns the index set of side chain atoms
 IndexSet getSideChainMinusHAtomIndices()
          Returns the index set of side chain atoms except of hydrogen atoms
 Atom getZetaCarbon()
          Returns the zeta carbon atom if any (otherwise returns null)
 int hashCode()
          Returns hash code
 boolean isAminoAcid()
          Returns true if it is an amino acid
 boolean isAtomIn(Atom atom)
          Returns true if the given atom is in the residue
 boolean isOther()
          Returns true if it is not an amino acid neither a water molecule
 boolean isWater()
          Returns true if it is a water molecule
 void removeAtom(Atom atom)
          Removes the given atom from the residue
 void removeAtom(int i)
          Removes atom #i from the residue
 void removeAtomByIndex(int index)
          Removes atom of given index from the residue
 void setAllAtomCoordinates(PointList pointlist)
          Sets the coordinates of all atoms
 void setAtom(int i, Atom atom)
          Replaces atom #i of the residue with the given atom
 void setAtomCoordinates(int i, Point3D coordinates)
          Sets the coordinates of atom #i
 void setAtomOfIndexCoordinates(int index, Point3D coordinates)
          Sets the coordinates of atom of given index
 void setChainID(java.lang.String chainID)
          Sets the chain ID of residue
 void setIndex(int index)
          Sets the index of residue
 void setName(java.lang.String name)
          Sets the name of residue
 void setResidueType(ResidueType residueType)
          Sets the type of residue
 java.lang.String toString()
          Returns the String representation of residue
 java.lang.String toStringInfo()
          Returns summary information about the residue
 
Methods inherited from class java.lang.Object
getClass, notify, notifyAll, wait, wait, wait
 

Constructor Detail

Residue

public Residue()
Constructs a new Residue object


Residue

public Residue(java.lang.String name)
Constructs a new Residue object of given name


Residue

public Residue(int index,
               java.lang.String name,
               ResidueType residueType)
Constructs a new Residue object of given index, name and residue type


Residue

public Residue(int index,
               java.lang.String name,
               java.lang.String chainID,
               ResidueType residueType)
Constructs a new Residue object of given index, name, chainID and residue type

Method Detail

getIndex

public int getIndex()
Returns the index of residue

Returns:
index of the residue

setIndex

public void setIndex(int index)
Sets the index of residue

Parameters:
index - index of the residue

getName

public java.lang.String getName()
Returns the name of residue

Returns:
name of the residue

setName

public void setName(java.lang.String name)
Sets the name of residue

Parameters:
name - name of the residue

getResidueType

public ResidueType getResidueType()
Returns the type of residue

Returns:
residue type

setResidueType

public void setResidueType(ResidueType residueType)
Sets the type of residue

Parameters:
residueType - residue type

getChainID

public java.lang.String getChainID()
Returns the chain ID of residue

Returns:
chain ID of residue

setChainID

public void setChainID(java.lang.String chainID)
Sets the chain ID of residue

Parameters:
chainID - chain ID

getNumberOfAtoms

public int getNumberOfAtoms()
Returns the number of atoms in the residue

Returns:
number of atoms

getAtomsAsArrayList

public java.util.ArrayList<Atom> getAtomsAsArrayList()
Returns the list of atoms as an ArrayList object

Returns:
list of atoms as an ArrayList

getAtom

public Atom getAtom(int i)
Returns atom #i of the residue

Returns:
atom #i

getAtomByIndex

public Atom getAtomByIndex(int index)
Returns the atom of given index

Parameters:
index - index of atom
Returns:
atom of given index

getAtomByName

public Atom getAtomByName(java.lang.String name)
Returns the atom of name

Parameters:
name - name of atom
Returns:
atom of given name

get1LetterCode

public java.lang.String get1LetterCode()
Returns the 1 letter code of residue type

Returns:
1 letter code

get3LetterCode

public java.lang.String get3LetterCode()
Returns the 3 letter code of residue type

Returns:
3 letter code

getCombinedCode

public java.lang.String getCombinedCode()
Returns residue index and the 3 letter code of residue type

Returns:
combined code

addAtom

public void addAtom(Atom atom)
Adds a new atom to the residue

Parameters:
atom - new atom

setAtom

public void setAtom(int i,
                    Atom atom)
Replaces atom #i of the residue with the given atom

Parameters:
atom - new atom

removeAtom

public void removeAtom(int i)
Removes atom #i from the residue


removeAtom

public void removeAtom(Atom atom)
Removes the given atom from the residue

Parameters:
atom - the atom to be removed

removeAtomByIndex

public void removeAtomByIndex(int index)
Removes atom of given index from the residue

Parameters:
index - index of atom to be removed

getNTerminalNitrogen

public Atom getNTerminalNitrogen()
Returns the N-terminal nitrogen atom of amino acid

Returns:
N-terminal nitrogen atom

getAlphaCarbon

public Atom getAlphaCarbon()
Returns the alpha carbon atom if any (otherwise returns null)

Returns:
alpha carbon atom

getBetaCarbon

public Atom getBetaCarbon()
Returns the beta carbon atom if any (otherwise returns null)

Returns:
beta carbon atom

getGammaCarbon

public Atom getGammaCarbon()
Returns the gamma carbon atom if any (otherwise returns null)

Returns:
gamma carbon atom

getDeltaCarbon

public Atom getDeltaCarbon()
Returns the delta carbon atom if any (otherwise returns null)

Returns:
delta carbon atom

getEpsilonCarbon

public Atom getEpsilonCarbon()
Returns the epsilon carbon atom if any (otherwise returns null)

Returns:
epsilon carbon atom

getZetaCarbon

public Atom getZetaCarbon()
Returns the zeta carbon atom if any (otherwise returns null)

Returns:
zeta carbon atom

getCTerminalCarbon

public Atom getCTerminalCarbon()
Returns the C-terminal carbon atom of amino acid

Returns:
C-terminal carbon atom

getCarbonylOxygen

public Atom getCarbonylOxygen()
Returns the carbonyl oxygen atom of amino acid

Returns:
carbonyl oxygen atom

getAtomIndices

public IndexSet getAtomIndices()
Returns the index set of all atoms

Returns:
index set of all atoms

getHydrogenAtomIndices

public IndexSet getHydrogenAtomIndices()
Returns the index set of hydrogen atoms

Returns:
index set of hydrogen atoms

getHeavyAtomIndices

public IndexSet getHeavyAtomIndices()
Returns the index set of heavy atoms

Returns:
index set of heavy atoms

getBackBoneAtomIndices

public IndexSet getBackBoneAtomIndices()
Returns the index set of backbone atoms

Returns:
index set of backbone atoms

getMainChainAtomIndices

public IndexSet getMainChainAtomIndices()
Returns the index set of main chain atoms

Returns:
index set of main chain atoms

getMainChainPlusCbAtomIndices

public IndexSet getMainChainPlusCbAtomIndices()
Returns the index set of main chain atoms and beta carbon atom

Returns:
index set of main chain atoms and beta carbon atom

getMainChainPlusHAtomIndices

public IndexSet getMainChainPlusHAtomIndices()
Returns the index set of main chain atoms and hydrogen atoms

Returns:
index set of main chain atoms and hydrogen atoms

getSideChainAtomIndices

public IndexSet getSideChainAtomIndices()
Returns the index set of side chain atoms

Returns:
index set of side chain atoms

getSideChainMinusHAtomIndices

public IndexSet getSideChainMinusHAtomIndices()
Returns the index set of side chain atoms except of hydrogen atoms

Returns:
index set of side chain atoms except of hydrogens

setAtomCoordinates

public void setAtomCoordinates(int i,
                               Point3D coordinates)
Sets the coordinates of atom #i

Parameters:
coordinates - new atomic coordinates

setAtomOfIndexCoordinates

public void setAtomOfIndexCoordinates(int index,
                                      Point3D coordinates)
Sets the coordinates of atom of given index

Parameters:
coordinates - new atomic coordinates

setAllAtomCoordinates

public void setAllAtomCoordinates(PointList pointlist)
Sets the coordinates of all atoms

Parameters:
pointlist - list of new atomic coordinates

getAllAtomCoordinates

public PointList getAllAtomCoordinates()
Returns the coordinates of all atoms

Returns:
list of atomic coordinates

getAlphaCarbonCoordinates

public Point3D getAlphaCarbonCoordinates()
Returns the position of alpha-carbon atom

Returns:
position of alpha-carbon atom

distanceAlphaCarbons

public double distanceAlphaCarbons(Residue other)
Returns the Euclidean distance of alpha-carbon atoms of this amino acid an another

Parameters:
other - other amino acid
Returns:
Euclidean distance of alpha-carbon atoms

distanceClosest

public double distanceClosest(Residue other)
Returns the Euclidean distance of the closest atoms of this amino acid and another

Parameters:
other - other amino acid
Returns:
Euclidean distance of closest atoms

distanceClosestHeavy

public double distanceClosestHeavy(Residue other)
Returns the Euclidean distance of the closest heavy atoms of this amino acid and another

Parameters:
other - other amino acid
Returns:
Euclidean distance of closest heavy atoms

isAminoAcid

public boolean isAminoAcid()
Returns true if it is an amino acid


isWater

public boolean isWater()
Returns true if it is a water molecule


isOther

public boolean isOther()
Returns true if it is not an amino acid neither a water molecule


isAtomIn

public boolean isAtomIn(Atom atom)
Returns true if the given atom is in the residue

Parameters:
atom - the atom to search for

clone

public java.lang.Object clone()
Returns an identical Residue object

Overrides:
clone in class java.lang.Object
Returns:
clone of the residue

toString

public java.lang.String toString()
Returns the String representation of residue

Overrides:
toString in class java.lang.Object
Returns:
String representation

toStringInfo

public java.lang.String toStringInfo()
Returns summary information about the residue

Returns:
summary information

equals

public boolean equals(java.lang.Object other)
Returns true if the two residues are identical

Overrides:
equals in class java.lang.Object
Parameters:
other - the other residue

hashCode

public int hashCode()
Returns hash code

Overrides:
hashCode in class java.lang.Object