Amber masthead
Filler image AmberTools23 Amber22 Manuals Tutorials Force Fields Contacts History
Filler image

Useful links:

Amber Home
Download Amber
Installation
Amber Citations
GPU Support
Updates
Mailing Lists
For Educators
File Formats
Contributors

AMBER file formats


``PARM'' parameter/topology file specification

The following section provides a good overview of the general topology file format. A more complete prmtop specification, complete with chamber-style sections and expanded explanations of some sections and how they are used to calculate forces in Amber programs, is available in PDF-form by clicking here. (This document is taken from the appendix of Jason Swails’ doctoral dissertation.) This file is generated by the LEaP programs, either sleap or tleap/xleap.

%FLAG TITLE
%FORMAT(20a4)  (ITITL(i), i=1,20)
  ITITL  : title

%FLAG POINTERS
%FORMAT(10i8)  NATOM,    NTYPES, NBONH,  MBONA,  NTHETH, MTHETA,
               NPHIH,    MPHIA,  NHPARM, NPARM,  NNB,    NRES,
               NBONA,    NTHETA, NPHIA,  NUMBND, NUMANG, NPTRA,
               NATYP,    NPHB,   IFPERT, NBPER,  NGPER,  NDPER,
               MBPER,    MGPER,  MDPER,  IFBOX,  NMXRS,  IFCAP,
               NUMEXTRA, NCOPY
  NATOM    : total number of atoms 
  NTYPES   : total number of distinct atom types
  NBONH    : number of bonds containing hydrogen
  MBONA    : number of bonds not containing hydrogen
  NTHETH   : number of angles containing hydrogen
  MTHETA   : number of angles not containing hydrogen
  NPHIH    : number of dihedrals containing hydrogen
  MPHIA    : number of dihedrals not containing hydrogen
  NHPARM   : currently not used
  NPARM    : used to determine if addles created prmtop
  NNB      : number of excluded atoms
  NRES     : number of residues
  NBONA    : MBONA + number of constraint bonds
  NTHETA   : MTHETA + number of constraint angles
  NPHIA    : MPHIA + number of constraint dihedrals
  NUMBND   : number of unique bond types
  NUMANG   : number of unique angle types
  NPTRA    : number of unique dihedral types
  NATYP    : number of atom types in parameter file, see SOLTY below
  NPHB     : number of distinct 10-12 hydrogen bond pair types
  IFPERT   : set to 1 if perturbation info is to be read in
  NBPER    : number of bonds to be perturbed
  NGPER    : number of angles to be perturbed
  NDPER    : number of dihedrals to be perturbed
  MBPER    : number of bonds with atoms completely in perturbed group
  MGPER    : number of angles with atoms completely in perturbed group
  MDPER    : number of dihedrals with atoms completely in perturbed groups
  IFBOX    : set to 1 if standard periodic box, 2 when truncated octahedral
  NMXRS    : number of atoms in the largest residue
  IFCAP    : set to 1 if the CAP option from edit was specified
  NUMEXTRA : number of extra points found in topology
  NCOPY    : number of PIMD slices / number of beads

%FLAG ATOM_NAME
%FORMAT(20a4)  (IGRAPH(i), i=1,NATOM)
  IGRAPH : the user-specified atoms names

%FLAG CHARGE
%FORMAT(5E16.8)  (CHARGE(i), i=1,NATOM)
  CHARGE : the atom charges.  Amber internally uses units of charge such
           that E = q1*q2/r, where E is in kcal/mol, r is in Angstrom, 
           and q1,q2 are the values found in this section of the prmtop file.

           Codes that write prmtop files need to ensure that the values
           entered in this section satisfy the rule given above; similarly,
           codes that read prmtop files need to interpret these values in
           the same way that Amber does.  (This is, after all, an "Amber"
           file format!)

           [Aside: Conversion of these internal values to units of electron 
           charge depends of what values are chosen for various fundamental
           constants.  For reasons lost to history, Amber force field 
           development has been based on the formula: q(internal units) =
           q(electron charge units)*18.2223.  Gromacs and some other programs
           currently [2015] use a conversion factor of 18.222615, which
           more closely agrees with currently-accepted values for fundamental
           constants.  Conversions to or from the prmtop format may need 
           to take account of such small differences, if they are important.]

%FLAG ATOMIC_NUMBER
%FORMAT(10I8)  (ATNUM(i), i=1,NATOM)
  ATNUM : the atomic number of each atom.

%FLAG MASS
%FORMAT(5E16.8)  (AMASS(i), i=1,NATOM)
  AMASS  : the atom masses

%FLAG ATOM_TYPE_INDEX
%FORMAT(1OI8)  (IAC(i), i=1,NATOM)
  IAC    : index for the atom types involved in Lennard Jones (6-12) 
           interactions.  See ICO below.

%FLAG NUMBER_EXCLUDED_ATOMS
%FORMAT(10I8)  (NUMEX(i), i=1,NATOM)
  NUMEX  : total number of excluded atoms for atom "i".  Also called IBLO. See
           NATEX below.

%FLAG NONBONDED_PARM_INDEX
%FORMAT(10I8)  (ICO(i), i=1,NTYPES*NTYPES)
  ICO    : provides the index to the nonbon parameter
           arrays CN1, CN2 and ASOL, BSOL.  All possible 6-12
           or 10-12 atoms type interactions are represented.
           NOTE: A particular atom type can have either a 10-12
           or a 6-12 interaction, but not both.  The index is
           calculated as follows:
             index = ICO(NTYPES*(IAC(i)-1)+IAC(j))
           If index is positive, this is an index into the
           6-12 parameter arrays (CN1 and CN2) otherwise it
           is an index into the 10-12 parameter arrays (ASOL
           and BSOL).

%FLAG RESIDUE_LABEL
%FORMAT(20A4)  (LBRES(i), i=1,NRES)
  LBRES : names of each of the residues

%FLAG RESIDUE_POINTER
%FORMAT(10I8)  (IPRES(i), i=1,NRES)
  IPRES  : atoms in each residue are listed for atom "i" in
           IPRES(i) to IPRES(i+1)-1

%FLAG BOND_FORCE_CONSTANT
%FORMAT(5E16.8)  (RK(i), i=1,NUMBND)
  RK     : force constant for the bonds of each type, kcal/(mol Angstrom**2)

%FLAG BOND_EQUIL_VALUE
%FORMAT(5E16.8)  (REQ(i), i=1,NUMBND)
  REQ    : the equilibrium bond length for the bonds of each type, Angstroms

%FLAG ANGLE_FORCE_CONSTANT
%FORMAT(5E16.8)  (TK(i), i=1,NUMANG)
  TK     : force constant for the angles of each type, kcal/(mol radian**2)

%FLAG ANGLE_EQUIL_VALUE
%FORMAT(5E16.8)  (TEQ(i), i=1,NUMANG)
  TEQ    : the equilibrium angle for the angles of each type, radians

%FLAG DIHEDRAL_FORCE_CONSTANT
%FORMAT(5E16.8)  (PK(i), i=1,NPTRA)
  PK     : force constant for the dihedrals of each type, kcal/mol

%FLAG DIHEDRAL_PERIODICITY
%FORMAT(5E16.8)  (PN(i), i=1,NPTRA)
  PN     : periodicity of the dihedral of a given type

%FLAG DIHEDRAL_PHASE
%FORMAT(5E16.8)  (PHASE(i), i=1,NPTRA)
  PHASE  : phase of the dihedral of a given type, radians

%FLAG SCEE_SCALE_FACTOR
%FORMAT(5E16.8) (ONE_SCEE(i), i=1,NPTRA)
  ONE_SCEE : 1-4 electrostatic scaling constant. It is inverted right after
             it's read in for performance reasons. This allows variable 
             1-4 scaling. If not present, it defaults to 1.2 for all dihedrals.
             Therefore, the default ONE_SCEE value in the code is 1.0/1.2

%FLAG SCNB_SCALE_FACTOR
%FORMAT(5E16.8) (ONE_SCNB(i), i=1,NPTRA)
  ONE_SCNB : 1-4 VDW scaling constant. It is inverted right after
             it's read in. This allows variable 1-4 scaling. If not present,
             it defaults to 2.0 for all dihedrals. Therefore, the default
             ONE_SCNB value in the code is 1.0/2.0

%FLAG SOLTY
%FORMAT(5E16.8)  (SOLTY(i), i=1,NATYP)
  SOLTY  : currently unused (reserved for future use)

%FLAG LENNARD_JONES_ACOEF
%FORMAT(5E16.8)  (CN1(i), i=1,NTYPES*(NTYPES+1)/2)
  CN1  : Lennard Jones r**12 terms for all possible atom type interactions,
         indexed by ICO and IAC; for atom i and j where i < j, the index
         into this array is as follows (assuming the value of ICO(index) is
         positive): CN1(ICO(NTYPES*(IAC(i)-1)+IAC(j))).

%FLAG LENNARD_JONES_BCOEF
%FORMAT(5E16.8)  (CN2(i), i=1,NTYPES*(NTYPES+1)/2)
  CN2  : Lennard Jones r**6 terms for all possible
         atom type interactions.  Indexed like CN1 above.

NOTE: the atom numbers in the following arrays that describe bonds, angles, and dihedrals are coordinate array indexes for runtime speed. The true atom number equals the absolute value of the number divided by three, plus one. In the case of the dihedrals, if the fourth atom is negative, this implies that the dihedral is an improper. If the third atom is negative, this implies that the end group interations are to be ignored. End group interactions are ignored, for example, in dihedrals of various ring systems (to prevent double counting of 1-4 interactions) and in multiterm dihedrals.

%FLAG BONDS_INC_HYDROGEN
%FORMAT(10I8)  (IBH(i),JBH(i),ICBH(i), i=1,NBONH)
  IBH    : atom involved in bond "i", bond contains hydrogen
  JBH    : atom involved in bond "i", bond contains hydrogen
  ICBH   : index into parameter arrays RK and REQ

%FLAG BONDS_WITHOUT_HYDROGEN
%FORMAT(10I8)  (IB(i),JB(i),ICB(i), i=1,NBONA)
  IB     : atom involved in bond "i", bond does not contain hydrogen
  JB     : atom involved in bond "i", bond does not contain hydrogen
  ICB    : index into parameter arrays RK and REQ

%FLAG ANGLES_INC_HYDROGEN
%FORMAT(10I8)  (ITH(i),JTH(i),KTH(i),ICTH(i), i=1,NTHETH)
  ITH    : atom involved in angle "i", angle contains hydrogen
  JTH    : atom involved in angle "i", angle contains hydrogen
  KTH    : atom involved in angle "i", angle contains hydrogen
  ICTH   : index into parameter arrays TK and TEQ for angle
           ITH(i)-JTH(i)-KTH(i)

%FLAG ANGLES_WITHOUT_HYDROGEN
%FORMAT(10I8)  (IT(i),JT(i),KT(i),ICT(i), i=1,NTHETA)
  IT     : atom involved in angle "i", angle does not contain hydrogen
  JT     : atom involved in angle "i", angle does not contain hydrogen
  KT     : atom involved in angle "i", angle does not contain hydrogen
  ICT    : index into parameter arrays TK and TEQ for angle
           IT(i)-JT(i)-KT(i)

%FLAG DIHEDRALS_INC_HYDROGEN
%FORMAT(10I8)  (IPH(i),JPH(i),KPH(i),LPH(i),ICPH(i), i=1,NPHIH)
  IPH    : atom involved in dihedral "i", dihedral contains hydrogen
  JPH    : atom involved in dihedral "i", dihedral contains hydrogen
  KPH    : atom involved in dihedral "i", dihedral contains hydrogen
  LPH    : atom involved in dihedral "i", dihedral contains hydrogen
  ICPH   : index into parameter arrays PK, PN, PHASE, ONE_SCEE, and ONE_SCNB
           for dihedral IPH(i)-JPH(i)-KPH(i)-LPH(i)

%FLAG DIHEDRALS_WITHOUT_HYDROGEN
%FORMAT(10I8)  (IP(i),JP(i),KP(i),LP(i),ICP(i), i=1,NPHIA)
  IP     : atom involved in dihedral "i", dihedral does not contain hydrogen
  JP     : atom involved in dihedral "i", dihedral does not contain hydrogen
  KP     : atom involved in dihedral "i", dihedral does not contain hydrogen
  LP     : atom involved in dihedral "i", dihedral does not contain hydrogen
  ICP    : index into parameter arrays PK, PN, PHASE, ONE_SCEE, and ONE_SCNB
           for dihedral IPH(i)-JPH(i)-KPH(i)-LPH(i).  Note, if the
           periodicity is negative, this implies the following entry
           in the PK, PN, and PHASE arrays is another term in a
           multitermed dihedral.

%FLAG EXCLUDED_ATOMS_LIST
%FORMAT(10I8)  (INB(i), i=1,NNB)
  INB    : the excluded atom list.  To get the excluded list for atom 
           "i" you need to traverse the NUMEX list, adding up all
           the previous NUMEX values, since NUMEX(i) holds the number
           of excluded atoms for atom "i", not the index into the 
           NATEX list.  Let IEXCL = SUM(NUMEX(j), j=1,i-1), then
           excluded atoms are INB(IEXCL) to INB(IEXCL+NUMEX(i)). Note,
           this array was called NATEX at one point, and while in most
           places it is now INB, it is still called NATEX in some places
           (especially in pmemd)

%FLAG HBOND_ACOEF
%FORMAT(5E16.8)  (ASOL(i), i=1,NPHB)
  ASOL   : the value for the r**12 term for hydrogen bonds of all
           possible types.  Index into these arrays is equivalent
           to the CN1 and CN2 arrays, however the index is negative.
           For example, for atoms i and j, with i < j, the index is
           -ICO(NTYPES*(IAC(i)-1+IAC(j)). Note: Amber must be compiled
           with -DHAS_10_12 in order to make use of this term!

%FLAG HBOND_BCOEF
%FORMAT(5E16.8)  (BSOL(i), i=1,NPHB)
  BSOL   : the value for the r**10 term for hydrogen bonds of all
           possible types.  Indexed like ASOL. Note: same restriction
           applies for use.

%FLAG HBCUT
%FORMAT(5E16.8)  (HBCUT(i), i=1,NPHB)
  HBCUT  : no longer in use

%FLAG AMBER_ATOM_TYPE
%FORMAT(20A4)  (ISYMBL(i), i=1,NATOM)
  ISYMBL : the AMBER atom types for each atom

%FLAG TREE_CHAIN_CLASSIFICATION
%FORMAT(20A4)  (ITREE(i), i=1,NATOM)
  ITREE  : the list of tree joining information, classified into five
           types.  M -- main chain, S -- side chain, B -- branch point, 
           3 -- branch into three chains, E -- end of the chain

%FLAG JOIN_ARRAY
%FORMAT(10I8)  (JOIN(i), i=1,NATOM)
  JOIN   : tree joining information, potentially used in ancient
           analysis programs.  Currently unused in sander or gibbs.

%FLAG IROTAT
%FORMAT(10I8)  (IROTAT(i), i = 1, NATOM)
  IROTAT : apparently the last atom that would move if atom i was
           rotated, however the meaning has been lost over time.
           Currently unused in sander or gibbs.

%FLAG RADIUS_SET
%FORMAT(1a80) TYPE
  TYPE   : the radius set chosen inside LEaP

%FLAG RADII
%FORMAT(5E16.8) (RBORN(i), i=1,NATOM)
  RBORN  : Generalized Born intrinsic dielectric radii

%FLAG SCREEN
%FORMAT(5E16.8) (FS(i), i=1,NATOM)
  FS     : Screening parameters used in Generalized Born

The following are only present if IFBOX .gt. 0

%FLAG SOLVENT_POINTERS
%FORMAT(3I8)  IPTRES, NSPM, NSPSOL
  IPTRES : final residue that is considered part of the solute,
           reset in sander and gibbs
  NSPM   : total number of molecules
  NSPSOL : the first solvent "molecule"

%FLAG ATOMS_PER_MOLECULE
%FORMAT(10I8)  (NSP(i), i=1,NSPM)
  NSP    : the total number of atoms in each molecule,
           necessary to correctly perform the pressure
           scaling.

%FLAG BOX_DIMENSIONS
%FORMAT(5E16.8)  OLDBETA, BOX(1), BOX(2), BOX(3)
  OLDBETA   : periodic box, angle between the XY and YZ planes in
              degrees. This is now redundant, as it is present in the
              inpcrd files. It is ignored here.
  BOX       : the periodic box lengths in the X, Y, and Z directions.
              This is now redundant, as it is present in the inpcrd files.
              It is ignored here.

The following are only present if IFCAP .gt. 0

%FLAG CAP_INFO
%FORMAT(10I8)  NATCAP
  NATCAP : last atom before the start of the cap of waters
           placed by edit

%FLAG CAP_INFO2
%FORMAT(5E16.8)  CUTCAP, XCAP, YCAP, ZCAP
  CUTCAP : the distance from the center of the cap to the outside
  XCAP   : X coordinate for the center of the cap
  YCAP   : Y coordinate for the center of the cap
  ZCAP   : Z coordinate for the center of the cap

The following is only present if IFPERT .gt. 0 Note that the initial state, or equivalently the prep/link/edit state, is represented by lambda=1 and the perturbed state, or final state specified in parm, is the lambda=0 state. This information is only used for GIBBS and is unused in either SANDER or PMEMD.

FORMAT(12I6)  (IBPER(i), JBPER(i), i=1,NBPER)
  IBPER  : atoms involved in perturbed bonds
  JBPER  : atoms involved in perturbed bonds

FORMAT(12I6)  (ICBPER(i), i=1,2*NBPER)
  ICBPER : pointer into the bond parameter arrays RK and REQ for the
           perturbed bonds.  ICBPER(i) represents lambda=1 and 
           ICBPER(i+NBPER) represents lambda=0.

FORMAT(12I6)  (ITPER(i), JTPER(i), KTPER(i), i=1,NGPER)
  IPTER  : atoms involved in perturbed angles
  JTPER  : atoms involved in perturbed angles
  KTPER  : atoms involved in perturbed angles

FORMAT(12I6)  (ICTPER(i), i=1,2*NGPER)
  ICTPER : pointer into the angle parameter arrays TK and TEQ for 
           the perturbed angles.  ICTPER(i) represents lambda=0 and 
           ICTPER(i+NGPER) represents lambda=1.

FORMAT(12I6)  (IPPER(i), JPPER(i), KPPER(i), LPPER(i), i=1,NDPER)
  IPPER  : atoms involved in perturbed dihedrals
  JPPER  : atoms involved in perturbed dihedrals
  KPPER  : atoms involved in perturbed dihedrals
  LPPER  : atoms involved in pertrubed dihedrals

FORMAT(12I6)  (ICPPER(i), i=1,2*NDPER)
  ICPPER : pointer into the dihedral parameter arrays PK, PN and
           PHASE for the perturbed dihedrals.  ICPPER(i) represents 
           lambda=1 and ICPPER(i+NGPER) represents lambda=0.

FORMAT(20A4)  (LABRES(i), i=1,NRES)
  LABRES : residue names at lambda=0

FORMAT(20A4)  (IGRPER(i), i=1,NATOM)
  IGRPER : atomic names at lambda=0

FORMAT(20A4)  (ISMPER(i), i=1,NATOM)
  ISMPER : atomic symbols at lambda=0

FORMAT(5E16.8)  (ALMPER(i), i=1,NATOM)
  ALMPER : unused currently in gibbs

FORMAT(12I6)  (IAPER(i), i=1,NATOM)
  IAPER  : IAPER(i) = 1 if the atom is being perturbed

FORMAT(12I6)  (IACPER(i), i=1,NATOM)
  IACPER : index for the atom types involved in Lennard Jones
           interactions at lambda=0.  Similar to IAC above.  
           See ICO above.

FORMAT(5E16.8)  (CGPER(i), i=1,NATOM)
  CGPER  : atomic charges at lambda=0

The following is only present if IPOL .eq. 1

%FLAG POLARIZABILITY
%FORMAT(5E18.8) (ATPOL(i), i=1,NATOM)
  ATPOL  : atomic polarizabilities

The following is only present if IPOL == 1 .and. IFPERT == 1

FORMAT(5E18.8) (ATPOL1(i), i=1,NATOM)
  ATPOL1 : atomic polarizabilities at lambda = 1 (above is at lambda = 0)

AMBER coordinate/restart file specification

The 'coord' version of this file is generated by the PARM or LEaP programs. The 'restrt' version is the result of energy minimization or molecular dynamics in SANDER or GIBBS and may contain velocity and periodic box information.

FORMAT(20A4) ITITL
  ITITL  : the title of the current run, from the AMBER
           parameter/topology file

FORMAT(I5,5E15.7) NATOM,TIME
  NATOM  : total number of atoms in coordinate file
  TIME   : option, current time in the simulation (picoseconds)

FORMAT(6F12.7) (X(i), Y(i), Z(i), i = 1,NATOM)
  X,Y,Z  : coordinates

IF dynamics

FORMAT(6F12.7) (VX(i), VY(i), VZ(i), i = 1,NATOM)
  VX,VY,VZ : velocities (units: Angstroms per 1/20.455 ps)

IF constant pressure (in 4.1, also constant volume)

FORMAT(6F12.7) BOX(1), BOX(2), BOX(3)
  BOX    : size of the periodic box

Note: in AMBER 4.1 if the EWALD option is turned on, the box angles will also be written out in the same format.


AMBER trajectory (coordinate or velocity) file specification

This file is optionally written during dynamics in SANDER or GIBBS.

FORMAT(20A4) ITITL
  ITITL  : the title of the current run, from the AMBER
           parameter/topology file

The following snapshot is written every NTWX steps in the trajectory (specified in the control input file):

FORMAT(10F8.3)  (X(i), Y(i), Z(i), i=1,NATOM)
  X,Y,Z  : coordinates or velocities (velocity units: Angstroms per 1/20.455 ps)

IF constant pressure (in 4.1, also constant volume), for each snapshot:

FORMAT(10F8.3)  BOX(1), BOX(2), BOX(3)
  BOX    : size of periodic box

AMBER Trajectory NetCDF Convention

Beginning with AMBER 9, a binary file format for trajectory data based on NetCDF is supported by sander, pmemd and ptraj. This format is supported by VMD beginning with version 1.8.4.

The AMBER Trajectory NetCDF Convention, describes the layout of dimensions, variables and attributes within the file. It is also available as a PDF.


Force field parameter file specification

      Force field information on the file frcfld: The following  sec-
      tion  of  this document describes the format of the AMBER Force
      Field Parameter File.  It is not expected that  the  user  will
      ordinarily  modify this file; rather modifications should ordi-
      narily be entered through the  frcmod  file  described  further
      below.   Of course, major changes, such as using the AMBER/OPLS
      force field rather than the AMBER one, would best  be  made  by
      changing this file. WARNING: multiple entries for the same atom
      symbols within a single frcfld or frcmod file can lead to unde-
      fined results, e.g. if there are two definitions of angle ener-
      gies between atom types A, B and C one of them is picked  arbi-
      trarily.

         - 1  -     ITITL

                        FORMAT(20A4)

         ITITL      A title for identification of the parameter set.

      ------------------------------------------------------------------------

         - 2 -      ***** INPUT FOR ATOM SYMBOLS AND MASSES *****

                    KNDSYM , AMASS, ATPOL

                        FORMAT(A2,2X,F10.2x,f10.2)

         KNDSYM     The unique atom symbol used in the system.

         AMASS      Atomic mass of the center having the symbol "KNDSYM".

         ATPOL      The atomic polarizability for each atom (in A**3)
                    This is the type of polarizability used in sander
                    and gibbs. No parameters are supplied for this since
                    the feature is still in development (Amber 4.1).

              NOTE: All the unique atomic symbols and their masses must
                    be read.  The input is terminated by a blank card.


      ------------------------------------------------------------------------

         - 3 -      ***** INPUT FOR ATOM SYMBOLS THAT ARE HYDROPHILIC *****

                    JSOLTY(I)

                         FORMAT(20(A2,2X))

         JSOLTY(I)  The atom symbols which are hydrophilic in solution.
                    This information is read but not used.

                    The input is terminated when a blank value is read for
                    the atom symbol.

      ------------------------------------------------------------------------

         - 4 -      ***** INPUT FOR BOND LENGTH PARAMETERS *****

                    IBT , JBT , RK , REQ

                        FORMAT(A2,1X,A2,2F10.2)

         IBT,JBT    Atom symbols for the two bonded atoms.

         RK         The harmonic force constant for the bond "IBT"-"JBT".
                    The unit is kcal/mol/(A**2).

         REQ        The equilibrium bond length for the above bond in Angstroms

                    The input is terminated by a blank card.

      ------------------------------------------------------------------------

         - 5 -      ***** INPUT FOR BOND ANGLE PARAMETERS *****

                    ITT , JTT , KTT , TK , TEQ

                        FORMAT(A2,1X,A2,1X,A2,2F10.2)

         ITT,...    The atom symbols for the atoms making an angle.

         TK         The harmonic force constants for the angle "ITT"-"JTT"-
                    "KTT" in units of kcal/mol/(rad**2) (radians are the
                    traditional unit for angle parameters in force fields).

         TEQ        The equilibrium bond angle for the above angle in degrees.

                    The input is terminated by a blank card.

      ------------------------------------------------------------------------

         - 6 -      ***** INPUT FOR DIHEDRAL PARAMETERS *****

                    IPT , JPT , KPT , LPT , IDIVF , PK , PHASE , PN

                        FORMAT(A2,1X,A2,1X,A2,1X,A2,I4,3F15.2)

         IPT, ...   The atom symbols for the atoms forming a dihedral
                    angle.  If IPT .eq. 'X ' .and. LPT .eq. 'X ' then
                    any dihedrals in the system involving the atoms "JPT" and
                    and "KPT" are assigned the same parameters.  This is
                    called the general dihedral type and is of the form
                    "X "-"JPT"-"KPT"-"X ".

         IDIVF      The factor by which the torsional barrier is divided.
                    Consult Weiner, et al., JACS 106:765 (1984) p. 769 for
                    details. Basically, the actual torsional potential is

                           (PK/IDIVF) * (1 + cos(PN*phi - PHASE))

         PK         The barrier height divided by a factor of 2.

         PHASE      The phase shift angle in the torsional function.

                    The unit is degrees.

         PN         The periodicity of the torsional barrier.
                    NOTE: If PN .lt. 0.0 then the torsional potential
                          is assumed to have more than one term, and the
                          values of the rest of the terms are read from the
                          next cards until a positive PN is encountered.  The
                          negative value of pn is used only for identifying
                          the existence of the next term and only the
                          absolute value of PN is kept.

            The input is terminated by a blank card.

      ------------------------------------------------------------------------

         - 7 -      ***** INPUT FOR IMPROPER DIHEDRAL PARAMETERS *****

                    IPT , JPT , KPT , LPT , IDIVF , PK , PHASE , PN

                        FORMAT(A2,1X,A2,1X,A2,1X,A2,I4,3F15.2)

                    The input is the same as in for the dihedrals except that
                    the torsional barrier height is NOT divided by the factor
                    idivf.  The improper torsions are defined between any four
                    atoms not bonded (in a successive fashion) with each other
                    as in the case of "regular" or "proper" dihedrals.  Improper
                    dihedrals are used to keep certain groups planar and to
                    prevent the racemization of certain centers in the united
                    atom model.  Consult the above reference for details.

             Important note: all general type improper dihedrals
                             (e.g. x -x -ct-hc) should appear before all
                             specifics (ct-ct-ct-hc) in the parm list.
                             Otherwise the generals will override the
                             specific with no warning.

             The input is terminated by a blank card.

      ------------------------------------------------------------------------

         - 8 -      ***** INPUT FOR H-BOND 10-12 POTENTIAL PARAMETERS *****

                    KT1 , KT2 , A , B , ASOLN , BSOLN , HCUT , IC

                        FORMAT(2X,A2,2X,A2,2x,5F10.2,I2)

         KT1,KT2    The atom symbols for the atom pairs for which the
                    parameters are defined.

         A          The coefficient of the 12th power term (A/(r**12)).

         B          The coefficient of the 10th power term (-B/(r**10)).

         ASOLN      Not used

         BSOLN      Not used

         HCUT       Not used

         IC         Not used

      ------------------------------------------------------------------------

         - 9 -      ***** INPUT FOR EQUIVALENCING ATOM SYMBOLS FOR
                          THE NON-BONDED 6-12 POTENTIAL PARAMETERS *****

                          IORG , IEQV(I) , I = 1 , 19

                              FORMAT(20(A2,2X))

         IORG        The atom symbols to which other atom symbols are to be
                     equivalenced in generating the 6-12 potential parameters.

         IEQV(I)     The atoms symbols which are to be equivalenced to the
                     atom symbol "IORG".  If more than 19 atom symbols have
                     to be equivalenced to a given atom symbol they can be
                     included as extra cards.

                     It is advisable not to equivalence any hydrogen bond
                     atom type atoms with any other atom types.

          NOTE: The input is terminated by a blank card.

      ------------------------------------------------------------------------

         - 10 -      ***** INPUT FOR THE 6-12 POTENTIAL PARAMETERS *****

                     LABEL , KINDNB

                         FORMAT(A4,6X,A2)

         LABEL       The name of the non-bonded input parameter to be
                     used.  It has to be matched with "NAMNB" read through
                     unit 5.  The program searches the file to load the
                     the required non-bonded parameters.  If that name is
                     not found the run will be terminated.

         KINDNB      Flag for the type of 6-12 parameters.

          'SK'       Slater-Kirkwood parameters are input.
                     see "caution" below.

          'RE'       van der Waals radius and the potential well depth
                     parameters are read.

          'AC'       The 6-12 potential coefficients are read.

		NOTE: All the non equivalenced atoms' parameters have to
                   be given.

           The input is terminated when label .eq. 'END'

      ------------------------------------------------------------------------
      CAUTION: the polarizabilities mentioned below are NOT the
               polarizabilities used in the sander (min/md) code.
               KINDNB 'SK' parameters are not currently part of
               the AMBER force field. See card 2, ATPOL for sander
               polarizability.


         - 10A -     ***** ONLY IF KINDNB .EQ. 'SK' *****

                     LTYNB , POL , XNEFF , RMIN

                         FORMAT(2X,A2,6X,3F10.6)

         LTYNB       Atom symbol.

         POL         Atomic polarizability for the atom centers having the
                     the above symbol.

         XNEFF       Effective number of electrons on the atom centers having
                     the above symbol.

         RMIN        van der Waals radius of the atom center having the above
                     symbol.

      ------------------------------------------------------------------------

         - 10B -     ***** ONLY IF KINDNB .EQ. 'RE' *****

                     LTYNB , R , EDEP

         LTYNB       Atom symbol.

         R           The van der Waals radius of the atoms having the symbol
                     "LTYNB"  (Angstoms)

         EDEP        The 6-12 potential well depth. (kcal/mol)

      ------------------------------------------------------------------------

         - 10C -     ***** ONLY IF KINDNB .EQ. 'AC' *****

                     LTYNB , A , C

         LTYNB       Atom symbol.

         A           The coefficient of the 12th power term (A/r**12).

         C           The coefficient of the 6th power term (-C/r**6).


Force field parameter modification file specification


      Modified  force  field  parameters in file frcmod: This file is
      normally the one that will be changed by the user.  It consists
      of  a 1-card title, followed by a blank line, then keyword sec-
      tions.  The allowed keywords (appearing in columns 1-4) are:


      MASS          follow this card by card of type - 2 - listed in the
                    unit 10 instructions above.  End with a blank line.

      BOND          follow this card by card of type - 4 - listed in the
                    unit 10 instructions above.  End with a blank line.

      ANGL          follow this card by card of type - 5 - listed in the
                    unit 10 instructions above.  End with a blank line.

      DIHE          follow this card by card of type - 6 - listed in the
                    unit 10 instructions above.  End with a blank line.

      IMPR          follow this card by card of type - 7 - listed in the
                    unit 10 instructions above.  End with a blank line.

      HBON          follow this card by card of type - 8 - listed in the
                    unit 10 instructions above.  End with a blank line.

      NONB          follow this card by card of type - 10A, B or C - listed
                    in the unit 10 instructions above.  E.g. if you specify
                    STDA in parm.in for the "regular" parm.dat file, this is
                    the convention that will be used when reading frcmod.
                    End with a blank line.


           Any or all of the keywords may be missing, if you have  no
      changes  for  that  section.  The entire file can be missing if
      you have no changes at all to make to the standard force field.
      Restrictions:  note that you cannot modify the equivalence pat-
      tern set up in the standard force field.

           If you have parameters in the frcmod file that modify val-
      ues  in the standard parameter file, and if iparml is set to 2,
      then both the original and  the  modified  parameters  will  be
      printed  to  the  output file.  The modified parameters will be
      marked with asterisks, and it is these values that will be used
      in subsequent calculations.

"How's that for maxed out?"

Last modified: May 3, 2020