Main Page | Namespace List | Class Hierarchy | Alphabetical List | Class List | File List | Class Members | File Members

GromacsTopFile.h

Go to the documentation of this file.
00001 #ifndef GROMACSTOPFILE_H
00002 #define GROMACSTOPFILE_H
00003 
00004 typedef float Real;
00005 
00006 #define NAMESIZE 5
00007 #define LONGNAMESIZE 20
00008 
00009 /* A GenericBond represents a bond between two atoms in a template
00010    that can be used one or more times within a simulation.  This is an
00011    immutable type.
00012    Here, atomi and atomj start at 0. */
00013 class GenericBond {
00014  private:
00015   int atomi; /* the first atom of the bond */
00016   int atomj; /* the second atom of the bond */
00017   int type; /* an index into the bond's parameters, stored elsewhere */
00018 
00019  public:
00020   /* initializes this to be a bond between atoms <i> and <j> of type
00021      <type> */
00022   GenericBond(int i, int j, int theType);
00023 
00024   int getAtomi() const {return atomi;}
00025   int getAtomj() const {return atomj;}
00026   int getType() const {return type;}
00027 
00028 };
00029 
00030 /* A GenericAngle represents a angle between two atoms in a template
00031    that can be used one or more times within a simulation.  This is an
00032    immutable type.
00033    Here, atomi and atomj start at 0. */
00034 class GenericAngle {
00035  private:
00036   int atomi; /* the first atom of the angle */
00037   int atomj; /* the second atom of the angle */
00038   int atomk; /* the third atom of the angle */
00039   int type; /* an index into the angle's parameters, stored elsewhere */
00040 
00041  public:
00042   /* initializes this to be a angle between atoms <i>, <j>, and <k> of
00043      type <type> */
00044   GenericAngle(int i, int j, int k, int theType);
00045 
00046   int getAtomi() const {return atomi;}
00047   int getAtomj() const {return atomj;}
00048   int getAtomk() const {return atomk;}
00049   int getType() const {return type;}
00050 
00051 };
00052 
00053 /* A GenericDihedral represents a dihedral angle between four atoms in
00054    a template that can be used one or more times within a simulation.
00055    This is an immutable type.
00056    Here, atomi and atomj start at 0. */
00057 class GenericDihedral {
00058  private:
00059   int atomi; /* the first atom of the angle */
00060   int atomj; /* the second atom of the angle */
00061   int atomk; /* the third atom of the angle */
00062   int atoml; /* the fourth atom of the angle */
00063   int type; /* an index into the parameters, stored elsewhere */
00064 
00065  public:
00066   /* initializes this to be a angle between atoms <i>, <j>, <k>, and
00067      <l> of type <type> */
00068   GenericDihedral(int i, int j, int k, int l, int theType);
00069 
00070   int getAtomi() const {return atomi;}
00071   int getAtomj() const {return atomj;}
00072   int getAtomk() const {return atomk;}
00073   int getAtoml() const {return atoml;}
00074   int getType() const {return type;}
00075 
00076 };
00077 
00078 /* A GenericAtom represents an atom in a template that can be used one
00079    or more times within a simulation.  This is an immutable type. */
00080 class GenericAtom {
00081  private:
00082   char type[NAMESIZE+1];
00083   int typeNum;
00084   int resNum;
00085   char resType[NAMESIZE+1];
00086   char atomName[NAMESIZE+1];
00087   Real charge;
00088   Real mass;
00089 
00090  public: 
00091 
00092   /* Initializes this to be the atom specified by the given
00093      parameters.  All the strings will be copied, and must be at most
00094      four characters long. */
00095   GenericAtom(const char *theType, int theTypeNum, int theResNum,
00096               const char *theResType,
00097               const char *theAtomName, Real theCharge, Real theMass);
00098 
00099   /* The following just return the parameters of the atom. */
00100   const char *getType() const { return type; }
00101   int getTypeNum() const { return typeNum; }
00102   int getResNum() const {return resNum;}
00103   const char *getResType() const {return resType;}
00104   const char *getAtomName() const {return atomName;}
00105   Real getCharge() const {return charge;}
00106   Real getMass() const {return mass;}
00107 
00108 };
00109 
00110 /* A GenericMol represents a complete molecule template.  This will be
00111    a mutable type, for convenience's sake.  The atoms here are
00112    numbered starting at ZERO */
00113 class GenericMol {
00114  private:
00115   /* the name of this Molecule, so it can be identified later */
00116   const char *name;
00117 
00118   /* the rest of the list */
00119   const GenericMol *next;
00120 
00121   /* the list of GenericAtoms */
00122   ResizeArray <const GenericAtom *> atomList;
00123 
00124   /* the list of GenericBonds */
00125   ResizeArray <const GenericBond *> bondList;
00126   
00127   /* the list of GenericAngles */
00128   ResizeArray <const GenericAngle *> angleList;
00129   
00130   /* the list of GenericDihedrals */
00131   ResizeArray <const GenericDihedral *> dihedralList;
00132   
00133  public:
00134   /* initializes this to be the molecule with name <name> */
00135   GenericMol(const char *theName);
00136 
00137   /* getNumAtoms returns the number of atoms stored here. */
00138   int getNumAtoms() const { return atomList.size(); }
00139   int getNumBonds() const { return bondList.size(); }
00140   int getNumAngles() const { return angleList.size(); }
00141   int getNumDihedrals() const { return dihedralList.size(); }
00142   int getNumRes() const { return atomList[atomList.size()-1]->getResNum(); }
00143 
00144   /* adds an atom to the end of list.  The residue numbers must start
00145      at zero and count up in unit steps */
00146   void addAtom(const char *theType, int theTypeNum, int theResNum,
00147                const char *theResType,
00148                const char *theAtomName, Real theCharge, Real theMass);
00149 
00150   /* adds a bond/angle/dihedral to the molecule */
00151   void addBond(int atomi, int atomj, int type);
00152   void addAngle(int atomi, int atomj, int atomk, int type);
00153   void addDihedral(int atomi, int atomj, int atomk, int atoml, int type);
00154 
00155   /* gets the name */
00156   const char *getName() const { return name; }
00157 
00158   /* gets atom <n>, where atom 0 is the first atom in the list. */
00159   const GenericAtom *getAtom(int n) const;
00160 
00161   /* gets bond <n>, where bond 0 is the first bond in the list. */
00162   const GenericBond *getBond(int n) const;
00163 
00164   /* gets angle <n>, where angle 0 is the first angle in the list. */
00165   const GenericAngle *getAngle(int n) const;
00166 
00167   /* gets dihedral <n>, where angle 0 is the first one in the list. */
00168   const GenericDihedral *getDihedral(int n) const;
00169 
00170 };
00171 
00172 /* A MolInst represents a number of instances of a molecule.  These
00173    objects are meant to be assembled (in order) into a list, to
00174    represent all the instances of molecules in a system. */
00175 class MolInst {
00176   int num;
00177   const GenericMol *mol;
00178  public:
00179   /* initializes this to represent <theNum> copies of <theMol> */
00180   MolInst(const GenericMol *theMol, int theNum);
00181 
00182   /* get the number */
00183   int getNum() const { return num; }
00184 
00185   /* get the total number of X here */
00186   int getNumAtoms() const;
00187   int getNumRes() const;
00188   int getNumBonds() const;
00189   int getNumAngles() const;
00190   int getNumDihedrals() const;
00191   
00192   /* get the molecule */
00193   const GenericMol *getMol() const { return mol; }
00194 };
00195 
00196 /* An AtomTable stores all the known types of atoms.
00197    mass: the mass in amu
00198    charge: the charge in e
00199    c6: the LJ parameter c6 in kcal/mol A^6
00200    c12: the LJ parameter c12 in kcal/mol A^12 */
00201 class AtomTable {
00202  private:
00203   /* implementation of the extensible array */
00204   ResizeArray<Real> mArray;
00205   ResizeArray<Real> qArray;
00206   ResizeArray<Real> c6Array;
00207   ResizeArray<Real> c12Array;
00208   ResizeArray<const char *> typeArray; /* should never be NULL */
00209 
00210  public:
00211   /* this adds a entry for an atom type to the table */
00212   void addType(const char *type, Real m, Real q, Real c6, Real c12);
00213 
00214   /* This looks up the atom by type - if it cannot be found, this
00215      returns -1, otherwise this returns the index of the atomtype */
00216   int getParams(const char *type, Real *m, Real *q,
00217                      Real *c6, Real *c12) const;
00218 
00219   /* This looks up the atom type by number, returning it in the string
00220      <type>, which must be at least 11 characters long. */
00221   void getType(int num, char *type) const;
00222 
00223 
00224   /* returns the number of types in the table */
00225   int size() const { return mArray.size(); }
00226 };
00227 
00228 /* A BondTable stores all the different types of bonds, so that we
00229    only have to have one copy of each.
00230    b0: the natural bond length in A.
00231    kB: the spring constant in kcal/mol/A^2, where E=kx^2 to 1st order.
00232    funct: 1 for harmonic potential, 2 for fourth-order GROMACS96
00233    potential. */
00234 class BondTable {
00235  private:
00236   /* implementation of the extensible array */
00237   ResizeArray<Real> b0Array;
00238   ResizeArray<Real> kBArray;
00239   ResizeArray<int> functArray;
00240   ResizeArray<const char *> typeaArray; /* NULL if no type is set */
00241   ResizeArray<const char *> typebArray; /* NULL if no type is set */
00242 
00243  public:
00244   /* returns the number of types in the table */
00245   int size() const { return b0Array.size(); }
00246 
00247   /* this gets the index of a bond in the table (adding an entry if
00248      none exists). */
00249   int getIndex(Real b0, Real kB, int funct);
00250 
00251   /* this adds a entry for a bond type to the table, including the two
00252      atoms involved in the bond */
00253   void addType(const char *typea, const char *typeb, Real b0,
00254                    Real kB, int funct);
00255 
00256   /* getParams puts the parameters for bond-type <num> into the
00257      spaces pointed to by the other arguments. */
00258   void getParams(int num, Real *b0, Real *kB, int *funct) const;
00259 
00260   /* This version looks up the bond by atom type - the order of the
00261      types doesn't matter!  If the specified bond/function combination
00262      cannot be found, it returns -1, otherwise it returns the index of
00263      the bondtype. */
00264   int getParams(const char *typea, const char *typeb, int funct,
00265                      Real *b0, Real *kB) const;
00266 };
00267 
00268 /* A AngleTable stores all the different types of angles, so that we
00269    only have to have one copy of each.  Units:
00270    th0 - degrees
00271    kth - kcal/mol/rad^2 (funct 1)  OR  kcal/mol (funct 2)
00272    funct 1: harmonic angle potential
00273    funct 2: harmonic G96 cosine potential
00274    The atoms are in standard geometry order - A--B--C for angle /ABC
00275 */
00276 class AngleTable {
00277  private:
00278   /* implementation of the extensible array */
00279   ResizeArray<Real> th0Array;
00280   ResizeArray<Real> kthArray;
00281   ResizeArray<int> functArray;
00282   ResizeArray<const char *> typeaArray; /* NULL if no type is set */
00283   ResizeArray<const char *> typebArray; /* NULL if no type is set */
00284   ResizeArray<const char *> typecArray; /* NULL if no type is set */
00285 
00286  public:
00287   /* returns the number of types in the table */
00288   int size() const { return th0Array.size(); }
00289 
00290   /* this gets the index of a angle in the table (adding an entry if
00291      none exists).
00292   */
00293   int getIndex(Real th0, Real kth, int funct);
00294 
00295   /* this adds a entry for a angle type to the table, including the
00296      three atoms involved in the angle */
00297   void addType(const char *typea, const char *typeb,
00298                     const char *typec, Real th0,
00299                     Real kth, int funct);
00300 
00301   /* getAngleParams puts the parameters for angle-type <num> into the
00302      spaces pointed to by the other arguments. 
00303   */
00304   void getParams(int num, Real *th0, Real *kth, int *funct) const;
00305 
00306   /* This version looks up the angle by atom type - the direction of
00307      the types doesn't matter (it can be A--B--C or C--B--A)!  If the
00308      specified angle/function combination cannot be found, it returns
00309      -1, otherwise it returns the index of the angletype. */
00310   int getParams(const char *typea, const char *typeb,
00311                      const char *typec, int funct,
00312                      Real *th0, Real *kth) const;
00313 };
00314 
00315 /* A DihedralTable stores all the different types of angles, so that we
00316    only have to have one copy of each.  Units:
00317    funct 1: proper dihedral
00318      c0 - thmax (deg)
00319      c1 - fc (kcal/mol)
00320      a,b are the two inner atoms
00321    funct 2: improper dihedral
00322      c0 - th0 (deg)
00323      c1 - fc (kcal/mol)
00324      a,b are the two outer atoms
00325    funct 3: RB dihedral
00326      c0-c5 - C0-C5 (kcal/mol)
00327      a,b are the two inner atoms
00328 
00329    The atoms are in standard geometry order - A--B--C--D for the
00330    dihedral between the planes containing ABC and BCD.
00331 */
00332 class DihedralTable {
00333  private:
00334   /* implementation of the extensible array */
00335   ResizeArray<int> functArray;
00336   ResizeArray<const char *> typeaArray; /* NULL if no type is set */
00337   ResizeArray<const char *> typebArray; /* NULL if no type is set */
00338   ResizeArray<int> multArray; /* 0 if funct is not 1 */
00339 
00340   /* cArray is where the actual paramters live - each Real * is a
00341      pointer to a list of Reals that is at least as long as the
00342      number of parameters required for the given type. */
00343   ResizeArray<const Real *> cArray;
00344 
00345  public:
00346   /* returns the number of types in the table */
00347   int size() const { return functArray.size(); }
00348 
00349   /* this function gets the index of a dihedral angle in the table
00350      (adding an entry if none exists).  The required number of
00351      parameters (see notes on class DihedralTable) must be stored in
00352      the array <c> */
00353   int getIndex(const Real *c, int mult, int funct);
00354 
00355   /* this adds a entry for a angle type to the table, including the
00356      two atoms involved in the angle.  The required number of
00357      parameters (see notes on class DihedralTable) must be stored in
00358      the array <c>.  Note that these two angles really refer to either
00359      atoms A and D or B and C depending on the dihedral type.  */
00360   void addType(const char *typea, const char *typeb,
00361                const Real *c, int mult, int funct);
00362 
00363   /* getParams puts the parameters for angle-type <num> into the
00364      spaces pointed to by the other arguments.  The required number of
00365      parameters (see notes on class DihedralTable) will be stored in
00366      the array <c>, so make sure that c has size >= 6 */
00367   void getParams(int num, Real *c, int *mult, int *funct) const;
00368 
00369   /* This version looks up the angle by atom type - the direction of
00370      the types doesn't matter (it can be ABCD or DCBA)!  If the
00371      specified angle/function combination cannot be found, it returns
00372      -1, otherwise it returns the index of the dihedraltye.  This is a
00373      little strange because the different dihedral types use different
00374      atoms - but since we always will know all four AND the funct when
00375      we need to look up a dihedral, this will work.  The required
00376      number of parameters (see notes on class DihedralTable) will be
00377      stored in the array <c>, so make sure that c is big enough for
00378      the specified function. */
00379   int getParams(const char *typea, const char *typeb,
00380                 const char *typec, const char *typed, int funct,
00381                 Real *c, int *mult) const;
00382 };
00383 
00384 /* this stores the VDW interaction parameters as a function of the
00385    atom types. */
00386 class VDWTable {
00387  private:
00388   /* atom types a and b */
00389   ResizeArray<const char *> typeAArray; /* NULL if no type is set */
00390   ResizeArray<const char *> typeBArray; /* NULL if no type is set */
00391 
00392   /* the VDW parameters */
00393   ResizeArray<Real> c6Array;
00394   ResizeArray<Real> c12Array;
00395   ResizeArray<Real> c6PairArray;  // these two are
00396   ResizeArray<Real> c12PairArray; // for 1-4 pairs
00397   
00398   /* finds the index from the two interacting types, or returns -1 */
00399   int getIndex(const char *typea, const char *typeb) const;
00400   
00401  public:
00402   /* these add VDW parameters to the list */
00403   void addType(const char *typea, const char *typeb, Real c6, Real c12);
00404   void add14Type(const char *typea, const char *typeb,
00405                  Real c6pair, Real c12pair);
00406 
00407   /* this returns the VDW interaction parameters that were added to
00408      the list earlier, and returns the index or the parameters (a
00409      useless number) or -1 if it can't find the specified types. */
00410   int getParams(const char *typea, const char *typeb,
00411                  Real *c6, Real *c12, Real *c6pair, Real *c12pair) const;
00412 };
00413 
00414 /* A GromacsTopFile represents the information stored in a GROMACS
00415    topology file.  This is an immutable type. */
00416 class GromacsTopFile {
00417  private:
00418   /* The system name */
00419   char *systemName;
00420 
00421   /* 1-4 parameter scaling values */
00422   Real fudgeLJ, fudgeQQ;
00423 
00424   /* whether or not to generate the 1-4 interactions automatically */
00425   int genPairs;
00426 
00427   /* a list of molecule templates */
00428   ResizeArray <GenericMol *> genericMols;
00429 
00430   /* the list of molecule instances */
00431   ResizeArray <const MolInst *> molInsts;
00432  
00433   /* the table of bonds/angles/atoms */
00434   AtomTable atomTable;
00435   BondTable bondTable;
00436   AngleTable angleTable;
00437   DihedralTable dihedralTable;
00438   VDWTable vdwTable;
00439  
00440  public:
00441 
00442   /* initializes this to represent the data stored in the topology
00443      file <filename>, or exits on error. */
00444   GromacsTopFile(char *filename);
00445 
00446   /* getSystemName returns the name of the system */
00447   char *getSystemName() const { return systemName; }
00448 
00449   /* returns the number of atoms/bonds/angles/dihedrals in the file */
00450   int getNumAtoms() const;
00451   int getNumBonds() const;
00452   int getNumAngles() const;
00453   int getNumDihedrals() const;
00454 
00455   /* returns the number of atoms/bonds/angles/dihedrals params in the file */
00456   int getNumAtomParams() const { return atomTable.size(); }
00457   int getNumBondParams() const { return bondTable.size(); }
00458   int getNumAngleParams() const { return angleTable.size(); }
00459   int getNumDihedralParams() const { return dihedralTable.size(); }
00460 
00461   /* getAtom puts the information about atom number <num> into the
00462      spaces pointed to by the other arguments.  The string buffers
00463      must be at least 11 characters long.  Atom number 0 is the first
00464      atom in the list, so that it corresponds with our normal idea of
00465      atom numbers.
00466      charge - e
00467      mass   - amu
00468   */
00469   void getAtom(int num, int *residue_number, char *residue_name,
00470                char *atom_name, char *atom_type, int *atom_typenum,
00471                Real *charge, Real *mass)
00472     const;
00473 
00474   /* getAtomParams puts the parameters for atom-type <num> into the
00475      spaces pointed to by the other arguements.  Currently, all anyone
00476      cares about is the atom type string, so that's all we're
00477      returning!  This may change later, but for now we get everything
00478      else with the getAtom() method.  The string buffer must be at
00479      least 11 characters long. */
00480   void getAtomParams(int num, char *type) const {
00481     atomTable.getType(num,type);
00482   }
00483 
00484   /* getBond puts the information about bond number <num> into the
00485      spaces pointed to by the other arguments.  Bond number 0 is the
00486      first bond in the list.
00487   */
00488   void getBond(int num, int *atomi, int *atomj, int *bondtype) const;
00489 
00490   /* getBondParams puts the parameters for bond-type <num> into the
00491      spaces pointed to by the other arguments. 
00492      b0 - natural length in A
00493      kB - spring constant in kcal/mol/A^2 - E=kx^2
00494      funct - 1 for harmonic, 2 for special fourth-order GROMOS96 thing */
00495   void getBondParams(int num, Real *b0, Real *kB, int *funct) const;
00496 
00497   /* getAngle puts the information about angle number <num> into the
00498      spaces pointed to by the other arguments.  Angle number 0 is the
00499      first angle in the list.
00500   */
00501   void getAngle(int num, int *atomi, int *atomj, int *atomk,
00502                 int *angletype) const;
00503 
00504   /* getAngleParams puts the parameters for angle-type <num> into the
00505      spaces pointed to by the other arguments.
00506      th0 - natural angle in degrees
00507      kth - spring constant in kcal/rad - E=kx^2
00508      funct - 1 for harmonic, 2 for special GROMOS96 thing */
00509   void getAngleParams(int num, Real *th0, Real *kth, int *funct) const;
00510 
00511   /* getDihedral puts the information about dihedral number <num> into
00512      the spaces pointed to by the other arguments.  Dihedral number 0
00513      is the first angle in the list. */
00514   void getDihedral(int num, int *atomi, int *atomj, int *atomk,
00515                    int *atoml, int *type) const;
00516 
00517   /* getDihedralParams puts the parameters for dihedral-type <num>
00518      into the spaces pointed to by the other arguments.  The array of
00519      Reals <c> must be of size >= 6, since up to 6 parameters will be
00520      set there.  These parameters are as follows:
00521    funct 1: proper dihedral
00522      c0 - thmax (deg)
00523      c1 - fc (kcal/mol)
00524      mult is the multiplicity
00525    funct 2: improper dihedral
00526      c0 - th0 (deg)
00527      c1 - fc (kcal/mol/rad2)
00528    funct 3: RB dihedral
00529      c0-c5 - C0-C5 (kcal/mol) 
00530   */
00531   void getDihedralParams(int num, Real *c, int *mult, int *funct) const;
00532 
00533   /* getVDWParams returs the Lennard-Jones bonding parameters for the
00534      specified two atom types, and the modified bonding parameters
00535      for 1-4 L-J interactons (c6pair, c12pair).
00536      Note that unlike the other types of parameters, you don't refer
00537      to this one by number - any two atom types define a VDWParam */
00538   void getVDWParams(int typea, int typeb,
00539                     Real *c6, Real *c12, Real *c6pair, Real *c7) const;
00540     
00541 };
00542 
00543 #endif // GROMACSTOPFILE_H

Generated on Mon Nov 23 04:59:21 2009 for NAMD by  doxygen 1.3.9.1