jgromacs.data
Class Structure

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

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

Objects of this class represent a single structure


Constructor Summary
Structure()
          Constructs a new Structure object
Structure(java.lang.String name)
          Constructs a new Structure object of given name
 
Method Summary
 void addAtomToResidue(int i, Atom atom)
          Adds the given atom to residue #i
 void addAtomToResidueOfIndex(int index, Atom atom)
          Adds the given atom to the residue of given index
 void addResidue(Residue residue)
          Adds a new residue to the structure
 java.lang.Object clone()
          Returns an identical Structure object
 java.util.ArrayList<java.lang.Integer> convertIndicesToArrayListIndices(IndexSet indices)
          Converts atomic indices into list indices
 boolean equals(java.lang.Object other)
          Returns true if the two structures are identical
 PointList getAllAtomCoordinates()
          Returns the coordinates of all atoms
 IndexSet getAlphaCarbonIndexSet()
          Returns the index set of alpha carbon atoms
 Atom getAtom(int i)
          Returns atom #i of the structure
 Atom getAtomByIndex(int index)
          Returns the atom of given index
 Atom getAtomInResidueOfIndex(int residueindex, int i)
          Returns atom #i of the residue of given index
 Atom getAtomInResidueOfIndex(int residueindex, java.lang.String chainID, int i)
          Returns atom #i of the residue of given index and chain ID
 IndexSet getBackboneIndexSet()
          Returns the index set of backbone atoms
 java.util.ArrayList<java.lang.String> getChainIDs()
          Returns the list of chain IDs in the structure
 Structure[] getChains()
          Returns all chains in an array of Structure objects
 IndexSetList getDefaultIndexSetList()
          Returns the list of default index sets
 IndexSet getHeavyProteinIndexSet()
          Returns the index set of protein atoms except of hydrogen atoms
 IndexSet getIndexSetOfChainID(java.lang.String chainID)
          Returns the index set of atoms of a given chain ID
 IndexSet getMainChainIndexSet()
          Returns the index set of main chain atoms
 IndexSet getMainChainPlusCbIndexSet()
          Returns the index set of main chain plus beta carbon atoms
 IndexSet getMainChainPlusHIndexSet()
          Returns the index set of main chain plus hydrogen atoms
 java.lang.String getName()
          Returns the name of structure
 IndexSet getNonProteinIndexSet()
          Returns the index set of non-protein atoms
 int getNumberOfAtoms()
          Returns the number of atoms in the structure
 int getNumberOfChains()
          Returns the number of chains in the structure
 int getNumberOfResidues()
          Returns the number of residues in the structure
 IndexSet getProteinIndexSet()
          Returns the index set of protein atoms
 Residue getResidue(int i)
          Returns residue #i of the structure
 Residue getResidueByIndex(int index)
          Returns the residue of given index
 Residue getResidueByIndex(int index, java.lang.String chainID)
          Returns the residue of given index and given chain ID
 java.util.ArrayList<Residue> getResiduesAsArrayList()
          Returns residues as an ArrayList object
 Sequence getSequence()
          Returns the amino acid sequence of the structure
 IndexSet getSideChainIndexSet()
          Returns the index set of side chain atoms
 IndexSet getSideChainMinusHIndexSet()
          Returns the index set of side chain atoms except of hydrogen atoms
 Structure getSubStructure(IndexSet indices)
          Returns the substructure defined by the given index set
 IndexSet getSystemIndexSet()
          Returns the index set of all atoms in the system
 IndexSet getWaterIndexSet()
          Returns the index set of water atoms
 int hashCode()
          Returns hash code
 java.util.ArrayList<java.lang.Integer> mapAtomIndicesToResidueIndices(IndexSet indices)
          Maps atomic indices to residue indices (i.e.
 IndexSet mapResidueIndicesToAtomIndices(java.util.ArrayList<java.lang.Integer> indices)
          Maps residue indices to atomic indices (i.e.
 void removeAtom(Atom atom)
          Removes the given atom from the structure
 void removeAtom(int i)
          Removes atom #i from the structure
 void removeAtomByIndex(int index)
          Removes the atom of given index from the structure
 void removeResidue(int i)
          Removes residue #i from the structure
 void removeResidue(Residue residue)
          Removes a residue from the structure
 void setAllAtomCoordinates(PointList pointlist)
          Sets the coordinates of all atoms
 void setAtomCoordinates(int i, Point3D coordinates)
          Sets the coordinates of atom #i
 void setAtomInResidue(int i, int j, Atom atom)
          Replaces atom #j of residue #i with a new atom
 void setAtomInResidueOfIndex(int index, int i, Atom atom)
          Replaces atom #i of the residue of given index with a new atom
 void setAtomOfIndexCoordinates(int index, Point3D coordinates)
          Sets the coordinates of atom of given index
 void setName(java.lang.String name)
          Sets the name of structure
 void setResidue(int i, Residue residue)
          Replaces residue #i with a new residue
 java.lang.String toString()
          Returns the String representation of structure
 java.lang.String toStringAsGRO()
          Returns the String representation of structure in GRO format
 java.lang.String toStringAsPDB()
          Returns the String representation of structure in PDB format
 java.lang.String toStringInfo()
          Returns summary information about the structure
 
Methods inherited from class java.lang.Object
getClass, notify, notifyAll, wait, wait, wait
 

Constructor Detail

Structure

public Structure()
Constructs a new Structure object


Structure

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

Method Detail

getName

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

Returns:
name of structure

setName

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

Parameters:
name - of structure

getResiduesAsArrayList

public java.util.ArrayList<Residue> getResiduesAsArrayList()
Returns residues as an ArrayList object

Returns:
residues as an ArrayList

getNumberOfResidues

public int getNumberOfResidues()
Returns the number of residues in the structure

Returns:
number of residues

getNumberOfAtoms

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

Returns:
number of atoms

getNumberOfChains

public int getNumberOfChains()
Returns the number of chains in the structure

Returns:
number of chains

getChainIDs

public java.util.ArrayList<java.lang.String> getChainIDs()
Returns the list of chain IDs in the structure

Returns:
list of chain IDs

getResidue

public Residue getResidue(int i)
Returns residue #i of the structure

Returns:
residue #i

getResidueByIndex

public Residue getResidueByIndex(int index,
                                 java.lang.String chainID)
Returns the residue of given index and given chain ID

Parameters:
index - index of residue
chainID - chain ID of residue
Returns:
residue of given index and chain ID

getResidueByIndex

public Residue getResidueByIndex(int index)
Returns the residue of given index

Parameters:
index - index of residue
Returns:
residue of given index

getAtom

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

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

getAtomInResidueOfIndex

public Atom getAtomInResidueOfIndex(int residueindex,
                                    java.lang.String chainID,
                                    int i)
Returns atom #i of the residue of given index and chain ID

Parameters:
residueindex - index of residue
chainID - chain ID of residue
Returns:
atom #i in residue of given index and chain ID

getAtomInResidueOfIndex

public Atom getAtomInResidueOfIndex(int residueindex,
                                    int i)
Returns atom #i of the residue of given index

Parameters:
residueindex - index of residue
Returns:
atom #i in residue of given index

addResidue

public void addResidue(Residue residue)
Adds a new residue to the structure

Parameters:
residue - new residue

removeResidue

public void removeResidue(Residue residue)
Removes a residue from the structure

Parameters:
residue - the residue to be removed

removeResidue

public void removeResidue(int i)
Removes residue #i from the structure


setResidue

public void setResidue(int i,
                       Residue residue)
Replaces residue #i with a new residue

Parameters:
residue - new residue

removeAtom

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


removeAtomByIndex

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

Parameters:
index - index of atom

removeAtom

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

Parameters:
atom - the atom to be removed

addAtomToResidue

public void addAtomToResidue(int i,
                             Atom atom)
Adds the given atom to residue #i

Parameters:
atom - new atom

addAtomToResidueOfIndex

public void addAtomToResidueOfIndex(int index,
                                    Atom atom)
Adds the given atom to the residue of given index

Parameters:
index - index of residue
atom - new atom

setAtomInResidue

public void setAtomInResidue(int i,
                             int j,
                             Atom atom)
Replaces atom #j of residue #i with a new atom

Parameters:
atom - new atom

setAtomInResidueOfIndex

public void setAtomInResidueOfIndex(int index,
                                    int i,
                                    Atom atom)
Replaces atom #i of the residue of given index with a new atom

Parameters:
index - index of residue
atom - new atom

getSystemIndexSet

public IndexSet getSystemIndexSet()
Returns the index set of all atoms in the system

Returns:
index set of all atoms

getProteinIndexSet

public IndexSet getProteinIndexSet()
Returns the index set of protein atoms

Returns:
index set of protein atoms

getHeavyProteinIndexSet

public IndexSet getHeavyProteinIndexSet()
Returns the index set of protein atoms except of hydrogen atoms

Returns:
index set of protein atoms except of hydrogen atoms

getAlphaCarbonIndexSet

public IndexSet getAlphaCarbonIndexSet()
Returns the index set of alpha carbon atoms

Returns:
index set of alpha carbon atoms

getBackboneIndexSet

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

Returns:
index set of backbone atoms

getMainChainIndexSet

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

Returns:
index set of main chain atoms

getMainChainPlusCbIndexSet

public IndexSet getMainChainPlusCbIndexSet()
Returns the index set of main chain plus beta carbon atoms

Returns:
index set of main chain plus beta carbon atoms

getMainChainPlusHIndexSet

public IndexSet getMainChainPlusHIndexSet()
Returns the index set of main chain plus hydrogen atoms

Returns:
index set of main chain plus hydrogen atoms

getSideChainIndexSet

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

Returns:
index set of side chain atoms

getSideChainMinusHIndexSet

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

Returns:
index set of side chain atoms except of hydrogen atoms

getNonProteinIndexSet

public IndexSet getNonProteinIndexSet()
Returns the index set of non-protein atoms

Returns:
index set of non-protein atoms

getWaterIndexSet

public IndexSet getWaterIndexSet()
Returns the index set of water atoms

Returns:
index set of water atoms

getDefaultIndexSetList

public IndexSetList getDefaultIndexSetList()
Returns the list of default index sets

Returns:
list of default index sets

getIndexSetOfChainID

public IndexSet getIndexSetOfChainID(java.lang.String chainID)
Returns the index set of atoms of a given chain ID

Parameters:
chainID - selected chain ID
Returns:
index set defined by chain ID

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

getSequence

public Sequence getSequence()
Returns the amino acid sequence of the structure

Returns:
amino acid sequence

getChains

public Structure[] getChains()
Returns all chains in an array of Structure objects

Returns:
array of chains

getSubStructure

public Structure getSubStructure(IndexSet indices)
Returns the substructure defined by the given index set

Parameters:
indices - index set of substructure
Returns:
substructure

convertIndicesToArrayListIndices

public java.util.ArrayList<java.lang.Integer> convertIndicesToArrayListIndices(IndexSet indices)
Converts atomic indices into list indices

Parameters:
indices - set of atomic indices
Returns:
list indices

mapAtomIndicesToResidueIndices

public java.util.ArrayList<java.lang.Integer> mapAtomIndicesToResidueIndices(IndexSet indices)
Maps atomic indices to residue indices (i.e. returns which residues contain the given atoms). If an atom is not present in any residue it is mapped to residue index -1

Parameters:
indices - set of atomic indices
Returns:
residue indices

mapResidueIndicesToAtomIndices

public IndexSet mapResidueIndicesToAtomIndices(java.util.ArrayList<java.lang.Integer> indices)
Maps residue indices to atomic indices (i.e. returns which atoms are present in the given residues).

Parameters:
indices - ArrayList of residue indices
Returns:
IndexSet of atomic indices

toString

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

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

toStringAsGRO

public java.lang.String toStringAsGRO()
Returns the String representation of structure in GRO format

Returns:
String representation in GRO format

toStringAsPDB

public java.lang.String toStringAsPDB()
Returns the String representation of structure in PDB format

Returns:
String representation in PDB format

toStringInfo

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

Returns:
summary information

clone

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

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

equals

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

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

hashCode

public int hashCode()
Returns hash code

Overrides:
hashCode in class java.lang.Object